meshToMesh_fvMeshTopoChanger.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-2025 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 "surfaceInterpolate.H"
28 #include "pointFields.H"
31 #include "surfaceToVolVelocity.H"
33 #include "polyMeshMap.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace fvMeshTopoChangers
41 {
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
49 
50 bool Foam::fvMeshTopoChangers::meshToMesh::forward() const
51 {
52  return
53  cycle_ == 0
54  || int((mesh().time().userTimeValue() - begin_)/cycle_) % 2 == 0;
55 }
56 
57 
58 Foam::scalar Foam::fvMeshTopoChangers::meshToMesh::meshTime() const
59 {
60  const Time& time = mesh().time();
61 
62  if (repeat_ > 0)
63  {
64  return begin_ + fmod(time.userTimeValue() - begin_, repeat_);
65  }
66  else if (cycle_ > 0)
67  {
68  if (forward())
69  {
70  return
71  begin_
72  + fmod(time.userTimeValue() - begin_, cycle_);
73  }
74  else
75  {
76  return
77  begin_ + cycle_
78  - fmod(time.userTimeValue() - begin_, cycle_);
79  }
80  }
81  else
82  {
83  return time.userTimeValue();
84  }
85 }
86 
87 
88 void Foam::fvMeshTopoChangers::meshToMesh::interpolateUfs()
89 {
90  // Interpolate U to Uf
91  UPtrList<surfaceVectorField> Ufs(mesh().curFields<surfaceVectorField>());
92 
93  forAll(Ufs, i)
94  {
95  surfaceVectorField& Uf = Ufs[i];
96 
98 
99  if (!isNull(U))
100  {
101  Uf.reset(fvc::interpolate(U));
102  }
103  }
104 }
105 
106 
107 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
108 
110 (
111  fvMesh& mesh,
112  const dictionary& dict
113 )
114 :
116  times_(dict.lookup<scalarList>("times", unitNone)),
117  timeDelta_(dict.lookup<scalar>("timeDelta", unitNone)),
118  begin_
119  (
120  dict.lookupOrDefault<scalar>
121  (
122  "begin",
123  unitNone,
124  mesh().time().beginTime().value()
125  )
126  ),
127  repeat_(dict.lookupOrDefault<scalar>("repeat", unitNone, 0)),
128  cycle_(dict.lookupOrDefault<scalar>("cycle", unitNone, 0)),
129  timeIndex_(-1),
130  mapped_(false)
131 {
132  if (repeat_ > 0 && cycle_ > 0)
133  {
135  << "Both 'repeat' and 'cycle' options specified"
136  << exit(FatalIOError);
137  }
138 
139  forAll(times_, i)
140  {
141  timeIndices_.insert(int64_t((times_[i] + timeDelta_/2.0)/timeDelta_));
142  }
143 }
144 
145 
146 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
147 
149 {}
150 
151 
152 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
153 
155 {
156  const Time& time = mesh().time();
157 
158  if
159  (
160  repeat_ > 0
161  || cycle_ > 0
162  || time.userTimeValue() + timeDelta_ < times_.last()
163  )
164  {
165  const scalar meshTime = this->meshTime();
166 
167  if (cycle_ == 0 || int((time.userTimeValue() - begin_)/cycle_) % 2 == 0)
168  {
169  forAll(times_, i)
170  {
171  if (times_[i] > meshTime + timeDelta_)
172  {
173  return time.userTimeToTime(times_[i] - meshTime);
174  }
175  }
176  }
177  else
178  {
179  forAllReverse(times_, i)
180  {
181  if (times_[i] < meshTime - timeDelta_)
182  {
183  return time.userTimeToTime(meshTime - times_[i]);
184  }
185  }
186  }
187  }
188 
189  return vGreat;
190 }
191 
192 
194 {
195  return mapped_;
196 }
197 
198 
200 {
201  mapped_ = false;
202 
203  const Time& time = mesh().time();
204 
205  // Add the meshToMeshAdjustTimeStepFunctionObject functionObject
206  // if not already present
207  if
208  (
209  time.functionObjects().findObjectID("meshToMeshAdjustTimeStep")
210  == -1
211  )
212  {
213  const_cast<Time&>(time).functionObjects().append
214  (
216  (
217  "meshToMeshAdjustTimeStep",
218  mesh()
219  )
220  );
221  }
222 
223  // Only map on the first call in a time-step
224  if (timeIndex_ != time.timeIndex())
225  {
226  timeIndex_ = time.timeIndex();
227  }
228  else
229  {
230  return false;
231  }
232 
233  // Obtain the mesh time and index from the user time
234  // repeating or cycling the mesh sequence if required
235  scalar meshTime = this->meshTime();
236  const int64_t timeIndex = int64_t((meshTime + timeDelta_/2)/timeDelta_);
237 
238  if (timeIndices_.found(timeIndex))
239  {
240  // Reset the meshTime from the timeIndex
241  // to avoid round-off errors at start of cyclic or repeat sequences
242  meshTime = timeIndex*timeDelta_;
243 
244  const word meshTimeName(time.timeName(meshTime));
245 
246  Info<< "Mapping to mesh time " << meshTimeName << endl;
247 
248  fvMesh otherMesh
249  (
250  IOobject
251  (
252  mesh().dbDir(),
253  time.constant(),
254  "meshes"/meshTimeName,
255  time,
257  ),
258  false
259  );
260 
261  mesh().preChange();
262 
263  mesh().swap(otherMesh);
264 
265  fvMeshToFvMesh mapper
266  (
267  otherMesh,
268  mesh(),
269  cellsToCellss::intersection::typeName
270  );
271 
272  // Ensure the deltaCoeffs are available for constraint patch evaluation
273  mesh().deltaCoeffs();
274 
275  // Map all the volFields in the objectRegistry
276  #define mapVolFieldType(Type, nullArg) \
277  MeshToMeshMapVolFields<Type>(mesh(), mapper);
279 
280  // Map all the volFields in the objectRegistry
281  #define mapVolInternalFieldType(Type, nullArg) \
282  MeshToMeshMapVolInternalFields<Type>(mesh(), mapper);
284 
285  // Set all the surfaceFields in the objectRegistry to NaN
286  #define NaNSurfaceFieldType(Type, nullArg) \
287  NaNGeometricFields<Type, surfaceMesh>(mesh());
289 
290  // Set all the pointFields in the objectRegistry to NaN
291  #define NaNPointFieldType(Type, nullArg) \
292  NaNGeometricFields<Type, pointMesh>(mesh());
294 
295  // Interpolate U's to Uf's
296  interpolateUfs();
297 
298  polyMeshMap map(mesh(), mapper);
299 
300  mesh().mapMesh(map);
301 
302  mapped_ = true;
303 
304  return true;
305  }
306  else
307  {
308  return false;
309  }
310 }
311 
312 
314 (
315  const polyTopoChangeMap& map
316 )
317 {}
318 
319 
321 {}
322 
323 
325 (
326  const polyDistributionMap& map
327 )
328 {}
329 
330 
331 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: UList.H:445
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:833
scalar userTimeValue() const
Return current user time value.
Definition: Time.C:819
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:420
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
label findObjectID(const word &name) const
Find the ID of a given function object by name.
Adjusts time-step for meshToMesh mapping.
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.
bool mapped() const
Return true if the mesh has been mapped this time-step,.
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:96
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
void preChange()
Prepare for a mesh change.
Definition: fvMesh.C:1222
void swap(fvMesh &)
Swap mesh.
Definition: fvMesh.C:819
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: fvMesh.C:1429
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.
const surfaceScalarField & deltaCoeffs() const
Return reference to cell-centre difference coefficients.
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
U
Definition: pEqn.H:72
#define NaNPointFieldType(Type, nullArg)
#define mapVolInternalFieldType(Type, nullArg)
#define NaNSurfaceFieldType(Type, nullArg)
#define mapVolFieldType(Type, nullArg)
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:63
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
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.