All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 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  {
37  defineTypeNameAndDebug(interfaceTurbulenceDamping, 0);
38 
40  (
41  fvModel,
42  interfaceTurbulenceDamping,
43  dictionary
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 
74  const surfaceScalarField alphaf(fvc::interpolate(alpha));
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 
100  forAll(mesh.boundary(), patchi)
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 dictionary& dict,
183  const fvMesh& mesh
184 )
185 :
186  fvModel(sourceName, modelType, dict, mesh),
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.lookupObject<phaseCompressible::momentumTransportModel>
196  (
197  IOobject::groupName
198  (
199  momentumTransportModel::typeName,
200  phaseName_
201  )
202  )
203  ),
204  C2_("C2", dimless, 0),
205  betaStar_("betaStar", dimless, 0),
206  beta_("beta", dimless, 0)
207 {
208  const word epsilonName(IOobject::groupName("epsilon", phaseName_));
209  const word omegaName(IOobject::groupName("omega", phaseName_));
210 
211  if (mesh.foundObject<volScalarField>(epsilonName))
212  {
213  fieldName_ = epsilonName;
214  C2_.read(turbulence_.coeffDict());
215  }
216  else if (mesh.foundObject<volScalarField>(omegaName))
217  {
218  fieldName_ = omegaName;
219  betaStar_.read(turbulence_.coeffDict());
220 
221  // Read beta for k-omega models or beta1 for k-omega SST
222  if (turbulence_.coeffDict().found("beta"))
223  {
224  beta_.read(turbulence_.coeffDict());
225  }
226  else
227  {
228  beta_ =
229  dimensionedScalar("beta1", dimless, turbulence_.coeffDict());
230  }
231  }
232  else
233  {
235  << "Cannot find either " << epsilonName << " or " << omegaName
236  << " field for fvModel " << typeName << exit(FatalIOError);
237  }
238 }
239 
240 
241 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
242 
244 {
245  return wordList(1, fieldName_);
246 }
247 
248 
250 (
251  fvMatrix<scalar>& eqn,
252  const word& fieldName
253 ) const
254 {
255  addRhoSup(one(), eqn, fieldName);
256 }
257 
258 
260 (
261  const volScalarField& rho,
262  fvMatrix<scalar>& eqn,
263  const word& fieldName
264 ) const
265 {
266  addRhoSup(rho(), eqn, fieldName);
267 }
268 
269 
271 (
272  const volScalarField& alpha,
273  const volScalarField& rho,
274  fvMatrix<scalar>& eqn,
275  const word& fieldName
276 ) const
277 {
278  if (debug)
279  {
280  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
281  }
282 
284  (
285  alpha*sqr(phase_.thermo().nu()()())
286  );
287 
288  if (fieldName == IOobject::groupName("epsilon", phaseName_))
289  {
290  eqn += rho()*interfaceFraction(alpha)
291  *C2_*aSqrnu*turbulence_.k()()/pow4(delta_);
292  }
293  else if (fieldName == IOobject::groupName("omega", phaseName_))
294  {
295  eqn += rho()*interfaceFraction(alpha)
296  *beta_*aSqrnu/(sqr(betaStar_)*pow4(delta_));
297  }
298  else
299  {
301  << "Support for field " << fieldName << " is not implemented"
302  << exit(FatalError);
303  }
304 }
305 
306 
307 void Foam::fv::interfaceTurbulenceDamping::topoChange(const polyTopoChangeMap&)
308 {}
309 
310 
311 void Foam::fv::interfaceTurbulenceDamping::mapMesh(const polyMeshMap& map)
312 {}
313 
314 
316 (
317  const polyDistributionMap&
318 )
319 {}
320 
321 
323 {
324  return true;
325 }
326 
327 
328 // ************************************************************************* //
interfaceTurbulenceDamping(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
fvsPatchField< scalar > fvsPatchScalarField
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
defineTypeNameAndDebug(fixedTemperatureConstraint, 0)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
fluidReactionThermo & thermo
Definition: createFields.H:28
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &)
Return a temporary field constructed from name, mesh.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
label phasei
Definition: pEqn.H:27
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
error FatalError
virtual wordList addSupFields() const
Return the list of fields for which the option adds source term.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const dimensionSet dimless
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:59
fvMesh & mesh
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
virtual bool movePoints()
Update for mesh motion.
Macros for easy insertion into run-time selection tables.
UList< label > labelUList
Definition: UList.H:65
const dimensionSet dimLength
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
compressibleMomentumTransportModel momentumTransportModel
stressControl lookup("compactNormalStress") >> compactNormalStress
Calculate the gradient of the given field.
labelList fv(nPoints)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
virtual const word & type() const =0
Runtime type information.
UPtrList< phaseModel > phaseModelPartialList
Definition: phaseSystem.H:84
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
addToRunTimeSelectionTable(fvConstraint, fixedTemperatureConstraint, dictionary)
List< word > wordList
A List of words.
Definition: fileName.H:54
virtual void addSup(fvMatrix< scalar > &eqn, const word &fieldName) const
Add explicit contribution to mixture epsilon or omega equation.
label patchi
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
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
dimensionedScalar pow4(const dimensionedScalar &ds)
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Namespace for OpenFOAM.
IOerror FatalIOError