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-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 
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.field();
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  fvMatrix<scalar>& eqn,
135  const word& fieldName
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].thermo().nu()()())
149  );
150 
151  for (label phasei=1; phasei<movingPhases.size(); phasei++)
152  {
153  aSqrnu +=
154  movingPhases[phasei]*sqr(movingPhases[phasei].thermo().nu()()());
155  }
156 
157  if (fieldName == "epsilon")
158  {
159  eqn += rho*interfaceFraction(phase_)*C2_*aSqrnu*turbulence_.k()()
160  /pow4(delta_);
161  }
162  else if (fieldName == "omega")
163  {
164  eqn += rho*interfaceFraction(phase_)*beta_*aSqrnu
165  /(sqr(betaStar_)*pow4(delta_));
166  }
167  else
168  {
170  << "Support for field " << fieldName << " is not implemented"
171  << exit(FatalError);
172  }
173 }
174 
175 
176 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
177 
179 (
180  const word& sourceName,
181  const word& modelType,
182  const fvMesh& mesh,
183  const dictionary& dict
184 )
185 :
186  fvModel(sourceName, modelType, mesh, dict),
187  phaseName_(dict.lookup("phase")),
188  delta_("delta", dimLength, dict),
189  phase_
190  (
191  mesh.lookupObject<phaseModel>(IOobject::groupName("alpha", phaseName_))
192  ),
193  turbulence_
194  (
195  mesh.lookupType<phaseCompressible::momentumTransportModel>(phaseName_)
196  ),
197  C2_("C2", dimless, 0),
198  betaStar_("betaStar", dimless, 0),
199  beta_("beta", dimless, 0)
200 {
201  const word epsilonName(IOobject::groupName("epsilon", phaseName_));
202  const word omegaName(IOobject::groupName("omega", phaseName_));
203 
204  if (mesh.foundObject<volScalarField>(epsilonName))
205  {
206  fieldName_ = epsilonName;
207  C2_.read(turbulence_.coeffDict());
208  }
209  else if (mesh.foundObject<volScalarField>(omegaName))
210  {
211  fieldName_ = omegaName;
212  betaStar_.read(turbulence_.coeffDict());
213 
214  // Read beta for k-omega models or beta1 for k-omega SST
215  if (turbulence_.coeffDict().found("beta"))
216  {
217  beta_.read(turbulence_.coeffDict());
218  }
219  else
220  {
221  beta_ =
222  dimensionedScalar("beta1", dimless, turbulence_.coeffDict());
223  }
224  }
225  else
226  {
228  << "Cannot find either " << epsilonName << " or " << omegaName
229  << " field for fvModel " << typeName << exit(FatalIOError);
230  }
231 }
232 
233 
234 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
235 
237 {
238  return wordList(1, fieldName_);
239 }
240 
241 
243 (
244  fvMatrix<scalar>& eqn,
245  const word& fieldName
246 ) const
247 {
248  addRhoSup(one(), eqn, fieldName);
249 }
250 
251 
253 (
254  const volScalarField& rho,
255  fvMatrix<scalar>& eqn,
256  const word& fieldName
257 ) const
258 {
259  addRhoSup(rho(), eqn, fieldName);
260 }
261 
262 
264 (
265  const volScalarField& alpha,
266  const volScalarField& rho,
267  fvMatrix<scalar>& eqn,
268  const word& fieldName
269 ) const
270 {
271  if (debug)
272  {
273  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
274  }
275 
277  (
278  alpha*sqr(phase_.thermo().nu()()())
279  );
280 
281  if (fieldName == IOobject::groupName("epsilon", phaseName_))
282  {
283  eqn += rho()*interfaceFraction(alpha)
284  *C2_*aSqrnu*turbulence_.k()()/pow4(delta_);
285  }
286  else if (fieldName == IOobject::groupName("omega", phaseName_))
287  {
288  eqn += rho()*interfaceFraction(alpha)
289  *beta_*aSqrnu/(sqr(betaStar_)*pow4(delta_));
290  }
291  else
292  {
294  << "Support for field " << fieldName << " is not implemented"
295  << exit(FatalError);
296  }
297 }
298 
299 
301 {}
302 
303 
305 {}
306 
307 
309 (
310  const polyDistributionMap&
311 )
312 {}
313 
314 
316 {
317  return true;
318 }
319 
320 
321 // ************************************************************************* //
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
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
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:451
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:893
const surfaceVectorField & Sf() const
Return cell face area vectors.
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:457
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 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(fvMatrix< scalar > &eqn, const word &fieldName) const
Add explicit contribution 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:127
const dimensionedScalar & deltaN() const
Stabilisation for normalisation of the interface normal.
Definition: phaseSystemI.H:187
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:318
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
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
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
List< word > wordList
A List of words.
Definition: fileName.H:54
VolField< vector > volVectorField
Definition: volFieldsFwd.H:62
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:251
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:61
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
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
fluidMulticomponentThermo & thermo
Definition: createFields.H:31