26 #include "VoFTurbulenceDamping.H" 55 const word& sourceName,
56 const word& modelType,
57 const dictionary& dict,
61 fvModel(sourceName, modelType, dict, mesh),
62 phaseName_(dict.lookupOrDefault(
"phase", word::null)),
66 mesh.lookupObject<incompressibleTwoPhaseMixture>
71 interface_(
refCast<const interfaceProperties>(mixture_)),
74 mesh.lookupObject<incompressibleMomentumTransportModel>
78 betaStar_(
"betaStar",
dimless, 0),
86 fieldName_ = epsilonName;
87 C2_.read(turbulence_.coeffDict());
91 fieldName_ = omegaName;
92 betaStar_.read(turbulence_.coeffDict());
95 if (turbulence_.coeffDict().found(
"beta"))
97 beta_.read(turbulence_.coeffDict());
108 <<
"Cannot find either " << epsilonName <<
" or " << omegaName
124 fvMatrix<scalar>& eqn,
125 const word& fieldName
130 Info<<
type() <<
": applying source to " << eqn.psi().name() <<
endl;
135 mixture_.alpha1()()*
sqr(mixture_.nuModel1().nu()()())
136 + mixture_.alpha2()()*
sqr(mixture_.nuModel2().nu()()())
139 if (fieldName ==
"epsilon")
141 eqn += interface_.fraction()*C2_*aSqrnu*turbulence_.k()()/
pow4(delta_);
143 else if (fieldName ==
"omega")
145 eqn += interface_.fraction()*beta_*aSqrnu/(
sqr(betaStar_)*
pow4(delta_));
150 <<
"Support for field " << fieldName <<
" is not implemented" 160 fvMatrix<scalar>& eqn,
161 const word& fieldName
166 Info<<
type() <<
": applying source to " << eqn.psi().name() <<
endl;
169 tmp<volScalarField::Internal> taSqrnu;
171 if (mixture_.alpha1().name() == alpha.name())
173 taSqrnu = mixture_.alpha1()()*
sqr(mixture_.nuModel1().nu()()());
175 else if (mixture_.alpha2().name() == alpha.name())
177 taSqrnu = mixture_.alpha2()()*
sqr(mixture_.nuModel2().nu()()());
182 <<
"Unknown phase-fraction " << alpha.name()
188 eqn += interface_.fraction()*C2_*taSqrnu*turbulence_.k()()
193 eqn += interface_.fraction()*beta_*taSqrnu
194 /(
sqr(betaStar_)*
pow4(delta_));
199 <<
"Support for field " << fieldName <<
" is not implemented" virtual bool movePoints()
Update for mesh motion.
defineTypeNameAndDebug(fixedTemperatureConstraint, 0)
errorManipArg< error, int > exit(error &err, const int errNo=1)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
virtual void addSup(const volScalarField &rho, fvMatrix< scalar > &eqn, const word &fieldName) const
Add explicit contribution to epsilon or omega equation.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
To & refCast(From &r)
Reference type cast template function.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Ostream & endl(Ostream &os)
Add newline and flush stream.
VoFTurbulenceDamping(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
const dimensionSet dimless
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
Macros for easy insertion into run-time selection tables.
const dimensionSet dimLength
GeometricField< scalar, fvPatchField, volMesh > volScalarField
compressibleMomentumTransportModel momentumTransportModel
static word groupName(Name name, const word &group)
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
addToRunTimeSelectionTable(fvConstraint, fixedTemperatureConstraint, dictionary)
List< word > wordList
A List of words.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual wordList addSupFields() const
Return the list of fields for which the option adds source term.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
dimensionedScalar pow4(const dimensionedScalar &ds)