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...
Public Types | |
enum | pressureDropModel { pVolumetricFlowRateTable, pConstant, pDarcyForchheimer } |
Modes of pressure drop. More... | |
![]() | |
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 &) const |
Write the source properties. More... | |
virtual bool | read (const dictionary &dict) |
Read source dictionary. More... | |
![]() | |
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 () |
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 selectionModeType & | selectionMode () const |
Return const access to the cell selection mode. More... | |
const word & | cellSetName () 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 labelList & | cells () 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... | |
![]() | |
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< option > | clone () const |
Return clone. More... | |
virtual | ~option () |
Destructor. More... | |
const word & | name () const |
Return const access to the source name. More... | |
const fvMesh & | mesh () const |
Return const access to the mesh database. More... | |
const dictionary & | coeffs () 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... | |
Switch & | active () |
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 | writeHeader (Ostream &) const |
Write the source header information. More... | |
virtual void | writeFooter (Ostream &) const |
Write the source footer information. More... | |
Private Member Functions | |
void | initialise () |
Init. More... | |
void | writeProps (const vectorField &gradP) const |
Write the pressure gradient to file (for restarts etc) More... | |
void | update (fvMatrix< vector > &eqn) |
Correct driving force for a constant mass flow rate. More... | |
directionalPressureGradientExplicitSource (const directionalPressureGradientExplicitSource &) | |
Disallow default bitwise copy construct. More... | |
void | operator= (const directionalPressureGradientExplicitSource &) |
Disallow default bitwise assignment. More... | |
Private Attributes | |
pressureDropModel | model_ |
Pressure drop model. More... | |
vectorField | gradP0_ |
Pressure gradient before correction. More... | |
vectorField | dGradP_ |
Change in pressure gradient. More... | |
vectorField | gradPporous_ |
Pressure drop due to porous media. More... | |
vector | flowDir_ |
Flow direction. More... | |
autoPtr< volScalarField > | invAPtr_ |
Matrix 1/A coefficients field pointer. More... | |
scalar | D_ |
Darcy pressure loss coefficient. More... | |
scalar | I_ |
Inertia pressure lost coefficient. More... | |
scalar | length_ |
Porous media length. More... | |
scalar | pressureDrop_ |
Constant pressure drop. More... | |
interpolationTable< scalar > | flowRate_ |
Volumetric flow rate vs pressure drop table. More... | |
word | faceZoneName_ |
Name of the faceZone at the heat exchange inlet. More... | |
label | zoneID_ |
Id for the face zone. More... | |
labelList | faceId_ |
Local list of face IDs. More... | |
labelList | facePatchId_ |
Local list of patch ID per face. More... | |
scalar | relaxationFactor_ |
Relaxation factor. More... | |
labelList | cellFaceMap_ |
Cells faces mapping. More... | |
Static Private Attributes | |
static const NamedEnum< pressureDropModel, 3 > | PressureDropModelNames_ |
Additional Inherited Members | |
![]() | |
static autoPtr< option > | New (const word &name, const dictionary &dict, const fvMesh &mesh) |
Return a reference to the selected fvOption model. More... | |
![]() | |
static const NamedEnum< selectionModeType, 4 > | selectionModeTypeNames_ |
Word list of selection mode type names. More... | |
![]() | |
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... | |
![]() | |
scalar | timeStart_ |
Time start. More... | |
scalar | duration_ |
Duration. More... | |
selectionModeType | selectionMode_ |
Cell selection mode. More... | |
word | cellSetName_ |
Name of cell set for "cellSet" and "cellZone" selectionMode. More... | |
List< point > | points_ |
List of points for "points" selectionMode. More... | |
labelList | cells_ |
Set of cells to apply source to. More... | |
scalar | V_ |
Sum of cell volumes. More... | |
![]() | |
const word | name_ |
Source name. More... | |
const word | modelType_ |
Model type. More... | |
const fvMesh & | mesh_ |
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< bool > | applied_ |
Applied flag list - corresponds to each fieldNames_ entry. More... | |
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.
airDeflection { type directionalPressureGradientExplicitSource; active true; directionalPressureGradientExplicitSourceCoeffs { selectionMode cellZone; cellZone cZone; fieldNames (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.
Definition at line 104 of file directionalPressureGradientExplicitSource.H.
enum pressureDropModel |
Modes of pressure drop.
Enumerator | |
---|---|
pVolumetricFlowRateTable | |
pConstant | |
pDarcyForchheimer |
Definition at line 113 of file directionalPressureGradientExplicitSource.H.
|
private |
Disallow default bitwise copy construct.
directionalPressureGradientExplicitSource | ( | const word & | sourceName, |
const word & | modelType, | ||
const dictionary & | dict, | ||
const fvMesh & | mesh | ||
) |
Construct from explicit source name and mesh.
Definition at line 173 of file directionalPressureGradientExplicitSource.C.
References Foam::endl(), Foam::exit(), Foam::FatalError, FatalErrorInFunction, IOstream::good(), Foam::Info, Foam::mag(), Foam::name(), Foam::nl, dictionary::null, propsDict(), and Foam::type().
|
private |
Init.
Definition at line 86 of file directionalPressureGradientExplicitSource.C.
References polyMesh::boundaryMesh(), faceId(), directionalPressureGradientExplicitSource::faceId_, directionalPressureGradientExplicitSource::facePatchId_, polyMesh::faceZones(), forAll, primitiveMesh::isInternalFace(), option::mesh_, List::setSize(), List::size(), polyPatch::start(), polyPatch::whichFace(), polyBoundaryMesh::whichPatch(), and directionalPressureGradientExplicitSource::zoneID_.
|
private |
Write the pressure gradient to file (for restarts etc)
Definition at line 144 of file directionalPressureGradientExplicitSource.C.
References gradP(), IOobject::NO_READ, IOobject::NO_WRITE, and propsDict().
Correct driving force for a constant mass flow rate.
|
private |
Disallow default bitwise assignment.
TypeName | ( | "directionalPressureGradientExplicitSource" | ) |
Runtime type information.
|
virtual |
Correct the pressure gradient.
Reimplemented from option.
Definition at line 267 of file directionalPressureGradientExplicitSource.C.
References cellId, fvPatchField::coupled(), 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, TurbulenceModel< Alpha, Rho, BasicTurbulenceModel, TransportModel >::nu(), nu, polyBoundaryMesh::patchID(), fvPatchField::patchNeighbourField(), phi, turbulenceModel::propertiesName, rAU(), Foam::reduce(), rho, Foam::rotationTensor(), List::size(), PtrList::size(), faceZone::slaveCells(), U, w(), and Vector< scalar >::zero.
Add explicit contribution to momentum equation.
Reimplemented from option.
Definition at line 473 of file directionalPressureGradientExplicitSource.C.
References fvMatrix::dimensions(), Foam::dimVolume, IOobject::NO_READ, IOobject::NO_WRITE, Foam::fvc::Su(), and Vector< scalar >::zero.
|
virtual |
Add explicit contribution to compressible momentum equation.
Reimplemented from option.
Definition at line 499 of file directionalPressureGradientExplicitSource.C.
Set 1/A coefficient.
Reimplemented from option.
Definition at line 510 of file directionalPressureGradientExplicitSource.C.
References fvMatrix::A(), IOobject::NO_READ, IOobject::NO_WRITE, and Vector< scalar >::zero.
|
virtual |
Write the source properties.
Reimplemented from option.
Definition at line 31 of file directionalPressureGradientExplicitSourceIO.C.
References NotImplemented.
|
virtual |
Read source dictionary.
Reimplemented from cellSetOption.
Definition at line 40 of file directionalPressureGradientExplicitSourceIO.C.
References dict, dictionary::lookup(), dictionary::lookupOrDefault(), and Foam::mag().
|
staticprivate |
Definition at line 125 of file directionalPressureGradientExplicitSource.H.
|
private |
Pressure drop model.
Definition at line 128 of file directionalPressureGradientExplicitSource.H.
|
private |
Pressure gradient before correction.
Definition at line 131 of file directionalPressureGradientExplicitSource.H.
|
private |
Change in pressure gradient.
Definition at line 134 of file directionalPressureGradientExplicitSource.H.
|
private |
Pressure drop due to porous media.
Definition at line 137 of file directionalPressureGradientExplicitSource.H.
|
private |
Flow direction.
Definition at line 140 of file directionalPressureGradientExplicitSource.H.
|
private |
Matrix 1/A coefficients field pointer.
Definition at line 143 of file directionalPressureGradientExplicitSource.H.
|
private |
Darcy pressure loss coefficient.
Definition at line 146 of file directionalPressureGradientExplicitSource.H.
|
private |
Inertia pressure lost coefficient.
Definition at line 149 of file directionalPressureGradientExplicitSource.H.
|
private |
Porous media length.
Definition at line 152 of file directionalPressureGradientExplicitSource.H.
|
private |
Constant pressure drop.
Definition at line 155 of file directionalPressureGradientExplicitSource.H.
|
private |
Volumetric flow rate vs pressure drop table.
Definition at line 158 of file directionalPressureGradientExplicitSource.H.
|
private |
Name of the faceZone at the heat exchange inlet.
Definition at line 161 of file directionalPressureGradientExplicitSource.H.
|
private |
Id for the face zone.
Definition at line 164 of file directionalPressureGradientExplicitSource.H.
Referenced by directionalPressureGradientExplicitSource::initialise().
|
private |
Local list of face IDs.
Definition at line 167 of file directionalPressureGradientExplicitSource.H.
Referenced by directionalPressureGradientExplicitSource::initialise().
|
private |
Local list of patch ID per face.
Definition at line 170 of file directionalPressureGradientExplicitSource.H.
Referenced by directionalPressureGradientExplicitSource::initialise().
|
private |
Relaxation factor.
Definition at line 173 of file directionalPressureGradientExplicitSource.H.
|
private |
Cells faces mapping.
Definition at line 176 of file directionalPressureGradientExplicitSource.H.
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.