55 void Foam::fv::propellerDisk::readCoeffs(
const dictionary&
dict)
69 <<
"disk direction vector is approximately zero"
73 n_ =
dict.lookup<scalar>(
"n");
83 log_ =
dict.lookupOrDefault<Switch>(
"log",
false);
90 {
"n",
"J",
"Jcorr",
"Udisk",
"Ucorr",
"Kt",
"Kq",
"T/rho",
"Q/rho"}
100 const Field<vector> zoneCellCentres(mesh().cellCentres(), set_.cells());
101 const Field<scalar> zoneCellVolumes(mesh().cellVolumes(), set_.cells());
103 return gSum(zoneCellVolumes*zoneCellCentres)/set_.V();
109 const Field<vector> zoneCellCentres(mesh().cellCentres(), set_.cells());
110 const scalar r2 =
gMax(
magSqr(zoneCellCentres - centre));
132 const scalar Uref = VUn/set_.V();
134 return mag(Uref/(n_*dProp_));
143 const word& modelType,
149 set_(mesh, coeffs()),
150 phaseName_(
word::null),
165 readCoeffs(coeffs());
189 addActuationDiskAxialInertialResistance
205 addActuationDiskAxialInertialResistance
222 addActuationDiskAxialInertialResistance
241 set_.topoChange(map);
256 set_.distribute(map);
#define forAll(list, i)
Loop across all elements in list.
Macros for easy insertion into run-time selection tables.
static autoPtr< Function1< Type > > New(const word &name, const Function1s::unitConversions &units, const dictionary &dict)
Select from dictionary.
Generic GeometricField class.
static word groupName(Name name, const word &group)
static bool master(const label communicator=0)
Am I the master process.
A list of keyword definitions, which are a keyword followed by any number of values (e....
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Mesh data needed to do the Finite Volume discretisation.
Finite volume model abstract base class.
const dictionary & coeffs() const
Return dictionary.
virtual bool read(const dictionary &dict)
Read source dictionary.
const fvMesh & mesh() const
Return const access to the mesh database.
const word & name() const
Return const access to the source name.
Disk momentum source which approximates a propeller based on a given propeller curve.
scalar dProp_
Propeller diameter.
virtual bool movePoints()
Update for mesh motion.
virtual void addSup(const volVectorField &U, fvMatrix< vector > &eqn) const
Source term to momentum equation.
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
word UName_
Name of the velocity field.
word phaseName_
The name of the phase to which this fvModel applies.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
scalar dHub_
Hub diameter.
virtual bool read(const dictionary &dict)
Read dictionary.
scalar diskThickness(const vector ¢re) const
Computes the thickness of the disk in streamwise direction.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
autoPtr< Function1< vector2D > > propellerFunction_
Propeller function.
vector diskNormal_
Disk normal direction.
scalar rotationDir_
Rotation direction (obtained from the sign of n_)
scalar J(const vectorField &U, const vector &nHat) const
Return the normalised flow-rate through the disk.
propellerDisk(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
autoPtr< functionObjects::logFile > logFile_
Optional log file.
scalar n_
Rotation speed [1/s].
Switch log_
Optional switch to enable logging of integral properties.
vector diskCentre() const
Computes the disk centre.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for handling words, derived from string.
static const word null
An empty word.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
errorManipArg< error, int > exit(error &err, const int errNo=1)
List< word > wordList
A List of words.
Type gSum(const FieldField< Field, Type > &f)
dimensionedScalar sign(const dimensionedScalar &ds)
word name(const bool)
Return a word representation of a bool.
const dimensionSet dimless
Vector< scalar > vector
A scalar version of the templated Vector.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
dimensioned< scalar > mag(const dimensioned< Type > &)
makeFoamTableReaders(avType, nullArg)
makeFunction1s(avType, nullArg)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
Type gMax(const FieldField< Field, Type > &f)