80 const word& sourceName,
81 const word& modelType,
86 fv::cellSetOption(sourceName, modelType,
dict,
mesh),
87 Ubar_(coeffs_.get<
vector>(
"Ubar")),
90 flowDir_(Ubar_/
mag(Ubar_)),
91 relaxation_(coeffs_.getOrDefault<scalar>(
"relaxation", 1)),
99 <<
"settings are:" << fieldNames_ << exit(FatalError);
110 if (propsFile.good())
112 Info<<
" Reading pressure gradient from file" << endl;
113 dictionary propsDict(propsFile);
114 propsDict.readEntry(
"gradient", gradP0_);
128 scalar magUbarAve = 0.0;
133 label celli = cells_[i];
134 scalar volCell = cv[celli];
135 magUbarAve += (flowDir_ &
U[celli])*volCell;
155 label celli = cells_[i];
156 scalar volCell = cv[celli];
157 rAUave +=
rAU[celli]*volCell;
166 scalar magUbarAve = this->magUbarAve(
U);
170 dGradP_ = relaxation_*(
mag(Ubar_) - magUbarAve)/rAUave;
175 label celli = cells_[i];
176 U[celli] += flowDir_*
rAU[celli]*dGradP_;
179 U.correctBoundaryConditions();
181 scalar
gradP = gradP0_ + dGradP_;
184 <<
", pressure gradient = " <<
gradP <<
endl;
198 name_ + fieldNames_[fieldi] +
"Sup",
203 auto&
Su = tSu.ref();
205 scalar
gradP = gradP0_ + dGradP_;
220 this->
addSup(eqn, fieldi);
243 rAPtr_() = 1.0/eqn.
A();
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
IOdictionary propsDict(dictIO)
static tmp< GeometricField< vector, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< vector >::calculatedType())
Input from file stream as an ISstream, normally using std::ifstream for the actual input.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
bool writeTime() const noexcept
True if this is a write interval.
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
fileName timePath() const
Return current time path = path/timeName.
A List with indirect addressing. Like IndirectList but does not store addressing.
void size(const label n)
Older name for setAddressableSize.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect,...
A special matrix type and solver, designed for finite volume solutions of scalar equations....
tmp< volScalarField > A() const
Return the central coefficient.
const dimensionSet & dimensions() const noexcept
Mesh data needed to do the Finite Volume discretisation.
const Time & time() const
Return the top-level database.
labelList cells_
Set of cells to apply source to.
cellSetOption(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
scalar V_
Sum of cell volumes.
Applies the force within a specified region to maintain the specified mean velocity for incompressibl...
vector Ubar_
Desired mean velocity.
meanVelocityForce(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
virtual bool read(const dictionary &dict)
Read source dictionary.
scalar dGradP_
Change in pressure gradient.
autoPtr< volScalarField > rAPtr_
Matrix 1/A coefficients field pointer.
virtual void constrain(fvMatrix< vector > &eqn, const label fieldi)
Set 1/A coefficient.
scalar gradP0_
Pressure gradient before correction.
scalar relaxation_
Relaxation factor.
virtual void addSup(fvMatrix< vector > &eqn, const label fieldi)
Add explicit contribution to momentum equation.
void writeProps(const scalar gradP) const
Write the pressure gradient to file (for restarts etc).
virtual scalar magUbarAve(const volVectorField &U) const
Calculate and return the magnitude of the mean velocity averaged over the selected cellSet.
virtual void correct(volVectorField &U)
Correct the pressure gradient.
vector flowDir_
Flow direction.
Base abstract class for handling finite volume options (i.e. fvOption).
const fvMesh & mesh_
Reference to the mesh database.
wordList fieldNames_
Field names to apply source to - populated by derived models.
dictionary coeffs_
Dictionary containing source coefficients.
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
const fvMesh & mesh() const noexcept
Return const access to the mesh database.
const word name_
Source name.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
volVectorField gradP(fvc::grad(p))
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
A special matrix type and solver, designed for finite volume solutions of scalar equations.
tmp< volScalarField > rAU
Namespace for finite-volume.
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
static constexpr const zero Zero
Global zero (0).
const dimensionSet dimVolume(pow3(dimLength))
void Su(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
constexpr char nl
The newline '\n' character (0x0a).
#define forAll(list, i)
Loop across all elements in list.