fvMeshTopoChangersMeshToMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2022-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "polyTopoChangeMap.H"
28 #include "volFields.H"
29 #include "surfaceInterpolate.H"
30 #include "pointFields.H"
32 #include "fvMeshToFvMesh.H"
34 #include "surfaceToVolVelocity.H"
36 #include "polyMeshMap.H"
37 #include "processorPolyPatch.H"
41 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 namespace fvMeshTopoChangers
47 {
50 }
51 }
52 
53 
54 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
55 
56 bool Foam::fvMeshTopoChangers::meshToMesh::forward() const
57 {
58  return
59  cycle_ == 0
60  || int((mesh().time().userTimeValue() - begin_)/cycle_) % 2 == 0;
61 }
62 
63 
64 Foam::scalar Foam::fvMeshTopoChangers::meshToMesh::meshTime() const
65 {
66  const Time& time = mesh().time();
67 
68  if (repeat_ > 0)
69  {
70  return begin_ + fmod(time.userTimeValue() - begin_, repeat_);
71  }
72  else if (cycle_ > 0)
73  {
74  if (forward())
75  {
76  return
77  begin_
78  + fmod(time.userTimeValue() - begin_, cycle_);
79  }
80  else
81  {
82  return
83  begin_ + cycle_
84  - fmod(time.userTimeValue() - begin_, cycle_);
85  }
86  }
87  else
88  {
89  return time.userTimeValue();
90  }
91 }
92 
93 
94 void Foam::fvMeshTopoChangers::meshToMesh::interpolateUfs()
95 {
96  // Interpolate U to Uf
97  UPtrList<surfaceVectorField> Ufs(mesh().curFields<surfaceVectorField>());
98 
99  forAll(Ufs, i)
100  {
101  surfaceVectorField& Uf = Ufs[i];
102 
104 
105  if (!isNull(U))
106  {
107  Uf.reset(fvc::interpolate(U));
108  }
109  }
110 }
111 
112 
113 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
114 
116 (
117  fvMesh& mesh,
118  const dictionary& dict
119 )
120 :
121  fvMeshTopoChanger(mesh),
122  dict_(dict),
123  times_(dict.lookup("times")),
124  timeDelta_(dict.lookup<scalar>("timeDelta")),
125  begin_(dict.lookupOrDefault("begin", mesh().time().beginTime().value())),
126  repeat_(dict.lookupOrDefault("repeat", 0.0)),
127  cycle_(dict.lookupOrDefault("cycle", 0.0)),
128  timeIndex_(-1)
129 {
130  if (repeat_ > 0 && cycle_ > 0)
131  {
133  << "Both 'repeat' and 'cycle' options specified"
134  << exit(FatalIOError);
135  }
136 
137  forAll(times_, i)
138  {
139  timeIndices_.insert(int64_t((times_[i] + timeDelta_/2.0)/timeDelta_));
140  }
141 }
142 
143 
144 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
145 
147 {}
148 
149 
150 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
151 
153 {
154  const Time& time = mesh().time();
155 
156  if
157  (
158  repeat_ > 0
159  || cycle_ > 0
160  || time.userTimeValue() + timeDelta_ < times_.last()
161  )
162  {
163  const scalar meshTime = this->meshTime();
164 
165  if (cycle_ == 0 || int((time.userTimeValue() - begin_)/cycle_) % 2 == 0)
166  {
167  forAll(times_, i)
168  {
169  if (times_[i] > meshTime + timeDelta_)
170  {
171  return time.userTimeToTime(times_[i] - meshTime);
172  }
173  }
174  }
175  else
176  {
177  forAllReverse(times_, i)
178  {
179  if (times_[i] < meshTime - timeDelta_)
180  {
181  return time.userTimeToTime(meshTime - times_[i]);
182  }
183  }
184  }
185  }
186 
187  return vGreat;
188 }
189 
190 
192 {
193  const Time& time = mesh().time();
194 
195  // Add the meshToMeshAdjustTimeStepFunctionObject functionObject
196  // if not already present
197  if
198  (
199  time.functionObjects().findObjectID("meshToMeshAdjustTimeStep")
200  == -1
201  )
202  {
203  const_cast<Time&>(time).functionObjects().append
204  (
206  (
207  "meshToMeshAdjustTimeStep",
208  time,
209  dict_
210  )
211  );
212  }
213 
214  // Only refine on the first call in a time-step
215  if (timeIndex_ != time.timeIndex())
216  {
217  timeIndex_ = time.timeIndex();
218  }
219  else
220  {
221  return false;
222  }
223 
224  // Obtain the mesh time from the user time
225  // repeating or cycling the mesh sequence if required
226 
227  const scalar meshTime = this->meshTime();
228 
229  if (timeIndices_.found((meshTime + timeDelta_/2)/timeDelta_))
230  {
231  const word otherMeshDir =
232  "meshToMesh_" + time.timeName(meshTime);
233 
234  Info << "Mapping to mesh " << otherMeshDir << endl;
235 
236  fvMesh otherMesh
237  (
238  IOobject
239  (
240  otherMeshDir,
241  time.constant(),
242  time,
244  ),
245  false,
247  );
248 
249  mesh().swap(otherMesh);
250 
251  fvMeshToFvMesh mapper
252  (
253  otherMesh,
254  mesh(),
255  cellsToCellss::intersection::typeName
256  );
257 
258  // Ensure the deltaCoeffs are available for constraint patch evaluation
259  mesh().deltaCoeffs();
260 
261  // Map all the volFields in the objectRegistry
262  #define mapVolFieldType(Type, nullArg) \
263  MeshToMeshMapVolFields<Type>(mesh(), mapper);
265 
266  // Map all the volFields in the objectRegistry
267  #define mapVolInternalFieldType(Type, nullArg) \
268  MeshToMeshMapVolInternalFields<Type>(mesh(), mapper);
270 
271  // Set all the surfaceFields in the objectRegistry to NaN
272  #define NaNSurfaceFieldType(Type, nullArg) \
273  NaNGeometricFields \
274  <Type, fvsPatchField, surfaceMesh, setSizeFvPatchFieldMapper> \
275  (mesh(), mapper);
277 
278  // Set all the pointFields in the objectRegistry to NaN
279  #define NaNPointFieldType(Type, nullArg) \
280  NaNGeometricFields \
281  <Type, pointPatchField, pointMesh, setSizePointPatchFieldMapper> \
282  (mesh(), mapper);
284 
285  // Interpolate U's to Uf's
286  interpolateUfs();
287 
288  polyMeshMap map(mesh(), mapper);
289 
290  mesh().mapMesh(map);
291 
292  return true;
293  }
294  else
295  {
296  return false;
297  }
298 }
299 
300 
302 (
303  const polyTopoChangeMap& map
304 )
305 {}
306 
307 
309 {}
310 
311 
313 (
314  const polyDistributionMap& map
315 )
316 {}
317 
318 
319 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: UList.H:446
Macros for easy insertion into run-time selection tables.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
static const word & constant()
Return constant name.
Definition: TimePaths.H:122
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:28
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
scalar userTimeToTime(const scalar tau) const
Convert the user-time (e.g. CA deg) to real-time (s).
Definition: Time.C:824
scalar userTimeValue() const
Return current user time value.
Definition: Time.C:818
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Definition: Time.C:630
const functionObjectList & functionObjects() const
Return the list of function objects.
Definition: Time.H:422
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
label findObjectID(const word &name) const
Find the ID of a given function object by name.
Abstract base class for fvMesh topology changers.
fvMesh & mesh()
Return the fvMesh.
fvMeshTopoChanger which maps the fields to a new mesh or sequence of meshes
virtual void topoChange(const polyTopoChangeMap &)
Update corresponding to the given map.
virtual void distribute(const polyDistributionMap &)
Update corresponding to the given distribution map.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
virtual bool update()
Update the mesh for both mesh motion and topology change.
meshToMesh(fvMesh &, const dictionary &dict)
Construct from fvMesh and dictionary.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const Time & time() const
Return time.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for handling words, derived from string.
Definition: word.H:62
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
#define NaNPointFieldType(Type, nullArg)
#define mapVolInternalFieldType(Type, nullArg)
#define NaNSurfaceFieldType(Type, nullArg)
#define mapVolFieldType(Type, nullArg)
U
Definition: pEqn.H:72
addToRunTimeSelectionTable(fvMeshTopoChanger, list, fvMesh)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
VolField< vector > volVectorField
Definition: volFieldsFwd.H:62
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
FOR_ALL_FIELD_TYPES(DefineContiguousFvWallLocationDataType)
IOerror FatalIOError
bool isNull(const T &t)
Return true if t is a reference to the nullObject of type T.
Definition: nullObjectI.H:52
SurfaceField< vector > surfaceVectorField
const volVectorField & surfaceToVolVelocity(const surfaceVectorField &Uf)
Get the cell velocity field corresponding to a given face velocity, or a.
dictionary dict
Return the vol-field velocity corresponding to a given surface-field velocity.