interfaceTurbulenceDamping.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) 2022-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 
27 #include "surfaceInterpolate.H"
28 #include "fvcGrad.H"
30 
31 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  namespace fv
36  {
38 
40  (
41  fvModel,
44  );
45  }
46 }
47 
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
52 Foam::fv::interfaceTurbulenceDamping::interfaceFraction
53 (
54  const volScalarField& alpha
55 ) const
56 {
57  const fvMesh& mesh = this->mesh();
58 
59  tmp<volScalarField::Internal> tA
60  (
62  (
63  "A",
64  mesh,
66  )
67  );
68  volScalarField::Internal& A = tA.ref();
69 
70  const surfaceVectorField& Sf = mesh.Sf();
71  const labelUList& own = mesh.owner();
72  const labelUList& nei = mesh.neighbour();
73 
75 
76  const volVectorField gradAlpha(fvc::grad(alpha));
78  (
79  gradAlpha()/(mag(gradAlpha()) + phase_.fluid().deltaN())
80  );
81 
82  const scalarField& ialpha = alpha;
83  const scalarField& ialphaf = alphaf;
84  scalarField sumnSf(mesh.nCells(), 0);
85 
86  forAll(own, facei)
87  {
88  {
89  const scalar nSf(mag(n[own[facei]] & Sf[facei]));
90  A[own[facei]] += nSf*(ialphaf[facei] - ialpha[own[facei]]);
91  sumnSf[own[facei]] += nSf;
92  }
93  {
94  const scalar nSf(mag(n[nei[facei]] & Sf[facei]));
95  A[nei[facei]] += nSf*(ialphaf[facei] - ialpha[nei[facei]]);
96  sumnSf[nei[facei]] += nSf;
97  }
98  }
99 
101  {
102  const labelUList& own = mesh.boundary()[patchi].faceCells();
103  const fvsPatchScalarField& palphaf = alphaf.boundaryField()[patchi];
104 
105  forAll(mesh.boundary()[patchi], facei)
106  {
107  const scalar nSf(mag(n[own[facei]] & Sf[facei]));
108  A[own[facei]] += nSf*(palphaf[facei] - ialpha[own[facei]]);
109  sumnSf[own[facei]] += nSf;
110  }
111  }
112 
113  scalarField& a = A.primitiveFieldRef();
114  forAll(a, i)
115  {
116  if (sumnSf[i] > small)
117  {
118  a[i] = 2*mag(a[i])/sumnSf[i];
119  }
120  else
121  {
122  a[i] = 0;
123  }
124  }
125 
126  return tA;
127 }
128 
129 
130 template<class RhoType>
131 void Foam::fv::interfaceTurbulenceDamping::addRhoSup
132 (
133  const RhoType& rho,
134  const volScalarField& field,
135  fvMatrix<scalar>& eqn
136 ) const
137 {
138  if (debug)
139  {
140  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
141  }
142 
143  const phaseSystem::phaseModelPartialList& movingPhases =
144  phase_.fluid().movingPhases();
145 
147  (
148  movingPhases[0]*sqr(movingPhases[0].fluidThermo().nu()()())
149  );
150 
151  for (label phasei=1; phasei<movingPhases.size(); phasei++)
152  {
153  aSqrnu +=
154  movingPhases[phasei]
155  *sqr(movingPhases[phasei].fluidThermo().nu()()());
156  }
157 
158  if (field.name() == "epsilon")
159  {
160  eqn += rho*interfaceFraction(phase_)*C2_*aSqrnu*turbulence_.k()()
161  /pow4(delta_);
162  }
163  else if (field.name() == "omega")
164  {
165  eqn += rho*interfaceFraction(phase_)*beta_*aSqrnu
166  /(sqr(betaStar_)*pow4(delta_));
167  }
168  else
169  {
171  << "Support for field " << field.name() << " is not implemented"
172  << exit(FatalError);
173  }
174 }
175 
176 
177 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
178 
180 (
181  const word& sourceName,
182  const word& modelType,
183  const fvMesh& mesh,
184  const dictionary& dict
185 )
186 :
187  fvModel(sourceName, modelType, mesh, dict),
188  phaseName_(dict.lookup("phase")),
189  delta_("delta", dimLength, dict),
190  phase_
191  (
192  mesh.lookupObject<phaseModel>(IOobject::groupName("alpha", phaseName_))
193  ),
194  turbulence_
195  (
196  mesh.lookupType<phaseCompressible::momentumTransportModel>(phaseName_)
197  ),
198  C2_("C2", dimless, 0),
199  betaStar_("betaStar", dimless, 0),
200  beta_("beta", dimless, 0)
201 {
202  const word epsilonName(IOobject::groupName("epsilon", phaseName_));
203  const word omegaName(IOobject::groupName("omega", phaseName_));
204 
205  if (mesh.foundObject<volScalarField>(epsilonName))
206  {
207  fieldName_ = epsilonName;
208  C2_.read(turbulence_.coeffDict());
209  }
210  else if (mesh.foundObject<volScalarField>(omegaName))
211  {
212  fieldName_ = omegaName;
213  betaStar_.read(turbulence_.coeffDict());
214 
215  // Read beta for k-omega models or beta1 for k-omega SST
216  if (turbulence_.coeffDict().found("beta"))
217  {
218  beta_.read(turbulence_.coeffDict());
219  }
220  else
221  {
222  beta_ =
223  dimensionedScalar("beta1", dimless, turbulence_.coeffDict());
224  }
225  }
226  else
227  {
229  << "Cannot find either " << epsilonName << " or " << omegaName
230  << " field for fvModel " << typeName << exit(FatalIOError);
231  }
232 }
233 
234 
235 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
236 
238 {
239  return wordList(1, fieldName_);
240 }
241 
242 
244 (
245  const volScalarField& field,
246  fvMatrix<scalar>& eqn
247 ) const
248 {
249  addRhoSup(one(), field, eqn);
250 }
251 
252 
254 (
255  const volScalarField& rho,
256  const volScalarField& field,
257  fvMatrix<scalar>& eqn
258 ) const
259 {
260  addRhoSup(rho(), field, eqn);
261 }
262 
263 
265 (
266  const volScalarField& alpha,
267  const volScalarField& rho,
268  const volScalarField& field,
269  fvMatrix<scalar>& eqn
270 ) const
271 {
272  if (debug)
273  {
274  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
275  }
276 
277  const volScalarField::Internal aSqrnu
278  (
279  alpha*sqr(phase_.fluidThermo().nu()()())
280  );
281 
282  if (field.name() == IOobject::groupName("epsilon", phaseName_))
283  {
284  eqn += rho()*interfaceFraction(alpha)
285  *C2_*aSqrnu*turbulence_.k()()/pow4(delta_);
286  }
287  else if (field.name() == IOobject::groupName("omega", phaseName_))
288  {
289  eqn += rho()*interfaceFraction(alpha)
290  *beta_*aSqrnu/(sqr(betaStar_)*pow4(delta_));
291  }
292  else
293  {
295  << "Support for field " << field.name() << " is not implemented"
296  << exit(FatalError);
297  }
298 }
299 
300 
302 {}
303 
304 
306 {}
307 
308 
310 (
311  const polyDistributionMap&
312 )
313 {}
314 
315 
317 {
318  return true;
319 }
320 
321 
322 // ************************************************************************* //
static const Foam::dimensionedScalar A("A", Foam::dimPressure, 611.21)
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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...
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &, const Field< Type > &)
Return a temporary field constructed from name, mesh,.
Generic GeometricField class.
DimensionedField< Type, GeoMesh > Internal
Type of the internal field from which this GeometricField is derived.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const word & name() const
Return name.
Definition: IOobject.H:310
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:162
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:509
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:99
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:469
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:857
const surfaceVectorField & Sf() const
Return cell face area vectors.
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:475
Finite volume model abstract base class.
Definition: fvModel.H:59
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:53
Free-surface phase 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.
virtual void addSup(const volScalarField &field, fvMatrix< scalar > &eqn) const
Add source to mixture 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.
interfaceTurbulenceDamping(const word &sourceName, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
Abstract base class for turbulence models (RAS, LES and laminar).
virtual const dictionary & coeffDict() const =0
Const access to the coefficients dictionary.
bool foundObject(const word &name) const
Is the named Type in registry.
A class representing the concept of 1 (scalar(1)) used to avoid unnecessary manipulations for objects...
Definition: one.H:51
const phaseSystem & fluid() const
Return the system to which this phase belongs.
Definition: phaseModel.C:163
const dimensionedScalar & deltaN() const
Stabilisation for normalisation of the interface normal.
Definition: phaseSystemI.H:248
UPtrList< phaseModel > phaseModelPartialList
Definition: phaseSystem.H:86
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.
label nCells() const
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:346
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
Calculate the gradient of the given field.
label patchi
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)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
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
VolField< vector > volVectorField
Definition: volFieldsFwd.H:65
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
const dimensionSet dimless
SurfaceField< scalar > surfaceScalarField
messageStream Info
const dimensionSet dimLength
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
dimensioned< scalar > mag(const dimensioned< Type > &)
dimensionedScalar pow4(const dimensionedScalar &ds)
IOerror FatalIOError
error FatalError
SurfaceField< vector > surfaceVectorField
UList< label > labelUList
Definition: UList.H:65
fvsPatchField< scalar > fvsPatchScalarField
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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