Loading...
Searching...
No Matches
kEpsilonLopesdaCosta.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2018 OpenFOAM Foundation
9 Copyright (C) 2018-2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
30#include "fvOptions.H"
32#include "bound.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
38namespace RASModels
39{
40
41// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
42
43template<class BasicTurbulenceModel>
45(
48)
49{
50 if (pm.dict().found(C.name()))
51 {
52 const labelList& cellZoneIDs = pm.cellZoneIDs();
53
54 const scalar Cpm = pm.dict().get<scalar>(C.name());
55
56 for (const label zonei : cellZoneIDs)
57 {
58 const labelList& cells = this->mesh_.cellZones()[zonei];
59
60 for (const label celli : cells)
61 {
62 C[celli] = Cpm;
63 }
64 }
65 }
66}
67
68
69template<class BasicTurbulenceModel>
71(
73 const porosityModels::powerLawLopesdaCosta& pm
74)
75{
76 if (pm.dict().found(C.name()))
77 {
78 const labelList& cellZoneIDs = pm.cellZoneIDs();
79 const scalarField& Sigma = pm.Sigma();
80
81 const scalar Cpm = pm.dict().get<scalar>(C.name());
82
83 for (const label zonei : cellZoneIDs)
84 {
85 const labelList& cells = this->mesh_.cellZones()[zonei];
86
87 forAll(cells, i)
88 {
89 const label celli = cells[i];
90 C[celli] = Cpm*Sigma[i];
91 }
92 }
93 }
94}
95
96
97template<class BasicTurbulenceModel>
99{
101
102 forAll(fvOptions, i)
103 {
105 {
106 const fv::explicitPorositySource& eps =
108
110 {
111 const porosityModels::powerLawLopesdaCosta& pm =
113 (
114 eps.model()
115 );
116
117 setPorosityCoefficient(Cmu_, pm);
118 Cmu_.correctBoundaryConditions();
119 setPorosityCoefficient(C1_, pm);
120 setPorosityCoefficient(C2_, pm);
121 setPorosityCoefficient(sigmak_, pm);
122 setPorosityCoefficient(sigmaEps_, pm);
123 sigmak_.correctBoundaryConditions();
124 sigmaEps_.correctBoundaryConditions();
125
126 setCdSigma(CdSigma_, pm);
127 setPorosityCoefficient(betap_, pm);
128 setPorosityCoefficient(betad_, pm);
129 setPorosityCoefficient(C4_, pm);
130 setPorosityCoefficient(C5_, pm);
132 }
133 }
134}
135
136
137template<class BasicTurbulenceModel>
139{
140 this->nut_ = Cmu_*sqr(k_)/epsilon_;
141 this->nut_.correctBoundaryConditions();
142 fv::options::New(this->mesh_).correct(this->nut_);
143
144 BasicTurbulenceModel::correctNut();
145}
146
147
148template<class BasicTurbulenceModel>
150(
151 const volScalarField::Internal& magU,
152 const volScalarField::Internal& magU3
153) const
155 return fvm::Su(CdSigma_*(betap_*magU3 - betad_*magU*k_()), k_);
156}
157
158
159template<class BasicTurbulenceModel>
162(
163 const volScalarField::Internal& magU,
164 const volScalarField::Internal& magU3
165) const
166{
167 return fvm::Su
168 (
169 CdSigma_
170 *(C4_*betap_*epsilon_()/k_()*magU3 - C5_*betad_*magU*epsilon_()),
171 epsilon_
172 );
173}
174
175
176// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
177
178template<class BasicTurbulenceModel>
179kEpsilonLopesdaCosta<BasicTurbulenceModel>::kEpsilonLopesdaCosta
180(
181 const alphaField& alpha,
182 const rhoField& rho,
183 const volVectorField& U,
184 const surfaceScalarField& alphaRhoPhi,
185 const surfaceScalarField& phi,
186 const transportModel& transport,
187 const word& propertiesName,
188 const word& type
189)
190:
191 eddyViscosity<RASModel<BasicTurbulenceModel>>
192 (
193 type,
194 alpha,
195 rho,
196 U,
197 alphaRhoPhi,
198 phi,
199 transport,
200 propertiesName
201 ),
202
203 Cmu_
204 (
206 (
207 "Cmu",
208 this->runTime_.timeName(),
209 this->mesh_
210 ),
211 this->mesh_,
212 dimensioned<scalar>::getOrAddToDict
213 (
214 "Cmu",
215 this->coeffDict_,
216 0.09
217 )
218 ),
219 C1_
220 (
222 (
223 "C1",
224 this->runTime_.timeName(),
225 this->mesh_
226 ),
227 this->mesh_,
228 dimensioned<scalar>::getOrAddToDict
229 (
230 "C1",
231 this->coeffDict_,
232 1.44
233 )
234 ),
235 C2_
236 (
238 (
239 "C2",
240 this->runTime_.timeName(),
241 this->mesh_
242 ),
243 this->mesh_,
244 dimensioned<scalar>::getOrAddToDict
245 (
246 "C2",
247 this->coeffDict_,
248 1.92
249 )
250 ),
251 sigmak_
252 (
254 (
255 "sigmak",
256 this->runTime_.timeName(),
257 this->mesh_
258 ),
259 this->mesh_,
260 dimensioned<scalar>::getOrAddToDict
261 (
262 "sigmak",
263 this->coeffDict_,
264 1.0
265 )
266 ),
267 sigmaEps_
268 (
270 (
271 "sigmaEps",
272 this->runTime_.timeName(),
273 this->mesh_
274 ),
275 this->mesh_,
276 dimensioned<scalar>::getOrAddToDict
277 (
278 "sigmaEps",
279 this->coeffDict_,
280 1.3
281 )
282 ),
283
284 CdSigma_
285 (
287 (
288 "CdSigma",
289 this->runTime_.timeName(),
290 this->mesh_
291 ),
292 this->mesh_,
293 dimensionedScalar("CdSigma", dimless/dimLength, 0)
294 ),
295 betap_
296 (
298 (
299 "betap",
300 this->runTime_.timeName(),
301 this->mesh_
302 ),
303 this->mesh_,
304 dimensionedScalar("betap", dimless, 0)
305 ),
306 betad_
307 (
309 (
310 "betad",
311 this->runTime_.timeName(),
312 this->mesh_
313 ),
314 this->mesh_,
315 dimensionedScalar("betad", dimless, 0)
316 ),
317 C4_
318 (
320 (
321 "C4",
322 this->runTime_.timeName(),
323 this->mesh_
324 ),
325 this->mesh_,
326 dimensionedScalar("C4", dimless, 0)
327 ),
328 C5_
329 (
331 (
332 "C5",
333 this->runTime_.timeName(),
334 this->mesh_
335 ),
336 this->mesh_,
337 dimensionedScalar("C5", dimless, 0)
338 ),
339
340 k_
341 (
343 (
344 IOobject::groupName("k", alphaRhoPhi.group()),
345 this->runTime_.timeName(),
346 this->mesh_,
347 IOobject::MUST_READ,
348 IOobject::AUTO_WRITE
349 ),
350 this->mesh_
351 ),
352 epsilon_
353 (
355 (
356 IOobject::groupName("epsilon", alphaRhoPhi.group()),
357 this->runTime_.timeName(),
358 this->mesh_,
359 IOobject::MUST_READ,
360 IOobject::AUTO_WRITE
361 ),
362 this->mesh_
363 )
364{
365 bound(k_, this->kMin_);
366 bound(epsilon_, this->epsilonMin_);
367
368 if (type == typeName)
369 {
370 this->printCoeffs(type);
371 }
372
374}
375
376
377// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
378
379template<class BasicTurbulenceModel>
381{
383 {
384 return true;
386
387 return false;
388}
389
390
391template<class BasicTurbulenceModel>
393{
394 if (!this->turbulence_)
395 {
396 return;
397 }
398
399 // Local references
400 const alphaField& alpha = this->alpha_;
401 const rhoField& rho = this->rho_;
402 const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
403 const volVectorField& U = this->U_;
404 volScalarField& nut = this->nut_;
405 fv::options& fvOptions(fv::options::New(this->mesh_));
406
407 eddyViscosity<RASModel<BasicTurbulenceModel>>::correct();
408
410 (
411 fvc::div(fvc::absolute(this->phi(), U))().v()
412 );
413
414 tmp<volTensorField> tgradU = fvc::grad(U);
416 (
417 this->GName(),
418 nut.v()*(devTwoSymm(tgradU().v()) && tgradU().v())
419 );
420 tgradU.clear();
421
422 // Update epsilon and G at the wall
423 epsilon_.boundaryFieldRef().updateCoeffs();
424 // Push any changed cell values to coupled neighbours
425 epsilon_.boundaryFieldRef().template evaluateCoupled<coupledFvPatch>();
426
428 volScalarField::Internal magU3(pow3(magU));
429
430 // Dissipation equation
431 tmp<fvScalarMatrix> epsEqn
432 (
433 fvm::ddt(alpha, rho, epsilon_)
434 + fvm::div(alphaRhoPhi, epsilon_)
435 - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
436 ==
437 C1_*alpha()*rho()*G*epsilon_()/k_()
438 - fvm::SuSp(((2.0/3.0)*C1_)*alpha()*rho()*divU, epsilon_)
439 - fvm::Sp(C2_*alpha()*rho()*epsilon_()/k_(), epsilon_)
440 + epsilonSource(magU, magU3)
441 + fvOptions(alpha, rho, epsilon_)
442 );
443
444 epsEqn.ref().relax();
445 fvOptions.constrain(epsEqn.ref());
446 epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
447 solve(epsEqn);
448 fvOptions.correct(epsilon_);
449 bound(epsilon_, this->epsilonMin_);
450
451 // Turbulent kinetic energy equation
452 tmp<fvScalarMatrix> kEqn
453 (
454 fvm::ddt(alpha, rho, k_)
455 + fvm::div(alphaRhoPhi, k_)
456 - fvm::laplacian(alpha*rho*DkEff(), k_)
457 ==
458 alpha()*rho()*G
459 - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
460 - fvm::Sp(alpha()*rho()*epsilon_()/k_(), k_)
461 + kSource(magU, magU3)
462 + fvOptions(alpha, rho, k_)
463 );
464
465 kEqn.ref().relax();
466 fvOptions.constrain(kEqn.ref());
467 solve(kEqn);
468 fvOptions.correct(k_);
469 bound(k_, this->kMin_);
470
471 correctNut();
472}
473
474
475// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
476
477} // End namespace RASModels
478} // End namespace Foam
479
480// ************************************************************************* //
static const Foam::dimensionedScalar C("", Foam::dimTemperature, 234.5)
Bound the given scalar field if it has gone unbounded.
fv::options & fvOptions
Graphite solid properties.
Definition C.H:49
DimensionedField< scalar, volMesh > Internal
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
Templated abstract base class for RAS turbulence models.
Definition RASModel.H:81
BasicTurbulenceModel::alphaField alphaField
BasicTurbulenceModel::rhoField rhoField
virtual tmp< fvScalarMatrix > kSource(const volScalarField::Internal &magU, const volScalarField::Internal &magU3) const
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
void setCdSigma(volScalarField::Internal &C, const porosityModels::powerLawLopesdaCosta &pm)
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
BasicTurbulenceModel::transportModel transportModel
virtual tmp< fvScalarMatrix > epsilonSource(const volScalarField::Internal &magU, const volScalarField::Internal &magU3) const
void setPorosityCoefficient(volScalarField::Internal &C, const porosityModels::powerLawLopesdaCosta &pm)
virtual bool read()
Re-read model coefficients if they have changed.
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
bool get(const label i) const
Return bool value at specified position, always false for out-of-range access.
Definition UList.H:868
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
Generic dimensioned Type class.
const word & name() const noexcept
Return const reference to name.
Eddy viscosity turbulence model base class.
eddyViscosity(const word &modelName, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName)
Applies an explicit porosity source to the momentum equation within a specified region.
const porosityModel & model() const noexcept
Access to the porosityModel.
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
optionList(const optionList &)=delete
No copy construct.
Finite-volume options, which is an IOdictionary of values and a fv::optionList.
Definition fvOptions.H:69
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present otherwise lookup and return.
Definition fvOptions.C:116
const labelList & cellZoneIDs() const
Return const access to the cell zone IDs.
const dictionary & dict() const
Return dictionary used for model construction.
const scalarField & Sigma() const
Return the porosity surface area per unit volume zone field.
Variant of the power law porosity model with spatially varying drag coefficient.
A class for managing temporary objects.
Definition tmp.H:75
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
Definition tmpI.H:289
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
Base-class for all transport models used by the incompressible turbulence models.
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
thermo correct()
scalar nut
const cellShapeList & cells
zeroField divU
Definition alphaSuSp.H:3
word timeName
Definition getTimeIndex.H:3
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcDiv.C:42
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition fvcMeshPhi.C:183
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition fvmDiv.C:41
zeroField Su(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition fvmDdt.C:41
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition int32.H:127
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition bound.C:29
GeometricField< scalar, fvPatchField, volMesh > volScalarField
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
volScalarField & alpha
CEqn solve()
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299