34 #ifndef MapGeometricFields_H
35 #define MapGeometricFields_H
44 template<
class Type,
class MeshMapper,
class GeoMesh>
55 const MeshMapper& mapper
63 template<
class Type,
class MeshMapper,
class GeoMesh>
66 const MeshMapper& mapper
71 mapper.thisDb().objectRegistry::template
78 iterator fieldIter =
fields.begin();
83 if (GeoMesh::Mesh::curGeometryFields.
found(fieldIter()->
name()))
89 if (&field.
mesh() == &mapper.mesh())
93 Info<<
"Mapping " << field.typeName <<
' ' << field.
name()
117 mapper.boundaryMap()[
patchi]
123 else if (polyMesh::debug)
125 Info<<
"Not mapping " << field.typeName <<
' ' << field.
name()
126 <<
" since originating mesh differs from that of mapper."
133 template<
class GeoField>
137 const label insertPatchi,
138 const word& defaultPatchFieldType
147 GeoField&
fld = *flds[fldNames[i]];
149 typename GeoField::Boundary& bfld =
fld.boundaryFieldRef();
151 if (bfld.size() !=
fld.mesh().boundary().size())
154 <<
" mesh size:" <<
fld.mesh().boundary().size()
163 defaultPatchFieldType,
164 fld.mesh().boundary()[insertPatchi],
172 template<
class GeoField>
185 GeoField&
fld = *flds[fldNames[i]];
187 typename GeoField::Boundary& bfld =
fld.boundaryFieldRef();
189 bfld.shuffle(newToOld);
#define forAll(list, i)
Loop across all elements in list.
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.
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...
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(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
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 LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
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)