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  {
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  (
65  mesh.lookupObject<incompressibleTwoPhaseVoFMixture>
66  (
67  "phaseProperties"
68  )
69  ),
70  turbulence_
71  (
72  mesh.lookupType<incompressibleMomentumTransportModel>(phaseName_)
73  ),
74  C2_("C2", dimless, 0),
75  betaStar_("betaStar", dimless, 0),
76  beta_("beta", dimless, 0)
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  C2_.read(turbulence_.coeffDict());
85  }
86  else if (mesh.foundObject<volScalarField>(omegaName))
87  {
88  fieldName_ = omegaName;
89  betaStar_.read(turbulence_.coeffDict());
90 
91  // Read beta for k-omega models or beta1 for k-omega SST
92  if (turbulence_.coeffDict().found("beta"))
93  {
94  beta_.read(turbulence_.coeffDict());
95  }
96  else
97  {
98  beta_ =
99  dimensionedScalar("beta1", dimless, turbulence_.coeffDict());
100  }
101  }
102  else
103  {
105  << "Cannot find either " << epsilonName << " or " << omegaName
106  << " field for fvModel " << typeName << exit(FatalIOError);
107  }
108 }
109 
110 
111 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
112 
114 {
115  return wordList(1, fieldName_);
116 }
117 
118 
120 (
121  fvMatrix<scalar>& eqn,
122  const word& fieldName
123 ) const
124 {
125  if (debug)
126  {
127  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
128  }
129 
130  const volScalarField::Internal aSqrnu
131  (
132  mixture_.alpha1()()*sqr(mixture_.nuModel1().nu()()())
133  + mixture_.alpha2()()*sqr(mixture_.nuModel2().nu()()())
134  );
135 
136  if (fieldName == "epsilon")
137  {
138  eqn += mixture_.interfaceFraction()*C2_*aSqrnu*turbulence_.k()()
139  /pow4(delta_);
140  }
141  else if (fieldName == "omega")
142  {
143  eqn += mixture_.interfaceFraction()*beta_*aSqrnu/(sqr(betaStar_)
144  *pow4(delta_));
145  }
146  else
147  {
149  << "Support for field " << fieldName << " is not implemented"
150  << exit(FatalError);
151  }
152 }
153 
154 
156 (
157  const volScalarField& alpha,
158  const volScalarField&,
159  fvMatrix<scalar>& eqn,
160  const word& fieldName
161 ) const
162 {
163  if (debug)
164  {
165  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
166  }
167 
169 
170  if (mixture_.alpha1().name() == alpha.name())
171  {
172  taSqrnu = mixture_.alpha1()()*sqr(mixture_.nuModel1().nu()()());
173  }
174  else if (mixture_.alpha2().name() == alpha.name())
175  {
176  taSqrnu = mixture_.alpha2()()*sqr(mixture_.nuModel2().nu()()());
177  }
178  else
179  {
181  << "Unknown phase-fraction " << alpha.name()
182  << exit(FatalError);
183  }
184 
185  if (fieldName == IOobject::groupName("epsilon", phaseName_))
186  {
187  eqn += mixture_.interfaceFraction()*C2_*taSqrnu*turbulence_.k()()
188  /pow4(delta_);
189  }
190  else if (fieldName == IOobject::groupName("omega", phaseName_))
191  {
192  eqn += mixture_.interfaceFraction()*beta_*taSqrnu
193  /(sqr(betaStar_)*pow4(delta_));
194  }
195  else
196  {
198  << "Support for field " << fieldName << " is not implemented"
199  << exit(FatalError);
200  }
201 }
202 
203 
205 {}
206 
207 
209 {}
210 
211 
213 {}
214 
215 
217 {
218  return true;
219 }
220 
221 
222 // ************************************************************************* //
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)
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 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(fvMatrix< scalar > &eqn, const word &fieldName) 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.
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))
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
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