Go to the documentation of this file.
41 namespace functionObjects
57 const auto&
U = lookupObject<volVectorField>(
UName_);
61 auto*
turb = findObject<cmpTurbModel>
87 "div(((rho*nuEff)*dev2(T(grad(U)))))"
100 const auto*
turb = findObject<icoTurbModel>
123 "div((nuEff*dev2(T(grad(U)))))"
150 const auto&
phi =lookupObject<surfaceScalarField>(phiName_);
160 word momName(scopedName(
"momentError"));
164 const fvMesh& subMesh = zoneSubSetPtr_->subsetter().subMesh();
172 scopedName(
"momentErrorMap"),
182 subMesh.objectRegistry::store(momentPtr);
184 momName = scopedName(
"momentErrorZone");
201 mesh_.objectRegistry::store(momentPtr);
214 if (
dict.readIfPresent<
word>(
"p", pName_))
218 if (
dict.readIfPresent<
word>(
"U", UName_))
223 if (
dict.readIfPresent<
word>(
"phi", phiName_))
225 Info<<
" phi: " << phiName_ <<
endl;
228 if (
dict.found(
"cellZones"))
234 zoneSubSetPtr_.reset(
nullptr);
246 const auto&
p = lookupObject<volScalarField>(pName_);
247 const auto&
U = lookupObject<volVectorField>(UName_);
248 const auto&
phi = lookupObject<surfaceScalarField>(phiName_);
252 const fvMeshSubset& subsetter = zoneSubSetPtr_->subsetter();
254 fvMesh& subMesh = zoneSubSetPtr_->subsetter().subMesh();
256 subMesh.fvSchemes::readOpt() = mesh_.fvSchemes::readOpt();
257 subMesh.fvSchemes::read();
262 scopedName(
"momentErrorMap")
279 lookupObjectRef<volVectorField>(scopedName(
"momentError"));
296 Log <<
" functionObjects::" <<
type() <<
" " <<
name();
300 Log <<
" writing field: " << scopedName(
"momentError") <<
endl;
302 const auto& momentErr =
303 lookupObjectRef<volVectorField>(scopedName(
"momentError"));
309 Log <<
" writing field: " << scopedName(
"momentErrorMap") <<
endl;
311 const fvMeshSubset& subsetter = zoneSubSetPtr_->subsetter();
314 const auto& momentErrMap =
317 scopedName(
"momentErrorMap")
321 zoneSubSetPtr_->mapToZone<
vector>(momentErrMap);
Defines the attributes of an object for which implicit objectRegistry management is supported,...
tmp< volVectorField > divDevRhoReff()
Return the effective viscous stress (laminar + turbulent).
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
A class for handling words, derived from Foam::string.
A class for managing temporary objects.
Given the original mesh and the list of selected cells, it creates the mesh consisting only of the de...
const dimensionSet dimVelocity
bool read(const char *buf, int32_t &val)
Same as readInt32.
static word timeName(const scalar t, const int precision=precision_)
static const word propertiesName
Default name of the turbulence properties dictionary.
Calculate the divergence of the given field.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Ostream & endl(Ostream &os)
Add newline and flush stream.
void calcMomentError()
Calculate the momentum error.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
virtual bool execute()
Execute.
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< volScalarField > trho
const fvMesh & subMesh() const
Return reference to subset mesh.
autoPtr< Detail::zoneSubSet > zoneSubSetPtr_
Sub-set mesh.
messageStream Info
Information stream (stdout output on master, null elsewhere)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
virtual bool read(const dictionary &)
Read the forces data.
const Type & lookupObject(const word &name, const bool recursive=false) const
word UName_
Name of velocity field.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
virtual bool read(const dictionary &dict)
Read optional controls.
Macros for easy insertion into run-time selection tables.
Mesh data needed to do the Finite Volume discretisation.
Type & lookupObjectRef(const word &name, const bool recursive=false) const
GeometricField< vector, fvPatchField, volMesh > volVectorField
static tmp< GeometricField< Type, fvPatchField, volMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &, const fvMesh &sMesh, const labelUList &patchMap, const labelUList &cellMap, const labelUList &faceMap, const bool allowUnmapped=false)
Map volume field. Optionally allow unmapped faces not to produce.
Calculate the laplacian of the given field.
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
Calculate the gradient of the given field.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
defineTypeNameAndDebug(ObukhovLength, 0)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
const Time & time() const
Return the top-level database.
Templated abstract base class for single-phase incompressible turbulence models.
static const GeometricField< vector, fvPatchField, volMesh > & null()
Return a null geometric field.
const dimensionSet dimVolume(pow3(dimLength))
compressible::turbulenceModel & turb
virtual bool write()
Write.
momentumError(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.