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-2025 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<class Type, class MeshMapper, class GeoMesh>
65 (
66  const MeshMapper& mapper
67 )
68 {
70  (
71  mapper.thisDb().objectRegistry::template
72  lookupClass<GeometricField<Type, GeoMesh>>()
73  );
74 
75  for
76  (
77  typename HashTable<const GeometricField<Type, GeoMesh>*>::
78  iterator fieldIter = fields.begin();
79  fieldIter != fields.end();
80  ++fieldIter
81  )
82  {
83  if (GeoMesh::Mesh::curGeometryFields.found(fieldIter()->name()))
84  continue;
85 
87  const_cast<GeometricField<Type, GeoMesh>&>(*fieldIter());
88 
89  if (&field.mesh() == &mapper.mesh())
90  {
91  if (polyMesh::debug)
92  {
93  Info<< "Mapping " << field.typeName << ' ' << field.name()
94  << endl;
95  }
96 
97  // Map the internal field
99  (
100  field.internalFieldRef(),
101  mapper
102  );
103 
104  // Map the patch fields
105  typename GeometricField<Type, GeoMesh>::Boundary& bfield =
106  field.boundaryFieldRef();
107  forAll(bfield, patchi)
108  {
109  // Cannot check sizes for patch fields because of
110  // empty fields in FV and because point fields get their size
111  // from the patch which has already been resized
112  //
113 
114  bfield[patchi].map
115  (
116  bfield[patchi],
117  mapper.boundaryMap()[patchi]
118  );
119  }
120 
121  field.instance() = field.time().name();
122  }
123  else if (polyMesh::debug)
124  {
125  Info<< "Not mapping " << field.typeName << ' ' << field.name()
126  << " since originating mesh differs from that of mapper."
127  << endl;
128  }
129  }
130 }
131 
132 
133 template<class GeoField>
134 void AddPatchFields
135 (
136  objectRegistry& obr,
137  const label insertPatchi,
138  const word& defaultPatchFieldType
139 )
140 {
141  HashTable<GeoField*> flds(obr.lookupClass<GeoField>());
142 
143  const wordList fldNames(flds.sortedToc());
144 
145  forAll(fldNames, i)
146  {
147  GeoField& fld = *flds[fldNames[i]];
148 
149  typename GeoField::Boundary& bfld = fld.boundaryFieldRef();
150 
151  if (bfld.size() != fld.mesh().boundary().size())
152  {
153  FatalErrorInFunction << "bfld size:" << bfld.size()
154  << " mesh size:" << fld.mesh().boundary().size()
155  << exit(FatalError);
156  }
157 
158  bfld.set
159  (
160  insertPatchi,
162  (
163  defaultPatchFieldType,
164  fld.mesh().boundary()[insertPatchi],
165  fld()
166  )
167  );
168  }
169 }
170 
171 
172 template<class GeoField>
174 (
175  objectRegistry& obr,
176  const labelUList& newToOld
177 )
178 {
179  HashTable<GeoField*> flds(obr.lookupClass<GeoField>());
180 
181  const wordList fldNames(flds.sortedToc());
182 
183  forAll(fldNames, i)
184  {
185  GeoField& fld = *flds[fldNames[i]];
186 
187  typename GeoField::Boundary& bfld = fld.boundaryFieldRef();
188 
189  bfld.shuffle(newToOld);
190  }
191 }
192 
193 
194 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
195 
196 } // End namespace Foam
197 
198 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
199 
200 #endif
201 
202 // ************************************************************************* //
bool found
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
const Mesh & mesh() const
Return mesh.
Internal & internalFieldRef()
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:242
const Time & time() const
Return time.
Definition: IOobject.C:315
fileName & instance() const
Return the instance directory, constant, system, <time> etc.
Definition: IOobject.C:352
const word & name() const
Return name.
Definition: IOobject.H:307
Generic internal field mapper. For "real" mapping, add template specialisations for mapping of intern...
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:334
label patchi
gmvFile<< "tracers "<< particles.size()<< nl;{ pointField positions(particles.size());label particlei=0;forAllConstIter(lagrangian::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(), lagrangian::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:234
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:258
messageStream Info
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
error FatalError
void MapGeometricFields(const MeshMapper &mapper)
Generic Geometric field mapper.
void AddPatchFields(objectRegistry &obr, const label insertPatchi, const word &defaultPatchFieldType)
void ReorderPatchFields(objectRegistry &obr, const labelUList &newToOld)