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...
Public Types | |
enum | pressureDropModel { pVolumetricFlowRateTable, pConstant, pDarcyForchheimer } |
![]() | |
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) |
![]() | |
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 word & | cellSetName () const noexcept |
scalar | V () const noexcept |
const labelList & | cells () const noexcept |
scalar | timeStart (scalar val) noexcept |
scalar | duration (scalar val) noexcept |
virtual bool | isActive () |
![]() | |
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< option > | clone () const |
virtual | ~option ()=default |
const word & | name () const noexcept |
const fvMesh & | mesh () const noexcept |
const dictionary & | coeffs () 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 autoPtr< option > | New (const word &name, const dictionary &dict, const fvMesh &mesh) |
![]() | |
bool | log |
![]() | |
static const Enum< selectionModeType > | selectionModeTypeNames_ |
![]() | |
void | setSelection (const dictionary &dict) |
void | setCellSelection () |
void | setVol () |
![]() | |
void | resetApplied () |
![]() | |
scalar | timeStart_ |
scalar | duration_ |
selectionModeType | selectionMode_ |
word | cellSetName_ |
List< point > | points_ |
labelList | cells_ |
scalar | V_ |
![]() | |
const word | name_ |
const word | modelType_ |
const fvMesh & | mesh_ |
dictionary | dict_ |
dictionary | coeffs_ |
wordList | fieldNames_ |
List< bool > | applied_ |
bool | active_ |
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.
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:
Definition at line 218 of file directionalPressureGradientExplicitSource.H.
enum pressureDropModel |
Enumerator | |
---|---|
pVolumetricFlowRateTable | |
pConstant | |
pDarcyForchheimer |
Definition at line 227 of file directionalPressureGradientExplicitSource.H.
directionalPressureGradientExplicitSource | ( | const word & | sourceName, |
const word & | modelType, | ||
const dictionary & | dict, | ||
const fvMesh & | mesh | ||
) |
Definition at line 150 of file directionalPressureGradientExplicitSource.C.
References Foam::endl(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, IOstream::good(), Foam::Info, Foam::name(), Foam::nl, propsDict(), option::resetApplied(), and Foam::type().
|
delete |
|
virtualdefault |
TypeName | ( | "directionalPressureGradientExplicitSource" | ) |
|
delete |
|
virtual |
Reimplemented from option.
Definition at line 243 of file directionalPressureGradientExplicitSource.C.
References cellId, fvPatchField::coupled(), D, Foam::expressions::patchExpr::debug, Foam::dimArea, Foam::dimVelocity, Foam::endl(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, faceZone::flipMap(), forAll, Foam::gSum(), Foam::gSumProd(), Foam::Info, Foam::mag(), faceZone::masterCells(), mesh, Foam::constant::physicoChemical::mu, nu, polyBoundaryMesh::patchID(), fvPatchField::patchNeighbourField(), phi, turbulenceModel::propertiesName, rAU, Foam::reduce(), rho, Foam::rotationTensor(), faceZone::slaveCells(), U, and Foam::Zero.
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.
|
virtual |
Reimplemented from option.
Definition at line 473 of file directionalPressureGradientExplicitSource.C.
Reimplemented from option.
Definition at line 484 of file directionalPressureGradientExplicitSource.C.
References fvMatrix::A(), IOobject::NO_READ, IOobject::NO_WRITE, and Foam::Zero.
|
virtual |
Reimplemented from option.
Definition at line 518 of file directionalPressureGradientExplicitSource.C.
References NotImplemented.
|
virtual |
Reimplemented from cellSetOption.
Definition at line 527 of file directionalPressureGradientExplicitSource.C.
References dict, dictionary::getOrDefault(), and dictionary::readEntry().
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.