Go to the documentation of this file.
94 #ifndef meanVelocityForce_H
95 #define meanVelocityForce_H
113 class meanVelocityForce
115 public fv::cellSetOption
137 autoPtr<volScalarField>
rAPtr_;
164 const word& sourceName,
165 const word& modelType,
volVectorField gradP(fvc::grad(p))
A class for handling words, derived from Foam::string.
Intermediate abstract class for handling cell-set options for the derived fvOptions.
meanVelocityForce(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
autoPtr< volScalarField > rAPtr_
Matrix 1/A coefficients field pointer.
scalar relaxation_
Relaxation factor.
void writeProps(const scalar gradP) const
Write the pressure gradient to file (for restarts etc)
void update(fvMatrix< vector > &eqn)
Correct driving force for a constant mass flow rate.
vector flowDir_
Flow direction.
vector Ubar_
Desired mean velocity.
virtual bool read(const dictionary &dict)
Read source dictionary.
scalar dGradP_
Change in pressure gradient.
virtual void addSup(fvMatrix< vector > &eqn, const label fieldi)
Add explicit contribution to momentum equation.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Mesh data needed to do the Finite Volume discretisation.
Vector< scalar > vector
A scalar version of the templated Vector.
virtual ~meanVelocityForce()=default
Destructor.
TypeName("meanVelocityForce")
Runtime type information.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Applies the force within a specified region to maintain the specified mean velocity for incompressibl...
const fvMesh & mesh() const noexcept
Return const access to the mesh database.
virtual void constrain(fvMatrix< vector > &eqn, const label fieldi)
Set 1/A coefficient.
virtual scalar magUbarAve(const volVectorField &U) const
scalar gradP0_
Pressure gradient before correction.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
void operator=(const meanVelocityForce &)=delete
No copy assignment.
virtual void correct(volVectorField &U)
Correct the pressure gradient.