Loading...
Searching...
No Matches
MRFZone Class Reference

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

#include <MRFZone.H>

Public Member Functions

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

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 64 of file MRFZone.H.

Constructor & Destructor Documentation

◆ MRFZone()

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

Construct from fvMesh.

Definition at line 230 of file MRFZone.C.

References dict, mesh, name(), read(), and Foam::Zero.

Here is the call graph for this function:

Member Function Documentation

◆ ClassName()

ClassName ( "MRFZone" )

References dict, mesh, name(), and word::null.

Here is the call graph for this function:

◆ clone()

autoPtr< MRFZone > clone ( ) const
inline

Return clone.

Definition at line 220 of file MRFZone.H.

References NotImplemented.

◆ name()

const Foam::word & name ( ) const
inline

Return const access to the MRF region name.

Definition at line 22 of file MRFZoneI.H.

Referenced by ClassName(), and MRFZone().

Here is the caller graph for this function:

◆ active()

bool active ( ) const
inline

Return const access to the MRF active flag.

Definition at line 28 of file MRFZoneI.H.

◆ origin()

const Foam::vector & origin ( ) const
inline

Return const access to the MRF origin.

Definition at line 34 of file MRFZoneI.H.

◆ axis()

const Foam::vector & axis ( ) const
inline

Return const access to the MRF axis.

Definition at line 40 of file MRFZoneI.H.

◆ Omega()

Foam::vector Omega ( ) const

Return the current Omega vector.

Definition at line 255 of file MRFZone.C.

Referenced by addCoriolis(), addCoriolis(), addCoriolis(), correctBoundaryVelocity(), makeAbsolute(), and makeRelative().

Here is the caller graph for this function:

◆ updateMesh()

void updateMesh ( const mapPolyMesh & mpm)
inline

Update the mesh corresponding to given map.

Definition at line 257 of file MRFZone.H.

◆ addCoriolis() [1/3]

void addCoriolis ( const volVectorField & U,
volVectorField & ddtU ) const

Add the Coriolis force contribution to the acceleration field.

Definition at line 261 of file MRFZone.C.

References cells, forAll, Omega(), GeometricField< Type, PatchField, GeoMesh >::primitiveFieldRef(), 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 286 of file MRFZone.C.

References cells, forAll, Omega(), Foam::rhs(), 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 319 of file MRFZone.C.

References cells, forAll, Omega(), rho, Foam::rhs(), U, and UEqn.

Here is the call graph for this function:

◆ makeRelative() [1/5]

void makeRelative ( volVectorField & U) const

Make the given absolute velocity relative within the MRF region.

Definition at line 357 of file MRFZone.C.

References C::C(), cells, forAll, Omega(), U, and Foam::Zero.

Here is the call graph for this function:

◆ makeRelative() [2/5]

void makeRelative ( surfaceScalarField & phi) const

Make the given absolute flux relative within the MRF region.

Definition at line 403 of file MRFZone.C.

References phi.

◆ makeRelative() [3/5]

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

Make the given absolute boundary flux relative.

within the MRF region

Definition at line 409 of file MRFZone.C.

References phi.

◆ makeRelative() [4/5]

void makeRelative ( Field< scalar > & phi,
const label patchi ) const

Make the given absolute patch flux relative.

within the MRF region

Definition at line 415 of file MRFZone.C.

References phi.

◆ makeRelative() [5/5]

void makeRelative ( const surfaceScalarField & rho,
surfaceScalarField & phi ) const

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

Definition at line 421 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 431 of file MRFZone.C.

References C::C(), cells, forAll, Omega(), 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.

Definition at line 476 of file MRFZone.C.

References phi.

◆ 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 482 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 492 of file MRFZone.C.

References forAll, Omega(), and U.

Here is the call graph for this function:

◆ zero()

template<class Type>
void zero ( GeometricField< Type, fvsPatchField, surfaceMesh > & phi) const

Zero the MRF region of the given field.

Definition at line 206 of file MRFZoneTemplates.C.

References forAll, phi, and Foam::Zero.

Referenced by MRFZoneList::zeroFilter().

Here is the caller graph for this function:

◆ update()

void update ( )

Update MRFZone faces if the mesh topology changes.

Definition at line 591 of file MRFZone.C.

◆ writeData()

void writeData ( Ostream & os) const

Write.

Definition at line 522 of file MRFZone.C.

References Foam::nl, and os().

Here is the call graph for this function:

◆ read()

bool read ( const dictionary & dict)

Read MRF dictionary.

Definition at line 542 of file MRFZone.C.

References dict, Foam::exit(), Foam::FatalError, FatalErrorInFunction, and Foam::returnReduceOr().

Referenced by MRFZone().

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

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