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-2023 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"
29 #include "fvMatrix.H"
31 
32 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  namespace fv
37  {
38  namespace compressible
39  {
41 
43  (
44  fvModel,
47  );
48  }
49  }
50 }
51 
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
56 (
57  const word& sourceName,
58  const word& modelType,
59  const fvMesh& mesh,
60  const dictionary& dict
61 )
62 :
63  fvModel(sourceName, modelType, mesh, dict),
64  phaseName_(dict.lookupOrDefault("phase", word::null)),
65  delta_("delta", dimLength, dict),
66  mixture_
67  (
68  mesh.lookupObject<compressibleTwoPhaseVoFMixture>
69  (
70  "phaseProperties"
71  )
72  ),
73  momentumTransport_
74  (
75  mesh.lookupType<compressibleMomentumTransportModel>(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(momentumTransport_.coeffDict());
88  }
89  else if (mesh.foundObject<volScalarField>(omegaName))
90  {
91  fieldName_ = omegaName;
92  betaStar_.read(momentumTransport_.coeffDict());
93 
94  // Read beta for k-omega models or beta1 for k-omega SST
95  if (momentumTransport_.coeffDict().found("beta"))
96  {
97  beta_.read(momentumTransport_.coeffDict());
98  }
99  else
100  {
101  beta_ = dimensionedScalar
102  (
103  "beta1",
104  dimless,
105  momentumTransport_.coeffDict()
106  );
107  }
108  }
109  else
110  {
112  << "Cannot find either " << epsilonName << " or " << omegaName
113  << " field for fvModel " << typeName << exit(FatalIOError);
114  }
115 }
116 
117 
118 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
119 
122 {
123  return wordList(1, fieldName_);
124 }
125 
126 
128 (
129  const volScalarField& rho,
130  fvMatrix<scalar>& eqn,
131  const word& fieldName
132 ) const
133 {
134  if (debug)
135  {
136  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
137  }
138 
139  const volScalarField::Internal aRhoSqrnu
140  (
141  mixture_.alpha1()()*mixture_.rho1()()*sqr(mixture_.thermo1().nu()()())
142  + mixture_.alpha2()()*mixture_.rho2()()*sqr(mixture_.thermo2().nu()()())
143  );
144 
145  if (fieldName == "epsilon")
146  {
147  eqn += mixture_.interfaceFraction()
148  *C2_*aRhoSqrnu*momentumTransport_.k()()/pow4(delta_);
149  }
150  else if (fieldName == "omega")
151  {
152  eqn += mixture_.interfaceFraction()
153  *beta_*aRhoSqrnu/(sqr(betaStar_)*pow4(delta_));
154  }
155  else
156  {
158  << "Support for field " << fieldName << " is not implemented"
159  << exit(FatalError);
160  }
161 }
162 
163 
165 (
166  const volScalarField& alpha,
167  const volScalarField&,
168  fvMatrix<scalar>& eqn,
169  const word& fieldName
170 ) const
171 {
172  if (debug)
173  {
174  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
175  }
176 
178 
179  if (mixture_.alpha1().name() == alpha.name())
180  {
181  taRhoSqrnu = mixture_.alpha1()()*mixture_.rho1()()
182  *sqr(mixture_.thermo1().nu()()());
183  }
184  else if (mixture_.alpha2().name() == alpha.name())
185  {
186  taRhoSqrnu = mixture_.alpha2()()*mixture_.rho2()()
187  *sqr(mixture_.thermo2().nu()()());
188  }
189  else
190  {
192  << "Unknown phase-fraction " << alpha.name()
193  << exit(FatalError);
194  }
195 
196  if (fieldName == IOobject::groupName("epsilon", phaseName_))
197  {
198  eqn += mixture_.interfaceFraction()
199  *C2_*taRhoSqrnu*momentumTransport_.k()()/pow4(delta_);
200  }
201  else if (fieldName == IOobject::groupName("omega", phaseName_))
202  {
203  eqn += mixture_.interfaceFraction()
204  *beta_*taRhoSqrnu/(sqr(betaStar_)*pow4(delta_));
205  }
206  else
207  {
209  << "Support for field " << fieldName << " is not implemented"
210  << exit(FatalError);
211  }
212 }
213 
214 
216 (
217  const polyTopoChangeMap&
218 )
219 {}
220 
221 
223 (
224  const polyMeshMap& map
225 )
226 {}
227 
228 
230 (
231  const polyDistributionMap&
232 )
233 {}
234 
235 
237 {
238  return true;
239 }
240 
241 
242 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
static word groupName(Name name, const word &group)
Base class for single-phase compressible turbulence models.
Class to represent a mixture of two rhoThermo-based phases.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:659
void read(const dictionary &)
Update the value of dimensioned<Type>
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
VolField< Type > & psi()
Definition: fvMatrix.H:289
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Finite volume model abstract base class.
Definition: fvModel.H:59
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:34
Free-surface momentumTransport damping function.
virtual bool movePoints()
Update for mesh motion.
virtual wordList addSupFields() const
Return the list of fields for which the option adds source term.
virtual void addSup(const volScalarField &rho, fvMatrix< scalar > &eqn, const word &fieldName) const
Add explicit contribution to epsilon or omega equation.
VoFTurbulenceDamping(const word &sourceName, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
virtual const dictionary & coeffDict() const =0
Const access to the coefficients dictionary.
bool foundObject(const word &name) const
Is the named Type in registry.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
defineTypeNameAndDebug(VoFCavitation, 0)
addToRunTimeSelectionTable(fvModel, VoFCavitation, dictionary)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
List< word > wordList
A List of words.
Definition: fileName.H:54
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const dimensionSet dimless
messageStream Info
const dimensionSet dimLength
dimensionedScalar pow4(const dimensionedScalar &ds)
IOerror FatalIOError
error FatalError
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
labelList fv(nPoints)
dictionary dict