MRFZone Class Reference

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

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...
 
const vectororigin () const
 Return const access to the MRF origin. More...
 
const vectoraxis () const
 Return const access to the MRF axis. 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 (Field< scalar > &phi, const label patchi) const
 Make the given absolute patch 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...
 
template<class Type >
void zero (GeometricField< Type, fvsPatchField, surfaceMesh > &phi) const
 Zero the MRF region of the given field. More...
 
void update ()
 Update MRFZone faces if the mesh topology changes. More...
 
void writeData (Ostream &os) const
 Write. More...
 
bool read (const dictionary &dict)
 Read MRF dictionary. 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 68 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 240 of file MRFZone.C.

References dict, and Foam::read().

Here is the call graph for this function:

Member Function Documentation

◆ ClassName()

ClassName ( "MRFZone"  )

◆ clone()

autoPtr<MRFZone> clone ( ) const
inline

Return clone.

Definition at line 178 of file MRFZone.H.

References NotImplemented.

◆ name()

const Foam::word & name ( ) const
inline

Return const access to the MRF region name.

Definition at line 29 of file MRFZoneI.H.

◆ active()

bool active ( ) const
inline

Return const access to the MRF active flag.

Definition at line 35 of file MRFZoneI.H.

◆ origin()

const Foam::vector & origin ( ) const
inline

Return const access to the MRF origin.

Definition at line 41 of file MRFZoneI.H.

◆ axis()

const Foam::vector & axis ( ) const
inline

Return const access to the MRF axis.

Definition at line 47 of file MRFZoneI.H.

◆ Omega()

Foam::vector Omega ( ) const

Return the current Omega vector.

Definition at line 264 of file MRFZone.C.

◆ updateMesh()

void updateMesh ( const mapPolyMesh mpm)
inline

Update the mesh corresponding to given map.

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

References cells, forAll, 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 295 of file MRFZone.C.

References cells, forAll, fvMatrix< Type >::psi(), fvMatrix< Type >::source(), 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 329 of file MRFZone.C.

References cells, forAll, fvMatrix< Type >::psi(), rho, fvMatrix< Type >::source(), 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 366 of file MRFZone.C.

References C::C(), cells, forAll, 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 412 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 418 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 424 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 431 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 440 of file MRFZone.C.

References C::C(), cells, forAll, 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 485 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 492 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 501 of file MRFZone.C.

References forAll, and U.

◆ zero()

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

Zero the MRF region of the given field.

Definition at line 214 of file MRFZoneTemplates.C.

References GeometricField< Type, PatchField, GeoMesh >::boundaryFieldRef(), forAll, phi, GeometricField< Type, PatchField, GeoMesh >::primitiveFieldRef(), and Foam::Zero.

Here is the call graph for this function:

◆ update()

void update ( )

Update MRFZone faces if the mesh topology changes.

Definition at line 614 of file MRFZone.C.

◆ writeData()

void writeData ( Ostream os) const

Write.

Definition at line 531 of file MRFZone.C.

References Ostream::beginBlock(), Ostream::endBlock(), Foam::nl, os(), and Ostream::writeEntry().

Here is the call graph for this function:

◆ read()

bool read ( const dictionary dict)

Read MRF dictionary.

Definition at line 551 of file MRFZone.C.

References dict, Foam::exit(), Foam::FatalError, FatalErrorInFunction, word::null, dictionary::readIfPresent(), and Foam::reduce().

Here is the call graph for this function:

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