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  {
39 
41  (
42  fvModel,
45  );
46  }
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
53 (
54  const word& sourceName,
55  const word& modelType,
56  const fvMesh& mesh,
57  const dictionary& dict
58 )
59 :
60  fvModel(sourceName, modelType, mesh, dict),
61  phaseName_(dict.lookupOrDefault("phase", word::null)),
62  delta_("delta", dimLength, dict),
63  mixture_
64  (
66  (
67  "phaseProperties"
68  )
69  ),
70  turbulence_
71  (
72  mesh.lookupType<incompressibleMomentumTransportModel>(phaseName_)
73  ),
74  C2_("C2", dimless, dict, 1.92),
75  betaStar_("betaStar", dimless, dict, 0.09),
76  beta_("beta", dimless, dict, 0.072)
77 {
78  const word epsilonName(IOobject::groupName("epsilon", phaseName_));
79  const word omegaName(IOobject::groupName("omega", phaseName_));
80 
81  if (mesh.foundObject<volScalarField>(epsilonName))
82  {
83  fieldName_ = epsilonName;
84  }
85  else if (mesh.foundObject<volScalarField>(omegaName))
86  {
87  fieldName_ = omegaName;
88  }
89  else
90  {
92  << "Cannot find either " << epsilonName << " or " << omegaName
93  << " field for fvModel " << typeName << exit(FatalIOError);
94  }
95 }
96 
97 
98 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
99 
101 {
102  return wordList(1, fieldName_);
103 }
104 
105 
107 (
108  const volScalarField& field,
109  fvMatrix<scalar>& eqn
110 ) const
111 {
112  if (debug)
113  {
114  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
115  }
116 
117  const volScalarField::Internal aSqrnu
118  (
119  mixture_.alpha1()()*sqr(mixture_.nuModel1().nu()()())
120  + mixture_.alpha2()()*sqr(mixture_.nuModel2().nu()()())
121  );
122 
123  if (field.name() == "epsilon")
124  {
125  eqn += mixture_.interfaceFraction()*C2_*aSqrnu*turbulence_.k()()
126  /pow4(delta_);
127  }
128  else if (field.name() == "omega")
129  {
130  eqn += mixture_.interfaceFraction()*beta_*aSqrnu/(sqr(betaStar_)
131  *pow4(delta_));
132  }
133  else
134  {
136  << "Support for field " << field.name() << " is not implemented"
137  << exit(FatalError);
138  }
139 }
140 
141 
143 (
144  const volScalarField& alpha,
145  const volScalarField& rho,
146  const volScalarField& field,
147  fvMatrix<scalar>& eqn
148 ) const
149 {
150  if (debug)
151  {
152  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
153  }
154 
156 
157  if (mixture_.alpha1().name() == alpha.name())
158  {
159  taSqrnu = mixture_.alpha1()()*sqr(mixture_.nuModel1().nu()()());
160  }
161  else if (mixture_.alpha2().name() == alpha.name())
162  {
163  taSqrnu = mixture_.alpha2()()*sqr(mixture_.nuModel2().nu()()());
164  }
165  else
166  {
168  << "Unknown phase-fraction " << alpha.name()
169  << exit(FatalError);
170  }
171 
172  if (field.name() == IOobject::groupName("epsilon", phaseName_))
173  {
174  eqn += mixture_.interfaceFraction()*C2_*taSqrnu*turbulence_.k()()
175  /pow4(delta_);
176  }
177  else if (field.name() == IOobject::groupName("omega", phaseName_))
178  {
179  eqn += mixture_.interfaceFraction()*beta_*taSqrnu
180  /(sqr(betaStar_)*pow4(delta_));
181  }
182  else
183  {
185  << "Support for field " << field.name() << " is not implemented"
186  << exit(FatalError);
187  }
188 }
189 
190 
192 {}
193 
194 
196 {}
197 
198 
200 {}
201 
202 
204 {
205  return true;
206 }
207 
208 
209 // ************************************************************************* //
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)
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 turbulence 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.
VoFTurbulenceDamping(const word &sourceName, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
virtual void addSup(const volScalarField &field, fvMatrix< scalar > &eqn) const
Add explicit contribution to epsilon or omega equation.
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.
Base class for single-phase incompressible turbulence models.
Class to represent a mixture of two constant density phases.
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))
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
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