MapGeometricFields.H
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) 2011-2021 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 Class
25  Foam::MapInternalField
26 
27 Description
28  Generic internal field mapper. For "real" mapping, add template
29  specialisations for mapping of internal fields depending on mesh
30  type.
31 
32 \*---------------------------------------------------------------------------*/
33 
34 #ifndef MapGeometricFields_H
35 #define MapGeometricFields_H
36 
37 #include "polyMesh.H"
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43 
44 template<class Type, class MeshMapper, class GeoMesh>
45 class MapInternalField
46 {
47 public:
48 
50  {}
51 
52  void operator()
53  (
55  const MeshMapper& mapper
56  ) const;
57 };
58 
59 
60 //- Generic Geometric field mapper.
61 // For "real" mapping, add template specialisations
62 // for mapping of internal fields depending on mesh type.
63 template
64 <
65  class Type,
66  template<class> class PatchField,
67  class MeshMapper,
68  class GeoMesh
69 >
71 (
72  const MeshMapper& mapper
73 )
74 {
76  (
77  mapper.thisDb().objectRegistry::template
78  lookupClass<GeometricField<Type, PatchField, GeoMesh>>()
79  );
80 
81  // It is necessary to enforce that all old-time fields are stored
82  // before the mapping is performed. Otherwise, if the
83  // old-time-level field is mapped before the field itself, sizes
84  // will not match.
85 
86  for
87  (
89  iterator fieldIter = fields.begin();
90  fieldIter != fields.end();
91  ++fieldIter
92  )
93  {
96  (*fieldIter());
97 
98  // Note: check can be removed once pointFields are actually stored on
99  // the pointMesh instead of now on the polyMesh!
100  if (&field.mesh() == &mapper.mesh())
101  {
102  field.storeOldTimes();
103  }
104  }
105 
106  for
107  (
109  iterator fieldIter = fields.begin();
110  fieldIter != fields.end();
111  ++fieldIter
112  )
113  {
116  (*fieldIter());
117 
118  if (&field.mesh() == &mapper.mesh())
119  {
120  if (polyMesh::debug)
121  {
122  Info<< "Mapping " << field.typeName << ' ' << field.name()
123  << endl;
124  }
125 
126  // Map the internal field
128  (
129  field.ref(),
130  mapper
131  );
132 
133  // Map the patch fields
135  ::Boundary& bfield = field.boundaryFieldRef();
136  forAll(bfield, patchi)
137  {
138  // Cannot check sizes for patch fields because of
139  // empty fields in FV and because point fields get their size
140  // from the patch which has already been resized
141  //
142 
143  bfield[patchi].autoMap(mapper.boundaryMap()[patchi]);
144  }
145 
146  field.instance() = field.time().timeName();
147  }
148  else if (polyMesh::debug)
149  {
150  Info<< "Not mapping " << field.typeName << ' ' << field.name()
151  << " since originating mesh differs from that of mapper."
152  << endl;
153  }
154  }
155 }
156 
157 
158 template<class GeoField>
160 (
161  objectRegistry& obr,
162  const label insertPatchi,
163  const dictionary& patchFieldDict,
164  const word& defaultPatchFieldType,
165  const typename GeoField::value_type& defaultPatchValue
166 )
167 {
168  HashTable<GeoField*> flds(obr.lookupClass<GeoField>());
169 
170  const wordList fldNames(flds.sortedToc());
171 
172  forAll(fldNames, i)
173  {
174  GeoField& fld = *flds[fldNames[i]];
175 
176  typename GeoField::Boundary& bfld = fld.boundaryFieldRef();
177 
178  if (bfld.size() != fld.mesh().boundary().size())
179  {
180  FatalErrorInFunction << "bfld size:" << bfld.size()
181  << " mesh size:" << fld.mesh().boundary().size()
182  << exit(FatalError);
183  }
184 
185  if (patchFieldDict.found(fld.name()))
186  {
187  bfld.set
188  (
189  insertPatchi,
191  (
192  fld.mesh().boundary()[insertPatchi],
193  fld(),
194  patchFieldDict.subDict(fld.name())
195  )
196  );
197  }
198  else
199  {
200  bfld.set
201  (
202  insertPatchi,
204  (
205  defaultPatchFieldType,
206  fld.mesh().boundary()[insertPatchi],
207  fld()
208  )
209  );
210  bfld[insertPatchi] == defaultPatchValue;
211  }
212  }
213 }
214 
215 
216 template<class GeoField>
218 (
219  objectRegistry& obr,
220  const labelUList& newToOld
221 )
222 {
223  HashTable<GeoField*> flds(obr.lookupClass<GeoField>());
224  const wordList fldNames(flds.sortedToc());
225  forAll(fldNames, i)
226  {
227  GeoField& fld = *flds[fldNames[i]];
228  typename GeoField::Boundary& bfld = fld.boundaryFieldRef();
229 
230  bfld.shuffle(newToOld);
231  }
232 }
233 
234 
235 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
236 
237 } // End namespace Foam
238 
239 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
240 
241 #endif
242 
243 // ************************************************************************* //
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:643
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const word & name() const
Return name.
Definition: IOobject.H:303
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:112
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
static const char *const typeName
Definition: Field.H:105
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Generic GeometricField class.
Info<< "Calculating turbulent flame speed field St\"<< endl;volScalarField St(IOobject("St", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:228
void MapGeometricFields(const MeshMapper &mapper)
Generic Geometric field mapper.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:982
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
void storeOldTimes() const
Store the old-time fields.
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
A class for handling words, derived from string.
Definition: word.H:59
void AddPatchFields(objectRegistry &obr, const label insertPatchi, const dictionary &patchFieldDict, const word &defaultPatchFieldType, const typename GeoField::value_type &defaultPatchValue)
iterator begin()
Iterator set to the beginning of the HashTable.
Definition: HashTableI.H:411
An STL-conforming hash table.
Definition: HashTable.H:61
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
const Mesh & mesh() const
Return mesh.
Internal & ref()
Return a reference to the dimensioned internal field.
Generic internal field mapper. For "real" mapping, add template specialisations for mapping of intern...
const fileName & instance() const
Definition: IOobject.H:390
label patchi
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
const Time & time() const
Return time.
Definition: IOobject.C:318
List< Key > sortedToc() const
Return the table of contents as a sorted list.
Definition: HashTable.C:217
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
messageStream Info
void ReorderPatchFields(objectRegistry &obr, const labelUList &newToOld)
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:46
rDeltaTY field()
Registry of regIOobjects.
HashTable< const Type * > lookupClass(const bool strict=false) const
Lookup and return all objects of the given Type.
Namespace for OpenFOAM.