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  {
91  (*fieldIter());
92 
93  if (&field.mesh() == &mapper.mesh())
94  {
95  if (polyMesh::debug)
96  {
97  Info<< "Mapping " << field.typeName << ' ' << field.name()
98  << endl;
99  }
100 
101  // Map the internal field
103  (
104  field.ref(),
105  mapper
106  );
107 
108  // Map the patch fields
110  ::Boundary& bfield = field.boundaryFieldRef();
111  forAll(bfield, patchi)
112  {
113  // Cannot check sizes for patch fields because of
114  // empty fields in FV and because point fields get their size
115  // from the patch which has already been resized
116  //
117 
118  bfield[patchi].map
119  (
120  bfield[patchi],
121  mapper.boundaryMap()[patchi]
122  );
123  }
124 
125  field.instance() = field.time().name();
126  }
127  else if (polyMesh::debug)
128  {
129  Info<< "Not mapping " << field.typeName << ' ' << field.name()
130  << " since originating mesh differs from that of mapper."
131  << endl;
132  }
133  }
134 }
135 
136 
137 template<class GeoField>
138 void AddPatchFields
139 (
140  objectRegistry& obr,
141  const label insertPatchi,
142  const dictionary& patchFieldDict,
143  const word& defaultPatchFieldType,
144  const typename GeoField::value_type& defaultPatchValue
145 )
146 {
147  HashTable<GeoField*> flds(obr.lookupClass<GeoField>());
148 
149  const wordList fldNames(flds.sortedToc());
150 
151  forAll(fldNames, i)
152  {
153  GeoField& fld = *flds[fldNames[i]];
154 
155  typename GeoField::Boundary& bfld = fld.boundaryFieldRef();
156 
157  if (bfld.size() != fld.mesh().boundary().size())
158  {
159  FatalErrorInFunction << "bfld size:" << bfld.size()
160  << " mesh size:" << fld.mesh().boundary().size()
161  << exit(FatalError);
162  }
163 
164  if (patchFieldDict.found(fld.name()))
165  {
166  bfld.set
167  (
168  insertPatchi,
170  (
171  fld.mesh().boundary()[insertPatchi],
172  fld(),
173  patchFieldDict.subDict(fld.name())
174  )
175  );
176  }
177  else
178  {
179  bfld.set
180  (
181  insertPatchi,
183  (
184  defaultPatchFieldType,
185  fld.mesh().boundary()[insertPatchi],
186  fld()
187  )
188  );
189  bfld[insertPatchi] == defaultPatchValue;
190  }
191  }
192 }
193 
194 
195 template<class GeoField>
197 (
198  objectRegistry& obr,
199  const labelUList& newToOld
200 )
201 {
202  HashTable<GeoField*> flds(obr.lookupClass<GeoField>());
203  const wordList fldNames(flds.sortedToc());
204  forAll(fldNames, i)
205  {
206  GeoField& fld = *flds[fldNames[i]];
207  typename GeoField::Boundary& bfld = fld.boundaryFieldRef();
208 
209  bfld.shuffle(newToOld);
210  }
211 }
212 
213 
214 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
215 
216 } // End namespace Foam
217 
218 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
219 
220 #endif
221 
222 // ************************************************************************* //
#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
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)