Go to the documentation of this file.
93 #ifndef functionObjects_nearWallFields_H
94 #define functionObjects_nearWallFields_H
105 namespace functionObjects
114 public fvMeshFunctionObject
152 PtrList<volScalarField>
vsf_;
153 PtrList<volVectorField>
vvf_;
155 PtrList<volSymmTensorField>
vSymmtf_;
156 PtrList<volTensorField>
vtf_;
222 virtual bool write();
TypeName("nearWallFields")
List< List< point > > cellToSamples_
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
List< Tuple2< word, word > > fieldSet_
A class for handling words, derived from Foam::string.
PtrList< volSphericalTensorField > vSpheretf_
autoPtr< mapDistribute > getPatchDataMapPtr_
Samples near-patch volume fields within an input distance range.
PtrList< volScalarField > vsf_
A HashTable with keys but without contents that is similar to std::unordered_set.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Given cell centre values and point (vertex) values decompose into tetrahedra and linear interpolate w...
void sampleBoundaryField(const interpolationCellPoint< Type > &interpolator, GeometricField< Type, fvPatchField, volMesh > &fld) const
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
PtrList< volSymmTensorField > vSymmtf_
void sampleFields(PtrList< GeometricField< Type, fvPatchField, volMesh >> &) const
A HashTable similar to std::unordered_map.
virtual bool read(const dictionary &dict)
void createFields(PtrList< GeometricField< Type, fvPatchField, volMesh >> &) const
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
List< labelList > labelListList
A List of labelList.
virtual ~nearWallFields()=default
const word & name() const noexcept
HashTable< word > fieldMap_
labelListList cellToWalls_
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
PtrList< volTensorField > vtf_
nearWallFields(const word &name, const Time &runTime, const dictionary &dict)
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys and label hasher.
Generic GeometricField class.
void operator=(const nearWallFields &)=delete
PtrList< volVectorField > vvf_
HashTable< word > reverseFieldMap_