31 template<
class GeoField>
36 mesh.objectRegistry::lookupClass<GeoField>()
41 const GeoField&
fld = *iter();
43 Pout<<
"Field:" << iter.key() <<
" internal size:" << fld.size()
49 <<
' ' << fld.boundaryField()[
patchi].patch().name()
50 <<
' ' << fld.boundaryField()[
patchi].type()
51 <<
' ' << fld.boundaryField()[
patchi].size()
58 template<
class T,
class Mesh>
59 void Foam::fvMeshDistribute::saveBoundaryFields
70 static_cast<const fvMesh&>(mesh_).objectRegistry::lookupClass<fldType>()
73 bflds.setSize(flds.
size());
79 const fldType&
fld = *iter();
81 bflds.set(i, fld.boundaryField().clone().ptr());
88 template<
class T,
class Mesh>
89 void Foam::fvMeshDistribute::mapBoundaryFields
104 mesh_.objectRegistry::lookupClass<fldType>()
107 if (flds.
size() != oldBflds.size())
117 fldType&
fld = *iter();
118 typename fldType::Boundary& bfld =
119 fld.boundaryFieldRef();
132 label oldFacei = faceMap[facei++];
135 forAll(oldPatchStarts, oldPatchi)
137 label oldLocalI = oldFacei - oldPatchStarts[oldPatchi];
139 if (oldLocalI >= 0 && oldLocalI < oldBfld[oldPatchi].size())
141 patchFld[i] = oldBfld[oldPatchi][oldLocalI];
151 void Foam::fvMeshDistribute::saveInternalFields
160 static_cast<const fvMesh&>(mesh_).objectRegistry::lookupClass<fldType>()
163 iflds.setSize(flds.
size());
169 const fldType&
fld = *iter();
171 iflds.set(i, fld.primitiveField().clone());
179 void Foam::fvMeshDistribute::mapExposedFaces
193 mesh_.objectRegistry::lookupClass<fldType>()
196 if (flds.
size() != oldFlds.size())
198 FatalErrorIn(
"fvMeshDistribute::mapExposedFaces(..)") <<
"problem" 207 fldType&
fld = *iter();
208 typename fldType::Boundary& bfld = fld.boundaryFieldRef();
209 const bool negateIfFlipped =
isFlux(fld);
211 const Field<T>& oldInternal = oldFlds[fieldI++];
221 const label facei = patchFld.
patch().start()+i;
222 const label oldFacei = faceMap[facei];
224 if (oldFacei < oldInternal.
size())
226 patchFld[i] = oldInternal[oldFacei];
230 patchFld[i] =
flipOp()(patchFld[i]);
239 template<
class GeoField>
240 void Foam::fvMeshDistribute::correctProcessorPatchFields()
243 processorPatchFieldType;
247 mesh_.objectRegistry::lookupClass<GeoField>()
252 GeoField&
fld = *iter();
254 typename GeoField::Boundary& bfld = fld.boundaryFieldRef();
266 if (isA<processorPatchFieldType>(bfld[
patchi]))
284 if (isA<processorPatchFieldType>(bfld[
patchi]))
295 forAll(patchSchedule, patchEvali)
297 if (isA<processorPatchFieldType>(bfld[patchEvali]))
299 if (patchSchedule[patchEvali].init)
301 bfld[patchSchedule[patchEvali].patch]
306 bfld[patchSchedule[patchEvali].patch]
316 template<
class GeoField>
317 void Foam::fvMeshDistribute::sendFields
349 Pout<<
"Subsetting field " << fieldNames[i]
350 <<
" for domain:" << domain <<
endl;
355 const GeoField&
fld =
369 template<
class GeoField>
370 void Foam::fvMeshDistribute::receiveFields
374 typename GeoField::Mesh&
mesh,
381 Pout<<
"Receiving fields " << fieldNames
382 <<
" from domain:" << domain <<
endl;
391 Pout<<
"Constructing field " << fieldNames[i]
392 <<
" from domain:" << domain <<
endl;
403 mesh.thisDb().time().timeName(),
409 fieldDicts.
subDict(fieldNames[i])
Class containing functor to negate primitives. Dummy for all other types.
#define forAll(list, i)
Loop across all elements in list.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
bool set(const label) const
Is element set.
A list of keyword definitions, which are a keyword followed by any number of values (e...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
const lduSchedule & patchSchedule() const
Order in which the patches should be initialised/evaluated.
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
void size(const label)
Override size to be inconsistent with allocated storage.
Ostream & endl(Ostream &os)
Add newline and flush stream.
static label nRequests()
Get number of outstanding requests.
Generic GeometricField class.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
label size() const
Return number of elements in table.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
const labelHashSet & flipFaceFlux() const
Map of flipped face flux faces.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
bool found(const Key &) const
Return true if hashedEntry is found in table.
const globalMeshData & globalData() const
Return parallel info.
An STL-conforming hash table.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
errorManip< error > abort(error &err)
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
static void printFieldInfo(const fvMesh &)
Print some field info.
Post-processing mesh subset tool. Given the original mesh and the list of selected cells...
static void waitRequests(const label start=0)
Wait until all requests (from start onwards) have finished.
const labelList & faceMap() const
Old face map.
const labelList & oldPatchStarts() const
Return list of the old patch start labels.
bool isFlux(const DimensionedField< Type, surfaceMesh > &df)
Check if surfaceField is a flux.
const fvPatch & patch() const
Return patch.
static commsTypes defaultCommsType
Default commsType.
static tmp< GeometricField< Type, fvPatchField, volMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &, const fvMesh &sMesh, const labelList &patchMap, const labelList &cellMap, const labelList &faceMap)
Map volume field.
static bool & parRun()
Is this a parallel run?
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
prefixOSstream Pout(cout, "Pout")
Mesh data needed to do the Finite Volume discretisation.
This boundary condition enables processor communication across patches.
A class for managing temporary objects.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
const fvMesh & baseMesh() const
Original mesh.