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 
98  HashTable<surfaceVectorField*> Ufs
99  (
100  mesh().lookupClass<surfaceVectorField>()
101  );
102 
103  forAllIter(HashTable<surfaceVectorField*>, Ufs, iter)
104  {
105  surfaceVectorField& Uf = *iter();
106 
108 
109  if (!isNull(U))
110  {
111  Uf.reset(fvc::interpolate(U));
112  }
113  }
114 }
115 
116 
117 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
118 
120 (
121  fvMesh& mesh,
122  const dictionary& dict
123 )
124 :
125  fvMeshTopoChanger(mesh),
126  dict_(dict),
127  times_(dict.lookup("times")),
128  timeDelta_(dict.lookup<scalar>("timeDelta")),
129  begin_(dict.lookupOrDefault("begin", mesh().time().beginTime().value())),
130  repeat_(dict.lookupOrDefault("repeat", 0.0)),
131  cycle_(dict.lookupOrDefault("cycle", 0.0)),
132  timeIndex_(-1)
133 {
134  if (repeat_ > 0 && cycle_ > 0)
135  {
137  << "Both 'repeat' and 'cycle' options specified"
138  << exit(FatalIOError);
139  }
140 
141  forAll(times_, i)
142  {
143  timeIndices_.insert(int64_t((times_[i] + timeDelta_/2.0)/timeDelta_));
144  }
145 }
146 
147 
148 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
149 
151 {}
152 
153 
154 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
155 
157 {
158  const Time& time = mesh().time();
159 
160  if
161  (
162  repeat_ > 0
163  || cycle_ > 0
164  || time.userTimeValue() + timeDelta_ < times_.last()
165  )
166  {
167  const scalar meshTime = this->meshTime();
168 
169  if (cycle_ == 0 || int((time.userTimeValue() - begin_)/cycle_) % 2 == 0)
170  {
171  forAll(times_, i)
172  {
173  if (times_[i] > meshTime + timeDelta_)
174  {
175  return time.userTimeToTime(times_[i] - meshTime);
176  }
177  }
178  }
179  else
180  {
181  forAllReverse(times_, i)
182  {
183  if (times_[i] < meshTime - timeDelta_)
184  {
185  return time.userTimeToTime(meshTime - times_[i]);
186  }
187  }
188  }
189  }
190 
191  return vGreat;
192 }
193 
194 
196 {
197  const Time& time = mesh().time();
198 
199  // Add the meshToMeshAdjustTimeStepFunctionObject functionObject
200  // if not already present
201  if
202  (
203  time.functionObjects().findObjectID("meshToMeshAdjustTimeStep")
204  == -1
205  )
206  {
207  const_cast<Time&>(time).functionObjects().append
208  (
210  (
211  "meshToMeshAdjustTimeStep",
212  time,
213  dict_
214  )
215  );
216  }
217 
218  // Only refine on the first call in a time-step
219  if (timeIndex_ != time.timeIndex())
220  {
221  timeIndex_ = time.timeIndex();
222  }
223  else
224  {
225  return false;
226  }
227 
228  // Obtain the mesh time from the user time
229  // repeating or cycling the mesh sequence if required
230 
231  const scalar meshTime = this->meshTime();
232 
233  if (timeIndices_.found((meshTime + timeDelta_/2)/timeDelta_))
234  {
235  const word otherMeshDir =
236  "meshToMesh_" + time.timeName(meshTime);
237 
238  Info << "Mapping to mesh " << otherMeshDir << endl;
239 
240  fvMesh otherMesh
241  (
242  IOobject
243  (
244  otherMeshDir,
245  time.constant(),
246  time,
248  ),
249  false,
251  );
252 
253  mesh().swap(otherMesh);
254 
255  fvMeshToFvMesh mapper
256  (
257  otherMesh,
258  mesh(),
259  cellsToCellss::intersection::typeName
260  );
261 
262  // Map all the volFields in the objectRegistry
263  #define mapVolFieldType(Type, nullArg) \
264  MeshToMeshMapVolFields<Type>(mesh(), mapper);
266 
267  // Map all the volFields in the objectRegistry
268  #define mapVolInternalFieldType(Type, nullArg) \
269  MeshToMeshMapVolInternalFields<Type>(mesh(), mapper);
271 
272  // Set all the surfaceFields in the objectRegistry to NaN
273  #define NaNSurfaceFieldType(Type, nullArg) \
274  NaNGeometricFields \
275  <Type, fvsPatchField, surfaceMesh, setSizeFvPatchFieldMapper> \
276  (mesh(), mapper);
278 
279  // Set all the pointFields in the objectRegistry to NaN
280  #define NaNPointFieldType(Type, nullArg) \
281  NaNGeometricFields \
282  <Type, pointPatchField, pointMesh, setSizePointPatchFieldMapper> \
283  (mesh(), mapper);
285 
286  // Interpolate U's to Uf's
287  interpolateUfs();
288 
289  polyMeshMap map(mesh(), mapper);
290 
291  mesh().mapMesh(map);
292 
293  return true;
294  }
295  else
296  {
297  return false;
298  }
299 }
300 
301 
303 (
304  const polyTopoChangeMap& map
305 )
306 {}
307 
308 
310 {}
311 
312 
314 (
315  const polyDistributionMap& map
316 )
317 {}
318 
319 
320 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
#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.