directionalPressureGradientExplicitSource Class Reference

Creates an explicit pressure gradient source in such a way to deflect the flow towards an specific direction (flowDir). Alternatively add an extra pressure drop in the flowDir direction using a model. More...

Inheritance diagram for directionalPressureGradientExplicitSource:
[legend]
Collaboration diagram for directionalPressureGradientExplicitSource:
[legend]

Public Types

enum  pressureDropModel { pVolumetricFlowRateTable, pConstant, pDarcyForchheimer }
 Modes of pressure drop. More...
 
- Public Types inherited from cellSetOption
enum  selectionModeType { smPoints, smCellSet, smCellZone, smAll }
 Enumeration for selection mode types. More...
 

Public Member Functions

 TypeName ("directionalPressureGradientExplicitSource")
 Runtime type information. More...
 
 directionalPressureGradientExplicitSource (const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
 Construct from explicit source name and mesh. More...
 
virtual void correct (volVectorField &U)
 Correct the pressure gradient. More...
 
virtual void addSup (fvMatrix< vector > &eqn, const label fieldI)
 Add explicit contribution to momentum equation. More...
 
virtual void addSup (const volScalarField &rho, fvMatrix< vector > &eqn, const label fieldI)
 Add explicit contribution to compressible momentum equation. More...
 
virtual void constrain (fvMatrix< vector > &eqn, const label fieldI)
 Set 1/A coefficient. More...
 
virtual void writeData (Ostream &os) const
 Write the source properties. More...
 
virtual bool read (const dictionary &dict)
 Read source dictionary. More...
 
- Public Member Functions inherited from cellSetOption
 TypeName ("cellSetOption")
 Runtime type information. More...
 
 cellSetOption (const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
 Construct from components. More...
 
virtual ~cellSetOption ()=default
 Destructor. More...
 
scalar timeStart () const
 Return const access to the time start. More...
 
scalar duration () const
 Return const access to the duration. More...
 
bool inTimeLimits (const scalar time) const
 Return true if within time limits. More...
 
const selectionModeTypeselectionMode () const
 Return const access to the cell selection mode. More...
 
const wordcellSetName () const
 Return const access to the name of cell set for "cellSet". More...
 
scalar V () const
 Return const access to the total cell volume. More...
 
const labelListcells () const
 Return const access to the cell set. More...
 
scalar & timeStart ()
 Return access to the time start. More...
 
scalar & duration ()
 Return access to the duration. More...
 
virtual bool isActive ()
 Is the source active? More...
 
- Public Member Functions inherited from option
 TypeName ("option")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, option, dictionary,(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh),(name, modelType, dict, mesh))
 
 option (const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
 Construct from components. More...
 
autoPtr< optionclone () const
 Return clone. More...
 
virtual ~option ()=default
 Destructor. More...
 
const wordname () const
 Return const access to the source name. More...
 
const fvMeshmesh () const
 Return const access to the mesh database. More...
 
const dictionarycoeffs () const
 Return dictionary. More...
 
bool active () const
 Return const access to the source active flag. More...
 
void setApplied (const label fieldi)
 Set the applied flag to true for field index fieldi. More...
 
Switchactive ()
 Return access to the source active flag. More...
 
virtual label applyToField (const word &fieldName) const
 Return index of field name if found in fieldNames list. More...
 
virtual void checkApplied () const
 Check that the source has been applied. More...
 
virtual void addSup (fvMatrix< scalar > &eqn, const label fieldi)
 
virtual void addSup (fvMatrix< symmTensor > &eqn, const label fieldi)
 
virtual void addSup (fvMatrix< sphericalTensor > &eqn, const label fieldi)
 
virtual void addSup (fvMatrix< tensor > &eqn, const label fieldi)
 
virtual void addSup (const volScalarField &rho, fvMatrix< scalar > &eqn, const label fieldi)
 
virtual void addSup (const volScalarField &rho, fvMatrix< symmTensor > &eqn, const label fieldi)
 
virtual void addSup (const volScalarField &rho, fvMatrix< sphericalTensor > &eqn, const label fieldi)
 
virtual void addSup (const volScalarField &rho, fvMatrix< tensor > &eqn, const label fieldi)
 
virtual void addSup (const volScalarField &alpha, const volScalarField &rho, fvMatrix< scalar > &eqn, const label fieldi)
 
virtual void addSup (const volScalarField &alpha, const volScalarField &rho, fvMatrix< vector > &eqn, const label fieldi)
 
virtual void addSup (const volScalarField &alpha, const volScalarField &rho, fvMatrix< symmTensor > &eqn, const label fieldi)
 
virtual void addSup (const volScalarField &alpha, const volScalarField &rho, fvMatrix< sphericalTensor > &eqn, const label fieldi)
 
virtual void addSup (const volScalarField &alpha, const volScalarField &rho, fvMatrix< tensor > &eqn, const label fieldi)
 
virtual void constrain (fvMatrix< scalar > &eqn, const label fieldi)
 
virtual void constrain (fvMatrix< sphericalTensor > &eqn, const label fieldi)
 
virtual void constrain (fvMatrix< symmTensor > &eqn, const label fieldi)
 
virtual void constrain (fvMatrix< tensor > &eqn, const label fieldi)
 
virtual void correct (volScalarField &field)
 
virtual void correct (volSphericalTensorField &field)
 
virtual void correct (volSymmTensorField &field)
 
virtual void correct (volTensorField &field)
 
virtual void writeHeader (Ostream &) const
 Write the source header information. More...
 
virtual void writeFooter (Ostream &) const
 Write the source footer information. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from option
static autoPtr< optionNew (const word &name, const dictionary &dict, const fvMesh &mesh)
 Return a reference to the selected fvOption model. More...
 
- Public Attributes inherited from option
bool log
 Switch write log to Info. More...
 
- Static Public Attributes inherited from cellSetOption
static const Enum< selectionModeTypeselectionModeTypeNames_
 List of selection mode type names. More...
 
- Protected Member Functions inherited from cellSetOption
void setSelection (const dictionary &dict)
 Set the cellSet or points selection. More...
 
void setCellSet ()
 Set the cell set based on the user input selection mode. More...
 
void setVol ()
 Recalculate the volume. More...
 
- Protected Attributes inherited from cellSetOption
scalar timeStart_
 Time start. More...
 
scalar duration_
 Duration. More...
 
selectionModeType selectionMode_
 Cell selection mode. More...
 
word cellSetName_
 Name of set/zone for "cellSet" and "cellZone" selectionMode. More...
 
List< pointpoints_
 List of points for "points" selectionMode. More...
 
labelList cells_
 Set of cells to apply source to. More...
 
scalar V_
 Sum of cell volumes. More...
 
- Protected Attributes inherited from option
const word name_
 Source name. More...
 
const word modelType_
 Model type. More...
 
const fvMeshmesh_
 Reference to the mesh database. More...
 
dictionary dict_
 Top level source dictionary. More...
 
dictionary coeffs_
 Dictionary containing source coefficients. More...
 
Switch active_
 Source active flag. More...
 
wordList fieldNames_
 Field names to apply source to - populated by derived models. More...
 
List< boolapplied_
 Applied flag list - corresponds to each fieldNames_ entry. More...
 

Detailed Description

Creates an explicit pressure gradient source in such a way to deflect the flow towards an specific direction (flowDir). Alternatively add an extra pressure drop in the flowDir direction using a model.


Source usage Example usage:

airDeflection
{
    type            directionalPressureGradientExplicitSource;
    active          true;

    directionalPressureGradientExplicitSourceCoeffs
    {
        selectionMode   cellZone;
        cellZone        cZone;

        fields      (U);            // Name of the field
        flowDir     (1 1 0);        // Desired flow direction
        faceZone    f0Zone;         // Face zone upstream cell zone
        relaxationFactor    0.3;    // Relaxation factor for flow
                                    // deflection (default 0.3)

        //Pressure drop model [Pa]
        model       volumetricFlowRateTable;//constant;//DarcyForchheimer;

        //DarcyForchheimer model
        // deltaP = (D*mu + 0.5*rho*magUn)*magUn*length_

        D           5e7;
        I           0;
        length      1e-3;

        //constant model
        pressureDrop    40;

        //volumetricFlowRateTable model
        outOfBounds     clamp;
        fileName        "volFlowRateTable";
    }
}

NOTE: In order to obtain the upwind velocities this function loops over the slaves cells of the faceZone specified in the dictionary, on the other hand, the cellZone to which this source term is applied should be composed of the master cells and they should be 'downwind' the faceZone.

Source files

Definition at line 107 of file directionalPressureGradientExplicitSource.H.

Member Enumeration Documentation

◆ pressureDropModel

Modes of pressure drop.

Enumerator
pVolumetricFlowRateTable 
pConstant 
pDarcyForchheimer 

Definition at line 116 of file directionalPressureGradientExplicitSource.H.

Constructor & Destructor Documentation

◆ directionalPressureGradientExplicitSource()

directionalPressureGradientExplicitSource ( const word sourceName,
const word modelType,
const dictionary dict,
const fvMesh mesh 
)

Construct from explicit source name and mesh.

Definition at line 163 of file directionalPressureGradientExplicitSource.C.

References Foam::endl(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, IOstream::good(), Foam::Info, Foam::name(), Foam::nl, dictionary::null, propsDict(), and Foam::type().

Here is the call graph for this function:

Member Function Documentation

◆ TypeName()

TypeName ( "directionalPressureGradientExplicitSource"  )

Runtime type information.

◆ correct()

◆ addSup() [1/2]

void addSup ( fvMatrix< vector > &  eqn,
const label  fieldI 
)
virtual

Add explicit contribution to momentum equation.

Reimplemented from option.

Definition at line 462 of file directionalPressureGradientExplicitSource.C.

References fvMatrix< Type >::dimensions(), Foam::dimVolume, IOobject::NO_READ, IOobject::NO_WRITE, Su, and Foam::Zero.

Here is the call graph for this function:

◆ addSup() [2/2]

void addSup ( const volScalarField rho,
fvMatrix< vector > &  eqn,
const label  fieldI 
)
virtual

Add explicit contribution to compressible momentum equation.

Reimplemented from option.

Definition at line 488 of file directionalPressureGradientExplicitSource.C.

◆ constrain()

void constrain ( fvMatrix< vector > &  eqn,
const label  fieldI 
)
virtual

Set 1/A coefficient.

Reimplemented from option.

Definition at line 499 of file directionalPressureGradientExplicitSource.C.

References fvMatrix< Type >::A(), IOobject::NO_READ, IOobject::NO_WRITE, and Foam::Zero.

Here is the call graph for this function:

◆ writeData()

void writeData ( Ostream os) const
virtual

Write the source properties.

Reimplemented from option.

Definition at line 533 of file directionalPressureGradientExplicitSource.C.

References NotImplemented.

◆ read()

bool read ( const dictionary dict)
virtual

Read source dictionary.

Reimplemented from cellSetOption.

Definition at line 542 of file directionalPressureGradientExplicitSource.C.

References dict, and dictionary::getOrDefault().

Here is the call graph for this function:

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