Public Member Functions | Private Member Functions | Private Attributes
MRFZone Class Reference

MRF zone definition based on cell zone and parameters obtained from a control dictionary constructed from the given stream. More...

Collaboration diagram for MRFZone:
Collaboration graph
[legend]

Public Member Functions

 ClassName ("MRFZone")
 
 MRFZone (const word &name, const fvMesh &mesh, const dictionary &dict, const word &cellZoneName=word::null)
 Construct from fvMesh. More...
 
autoPtr< MRFZoneclone () const
 Return clone. More...
 
const wordname () const
 Return const access to the MRF region name. More...
 
bool active () const
 Return const access to the MRF active flag. More...
 
vector Omega () const
 Return the current Omega vector. More...
 
void updateMesh (const mapPolyMesh &mpm)
 Update the mesh corresponding to given map. More...
 
void addCoriolis (const volVectorField &U, volVectorField &ddtU) const
 Add the Coriolis force contribution to the acceleration field. More...
 
void addCoriolis (fvVectorMatrix &UEqn, const bool rhs=false) const
 Add the Coriolis force contribution to the momentum equation. More...
 
void addCoriolis (const volScalarField &rho, fvVectorMatrix &UEqn, const bool rhs=false) const
 Add the Coriolis force contribution to the momentum equation. More...
 
void makeRelative (volVectorField &U) const
 Make the given absolute velocity relative within the MRF region. More...
 
void makeRelative (surfaceScalarField &phi) const
 Make the given absolute flux relative within the MRF region. More...
 
void makeRelative (FieldField< fvsPatchField, scalar > &phi) const
 Make the given absolute boundary flux relative. More...
 
void makeRelative (const surfaceScalarField &rho, surfaceScalarField &phi) const
 Make the given absolute mass-flux relative within the MRF region. More...
 
void makeAbsolute (volVectorField &U) const
 Make the given relative velocity absolute within the MRF region. More...
 
void makeAbsolute (surfaceScalarField &phi) const
 Make the given relative flux absolute within the MRF region. More...
 
void makeAbsolute (const surfaceScalarField &rho, surfaceScalarField &phi) const
 Make the given relative mass-flux absolute within the MRF region. More...
 
void correctBoundaryVelocity (volVectorField &U) const
 Correct the boundary velocity for the rotation of the MRF region. More...
 
void writeData (Ostream &os) const
 Write. More...
 
bool read (const dictionary &dict)
 Read MRF dictionary. More...
 

Private Member Functions

void setMRFFaces ()
 Divide faces in frame according to patch. More...
 
template<class RhoFieldType >
void makeRelativeRhoFlux (const RhoFieldType &rho, surfaceScalarField &phi) const
 Make the given absolute mass/vol flux relative within the MRF region. More...
 
template<class RhoFieldType >
void makeRelativeRhoFlux (const RhoFieldType &rho, FieldField< fvsPatchField, scalar > &phi) const
 Make the given absolute mass/vol flux relative within the MRF region. More...
 
template<class RhoFieldType >
void makeAbsoluteRhoFlux (const RhoFieldType &rho, surfaceScalarField &phi) const
 Make the given relative mass/vol flux absolute within the MRF region. More...
 
 MRFZone (const MRFZone &)
 Disallow default bitwise copy construct. More...
 
void operator= (const MRFZone &)
 Disallow default bitwise assignment. More...
 

Private Attributes

const fvMeshmesh_
 Reference to the mesh database. More...
 
const word name_
 Name of the MRF region. More...
 
dictionary coeffs_
 Coefficients dictionary. More...
 
bool active_
 MRF region active flag. More...
 
word cellZoneName_
 Name of cell zone. More...
 
label cellZoneID_
 Cell zone ID. More...
 
const wordReList excludedPatchNames_
 
labelList excludedPatchLabels_
 
labelList internalFaces_
 Internal faces that are part of MRF. More...
 
labelListList includedFaces_
 Outside faces (per patch) that move with the MRF. More...
 
labelListList excludedFaces_
 Excluded faces (per patch) that do not move with the MRF. More...
 
const vector origin_
 Origin of the axis. More...
 
vector axis_
 Axis vector. More...
 
autoPtr< DataEntry< scalar > > omega_
 Angular velocty (rad/sec) More...
 

Detailed Description

MRF zone definition based on cell zone and parameters obtained from a control dictionary constructed from the given stream.

The rotation of the MRF region is defined by an origin and axis of rotation and an angular speed.

Source files

Definition at line 65 of file MRFZone.H.

Constructor & Destructor Documentation

◆ MRFZone() [1/2]

MRFZone ( const MRFZone )
private

Disallow default bitwise copy construct.

◆ MRFZone() [2/2]

MRFZone ( const word name,
const fvMesh mesh,
const dictionary dict,
const word cellZoneName = word::null 
)

Construct from fvMesh.

Definition at line 237 of file MRFZone.C.

References Foam::exit(), Foam::FatalError, FatalErrorInFunction, forAllConstIter(), Foam::mag(), word::null, Foam::reduce(), and HashTable::size().

Here is the call graph for this function:

Member Function Documentation

◆ setMRFFaces()

void setMRFFaces ( )
private

◆ makeRelativeRhoFlux() [1/2]

void makeRelativeRhoFlux ( const RhoFieldType &  rho,
surfaceScalarField phi 
) const
private

Make the given absolute mass/vol flux relative within the MRF region.

Definition at line 36 of file MRFZoneTemplates.C.

References forAll, GeometricField::internalField(), phi, and rho.

Here is the call graph for this function:

◆ makeRelativeRhoFlux() [2/2]

void makeRelativeRhoFlux ( const RhoFieldType &  rho,
FieldField< fvsPatchField, scalar > &  phi 
) const
private

Make the given absolute mass/vol flux relative within the MRF region.

Definition at line 63 of file MRFZoneTemplates.C.

References GeometricField::boundaryField(), forAll, patchi, phi, and rho.

Here is the call graph for this function:

◆ makeAbsoluteRhoFlux()

void makeAbsoluteRhoFlux ( const RhoFieldType &  rho,
surfaceScalarField phi 
) const
private

Make the given relative mass/vol flux absolute within the MRF region.

Definition at line 102 of file MRFZoneTemplates.C.

References GeometricField::boundaryField(), forAll, GeometricField::internalField(), patchi, phi, and rho.

Here is the call graph for this function:

◆ operator=()

void operator= ( const MRFZone )
private

Disallow default bitwise assignment.

◆ ClassName()

ClassName ( "MRFZone"  )

◆ clone()

autoPtr<MRFZone> clone ( ) const
inline

Return clone.

Definition at line 164 of file MRFZone.H.

References NotImplemented.

◆ name()

const Foam::word & name ( ) const
inline

Return const access to the MRF region name.

Definition at line 26 of file MRFZoneI.H.

References MRFZone::name_.

Referenced by MRFZoneList::read().

Here is the caller graph for this function:

◆ active()

bool active ( ) const
inline

Return const access to the MRF active flag.

Definition at line 32 of file MRFZoneI.H.

◆ Omega()

Foam::vector Omega ( ) const

Return the current Omega vector.

Definition at line 304 of file MRFZone.C.

◆ updateMesh()

void updateMesh ( const mapPolyMesh mpm)
inline

Update the mesh corresponding to given map.

Definition at line 188 of file MRFZone.H.

References MRFZone::setMRFFaces().

Here is the call graph for this function:

◆ addCoriolis() [1/3]

void addCoriolis ( const volVectorField U,
volVectorField ddtU 
) const

Add the Coriolis force contribution to the acceleration field.

Definition at line 311 of file MRFZone.C.

References cells, forAll, GeometricField::internalField(), and U.

Here is the call graph for this function:

◆ addCoriolis() [2/3]

void addCoriolis ( fvVectorMatrix UEqn,
const bool  rhs = false 
) const

Add the Coriolis force contribution to the momentum equation.

Adds to the lhs of the equation; optionally add to rhs

Definition at line 335 of file MRFZone.C.

References cells, forAll, U, and UEqn().

Here is the call graph for this function:

◆ addCoriolis() [3/3]

void addCoriolis ( const volScalarField rho,
fvVectorMatrix UEqn,
const bool  rhs = false 
) const

Add the Coriolis force contribution to the momentum equation.

Adds to the lhs of the equation; optionally add to rhs

Definition at line 369 of file MRFZone.C.

References cells, forAll, rho, U, and UEqn().

Here is the call graph for this function:

◆ makeRelative() [1/4]

void makeRelative ( volVectorField U) const

Make the given absolute velocity relative within the MRF region.

Definition at line 406 of file MRFZone.C.

References C::C(), cells, forAll, patchi, U, and Vector< scalar >::zero.

Here is the call graph for this function:

◆ makeRelative() [2/4]

void makeRelative ( surfaceScalarField phi) const

Make the given absolute flux relative within the MRF region.

◆ makeRelative() [3/4]

void makeRelative ( FieldField< fvsPatchField, scalar > &  phi) const

Make the given absolute boundary flux relative.

within the MRF region

Definition at line 450 of file MRFZone.C.

References phi.

◆ makeRelative() [4/4]

void makeRelative ( const surfaceScalarField rho,
surfaceScalarField phi 
) const

Make the given absolute mass-flux relative within the MRF region.

Definition at line 457 of file MRFZone.C.

References phi, and rho.

◆ makeAbsolute() [1/3]

void makeAbsolute ( volVectorField U) const

Make the given relative velocity absolute within the MRF region.

Definition at line 466 of file MRFZone.C.

References C::C(), cells, forAll, patchi, and U.

Here is the call graph for this function:

◆ makeAbsolute() [2/3]

void makeAbsolute ( surfaceScalarField phi) const

Make the given relative flux absolute within the MRF region.

◆ makeAbsolute() [3/3]

void makeAbsolute ( const surfaceScalarField rho,
surfaceScalarField phi 
) const

Make the given relative mass-flux absolute within the MRF region.

Definition at line 511 of file MRFZone.C.

References phi, and rho.

◆ correctBoundaryVelocity()

void correctBoundaryVelocity ( volVectorField U) const

Correct the boundary velocity for the rotation of the MRF region.

Definition at line 520 of file MRFZone.C.

References forAll, patchi, and U.

◆ writeData()

void writeData ( Ostream os) const

Write.

Definition at line 543 of file MRFZone.C.

References token::BEGIN_BLOCK, Foam::decrIndent(), token::END_BLOCK, token::END_STATEMENT, Foam::incrIndent(), Foam::nl, Ostream::write(), and Ostream::writeKeyword().

Here is the call graph for this function:

◆ read()

bool read ( const dictionary dict)

Read MRF dictionary.

Definition at line 564 of file MRFZone.C.

References dict, and dictionary::lookupOrDefault().

Referenced by MRFZoneList::read().

Here is the call graph for this function:
Here is the caller graph for this function:

Field Documentation

◆ mesh_

const fvMesh& mesh_
private

Reference to the mesh database.

Definition at line 70 of file MRFZone.H.

Referenced by MRFZone::setMRFFaces().

◆ name_

const word name_
private

Name of the MRF region.

Definition at line 73 of file MRFZone.H.

Referenced by MRFZone::name().

◆ coeffs_

dictionary coeffs_
private

Coefficients dictionary.

Definition at line 76 of file MRFZone.H.

◆ active_

bool active_
private

MRF region active flag.

Definition at line 79 of file MRFZone.H.

◆ cellZoneName_

word cellZoneName_
private

Name of cell zone.

Definition at line 82 of file MRFZone.H.

◆ cellZoneID_

label cellZoneID_
private

Cell zone ID.

Definition at line 85 of file MRFZone.H.

Referenced by MRFZone::setMRFFaces().

◆ excludedPatchNames_

const wordReList excludedPatchNames_
private

Definition at line 87 of file MRFZone.H.

◆ excludedPatchLabels_

labelList excludedPatchLabels_
private

Definition at line 89 of file MRFZone.H.

Referenced by MRFZone::setMRFFaces().

◆ internalFaces_

labelList internalFaces_
private

Internal faces that are part of MRF.

Definition at line 92 of file MRFZone.H.

Referenced by MRFZone::setMRFFaces().

◆ includedFaces_

labelListList includedFaces_
private

Outside faces (per patch) that move with the MRF.

Definition at line 95 of file MRFZone.H.

Referenced by MRFZone::setMRFFaces().

◆ excludedFaces_

labelListList excludedFaces_
private

Excluded faces (per patch) that do not move with the MRF.

Definition at line 98 of file MRFZone.H.

Referenced by MRFZone::setMRFFaces().

◆ origin_

const vector origin_
private

Origin of the axis.

Definition at line 101 of file MRFZone.H.

◆ axis_

vector axis_
private

Axis vector.

Definition at line 104 of file MRFZone.H.

◆ omega_

autoPtr<DataEntry<scalar> > omega_
private

Angular velocty (rad/sec)

Definition at line 107 of file MRFZone.H.


The documentation for this class was generated from the following files: