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-2024 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"
38 #include "setSizeFieldMapper.H"
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 namespace fvMeshTopoChangers
46 {
49 }
50 }
51 
52 
53 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
54 
55 bool Foam::fvMeshTopoChangers::meshToMesh::forward() const
56 {
57  return
58  cycle_ == 0
59  || int((mesh().time().userTimeValue() - begin_)/cycle_) % 2 == 0;
60 }
61 
62 
63 Foam::scalar Foam::fvMeshTopoChangers::meshToMesh::meshTime() const
64 {
65  const Time& time = mesh().time();
66 
67  if (repeat_ > 0)
68  {
69  return begin_ + fmod(time.userTimeValue() - begin_, repeat_);
70  }
71  else if (cycle_ > 0)
72  {
73  if (forward())
74  {
75  return
76  begin_
77  + fmod(time.userTimeValue() - begin_, cycle_);
78  }
79  else
80  {
81  return
82  begin_ + cycle_
83  - fmod(time.userTimeValue() - begin_, cycle_);
84  }
85  }
86  else
87  {
88  return time.userTimeValue();
89  }
90 }
91 
92 
93 void Foam::fvMeshTopoChangers::meshToMesh::interpolateUfs()
94 {
95  // Interpolate U to Uf
96  UPtrList<surfaceVectorField> Ufs(mesh().curFields<surfaceVectorField>());
97 
98  forAll(Ufs, i)
99  {
100  surfaceVectorField& Uf = Ufs[i];
101 
103 
104  if (!isNull(U))
105  {
106  Uf.reset(fvc::interpolate(U));
107  }
108  }
109 }
110 
111 
112 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
113 
115 (
116  fvMesh& mesh,
117  const dictionary& dict
118 )
119 :
120  fvMeshTopoChanger(mesh),
121  dict_(dict),
122  times_(dict.lookup<scalarList>("times", unitNone)),
123  timeDelta_(dict.lookup<scalar>("timeDelta", unitNone)),
124  begin_
125  (
126  dict.lookupOrDefault<scalar>
127  (
128  "begin",
129  unitNone,
130  mesh().time().beginTime().value()
131  )
132  ),
133  repeat_(dict.lookupOrDefault<scalar>("repeat", unitNone, 0)),
134  cycle_(dict.lookupOrDefault<scalar>("cycle", unitNone, 0)),
135  timeIndex_(-1)
136 {
137  if (repeat_ > 0 && cycle_ > 0)
138  {
140  << "Both 'repeat' and 'cycle' options specified"
141  << exit(FatalIOError);
142  }
143 
144  forAll(times_, i)
145  {
146  timeIndices_.insert(int64_t((times_[i] + timeDelta_/2.0)/timeDelta_));
147  }
148 }
149 
150 
151 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
152 
154 {}
155 
156 
157 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
158 
160 {
161  const Time& time = mesh().time();
162 
163  if
164  (
165  repeat_ > 0
166  || cycle_ > 0
167  || time.userTimeValue() + timeDelta_ < times_.last()
168  )
169  {
170  const scalar meshTime = this->meshTime();
171 
172  if (cycle_ == 0 || int((time.userTimeValue() - begin_)/cycle_) % 2 == 0)
173  {
174  forAll(times_, i)
175  {
176  if (times_[i] > meshTime + timeDelta_)
177  {
178  return time.userTimeToTime(times_[i] - meshTime);
179  }
180  }
181  }
182  else
183  {
184  forAllReverse(times_, i)
185  {
186  if (times_[i] < meshTime - timeDelta_)
187  {
188  return time.userTimeToTime(meshTime - times_[i]);
189  }
190  }
191  }
192  }
193 
194  return vGreat;
195 }
196 
197 
199 {
200  const Time& time = mesh().time();
201 
202  // Add the meshToMeshAdjustTimeStepFunctionObject functionObject
203  // if not already present
204  if
205  (
206  time.functionObjects().findObjectID("meshToMeshAdjustTimeStep")
207  == -1
208  )
209  {
210  const_cast<Time&>(time).functionObjects().append
211  (
213  (
214  "meshToMeshAdjustTimeStep",
215  time,
216  dict_
217  )
218  );
219  }
220 
221  // Only refine on the first call in a time-step
222  if (timeIndex_ != time.timeIndex())
223  {
224  timeIndex_ = time.timeIndex();
225  }
226  else
227  {
228  return false;
229  }
230 
231  // Obtain the mesh time and index from the user time
232  // repeating or cycling the mesh sequence if required
233  scalar meshTime = this->meshTime();
234  const int64_t timeIndex = int64_t((meshTime + timeDelta_/2)/timeDelta_);
235 
236  if (timeIndices_.found(timeIndex))
237  {
238  // Reset the meshTime from the timeIndex
239  // to avoid round-off errors at start of cyclic or repeat sequences
240  meshTime = timeIndex*timeDelta_;
241 
242  const word otherMeshDir = "meshToMesh_" + time.timeName(meshTime);
243 
244  Info << "Mapping to mesh " << otherMeshDir << endl;
245 
246  fvMesh otherMesh
247  (
248  IOobject
249  (
250  otherMeshDir,
251  time.constant(),
252  mesh().dbDir(),
253  time,
255  ),
256  false
257  );
258 
259  mesh().swap(otherMesh);
260 
261  fvMeshToFvMesh mapper
262  (
263  otherMesh,
264  mesh(),
265  cellsToCellss::intersection::typeName
266  );
267 
268  // Ensure the deltaCoeffs are available for constraint patch evaluation
269  mesh().deltaCoeffs();
270 
271  // Map all the volFields in the objectRegistry
272  #define mapVolFieldType(Type, nullArg) \
273  MeshToMeshMapVolFields<Type>(mesh(), mapper);
275 
276  // Map all the volFields in the objectRegistry
277  #define mapVolInternalFieldType(Type, nullArg) \
278  MeshToMeshMapVolInternalFields<Type>(mesh(), mapper);
280 
281  // Set all the surfaceFields in the objectRegistry to NaN
282  #define NaNSurfaceFieldType(Type, nullArg) \
283  NaNGeometricFields \
284  <Type, fvsPatchField, surfaceMesh> \
285  (mesh(), mapper);
287 
288  // Set all the pointFields in the objectRegistry to NaN
289  #define NaNPointFieldType(Type, nullArg) \
290  NaNGeometricFields \
291  <Type, pointPatchField, pointMesh> \
292  (mesh(), mapper);
294 
295  // Interpolate U's to Uf's
296  interpolateUfs();
297 
298  polyMeshMap map(mesh(), mapper);
299 
300  mesh().mapMesh(map);
301 
302  return true;
303  }
304  else
305  {
306  return false;
307  }
308 }
309 
310 
312 (
313  const polyTopoChangeMap& map
314 )
315 {}
316 
317 
319 {}
320 
321 
323 (
324  const polyDistributionMap& map
325 )
326 {}
327 
328 
329 // ************************************************************************* //
#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:844
scalar userTimeValue() const
Return current user time value.
Definition: Time.C:830
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Definition: Time.C:641
const functionObjectList & functionObjects() const
Return the list of function objects.
Definition: Time.H:426
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
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:99
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:346
#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:65
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
messageStream Info
const unitConversion unitNone
IOerror FatalIOError
bool isNull(const T &t)
Return true if t is a reference to the nullObject of type T.
Definition: nullObjectI.H:58
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
SurfaceField< vector > surfaceVectorField
const volVectorField & surfaceToVolVelocity(const surfaceVectorField &Uf)
Get the cell velocity field corresponding to a given face velocity, or a.
label timeIndex
Definition: getTimeIndex.H:4
dictionary dict
Return the vol-field velocity corresponding to a given surface-field velocity.