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 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 "meshToMesh.H"
33 #include "cellVolumeWeightMethod.H"
34 #include "surfaceToVolVelocity.H"
36 #include "polyMeshMap.H"
37 #include "processorPolyPatch.H"
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 namespace Foam
43 {
44 namespace fvMeshTopoChangers
45 {
46  defineTypeNameAndDebug(meshToMesh, 0);
47  addToRunTimeSelectionTable(fvMeshTopoChanger, meshToMesh, fvMesh);
48 }
49 }
50 
51 
52 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
53 
54 void Foam::fvMeshTopoChangers::meshToMesh::interpolateUfs()
55 {
56  // Interpolate U to Uf
57 
58  HashTable<surfaceVectorField*> Ufs
59  (
60  mesh().lookupClass<surfaceVectorField>()
61  );
62 
63  forAllIter(HashTable<surfaceVectorField*>, Ufs, iter)
64  {
65  surfaceVectorField& Uf = *iter();
66 
67  const volVectorField& U = surfaceToVolVelocity(Uf);
68 
69  if (!isNull(U))
70  {
71  Uf.reset(fvc::interpolate(U));
72  }
73  }
74 }
75 
76 
77 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
78 
80 (
81  fvMesh& mesh,
82  const dictionary& dict
83 )
84 :
85  fvMeshTopoChanger(mesh),
86  dict_(dict),
87  cuttingPatches_(dict.lookupOrDefault("cuttingPatches", wordList::null())),
88  times_(dict.lookup("times")),
89  timeDelta_(dict.lookup<scalar>("timeDelta")),
90  timeIndex_(-1)
91 {
92  forAll(times_, i)
93  {
94  timeIndices_.insert(label(times_[i]/timeDelta_));
95  }
96 }
97 
98 
99 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
100 
102 {}
103 
104 
105 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
106 
108 {
109  if (timeIndex_ == -1)
110  {
111  const_cast<Time&>(mesh().time()).functionObjects().append
112  (
114  (
115  "meshToMeshAdjustTimeStep",
116  mesh().time(),
117  dict_
118  )
119  );
120  }
121 
122  bool hasChanged = false;
123 
124  // Only refine on the first call in a time-step
125  if (timeIndex_ != mesh().time().timeIndex())
126  {
127  timeIndex_ = mesh().time().timeIndex();
128  }
129  else
130  {
131  return hasChanged;
132  }
133 
134  const scalar userTime = mesh().time().userTimeValue();
135 
136  if (timeIndices_.found((userTime + timeDelta_/2)/timeDelta_))
137  {
138  const word meshDir = "meshToMesh_" + mesh().time().timeName(userTime);
139 
140  Info << "Mapping to mesh " << meshDir << endl;
141 
142  hasChanged = true;
143 
144  fvMesh newMesh
145  (
146  IOobject
147  (
148  meshDir,
149  mesh().time().constant(),
150  mesh().time(),
152  ),
153  false,
155  );
156 
158 
159  // Create mesh-to-mesh mapper with support for cuttingPatches
160  // if specified
161  if (cuttingPatches_.size())
162  {
163  HashSet<word> cuttingPatchTable;
164  forAll(cuttingPatches_, i)
165  {
166  cuttingPatchTable.insert(cuttingPatches_[i]);
167  }
168 
169  HashTable<word> patchMap(mesh().boundary().size());
170 
171  const polyBoundaryMesh& pbm = mesh().boundaryMesh();
172 
173  forAll(pbm, i)
174  {
175  if
176  (
177  !cuttingPatchTable.found(pbm[i].name())
178  && !isA<processorPolyPatch>(pbm[i])
179  )
180  {
181  patchMap.insert(pbm[i].name(), pbm[i].name());
182  }
183  }
184 
185  mapper = new Foam::meshToMesh
186  (
187  mesh(),
188  newMesh,
189  cellVolumeWeightMethod::typeName,
190  patchMap,
191  cuttingPatches_
192  );
193  }
194  else
195  {
196  mapper = new Foam::meshToMesh
197  (
198  mesh(),
199  newMesh,
200  cellVolumeWeightMethod::typeName
201  );
202  }
203 
204  mesh().reset(newMesh);
205 
206  mesh().deltaCoeffs();
207 
208  // Map all the volFields in the objectRegistry
209  #define mapVolFieldType(Type, nullArg) \
210  MeshToMeshMapVolFields<Type>(mapper);
212 
213  // Set all the surfaceFields in the objectRegistry to NaN
214  #define NaNSurfaceFieldType(Type, nullArg) \
215  NaNGeometricFields \
216  <Type, fvsPatchField, surfaceMesh, fvPatchFieldMapper>(mapper);
218 
219  // Set all the pointFields in the objectRegistry to NaN
220  #define NaNPointFieldType(Type, nullArg) \
221  NaNGeometricFields \
222  <Type, pointPatchField, pointMesh, pointPatchFieldMapper>(mapper);
224 
225  // Interpolate U's to Uf's
226  interpolateUfs();
227 
228  polyMeshMap map(mesh());
229  mesh().mapMesh(map);
230  }
231 
232  return hasChanged;
233 }
234 
235 
237 (
238  const polyTopoChangeMap& map
239 )
240 {}
241 
242 
244 {}
245 
246 
248 (
249  const polyDistributionMap& map
250 )
251 {}
252 
253 
254 // ************************************************************************* //
A HashTable with keys but without contents.
Definition: HashSet.H:59
Return the vol-field velocity corresponding to a given surface-field velocity.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const word & name() const
Return name.
Definition: IOobject.H:315
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define mapVolFieldType(Type, nullArg)
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const volVectorField & surfaceToVolVelocity(const surfaceVectorField &Uf)
Get the cell velocity field corresponding to a given face velocity, or a.
Class to calculate the cell-addressing between two overlapping meshes.
Definition: meshToMesh.H:60
static const List< T > & null()
Return a null List.
Definition: ListI.H:118
addToRunTimeSelectionTable(fvMeshTopoChanger, list, fvMesh)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:59
fvMesh & mesh
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:111
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Macros for easy insertion into run-time selection tables.
virtual bool update()
Update the mesh for both mesh motion and topology change.
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
bool isNull(const T &t)
Return true if t is a reference to the nullObject of type T.
Definition: nullObjectI.H:46
virtual void distribute(const polyDistributionMap &)
Update corresponding to the given distribution map.
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:113
A class for handling words, derived from string.
Definition: word.H:59
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
faceListList boundary(nPatches)
An STL-conforming hash table.
Definition: HashTable.H:61
Foam::polyBoundaryMesh.
virtual void topoChange(const polyTopoChangeMap &)
Update corresponding to the given map.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
meshToMesh(fvMesh &, const dictionary &dict)
Construct from fvMesh and dictionary.
#define NaNSurfaceFieldType(Type, nullArg)
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
label timeIndex
Definition: getTimeIndex.H:4
messageStream Info
FOR_ALL_FIELD_TYPES(DefineFvWallInfoType)
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:50
#define NaNPointFieldType(Type, nullArg)
An abstract class for the user time description.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
Abstract base class for fvMesh movers.
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864