34 #ifndef MapGeometricFields_H
35 #define MapGeometricFields_H
44 template<
class Type,
class MeshMapper,
class GeoMesh>
55 const MeshMapper& mapper
66 template<
class>
class PatchField,
72 const MeshMapper& mapper
77 mapper.thisDb().objectRegistry::template
84 iterator fieldIter =
fields.begin();
89 if (GeoMesh::Mesh::geometryFields.
found(fieldIter()->
name()))
continue;
95 if (&field.
mesh() == &mapper.mesh())
123 mapper.boundaryMap()[
patchi]
129 else if (polyMesh::debug)
132 <<
" since originating mesh differs from that of mapper."
139 template<
class GeoField>
143 const label insertPatchi,
145 const word& defaultPatchFieldType,
146 const typename GeoField::value_type& defaultPatchValue
155 GeoField&
fld = *flds[fldNames[i]];
157 typename GeoField::Boundary& bfld =
fld.boundaryFieldRef();
159 if (bfld.size() !=
fld.mesh().boundary().size())
162 <<
" mesh size:" <<
fld.mesh().boundary().size()
166 if (patchFieldDict.
found(
fld.name()))
173 fld.mesh().boundary()[insertPatchi],
186 defaultPatchFieldType,
187 fld.mesh().boundary()[insertPatchi],
191 bfld[insertPatchi] == defaultPatchValue;
197 template<
class GeoField>
208 GeoField&
fld = *flds[fldNames[i]];
209 typename GeoField::Boundary& bfld =
fld.boundaryFieldRef();
211 bfld.shuffle(newToOld);
#define forAll(list, i)
Loop across all elements in list.
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
static const char *const typeName
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Generic GeometricField class.
Internal & internalFieldRef()
Return a reference to the dimensioned internal field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
An STL-conforming hash table.
List< Key > sortedToc() const
Return the table of contents as a sorted list.
const Time & time() const
Return time.
fileName & instance() const
Return the instance directory, constant, system, <time> etc.
const word & name() const
Return name.
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....
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
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.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
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
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
errorManipArg< error, int > exit(error &err, const int errNo=1)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Ostream & endl(Ostream &os)
Add newline and flush stream.
word name(const bool)
Return a word representation of a bool.
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.
void ReorderPatchFields(objectRegistry &obr, const labelUList &newToOld)