Actuation disk source. More...


Public Member Functions | |
| TypeName ("actuationDisk") | |
| Runtime type information. More... | |
| actuationDisk (const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict) | |
| Construct from components. More... | |
| actuationDisk (const actuationDisk &)=delete | |
| Disallow default bitwise copy construction. More... | |
| virtual | ~actuationDisk () |
| Destructor. More... | |
| virtual wordList | addSupFields () const |
| Return the list of fields for which the fvModel adds source term. More... | |
| virtual void | addSup (const volVectorField &U, fvMatrix< vector > &eqn) const |
| Source term to momentum equation. More... | |
| virtual void | addSup (const volScalarField &rho, const volVectorField &U, fvMatrix< vector > &eqn) const |
| Source term to compressible momentum equation. More... | |
| virtual void | addSup (const volScalarField &alpha, const volScalarField &rho, const volVectorField &U, fvMatrix< vector > &eqn) const |
| Explicit and implicit sources for phase equations. More... | |
| virtual bool | movePoints () |
| Update for mesh motion. More... | |
| virtual void | topoChange (const polyTopoChangeMap &) |
| Update topology using the given map. More... | |
| virtual void | mapMesh (const polyMeshMap &) |
| Update from another mesh using the given map. More... | |
| virtual void | distribute (const polyDistributionMap &) |
| Redistribute or update using the given distribution map. More... | |
| virtual bool | read (const dictionary &dict) |
| Read dictionary. More... | |
| void | operator= (const actuationDisk &)=delete |
| Disallow default bitwise assignment. More... | |
Public Member Functions inherited from fvModel | |
| TypeName ("fvModel") | |
| Runtime type information. More... | |
| declareRunTimeSelectionTable (autoPtr, fvModel, dictionary,(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict),(name, modelType, mesh, dict)) | |
| fvModel (const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict) | |
| Construct from components. More... | |
| autoPtr< fvModel > | clone () const |
| Return clone. More... | |
| virtual | ~fvModel () |
| Destructor. More... | |
| const word & | name () const |
| Return const access to the source name. More... | |
| const word & | keyword () const |
| Return name as the keyword. More... | |
| const fvMesh & | mesh () const |
| Return const access to the mesh database. More... | |
| const dictionary & | coeffs (const dictionary &) const |
| Return the coefficients sub-dictionary. More... | |
| virtual bool | addsSupToField (const word &fieldName) const |
| Return true if the fvModel adds a source term to the given. More... | |
| virtual scalar | maxDeltaT () const |
| Return the maximum time-step for stable operation. More... | |
| virtual void | addSup (fvMatrix< scalar > &eqn) const |
| Add a source term to a field-less proxy equation. More... | |
| template<class Type > | |
| tmp< fvMatrix< Type > > | sourceProxy (const VolField< Type > &eqnField) const |
| Add a source term to an equation. More... | |
| template<class Type > | |
| tmp< fvMatrix< Type > > | source (const VolField< Type > &field) const |
| Return source for an equation. More... | |
| template<class Type > | |
| tmp< fvMatrix< Type > > | sourceProxy (const VolField< Type > &field, const VolField< Type > &eqnField) const |
| Return source for an equation. More... | |
| template<class Type > | |
| tmp< fvMatrix< Type > > | source (const volScalarField &rho, const VolField< Type > &field) const |
| Return source for a compressible equation. More... | |
| template<class Type > | |
| tmp< fvMatrix< Type > > | sourceProxy (const volScalarField &rho, const VolField< Type > &field, const VolField< Type > &eqnField) const |
| Return source for a compressible equation. More... | |
| template<class Type > | |
| tmp< fvMatrix< Type > > | source (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &field) const |
| Return source for a phase equation. More... | |
| template<class Type > | |
| tmp< fvMatrix< Type > > | sourceProxy (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &field, const VolField< Type > &eqnField) const |
| Return source for a phase equation. More... | |
| template<class Type > | |
| tmp< fvMatrix< Type > > | source (const volScalarField &alpha, const geometricOneField &rho, const VolField< Type > &field) const |
| Return source for a phase equation. More... | |
| template<class Type > | |
| tmp< fvMatrix< Type > > | source (const geometricOneField &alpha, const volScalarField &rho, const VolField< Type > &field) const |
| Return source for a phase equation. More... | |
| template<class Type > | |
| tmp< fvMatrix< Type > > | source (const geometricOneField &alpha, const geometricOneField &rho, const VolField< Type > &field) const |
| Return source for a phase equation. More... | |
| template<class Type > | |
| tmp< fvMatrix< Type > > | d2dt2 (const VolField< Type > &field) const |
| Return source for an equation with a second time derivative. More... | |
| virtual void | preUpdateMesh () |
| Prepare for mesh update. More... | |
| virtual void | correct () |
| Correct the fvModel. More... | |
| virtual bool | write (const bool write=true) const |
| Write fvModel data. More... | |
| template<class AlphaRhoFieldType , class ... AlphaRhoFieldTypes> | |
| Foam::dimensionSet | sourceDims (const dimensionSet &ds, const AlphaRhoFieldType &alphaRhoField, const AlphaRhoFieldTypes &... alphaRhoFields) |
| template<class AlphaRhoFieldType , class ... AlphaRhoFieldTypes> | |
| const Foam::word & | fieldName (const AlphaRhoFieldType &alphaRhoField, const AlphaRhoFieldTypes &... alphaRhoFields) |
| template<class AlphaRhoFieldType > | |
| const Foam::word & | fieldName (const AlphaRhoFieldType &alphaRhoField) |
| template<class Type , class ... AlphaRhoFieldTypes> | |
| Foam::tmp< Foam::fvMatrix< Type > > | sourceTerm (const VolField< Type > &eqnField, const dimensionSet &ds, const AlphaRhoFieldTypes &... alphaRhoFields) const |
| template<class Type > | |
| Foam::tmp< Foam::fvMatrix< Type > > | sourceProxy (const VolField< Type > &eqnField) const |
| template<class Type > | |
| Foam::tmp< Foam::fvMatrix< Type > > | source (const VolField< Type > &field) const |
| template<class Type > | |
| Foam::tmp< Foam::fvMatrix< Type > > | sourceProxy (const VolField< Type > &field, const VolField< Type > &eqnField) const |
| template<class Type > | |
| Foam::tmp< Foam::fvMatrix< Type > > | source (const volScalarField &rho, const VolField< Type > &field) const |
| template<class Type > | |
| Foam::tmp< Foam::fvMatrix< Type > > | sourceProxy (const volScalarField &rho, const VolField< Type > &field, const VolField< Type > &eqnField) const |
| template<class Type > | |
| Foam::tmp< Foam::fvMatrix< Type > > | source (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &field) const |
| template<class Type > | |
| Foam::tmp< Foam::fvMatrix< Type > > | sourceProxy (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &field, const VolField< Type > &eqnField) const |
| template<class Type > | |
| Foam::tmp< Foam::fvMatrix< Type > > | source (const geometricOneField &alpha, const geometricOneField &rho, const VolField< Type > &field) const |
| template<class Type > | |
| Foam::tmp< Foam::fvMatrix< Type > > | source (const volScalarField &alpha, const geometricOneField &rho, const VolField< Type > &field) const |
| template<class Type > | |
| Foam::tmp< Foam::fvMatrix< Type > > | source (const geometricOneField &alpha, const volScalarField &rho, const VolField< Type > &field) const |
| template<class Type > | |
| Foam::tmp< Foam::fvMatrix< Type > > | d2dt2 (const VolField< Type > &field) const |
Protected Attributes | |
| fvCellZone | zone_ |
| The cellZone the fvConstraint applies to. More... | |
| word | phaseName_ |
| The name of the phase to which this fvModel applies. More... | |
| word | UName_ |
| Name of the velocity field. More... | |
| vector | diskDir_ |
| Disk area normal. More... | |
| scalar | Cp_ |
| Power coefficient. More... | |
| scalar | Ct_ |
| Thrust coefficient. More... | |
| scalar | diskArea_ |
| Disk area. More... | |
| point | upstreamPoint_ |
| Upstream point sample. More... | |
| label | upstreamCellId_ |
| Upstream cell ID. More... | |
Additional Inherited Members | |
Static Public Member Functions inherited from fvModel | |
| static const dictionary & | coeffs (const word &modelType, const dictionary &) |
| Return the coefficients sub-dictionary for a given model type. More... | |
| template<class AlphaRhoFieldType , class ... AlphaRhoFieldTypes> | |
| static dimensionSet | sourceDims (const dimensionSet &ds, const AlphaRhoFieldType &alphaRhoField, const AlphaRhoFieldTypes &... alphaRhoFields) |
| Return the dimensions of the matrix of a source term. More... | |
| static const dimensionSet & | sourceDims (const dimensionSet &ds) |
| Return the dimensions of the matrix of a source term (base. More... | |
| template<class AlphaRhoFieldType , class ... AlphaRhoFieldTypes> | |
| static const word & | fieldName (const AlphaRhoFieldType &alphaRhoField, const AlphaRhoFieldTypes &... alphaRhoFields) |
| Return the name of the field associated with a source term. More... | |
| template<class AlphaRhoFieldType > | |
| static const word & | fieldName (const AlphaRhoFieldType &alphaRhoField) |
| Return the name of the field associated with a source term (base. More... | |
| static const word & | fieldName () |
| Return the name of the field associated with a source term. Special. More... | |
| static autoPtr< fvModel > | New (const word &name, const fvMesh &mesh, const dictionary &dict) |
| Return a reference to the selected fvModel. More... | |
Static Public Attributes inherited from fvModel | |
| static const wordHashSet | keywords |
| The keywords read by this class. More... | |
Protected Member Functions inherited from fvModel | |
| template<class Type > | |
| void | addSupType (const VolField< Type > &field, fvMatrix< Type > &eqn) const |
| Add a source term to an equation. More... | |
| template<class Type > | |
| void | addSupType (const volScalarField &rho, const VolField< Type > &field, fvMatrix< Type > &eqn) const |
| Add a source term to a compressible equation. More... | |
| template<class Type > | |
| void | addSupType (const volScalarField &alpha, const volScalarField &rho, const VolField< Type > &field, fvMatrix< Type > &eqn) const |
| Add a source term to a phase equation. More... | |
| template<class Type , class ... AlphaRhoFieldTypes> | |
| tmp< fvMatrix< Type > > | sourceTerm (const VolField< Type > &eqnField, const dimensionSet &ds, const AlphaRhoFieldTypes &... alphaRhoFields) const |
| Return a source for an equation. More... | |
Actuation disk source.
Constant values for momentum source for actuation disk
where:
| = | Disk area | |
| = | Unit disk direction | |
| = | Upstream velocity | |
| = | 1 - Cp/Ct | |
| = | Power coefficient | |
| = | Thrust coefficient |
actuationDisk1
{
type actuationDisk;
cellZone actuationDisk1;
diskDir (-1 0 0); // Disk direction
Cp 0.1; // Power coefficient
Ct 0.5; // Thrust coefficient
diskArea 5.0; // Disk area
upstreamPoint (0 0 0); // Upstream point
}
Definition at line 108 of file actuationDisk.H.
| actuationDisk | ( | const word & | name, |
| const word & | modelType, | ||
| const fvMesh & | mesh, | ||
| const dictionary & | dict | ||
| ) |
Construct from components.
Definition at line 134 of file actuationDisk.C.
References fvModel::coeffs(), and dict.

|
delete |
Disallow default bitwise copy construction.
|
inlinevirtual |
Destructor.
Definition at line 186 of file actuationDisk.H.
| TypeName | ( | "actuationDisk" | ) |
Runtime type information.
|
virtual |
Return the list of fields for which the fvModel adds source term.
to the transport equation
Reimplemented from fvModel.
Definition at line 159 of file actuationDisk.C.
|
virtual |
Source term to momentum equation.
Reimplemented in radialActuationDisk.
Definition at line 165 of file actuationDisk.C.
References mesh, fvMatrix< Type >::source(), and U.

|
virtual |
Source term to compressible momentum equation.
Reimplemented in radialActuationDisk.
Definition at line 183 of file actuationDisk.C.
References mesh, rho, fvMatrix< Type >::source(), and U.

|
virtual |
Explicit and implicit sources for phase equations.
Definition at line 202 of file actuationDisk.C.
References alpha(), mesh, rho, fvMatrix< Type >::source(), and U.

|
virtual |
|
virtual |
Update topology using the given map.
Implements fvModel.
Definition at line 229 of file actuationDisk.C.
|
virtual |
Update from another mesh using the given map.
Implements fvModel.
Definition at line 235 of file actuationDisk.C.
|
virtual |
Redistribute or update using the given distribution map.
Implements fvModel.
Definition at line 241 of file actuationDisk.C.
|
virtual |
Read dictionary.
Reimplemented from fvModel.
Reimplemented in radialActuationDisk.
Definition at line 247 of file actuationDisk.C.
References dict, and fvModel::read().
Referenced by radialActuationDisk::read().


|
delete |
Disallow default bitwise assignment.
|
protected |
The cellZone the fvConstraint applies to.
Definition at line 117 of file actuationDisk.H.
|
protected |
The name of the phase to which this fvModel applies.
Definition at line 120 of file actuationDisk.H.
|
protected |
Name of the velocity field.
Definition at line 123 of file actuationDisk.H.
|
protected |
Disk area normal.
Definition at line 126 of file actuationDisk.H.
|
protected |
Power coefficient.
Definition at line 129 of file actuationDisk.H.
|
protected |
Thrust coefficient.
Definition at line 132 of file actuationDisk.H.
|
protected |
Disk area.
Definition at line 135 of file actuationDisk.H.
|
protected |
Upstream point sample.
Definition at line 138 of file actuationDisk.H.
|
protected |
Upstream cell ID.
Definition at line 141 of file actuationDisk.H.