Go to the documentation of this file.
34 Foam::MRFZoneList::MRFZoneList
62 a = a || this->operator[](i).active();
65 if (warn && this->size() && !a)
67 Info<<
" No MRF zones active" <<
endl;
112 MRFZone& pm = this->operator[](i);
114 allOk = (allOk && ok);
125 this->operator[](i).writeData(os);
140 operator[](i).addCoriolis(
U, ddtU);
149 operator[](i).addCoriolis(
UEqn);
162 operator[](i).addCoriolis(
rho,
UEqn);
178 "MRFZoneList:acceleration",
179 U.mesh().time().timeName(),
190 operator[](i).addCoriolis(
U, acceleration);
193 return tacceleration;
211 operator[](i).makeRelative(
U);
220 operator[](i).makeRelative(
phi);
237 "relative(" + tphi().
name() +
')',
268 operator[](i).makeRelative(rphi.
ref());
295 operator[](i).makeRelative(rphi.
ref(), patchi);
317 operator[](i).makeRelative(
rho,
phi);
326 operator[](i).makeAbsolute(
U);
335 operator[](i).makeAbsolute(
phi);
352 "absolute(" + tphi().
name() +
')',
379 operator[](i).makeAbsolute(
rho,
phi);
388 operator[](i).correctBoundaryVelocity(
U);
401 relative(mesh_.Sf().boundaryField() &
U.boundaryField())
407 forAll(mesh_.boundary(), patchi)
411 isA<fixedValueFvsPatchScalarField>(phibf[patchi])
414 phibf[patchi] ==
Uf[patchi];
422 if (mesh_.topoChanging())
426 operator[](i).update();
440 models.writeData(os);
MRF zone definition based on cell zone and parameters obtained from a control dictionary constructed ...
A keyword and a list of tokens is an 'entry'.
void addAcceleration(const volVectorField &U, volVectorField &ddtU) const
Add the frame acceleration.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
A class for handling words, derived from Foam::string.
A field of fields is a PtrList of fields with reference counting.
void clear() const noexcept
bool writeData(Ostream &os) const
Write data to Ostream.
A class for managing temporary objects.
static constexpr const zero Zero
Global zero (0)
List container for MRF zomes.
void reset(const dictionary &dict)
Reset the source list.
autoPtr< surfaceVectorField > Uf
Ostream & endl(Ostream &os)
Add newline and flush stream.
void correctBoundaryVelocity(volVectorField &U) const
Correct the boundary velocity for the rotation of the MRF region.
bool active(const bool warn=false) const
Return active status.
#define forAll(list, i)
Loop across all elements in list.
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
bool read(const dictionary &dict)
Read MRF dictionary.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
patchWriters resize(patchIds.size())
messageStream Info
Information stream (uses stdout - output is on the master only)
void makeRelative(volVectorField &U) const
Make the given absolute velocity relative within the MRF region.
word name(const complex &c)
Return string representation of complex.
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &phi) const
Return the given relative flux absolute within the MRF region.
tmp< volVectorField > DDt(const volVectorField &U) const
Return the frame acceleration.
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Mesh data needed to do the Finite Volume discretisation.
tmp< surfaceScalarField > relative(const tmp< surfaceScalarField > &phi) const
Return the given absolute flux relative within the MRF region.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
void makeAbsolute(volVectorField &U) const
Make the given relative velocity absolute within the MRF region.
void correctBoundaryFlux(const volVectorField &U, surfaceScalarField &phi) const
Correct the boundary flux for the rotation of the MRF region.
void update()
Update MRFZone faces if the mesh topology changes.
tmp< surfaceScalarField > relative(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given absolute flux in relative form.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
bool read(const dictionary &dict)
Read dictionary.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
~MRFZoneList()
Destructor.
void makeAbsolute(surfaceScalarField &phi, const volVectorField &U)
Make the given flux absolute.
bool good() const
Return true if next operation might succeed.
const word & name() const
Return const access to the MRF region name.
tmp< GeometricField< Type, fvPatchField, volMesh > > DDt(const surfaceScalarField &phi, const GeometricField< Type, fvPatchField, volMesh > &psi)