Public Types | Public Member Functions | List of all members
directionalPressureGradientExplicitSource Class Reference

Applies 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:
Inheritance graph
[legend]
Collaboration diagram for directionalPressureGradientExplicitSource:
Collaboration graph
[legend]

Public Types

enum  pressureDropModel { pVolumetricFlowRateTable, pConstant, pDarcyForchheimer }
 
- Public Types inherited from cellSetOption
enum  selectionModeType { smAll, smCellSet, smCellZone, smPoints }
 

Public Member Functions

 TypeName ("directionalPressureGradientExplicitSource")
 
 directionalPressureGradientExplicitSource (const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
 
 directionalPressureGradientExplicitSource (const directionalPressureGradientExplicitSource &)=delete
 
void operator= (const directionalPressureGradientExplicitSource &)=delete
 
virtual ~directionalPressureGradientExplicitSource ()=default
 
virtual void correct (volVectorField &U)
 
virtual void addSup (fvMatrix< vector > &eqn, const label fieldI)
 
virtual void addSup (const volScalarField &rho, fvMatrix< vector > &eqn, const label fieldI)
 
virtual void constrain (fvMatrix< vector > &eqn, const label fieldI)
 
virtual void writeData (Ostream &os) const
 
virtual bool read (const dictionary &dict)
 
- Public Member Functions inherited from cellSetOption
 TypeName ("cellSetOption")
 
 cellSetOption (const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
 
virtual ~cellSetOption ()=default
 
scalar timeStart () const noexcept
 
scalar duration () const noexcept
 
bool inTimeLimits (const scalar timeValue) const
 
selectionModeType selectionMode () const noexcept
 
bool useSubMesh () const noexcept
 
const wordcellSetName () const noexcept
 
scalar V () const noexcept
 
const labelListcells () const noexcept
 
scalar timeStart (scalar val) noexcept
 
scalar duration (scalar val) noexcept
 
virtual bool isActive ()
 
- Public Member Functions inherited from option
 TypeName ("option")
 
 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)
 
autoPtr< optionclone () const
 
virtual ~option ()=default
 
const wordname () const noexcept
 
const fvMeshmesh () const noexcept
 
const dictionarycoeffs () const noexcept
 
bool active () const noexcept
 
void setApplied (const label fieldi)
 
bool active (const bool on) noexcept
 
virtual label applyToField (const word &fieldName) const
 
virtual void checkApplied () const
 
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 postProcessSens (scalarField &sensField, const word &fieldName=word::null, const word &designVariablesName=word::null)
 
virtual void postProcessSens (vectorField &sensField, const word &fieldName=word::null, const word &designVariablesName=word::null)
 
virtual void postProcessSens (tensorField &sensField, const word &fieldName=word::null, const word &designVariablesName=word::null)
 
virtual void writeHeader (Ostream &) const
 
virtual void writeFooter (Ostream &) const
 

Additional Inherited Members

- Static Public Member Functions inherited from option
static autoPtr< optionNew (const word &name, const dictionary &dict, const fvMesh &mesh)
 
- Public Attributes inherited from option
bool log
 
- Static Public Attributes inherited from cellSetOption
static const Enum< selectionModeTypeselectionModeTypeNames_
 
- Protected Member Functions inherited from cellSetOption
void setSelection (const dictionary &dict)
 
void setCellSelection ()
 
void setVol ()
 
- Protected Member Functions inherited from option
void resetApplied ()
 
- Protected Attributes inherited from cellSetOption
scalar timeStart_
 
scalar duration_
 
selectionModeType selectionMode_
 
word cellSetName_
 
List< pointpoints_
 
labelList cells_
 
scalar V_
 
- Protected Attributes inherited from option
const word name_
 
const word modelType_
 
const fvMeshmesh_
 
dictionary dict_
 
dictionary coeffs_
 
wordList fieldNames_
 
List< boolapplied_
 
bool active_
 

Detailed Description

Applies 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.

Usage
Minimal example by using constant/fvOptions:
directionalPressureGradientExplicitSource1
{
    // Mandatory entries (unmodifiable)
    type            directionalPressureGradientExplicitSource;

    // Mandatory entries (unmodifiable)
    model       <modelName>;
    fields      (<fieldName>);

    // Mandatory entries (runtime modifiable)
    flowDir     (1 1 0);
    faceZone    <faceZoneName>;

    // Conditional mandatory entries (unmodifiable)

        // when <timePath>/uniform/<name>Properties file exists
        gradient    <vectorField>; // reading from the aforementioned file

        // when model=DarcyForchheimer
        // deltaP = (D*mu + 0.5*rho*magUn)*magUn*length
        D           5e7;
        I           0;
        length      1e-3;

        // when model=constant
        pressureDrop    40;

        // when model=volumetricFlowRateTable
        outOfBounds     clamp;
        fileName        "volFlowRateTable";

    // Optional entries (runtime modifiable)
    relaxationFactor    0.3;

    // Mandatory/Optional (inherited) entries
    ...
}

where the entries mean:

Property Description Type Reqd Dflt
type Type name: directionalPressureGradientExplicitSource word yes -
model Pressure drop model [Pa] word yes -
fields Name of operand field word yes -
gradient Initial pressure gradient field vectorField cndtnl -
flowDir Deflection flow direction vector yes -
faceZone Name of upstream faceZone word yes -
relaxationFactor Relaxation factor for flow deflection scalar no 0.3
D Darcy pressure loss coefficient scalar cndtnl -
I Inertia pressure lost coefficient scalar cndtnl -
length Porous media length scalar cndtnl -
presssureDrop Constant pressure drop scalar cndtnl -
fileName Interpolation table for volumetric flow rate interpolationTable cndtnl -

Options for the model entry:

      volumetricFlowRateTable | Pressure-gradient file
      constant                | Constant pressure drop
      DarcyForchheimer        | Darcy-Forchheimer model

The inherited entries are elaborated in:

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 218 of file directionalPressureGradientExplicitSource.H.

Member Enumeration Documentation

◆ pressureDropModel

Enumerator
pVolumetricFlowRateTable 
pConstant 
pDarcyForchheimer 

Definition at line 227 of file directionalPressureGradientExplicitSource.H.

Constructor & Destructor Documentation

◆ directionalPressureGradientExplicitSource() [1/2]

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

◆ directionalPressureGradientExplicitSource() [2/2]

◆ ~directionalPressureGradientExplicitSource()

virtual ~directionalPressureGradientExplicitSource ( )
virtualdefault

Member Function Documentation

◆ TypeName()

TypeName ( "directionalPressureGradientExplicitSource"  )

◆ operator=()

void operator= ( const directionalPressureGradientExplicitSource )
delete

◆ correct()

void correct ( volVectorField U)
virtual

◆ addSup() [1/2]

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

Reimplemented from option.

Definition at line 447 of file directionalPressureGradientExplicitSource.C.

References fvMatrix::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

Reimplemented from option.

Definition at line 473 of file directionalPressureGradientExplicitSource.C.

◆ constrain()

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

Reimplemented from option.

Definition at line 484 of file directionalPressureGradientExplicitSource.C.

References fvMatrix::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

Reimplemented from option.

Definition at line 518 of file directionalPressureGradientExplicitSource.C.

References NotImplemented.

◆ read()

bool read ( const dictionary dict)
virtual

Reimplemented from cellSetOption.

Definition at line 527 of file directionalPressureGradientExplicitSource.C.

References dict, dictionary::getOrDefault(), and dictionary::readEntry().

Here is the call graph for this function:

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