mixtureKEpsilon.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) 2013-2017 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
29 #include "mixtureKEpsilon.H"
30 #include "fvOptions.H"
31 #include "bound.H"
32 #include "twoPhaseSystem.H"
33 #include "dragModel.H"
34 #include "virtualMassModel.H"
37 #include "fvmSup.H"
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43 namespace RASModels
44 {
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
48 template<class BasicTurbulenceModel>
49 mixtureKEpsilon<BasicTurbulenceModel>::mixtureKEpsilon
50 (
51  const alphaField& alpha,
52  const rhoField& rho,
53  const volVectorField& U,
54  const surfaceScalarField& alphaRhoPhi,
55  const surfaceScalarField& phi,
56  const transportModel& transport,
57  const word& propertiesName,
58  const word& type
59 )
60 :
62  (
63  type,
64  alpha,
65  rho,
66  U,
67  alphaRhoPhi,
68  phi,
69  transport,
70  propertiesName
71  ),
72 
73  liquidTurbulencePtr_(nullptr),
74 
75  Cmu_
76  (
78  (
79  "Cmu",
80  this->coeffDict_,
81  0.09
82  )
83  ),
84  C1_
85  (
87  (
88  "C1",
89  this->coeffDict_,
90  1.44
91  )
92  ),
93  C2_
94  (
96  (
97  "C2",
98  this->coeffDict_,
99  1.92
100  )
101  ),
102  C3_
103  (
105  (
106  "C3",
107  this->coeffDict_,
108  C2_.value()
109  )
110  ),
111  Cp_
112  (
114  (
115  "Cp",
116  this->coeffDict_,
117  0.25
118  )
119  ),
120  sigmak_
121  (
123  (
124  "sigmak",
125  this->coeffDict_,
126  1.0
127  )
128  ),
129  sigmaEps_
130  (
132  (
133  "sigmaEps",
134  this->coeffDict_,
135  1.3
136  )
137  ),
138 
139  k_
140  (
141  IOobject
142  (
143  IOobject::groupName("k", alphaRhoPhi.group()),
144  this->runTime_.timeName(),
145  this->mesh_,
148  ),
149  this->mesh_
150  ),
151  epsilon_
152  (
153  IOobject
154  (
155  IOobject::groupName("epsilon", alphaRhoPhi.group()),
156  this->runTime_.timeName(),
157  this->mesh_,
160  ),
161  this->mesh_
162  )
163 {
164  bound(k_, this->kMin_);
165  bound(epsilon_, this->epsilonMin_);
166 
167  if (type == typeName)
168  {
169  this->printCoeffs(type);
170  }
171 }
172 
173 
174 template<class BasicTurbulenceModel>
176 (
177  const volScalarField& epsilon
178 ) const
179 {
180  const volScalarField::Boundary& ebf = epsilon.boundaryField();
181 
182  wordList ebt = ebf.types();
183 
184  forAll(ebf, patchi)
185  {
186  if (isA<fixedValueFvPatchScalarField>(ebf[patchi]))
187  {
188  ebt[patchi] = fixedValueFvPatchScalarField::typeName;
189  }
190  }
191 
192  return ebt;
193 }
194 
195 
196 template<class BasicTurbulenceModel>
198 (
199  volScalarField& vsf,
200  const volScalarField& refVsf
201 ) const
202 {
203  volScalarField::Boundary& bf = vsf.boundaryFieldRef();
204  const volScalarField::Boundary& refBf =
205  refVsf.boundaryField();
206 
207  forAll(bf, patchi)
208  {
209  if
210  (
211  isA<inletOutletFvPatchScalarField>(bf[patchi])
212  && isA<inletOutletFvPatchScalarField>(refBf[patchi])
213  )
214  {
215  refCast<inletOutletFvPatchScalarField>
216  (bf[patchi]).refValue() =
217  refCast<const inletOutletFvPatchScalarField>
218  (refBf[patchi]).refValue();
219  }
220  }
221 }
222 
223 
224 template<class BasicTurbulenceModel>
226 {
227  if (rhom_) return;
228 
229  // Local references to gas-phase properties
230  const volScalarField& kg = this->k_;
231  const volScalarField& epsilong = this->epsilon_;
232 
233  // Local references to liquid-phase properties
234  mixtureKEpsilon<BasicTurbulenceModel>& turbc = this->liquidTurbulence();
235  const volScalarField& kl = turbc.k_;
236  const volScalarField& epsilonl = turbc.epsilon_;
237 
238  word startTimeName
239  (
240  this->runTime_.timeName(this->runTime_.startTime().value())
241  );
242 
243  Ct2_.reset
244  (
245  new volScalarField
246  (
247  IOobject
248  (
249  "Ct2",
250  startTimeName,
251  this->mesh_,
254  ),
255  Ct2()
256  )
257  );
258 
259  rhom_.reset
260  (
261  new volScalarField
262  (
263  IOobject
264  (
265  "rhom",
266  startTimeName,
267  this->mesh_,
270  ),
271  rhom()
272  )
273  );
274 
275  km_.reset
276  (
277  new volScalarField
278  (
279  IOobject
280  (
281  "km",
282  startTimeName,
283  this->mesh_,
286  ),
287  mix(kl, kg),
288  kl.boundaryField().types()
289  )
290  );
291  correctInletOutlet(km_(), kl);
292 
293  epsilonm_.reset
294  (
295  new volScalarField
296  (
297  IOobject
298  (
299  "epsilonm",
300  startTimeName,
301  this->mesh_,
304  ),
305  mix(epsilonl, epsilong),
306  epsilonBoundaryTypes(epsilonl)
307  )
308  );
309  correctInletOutlet(epsilonm_(), epsilonl);
310 }
311 
312 
313 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
314 
315 template<class BasicTurbulenceModel>
317 {
319  {
320  Cmu_.readIfPresent(this->coeffDict());
321  C1_.readIfPresent(this->coeffDict());
322  C2_.readIfPresent(this->coeffDict());
323  C3_.readIfPresent(this->coeffDict());
324  Cp_.readIfPresent(this->coeffDict());
325  sigmak_.readIfPresent(this->coeffDict());
326  sigmaEps_.readIfPresent(this->coeffDict());
327 
328  return true;
329  }
330 
331  return false;
332 }
333 
334 
335 template<class BasicTurbulenceModel>
337 {
338  this->nut_ = Cmu_*sqr(k_)/epsilon_;
339  this->nut_.correctBoundaryConditions();
340  fv::options::New(this->mesh_).correct(this->nut_);
341 
342  BasicTurbulenceModel::correctNut();
343 }
344 
345 
346 template<class BasicTurbulenceModel>
349 {
350  if (!liquidTurbulencePtr_)
351  {
352  const volVectorField& U = this->U_;
353 
354  const transportModel& gas = this->transport();
355  const twoPhaseSystem& fluid =
356  refCast<const twoPhaseSystem>(gas.fluid());
357  const transportModel& liquid = fluid.otherPhase(gas);
358 
359  liquidTurbulencePtr_ =
361  (
362  U.db().lookupObject<mixtureKEpsilon<BasicTurbulenceModel>>
363  (
365  (
367  liquid.name()
368  )
369  )
370  );
371  }
372 
373  return *liquidTurbulencePtr_;
374 }
375 
376 
377 template<class BasicTurbulenceModel>
379 {
380  const mixtureKEpsilon<BasicTurbulenceModel>& liquidTurbulence =
381  this->liquidTurbulence();
382 
383  const transportModel& gas = this->transport();
384  const twoPhaseSystem& fluid = refCast<const twoPhaseSystem>(gas.fluid());
385  const transportModel& liquid = fluid.otherPhase(gas);
386 
387  const volScalarField& alphag = this->alpha_;
388 
389  volScalarField magUr(mag(liquidTurbulence.U() - this->U()));
390 
392  (
393  (6*this->Cmu_/(4*sqrt(3.0/2.0)))
394  *fluid.Kd()/liquid.rho()
395  *(liquidTurbulence.k_/liquidTurbulence.epsilon_)
396  );
397  volScalarField Ct0((3 + beta)/(1 + beta + 2*gas.rho()/liquid.rho()));
398  volScalarField fAlphad((180 + (-4.71e3 + 4.26e4*alphag)*alphag)*alphag);
399 
400  return sqr(1 + (Ct0 - 1)*exp(-fAlphad));
401 }
402 
403 
404 template<class BasicTurbulenceModel>
406 {
407  const transportModel& gas = this->transport();
408  const twoPhaseSystem& fluid = refCast<const twoPhaseSystem>(gas.fluid());
409  return fluid.otherPhase(gas).rho();
410 }
411 
412 
413 template<class BasicTurbulenceModel>
415 {
416  const transportModel& gas = this->transport();
417  const twoPhaseSystem& fluid = refCast<const twoPhaseSystem>(gas.fluid());
418  const virtualMassModel& virtualMass =
419  fluid.lookupSubModel<virtualMassModel>(gas, fluid.otherPhase(gas));
420  return
421  gas.rho()
422  + virtualMass.Cvm()*fluid.otherPhase(gas).rho();
423 }
424 
425 
426 template<class BasicTurbulenceModel>
428 {
429  const volScalarField& alphag = this->alpha_;
430  const volScalarField& alphal = this->liquidTurbulence().alpha_;
431 
432  return alphal*rholEff() + alphag*rhogEff();
433 }
434 
435 
436 template<class BasicTurbulenceModel>
438 (
439  const volScalarField& fc,
440  const volScalarField& fd
441 ) const
442 {
443  const volScalarField& alphag = this->alpha_;
444  const volScalarField& alphal = this->liquidTurbulence().alpha_;
445 
446  return (alphal*rholEff()*fc + alphag*rhogEff()*fd)/rhom_();
447 }
448 
449 
450 template<class BasicTurbulenceModel>
452 (
453  const volScalarField& fc,
454  const volScalarField& fd
455 ) const
456 {
457  const volScalarField& alphag = this->alpha_;
458  const volScalarField& alphal = this->liquidTurbulence().alpha_;
459 
460  return
461  (alphal*rholEff()*fc + alphag*rhogEff()*Ct2_()*fd)
462  /(alphal*rholEff() + alphag*rhogEff()*Ct2_());
463 }
464 
465 
466 template<class BasicTurbulenceModel>
468 (
469  const surfaceScalarField& fc,
470  const surfaceScalarField& fd
471 ) const
472 {
473  const volScalarField& alphag = this->alpha_;
474  const volScalarField& alphal = this->liquidTurbulence().alpha_;
475 
477  surfaceScalarField alphagf(fvc::interpolate(alphag));
478 
479  surfaceScalarField rholEfff(fvc::interpolate(rholEff()));
480  surfaceScalarField rhogEfff(fvc::interpolate(rhogEff()));
481 
482  return
483  (alphalf*rholEfff*fc + alphagf*rhogEfff*fvc::interpolate(Ct2_())*fd)
484  /(alphalf*rholEfff + alphagf*rhogEfff*fvc::interpolate(Ct2_()));
485 }
486 
487 
488 template<class BasicTurbulenceModel>
490 {
491  const mixtureKEpsilon<BasicTurbulenceModel>& liquidTurbulence =
492  this->liquidTurbulence();
493 
494  const transportModel& gas = this->transport();
495  const twoPhaseSystem& fluid = refCast<const twoPhaseSystem>(gas.fluid());
496  const transportModel& liquid = fluid.otherPhase(gas);
497 
498  const dragModel& drag = fluid.lookupSubModel<dragModel>(gas, liquid);
499 
500  volScalarField magUr(mag(liquidTurbulence.U() - this->U()));
501 
502  // Lahey model
503  tmp<volScalarField> bubbleG
504  (
505  Cp_
506  *liquid*liquid.rho()
507  *(
508  pow3(magUr)
509  + pow(drag.CdRe()*liquid.nu()/gas.d(), 4.0/3.0)
510  *pow(magUr, 5.0/3.0)
511  )
512  *gas
513  /gas.d()
514  );
515 
516  // Simple model
517  // tmp<volScalarField> bubbleG
518  // (
519  // Cp_*liquid*drag.K()*sqr(magUr)
520  // );
521 
522  return bubbleG;
523 }
524 
525 
526 template<class BasicTurbulenceModel>
528 {
529  return fvm::Su(bubbleG()/rhom_(), km_());
530 }
531 
532 
533 template<class BasicTurbulenceModel>
535 {
536  return fvm::Su(C3_*epsilonm_()*bubbleG()/(rhom_()*km_()), epsilonm_());
537 }
538 
539 
540 template<class BasicTurbulenceModel>
542 {
543  const transportModel& gas = this->transport();
544  const twoPhaseSystem& fluid = refCast<const twoPhaseSystem>(gas.fluid());
545 
546  // Only solve the mixture turbulence for the gas-phase
547  if (&gas != &fluid.phase1())
548  {
549  // This is the liquid phase but check the model for the gas-phase
550  // is consistent
551  this->liquidTurbulence();
552 
553  return;
554  }
555 
556  if (!this->turbulence_)
557  {
558  return;
559  }
560 
561  // Initialise the mixture fields if they have not yet been constructed
562  initMixtureFields();
563 
564  // Local references to gas-phase properties
565  tmp<surfaceScalarField> phig = this->phi();
566  const volVectorField& Ug = this->U_;
567  const volScalarField& alphag = this->alpha_;
568  volScalarField& kg = this->k_;
569  volScalarField& epsilong = this->epsilon_;
570  volScalarField& nutg = this->nut_;
571 
572  // Local references to liquid-phase properties
573  mixtureKEpsilon<BasicTurbulenceModel>& liquidTurbulence =
574  this->liquidTurbulence();
575  tmp<surfaceScalarField> phil = liquidTurbulence.phi();
576  const volVectorField& Ul = liquidTurbulence.U_;
577  const volScalarField& alphal = liquidTurbulence.alpha_;
578  volScalarField& kl = liquidTurbulence.k_;
579  volScalarField& epsilonl = liquidTurbulence.epsilon_;
580  volScalarField& nutl = liquidTurbulence.nut_;
581 
582  // Local references to mixture properties
583  volScalarField& rhom = rhom_();
584  volScalarField& km = km_();
585  volScalarField& epsilonm = epsilonm_();
586 
587  fv::options& fvOptions(fv::options::New(this->mesh_));
588 
590 
591  // Update the effective mixture density
592  rhom = this->rhom();
593 
594  // Mixture flux
595  surfaceScalarField phim("phim", mixFlux(phil, phig));
596 
597  // Mixture velocity divergence
598  volScalarField divUm
599  (
600  mixU
601  (
602  fvc::div(fvc::absolute(phil, Ul)),
604  )
605  );
606 
608  {
609  tmp<volTensorField> tgradUl = fvc::grad(Ul);
611  (
612  new volScalarField
613  (
614  this->GName(),
615  nutl*(tgradUl() && dev(twoSymm(tgradUl())))
616  )
617  );
618  tgradUl.clear();
619 
620  // Update k, epsilon and G at the wall
621  kl.boundaryFieldRef().updateCoeffs();
622  epsilonl.boundaryFieldRef().updateCoeffs();
623 
624  Gc.ref().checkOut();
625  }
626 
628  {
629  tmp<volTensorField> tgradUg = fvc::grad(Ug);
631  (
632  new volScalarField
633  (
634  this->GName(),
635  nutg*(tgradUg() && dev(twoSymm(tgradUg())))
636  )
637  );
638  tgradUg.clear();
639 
640  // Update k, epsilon and G at the wall
641  kg.boundaryFieldRef().updateCoeffs();
642  epsilong.boundaryFieldRef().updateCoeffs();
643 
644  Gd.ref().checkOut();
645  }
646 
647  // Mixture turbulence generation
648  volScalarField Gm(mix(Gc, Gd));
649 
650  // Mixture turbulence viscosity
651  volScalarField nutm(mixU(nutl, nutg));
652 
653  // Update the mixture k and epsilon boundary conditions
654  km == mix(kl, kg);
655  bound(km, this->kMin_);
656  epsilonm == mix(epsilonl, epsilong);
657  bound(epsilonm, this->epsilonMin_);
658 
659  // Dissipation equation
660  tmp<fvScalarMatrix> epsEqn
661  (
662  fvm::ddt(epsilonm)
663  + fvm::div(phim, epsilonm)
664  - fvm::Sp(fvc::div(phim), epsilonm)
665  - fvm::laplacian(DepsilonEff(nutm), epsilonm)
666  ==
667  C1_*Gm*epsilonm/km
668  - fvm::SuSp(((2.0/3.0)*C1_)*divUm, epsilonm)
669  - fvm::Sp(C2_*epsilonm/km, epsilonm)
670  + epsilonSource()
671  + fvOptions(epsilonm)
672  );
673 
674  epsEqn.ref().relax();
675  fvOptions.constrain(epsEqn.ref());
676  epsEqn.ref().boundaryManipulate(epsilonm.boundaryFieldRef());
677  solve(epsEqn);
678  fvOptions.correct(epsilonm);
679  bound(epsilonm, this->epsilonMin_);
680 
681 
682  // Turbulent kinetic energy equation
683  tmp<fvScalarMatrix> kmEqn
684  (
685  fvm::ddt(km)
686  + fvm::div(phim, km)
687  - fvm::Sp(fvc::div(phim), km)
688  - fvm::laplacian(DkEff(nutm), km)
689  ==
690  Gm
691  - fvm::SuSp((2.0/3.0)*divUm, km)
692  - fvm::Sp(epsilonm/km, km)
693  + kSource()
694  + fvOptions(km)
695  );
696 
697  kmEqn.ref().relax();
698  fvOptions.constrain(kmEqn.ref());
699  solve(kmEqn);
700  fvOptions.correct(km);
701  bound(km, this->kMin_);
703 
704  volScalarField Cc2(rhom/(alphal*rholEff() + alphag*rhogEff()*Ct2_()));
705  kl = Cc2*km;
707  epsilonl = Cc2*epsilonm;
708  epsilonl.correctBoundaryConditions();
709  liquidTurbulence.correctNut();
710 
711  Ct2_() = Ct2();
712  kg = Ct2_()*kl;
714  epsilong = Ct2_()*epsilonl;
715  epsilong.correctBoundaryConditions();
716  nutg = Ct2_()*(liquidTurbulence.nu()/this->nu())*nutl;
717 }
718 
719 
720 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
721 
722 } // End namespace RASModels
723 } // End namespace Foam
724 
725 // ************************************************************************* //
Foam::RASModels::mixtureKEpsilon::mixFlux
tmp< surfaceScalarField > mixFlux(const surfaceScalarField &fc, const surfaceScalarField &fd) const
Definition: mixtureKEpsilon.C:468
Foam::virtualMassModel
Definition: virtualMassModel.H:55
Foam::RASModels::mixtureKEpsilon::mix
tmp< volScalarField > mix(const volScalarField &fc, const volScalarField &fd) const
Definition: mixtureKEpsilon.C:438
drag
Info<< "Reading strained laminar flame speed field Su\n"<< endl;volScalarField Su(IOobject("Su", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field betav\n"<< endl;volScalarField betav(IOobject("betav", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field Lobs\n"<< endl;volScalarField Lobs(IOobject("Lobs", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field CT\n"<< endl;volSymmTensorField CT(IOobject("CT", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field Nv\n"<< endl;volScalarField Nv(IOobject("Nv", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field nsv\n"<< endl;volSymmTensorField nsv(IOobject("nsv", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);IOdictionary PDRProperties(IOobject("PDRProperties", runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE));autoPtr< PDRDragModel > drag
Definition: createFields.H:165
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::RASModels::mixtureKEpsilon::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: mixtureKEpsilon.H:205
Foam::RASModels::mixtureKEpsilon::mixU
tmp< volScalarField > mixU(const volScalarField &fc, const volScalarField &fd) const
Definition: mixtureKEpsilon.C:452
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:194
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
fvOptions.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::twoPhaseSystem
Class which solves the volume fraction equations for two phases.
Definition: twoPhaseSystem.H:53
Foam::fv::optionList::correct
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
Definition: fvOptionListTemplates.C:355
Foam::tmp::clear
void clear() const noexcept
Definition: tmpI.H:287
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::RASModels::mixtureKEpsilon::initMixtureFields
void initMixtureFields()
Definition: mixtureKEpsilon.C:225
Foam::RASModel
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:52
Foam::fv::options::New
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:103
Foam::RASModels::mixtureKEpsilon::correctInletOutlet
void correctInletOutlet(volScalarField &vsf, const volScalarField &refVsf) const
Definition: mixtureKEpsilon.C:198
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::bound
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:35
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::RASModels::mixtureKEpsilon::epsilonBoundaryTypes
wordList epsilonBoundaryTypes(const volScalarField &epsilon) const
Definition: mixtureKEpsilon.C:176
Foam::RASModels::mixtureKEpsilon::bubbleG
tmp< volScalarField > bubbleG() const
Definition: mixtureKEpsilon.C:489
Foam::RASModels::mixtureKEpsilon::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: mixtureKEpsilon.C:541
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
Foam::liquid
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Definition: liquid.H:54
Foam::RASModels::mixtureKEpsilon::correctNut
virtual void correctNut()
Definition: mixtureKEpsilon.C:336
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
fluid
twoPhaseSystem & fluid
Definition: setRegionFluidFields.H:3
mixtureKEpsilon.H
Foam::RASModels::mixtureKEpsilon::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: mixtureKEpsilon.H:206
Foam::fvm::Su
tmp< fvMatrix< Type > > Su(const DimensionedField< Type, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
rho
rho
Definition: readInitialConditions.H:88
Foam::RASModels::mixtureKEpsilon::epsilonSource
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: mixtureKEpsilon.C:534
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::RASModels::mixtureKEpsilon::rhom
tmp< volScalarField > rhom() const
Definition: mixtureKEpsilon.C:427
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
Foam::fvm::SuSp
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:23
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::virtualMassModel::Cvm
virtual tmp< volScalarField > Cvm() const =0
Return the virtual mass coefficient.
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
bound.H
Bound the given scalar field if it has gone unbounded.
correct
fvOptions correct(rho)
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:187
Foam::RASModels::mixtureKEpsilon
Mixture k-epsilon turbulence model for two-phase gas-liquid systems.
Definition: mixtureKEpsilon.H:86
beta
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
Foam::solve
SolverPerformance< Type > solve(faMatrix< Type > &, Istream &)
Solve returning the solution statistics given convergence tolerance.
Foam::RASModels::mixtureKEpsilon::k_
volScalarField k_
Definition: mixtureKEpsilon.H:123
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::RASModels::mixtureKEpsilon::rholEff
tmp< volScalarField > rholEff() const
Definition: mixtureKEpsilon.C:405
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::fv::options
Finite-volume options.
Definition: fvOptions.H:55
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:42
Foam::RASModels::mixtureKEpsilon::Ct2
tmp< volScalarField > Ct2() const
Definition: mixtureKEpsilon.C:378
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::fvc::interpolate
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.
Foam::RASModels::mixtureKEpsilon::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: mixtureKEpsilon.C:316
fvmSup.H
Calculate the matrix for implicit and explicit sources.
Foam::transportModel
Base-class for all transport models used by the incompressible turbulence models.
Definition: transportModel.H:53
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:940
Foam::liquid::rho
scalar rho(scalar p, scalar T) const
Liquid density [kg/m^3].
Definition: liquidI.H:28
Foam::RASModels::mixtureKEpsilon::rhogEff
tmp< volScalarField > rhogEff() const
Definition: mixtureKEpsilon.C:414
U
U
Definition: pEqn.H:72
phig
surfaceScalarField phig("phig", -rhorAUf *ghf *fvc::snGrad(rho) *mesh.magSf())
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:783
Foam::List< word >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
fixedValueFvPatchFields.H
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::RASModels::mixtureKEpsilon::epsilon_
volScalarField epsilon_
Definition: mixtureKEpsilon.H:124
alphal
alphal
Definition: alphavPsi.H:12
Foam::RASModels::mixtureKEpsilon::kSource
virtual tmp< fvScalarMatrix > kSource() const
Definition: mixtureKEpsilon.C:527
Foam::eddyViscosity
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:55
Foam::eddyViscosity< RASModel< BasicTurbulenceModel > >::nut_
volScalarField nut_
Definition: eddyViscosity.H:66
Foam::dragModel
Definition: dragModel.H:51
epsilon
scalar epsilon
Definition: evaluateNearWall.H:7
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Foam::RASModels::mixtureKEpsilon::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: mixtureKEpsilon.H:204
Foam::GeometricField< vector, fvPatchField, volMesh >
inletOutletFvPatchFields.H
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
Foam::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:95
Foam::fvc::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:190
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::IOobject::MUST_READ
Definition: IOobject.H:185
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:106