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-2024 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  (
69  (
70  "phaseProperties"
71  )
72  ),
73  momentumTransport_
74  (
75  mesh.lookupType<compressibleMomentumTransportModel>(phaseName_)
76  ),
77  C2_("C2", dimless, dict, 1.92),
78  betaStar_("betaStar", dimless, dict, 0.09),
79  beta_("beta", dimless, dict, 0.072)
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  }
88  else if (mesh.foundObject<volScalarField>(omegaName))
89  {
90  fieldName_ = omegaName;
91  }
92  else
93  {
95  << "Cannot find either " << epsilonName << " or " << omegaName
96  << " field for fvModel " << typeName << exit(FatalIOError);
97  }
98 }
99 
100 
101 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
102 
105 {
106  return wordList(1, fieldName_);
107 }
108 
109 
111 (
112  const volScalarField& rho,
113  const volScalarField& field,
114  fvMatrix<scalar>& eqn
115 ) const
116 {
117  if (debug)
118  {
119  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
120  }
121 
122  const volScalarField::Internal aRhoSqrnu
123  (
124  mixture_.alpha1()()*mixture_.rho1()()*sqr(mixture_.thermo1().nu()()())
125  + mixture_.alpha2()()*mixture_.rho2()()*sqr(mixture_.thermo2().nu()()())
126  );
127 
128  if (field.name() == "epsilon")
129  {
130  eqn += mixture_.interfaceFraction()
131  *C2_*aRhoSqrnu*momentumTransport_.k()()/pow4(delta_);
132  }
133  else if (field.name() == "omega")
134  {
135  eqn += mixture_.interfaceFraction()
136  *beta_*aRhoSqrnu/(sqr(betaStar_)*pow4(delta_));
137  }
138  else
139  {
141  << "Support for field " << field.name() << " is not implemented"
142  << exit(FatalError);
143  }
144 }
145 
146 
148 (
149  const volScalarField& alpha,
150  const volScalarField& rho,
151  const volScalarField& field,
152  fvMatrix<scalar>& eqn
153 ) const
154 {
155  if (debug)
156  {
157  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
158  }
159 
161 
162  if (mixture_.alpha1().name() == alpha.name())
163  {
164  taRhoSqrnu = mixture_.alpha1()()*mixture_.rho1()()
165  *sqr(mixture_.thermo1().nu()()());
166  }
167  else if (mixture_.alpha2().name() == alpha.name())
168  {
169  taRhoSqrnu = mixture_.alpha2()()*mixture_.rho2()()
170  *sqr(mixture_.thermo2().nu()()());
171  }
172  else
173  {
175  << "Unknown phase-fraction " << alpha.name()
176  << exit(FatalError);
177  }
178 
179  if (field.name() == IOobject::groupName("epsilon", phaseName_))
180  {
181  eqn += mixture_.interfaceFraction()
182  *C2_*taRhoSqrnu*momentumTransport_.k()()/pow4(delta_);
183  }
184  else if (field.name() == IOobject::groupName("omega", phaseName_))
185  {
186  eqn += mixture_.interfaceFraction()
187  *beta_*taRhoSqrnu/(sqr(betaStar_)*pow4(delta_));
188  }
189  else
190  {
192  << "Support for field " << field.name() << " is not implemented"
193  << exit(FatalError);
194  }
195 }
196 
197 
199 (
200  const polyTopoChangeMap&
201 )
202 {}
203 
204 
206 (
207  const polyMeshMap& map
208 )
209 {}
210 
211 
213 (
214  const polyDistributionMap&
215 )
216 {}
217 
218 
220 {
221  return true;
222 }
223 
224 
225 // ************************************************************************* //
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.
const word & name() const
Return name.
Definition: IOobject.H:307
static word groupName(Name name, const word &group)
Base class for single-phase compressible turbulence models.
Class to represent a mixture of two rhoFluidThermo-based phases.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
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:96
Finite volume model abstract base class.
Definition: fvModel.H:60
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:69
Free-surface momentumTransport damping function.
virtual bool movePoints()
Update for mesh motion.
virtual void addSup(const volScalarField &rho, const volScalarField &field, fvMatrix< scalar > &eqn) const
Add explicit contribution to epsilon or omega equation.
virtual wordList addSupFields() const
Return the list of fields for which the option adds source term.
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.
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
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
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
List< word > wordList
A List of words.
Definition: fileName.H:54
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
const dimensionSet dimless
messageStream Info
const dimensionSet dimLength
void pow4(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
IOerror FatalIOError
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
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