VoFTurbulenceDamping.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2021-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "VoFTurbulenceDamping.H"
28 #include "interfaceProperties.H"
30 #include "fvMatrix.H"
32 
33 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  namespace fv
38  {
39  defineTypeNameAndDebug(VoFTurbulenceDamping, 0);
40 
42  (
43  fvModel,
44  VoFTurbulenceDamping,
45  dictionary
46  );
47  }
48 }
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
54 (
55  const word& sourceName,
56  const word& modelType,
57  const dictionary& dict,
58  const fvMesh& mesh
59 )
60 :
61  fvModel(sourceName, modelType, dict, mesh),
62  phaseName_(dict.lookupOrDefault("phase", word::null)),
63  delta_("delta", dimLength, dict),
64  mixture_
65  (
66  mesh.lookupObject<compressibleTwoPhaseMixture>
67  (
68  "phaseProperties"
69  )
70  ),
71  interface_(refCast<const interfaceProperties>(mixture_)),
72  turbulence_
73  (
74  mesh.lookupObject<compressibleMomentumTransportModel>
75  (IOobject::groupName(momentumTransportModel::typeName, phaseName_))
76  ),
77  C2_("C2", dimless, 0),
78  betaStar_("betaStar", dimless, 0),
79  beta_("beta", dimless, 0)
80 {
81  const word epsilonName(IOobject::groupName("epsilon", phaseName_));
82  const word omegaName(IOobject::groupName("omega", phaseName_));
83 
84  if (mesh.foundObject<volScalarField>(epsilonName))
85  {
86  fieldName_ = epsilonName;
87  C2_.read(turbulence_.coeffDict());
88  }
89  else if (mesh.foundObject<volScalarField>(omegaName))
90  {
91  fieldName_ = omegaName;
92  betaStar_.read(turbulence_.coeffDict());
93 
94  // Read beta for k-omega models or beta1 for k-omega SST
95  if (turbulence_.coeffDict().found("beta"))
96  {
97  beta_.read(turbulence_.coeffDict());
98  }
99  else
100  {
101  beta_ =
102  dimensionedScalar("beta1", dimless, turbulence_.coeffDict());
103  }
104  }
105  else
106  {
108  << "Cannot find either " << epsilonName << " or " << omegaName
109  << " field for fvModel " << typeName << exit(FatalIOError);
110  }
111 }
112 
113 
114 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
115 
117 {
118  return wordList(1, fieldName_);
119 }
120 
121 
123 (
124  const volScalarField& rho,
125  fvMatrix<scalar>& eqn,
126  const word& fieldName
127 ) const
128 {
129  if (debug)
130  {
131  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
132  }
133 
134  const volScalarField::Internal aRhoSqrnu
135  (
136  mixture_.alpha1()()*mixture_.rho1()()*sqr(mixture_.thermo1().nu()()())
137  + mixture_.alpha2()()*mixture_.rho2()()*sqr(mixture_.thermo2().nu()()())
138  );
139 
140  if (fieldName == "epsilon")
141  {
142  eqn += interface_.fraction()*C2_*aRhoSqrnu*turbulence_.k()()
143  /pow4(delta_);
144  }
145  else if (fieldName == "omega")
146  {
147  eqn += interface_.fraction()*beta_*aRhoSqrnu
148  /(sqr(betaStar_)*pow4(delta_));
149  }
150  else
151  {
153  << "Support for field " << fieldName << " is not implemented"
154  << exit(FatalError);
155  }
156 }
157 
158 
160 (
161  const volScalarField& alpha,
162  const volScalarField&,
163  fvMatrix<scalar>& eqn,
164  const word& fieldName
165 ) const
166 {
167  if (debug)
168  {
169  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
170  }
171 
172  tmp<volScalarField::Internal> taRhoSqrnu;
173 
174  if (mixture_.alpha1().name() == alpha.name())
175  {
176  taRhoSqrnu = mixture_.alpha1()()*mixture_.rho1()()
177  *sqr(mixture_.thermo1().nu()()());
178  }
179  else if (mixture_.alpha2().name() == alpha.name())
180  {
181  taRhoSqrnu = mixture_.alpha2()()*mixture_.rho2()()
182  *sqr(mixture_.thermo2().nu()()());
183  }
184  else
185  {
187  << "Unknown phase-fraction " << alpha.name()
188  << exit(FatalError);
189  }
190 
191  if (fieldName == IOobject::groupName("epsilon", phaseName_))
192  {
193  eqn += interface_.fraction()*C2_*taRhoSqrnu*turbulence_.k()()
194  /pow4(delta_);
195  }
196  else if (fieldName == IOobject::groupName("omega", phaseName_))
197  {
198  eqn += interface_.fraction()*beta_*taRhoSqrnu
199  /(sqr(betaStar_)*pow4(delta_));
200  }
201  else
202  {
204  << "Support for field " << fieldName << " is not implemented"
205  << exit(FatalError);
206  }
207 }
208 
209 
210 void Foam::fv::VoFTurbulenceDamping::topoChange(const polyTopoChangeMap&)
211 {}
212 
213 
214 void Foam::fv::VoFTurbulenceDamping::mapMesh(const polyMeshMap& map)
215 {}
216 
217 
218 void Foam::fv::VoFTurbulenceDamping::distribute(const polyDistributionMap&)
219 {}
220 
221 
223 {
224  return true;
225 }
226 
227 
228 // ************************************************************************* //
virtual bool movePoints()
Update for mesh motion.
defineTypeNameAndDebug(fixedTemperatureConstraint, 0)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
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.
Definition: typeInfo.H:106
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
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
Definition: volFieldsFwd.H:58
compressibleMomentumTransportModel momentumTransportModel
labelList fv(nPoints)
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.
Definition: fileName.H:54
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
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.
Definition: POSIX.C:488
dimensionedScalar pow4(const dimensionedScalar &ds)
messageStream Info
Namespace for OpenFOAM.
IOerror FatalIOError