Go to the documentation of this file.
30 template<
class GeoField>
35 mesh.objectRegistry::lookupClass<GeoField>()
40 const GeoField&
fld = *iter();
42 Pout<<
"Field:" << iter.key() <<
" internalsize:" <<
fld.size()
49 <<
' ' <<
fld.boundaryField()[patchI].patch().name()
50 <<
' ' <<
fld.boundaryField()[patchI].type()
51 <<
' ' <<
fld.boundaryField()[patchI].size()
59 template<
class T,
class Mesh>
69 static_cast<const fvMesh&
>(mesh_).objectRegistry::lookupClass<fldType>()
72 bflds.setSize(flds.
size());
78 const fldType&
fld = *iter();
80 bflds.set(i,
fld.boundaryField().clone().ptr());
88 template<
class T,
class Mesh>
102 mesh_.objectRegistry::lookupClass<fldType>()
105 if (flds.
size() != oldBflds.size())
115 fldType&
fld = *iter();
116 typename fldType::GeometricBoundaryField& bfld =
fld.boundaryField();
132 forAll(oldPatchStarts, oldPatchI)
134 label oldLocalI = oldFaceI - oldPatchStarts[oldPatchI];
136 if (oldLocalI >= 0 && oldLocalI < oldBfld[oldPatchI].size())
138 patchFld[i] = oldBfld[oldPatchI][oldLocalI];
157 static_cast<const fvMesh&
>(mesh_).objectRegistry::lookupClass<fldType>()
160 iflds.setSize(flds.
size());
166 const fldType&
fld = *iter();
168 iflds.set(i,
fld.internalField().clone());
189 mesh_.objectRegistry::lookupClass<fldType>()
192 if (flds.
size() != oldFlds.size())
204 fldType&
fld = *iter();
205 typename fldType::GeometricBoundaryField& bfld =
fld.boundaryField();
207 const Field<T>& oldInternal = oldFlds[fieldI++];
221 if (oldFaceI < oldInternal.size())
223 patchFld[i] = oldInternal[oldFaceI];
227 patchFld[i] =
flipOp()(patchFld[i]);
237 template<
class GeoField,
class PatchFieldType>
240 const typename GeoField::value_type& initVal
245 mesh_.objectRegistry::lookupClass<GeoField>()
250 GeoField&
fld = *iter();
252 typename GeoField::GeometricBoundaryField& bfld =
257 if (isA<PatchFieldType>(bfld[patchI]))
259 bfld[patchI] == initVal;
267 template<
class GeoField>
272 mesh_.objectRegistry::lookupClass<GeoField>()
277 const GeoField&
fld = *iter();
278 fld.correctBoundaryConditions();
301 template<
class GeoField>
316 <<
" for domain:" << domain <<
endl;
321 const GeoField&
fld =
336 template<
class GeoField>
349 <<
" from domain:" << domain <<
endl;
359 <<
" from domain:" << domain <<
endl;
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.
const labelHashSet & flipFaceFlux() const
Map of flipped face flux faces.
static void receiveFields(const label domain, const wordList &fieldNames, fvMesh &, PtrList< GeoField > &, const dictionary &fieldDicts)
Receive fields. Opposite of sendFields.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
void saveBoundaryFields(PtrList< FieldField< fvsPatchField, T > > &bflds) const
Save boundary fields.
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
void mapBoundaryFields(const mapPolyMesh &map, const PtrList< FieldField< fvsPatchField, T > > &oldBflds)
Map boundary fields.
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
Post-processing mesh subset tool. Given the original mesh and the list of selected cells,...
Info<< "Creating field dpdt\n"<< endl;volScalarField dpdt(IOobject("dpdt", runTime.timeName(), mesh), mesh, dimensionedScalar("dpdt", p.dimensions()/dimTime, 0));Info<< "Creating field kinetic energy K\n"<< endl;volScalarField K("K", 0.5 *magSqr(U));volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);p_rgh=p - rho *gh;mesh.setFluxRequired(p_rgh.name());multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
static void printFieldInfo(const fvMesh &)
Print some field info.
Class containing functor to negate primitives. Dummy for all other types.
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Ostream & endl(Ostream &os)
Add newline and flush stream.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
static void sendFields(const label domain, const wordList &fieldNames, const fvMeshSubset &, Ostream &toNbr)
Send subset of fields.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
A list of keyword definitions, which are a keyword followed by any number of values (e....
label size() const
Return number of elements in table.
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.
Mesh data needed to do the Finite Volume discretisation.
errorManip< error > abort(error &err)
void saveInternalFields(PtrList< Field< T > > &iflds) const
Save internal fields of surfaceFields.
An STL-conforming hash table.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
const fvPatch & patch() const
Return patch.
prefixOSstream Pout(cout, "Pout")
void correctBoundaryConditions()
Call correctBoundaryConditions on fields.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
wordList fieldNames(const IOobjectList &objects, const bool syncPar)
Get sorted names of fields of type. If syncPar and running in parallel.
const labelList & oldPatchStarts() const
Return list of the old patch start labels.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
void initPatchFields(const typename GeoField::value_type &initVal)
Init patch fields of certain type.
label start() const
Return start label of this patch in the polyMesh face list.
const labelList & faceMap() const
Old face map.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const Time & time() const
Return the top-level database.
void mapExposedFaces(const mapPolyMesh &map, const PtrList< Field< T > > &oldFlds)
Set value of patch faces resulting from internal faces.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
void size(const label)
Override size to be inconsistent with allocated storage.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
const fvMesh & baseMesh() const
Original mesh.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Generic GeometricField class.