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-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 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
79  );
80 
81  for
82  (
84  iterator fieldIter = fields.begin();
85  fieldIter != fields.end();
86  ++fieldIter
87  )
88  {
89  if (GeoMesh::Mesh::geometryFields.found(fieldIter()->name())) continue;
90 
93  (*fieldIter());
94 
95  if (&field.mesh() == &mapper.mesh())
96  {
97  if (polyMesh::debug)
98  {
99  Info<< "Mapping " << field.typeName << ' ' << field.name()
100  << endl;
101  }
102 
103  // Map the internal field
105  (
106  field.ref(),
107  mapper
108  );
109 
110  // Map the patch fields
112  ::Boundary& bfield = field.boundaryFieldRef();
113  forAll(bfield, patchi)
114  {
115  // Cannot check sizes for patch fields because of
116  // empty fields in FV and because point fields get their size
117  // from the patch which has already been resized
118  //
119 
120  bfield[patchi].map
121  (
122  bfield[patchi],
123  mapper.boundaryMap()[patchi]
124  );
125  }
126 
127  field.instance() = field.time().name();
128  }
129  else if (polyMesh::debug)
130  {
131  Info<< "Not mapping " << field.typeName << ' ' << field.name()
132  << " since originating mesh differs from that of mapper."
133  << endl;
134  }
135  }
136 }
137 
138 
139 template<class GeoField>
140 void AddPatchFields
141 (
142  objectRegistry& obr,
143  const label insertPatchi,
144  const dictionary& patchFieldDict,
145  const word& defaultPatchFieldType,
146  const typename GeoField::value_type& defaultPatchValue
147 )
148 {
149  HashTable<GeoField*> flds(obr.lookupClass<GeoField>());
150 
151  const wordList fldNames(flds.sortedToc());
152 
153  forAll(fldNames, i)
154  {
155  GeoField& fld = *flds[fldNames[i]];
156 
157  typename GeoField::Boundary& bfld = fld.boundaryFieldRef();
158 
159  if (bfld.size() != fld.mesh().boundary().size())
160  {
161  FatalErrorInFunction << "bfld size:" << bfld.size()
162  << " mesh size:" << fld.mesh().boundary().size()
163  << exit(FatalError);
164  }
165 
166  if (patchFieldDict.found(fld.name()))
167  {
168  bfld.set
169  (
170  insertPatchi,
172  (
173  fld.mesh().boundary()[insertPatchi],
174  fld(),
175  patchFieldDict.subDict(fld.name())
176  )
177  );
178  }
179  else
180  {
181  bfld.set
182  (
183  insertPatchi,
185  (
186  defaultPatchFieldType,
187  fld.mesh().boundary()[insertPatchi],
188  fld()
189  )
190  );
191  bfld[insertPatchi] == defaultPatchValue;
192  }
193  }
194 }
195 
196 
197 template<class GeoField>
199 (
200  objectRegistry& obr,
201  const labelUList& newToOld
202 )
203 {
204  HashTable<GeoField*> flds(obr.lookupClass<GeoField>());
205  const wordList fldNames(flds.sortedToc());
206  forAll(fldNames, i)
207  {
208  GeoField& fld = *flds[fldNames[i]];
209  typename GeoField::Boundary& bfld = fld.boundaryFieldRef();
210 
211  bfld.shuffle(newToOld);
212  }
213 }
214 
215 
216 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
217 
218 } // End namespace Foam
219 
220 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
221 
222 #endif
223 
224 // ************************************************************************* //
bool found
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const Mesh & mesh() const
Return mesh.
void map(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 map from the given field
Definition: Field.C:265
static const char *const typeName
Definition: Field.H:105
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:47
Generic GeometricField class.
Internal & ref()
Return a reference to the dimensioned internal field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
An STL-conforming hash table.
Definition: HashTable.H:127
List< Key > sortedToc() const
Return the table of contents as a sorted list.
Definition: HashTable.C:217
const Time & time() const
Return time.
Definition: IOobject.C:318
fileName & instance() const
Return the instance directory, constant, system, <time> etc.
Definition: IOobject.C:355
const word & name() const
Return name.
Definition: IOobject.H:310
Generic internal field mapper. For "real" mapping, add template specialisations for mapping of intern...
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:998
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:659
const word & name() const
Return const reference to name.
Registry of regIOobjects.
HashTable< const Type * > lookupClass(const bool strict=false) const
Lookup and return all objects of the given Type.
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
label patchi
gmvFile<< "tracers "<< particles.size()<< nl;{ pointField positions(particles.size());label particlei=0;forAllConstIter(Cloud< passiveParticle >, particles, iter) { positions[particlei++]=iter().position(mesh);} for(i=0;i< pTraits< point >::nComponents;i++) { forAll(positions, particlei) { gmvFile<< component(positions[particlei], i)<< ' ';} gmvFile<< nl;}}forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:230
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
messageStream Info
void AddPatchFields(objectRegistry &obr, const label insertPatchi, const dictionary &patchFieldDict, const word &defaultPatchFieldType, const typename GeoField::value_type &defaultPatchValue)
void MapGeometricFields(const MeshMapper &mapper)
Generic Geometric field mapper.
error FatalError
void ReorderPatchFields(objectRegistry &obr, const labelUList &newToOld)