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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "mixtureKEpsilon.H"
29 #include "fvOptions.H"
30 #include "bound.H"
31 #include "twoPhaseSystem.H"
32 #include "dragModel.H"
33 #include "virtualMassModel.H"
36 #include "fvmSup.H"
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42 namespace RASModels
43 {
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 
47 template<class BasicTurbulenceModel>
48 mixtureKEpsilon<BasicTurbulenceModel>::mixtureKEpsilon
49 (
50  const alphaField& alpha,
51  const rhoField& rho,
52  const volVectorField& U,
53  const surfaceScalarField& alphaRhoPhi,
54  const surfaceScalarField& phi,
55  const transportModel& transport,
56  const word& propertiesName,
57  const word& type
58 )
59 :
61  (
62  type,
63  alpha,
64  rho,
65  U,
66  alphaRhoPhi,
67  phi,
68  transport,
69  propertiesName
70  ),
71 
72  liquidTurbulencePtr_(nullptr),
73 
74  Cmu_
75  (
77  (
78  "Cmu",
79  this->coeffDict_,
80  0.09
81  )
82  ),
83  C1_
84  (
86  (
87  "C1",
88  this->coeffDict_,
89  1.44
90  )
91  ),
92  C2_
93  (
95  (
96  "C2",
97  this->coeffDict_,
98  1.92
99  )
100  ),
101  C3_
102  (
104  (
105  "C3",
106  this->coeffDict_,
107  C2_.value()
108  )
109  ),
110  Cp_
111  (
113  (
114  "Cp",
115  this->coeffDict_,
116  0.25
117  )
118  ),
119  sigmak_
120  (
122  (
123  "sigmak",
124  this->coeffDict_,
125  1.0
126  )
127  ),
128  sigmaEps_
129  (
131  (
132  "sigmaEps",
133  this->coeffDict_,
134  1.3
135  )
136  ),
137 
138  k_
139  (
140  IOobject
141  (
142  IOobject::groupName("k", alphaRhoPhi.group()),
143  this->runTime_.timeName(),
144  this->mesh_,
147  ),
148  this->mesh_
149  ),
150  epsilon_
151  (
152  IOobject
153  (
154  IOobject::groupName("epsilon", alphaRhoPhi.group()),
155  this->runTime_.timeName(),
156  this->mesh_,
159  ),
160  this->mesh_
161  )
162 {
163  bound(k_, this->kMin_);
164  bound(epsilon_, this->epsilonMin_);
165 
166  if (type == typeName)
167  {
168  this->printCoeffs(type);
169  }
170 }
171 
172 
173 template<class BasicTurbulenceModel>
175 (
176  const volScalarField& epsilon
177 ) const
178 {
179  const volScalarField::Boundary& ebf = epsilon.boundaryField();
180 
181  wordList ebt = ebf.types();
182 
183  forAll(ebf, patchi)
184  {
185  if (isA<fixedValueFvPatchScalarField>(ebf[patchi]))
186  {
187  ebt[patchi] = fixedValueFvPatchScalarField::typeName;
188  }
189  }
190 
191  return ebt;
192 }
193 
194 
195 template<class BasicTurbulenceModel>
197 (
198  volScalarField& vsf,
199  const volScalarField& refVsf
200 ) const
201 {
202  volScalarField::Boundary& bf = vsf.boundaryFieldRef();
203  const volScalarField::Boundary& refBf =
204  refVsf.boundaryField();
205 
206  forAll(bf, patchi)
207  {
208  if
209  (
210  isA<inletOutletFvPatchScalarField>(bf[patchi])
211  && isA<inletOutletFvPatchScalarField>(refBf[patchi])
212  )
213  {
214  refCast<inletOutletFvPatchScalarField>
215  (bf[patchi]).refValue() =
216  refCast<const inletOutletFvPatchScalarField>
217  (refBf[patchi]).refValue();
218  }
219  }
220 }
221 
222 
223 template<class BasicTurbulenceModel>
225 {
226  if (rhom_.valid()) return;
227 
228  // Local references to gas-phase properties
229  const volScalarField& kg = this->k_;
230  const volScalarField& epsilong = this->epsilon_;
231 
232  // Local references to liquid-phase properties
233  mixtureKEpsilon<BasicTurbulenceModel>& turbc = this->liquidTurbulence();
234  const volScalarField& kl = turbc.k_;
235  const volScalarField& epsilonl = turbc.epsilon_;
236 
237  word startTimeName
238  (
239  this->runTime_.timeName(this->runTime_.startTime().value())
240  );
241 
242  Ct2_.reset
243  (
244  new volScalarField
245  (
246  IOobject
247  (
248  "Ct2",
249  startTimeName,
250  this->mesh_,
253  ),
254  Ct2()
255  )
256  );
257 
258  rhom_.reset
259  (
260  new volScalarField
261  (
262  IOobject
263  (
264  "rhom",
265  startTimeName,
266  this->mesh_,
269  ),
270  rhom()
271  )
272  );
273 
274  km_.reset
275  (
276  new volScalarField
277  (
278  IOobject
279  (
280  "km",
281  startTimeName,
282  this->mesh_,
285  ),
286  mix(kl, kg),
287  kl.boundaryField().types()
288  )
289  );
290  correctInletOutlet(km_(), kl);
291 
292  epsilonm_.reset
293  (
294  new volScalarField
295  (
296  IOobject
297  (
298  "epsilonm",
299  startTimeName,
300  this->mesh_,
303  ),
304  mix(epsilonl, epsilong),
305  epsilonBoundaryTypes(epsilonl)
306  )
307  );
308  correctInletOutlet(epsilonm_(), epsilonl);
309 }
310 
311 
312 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
313 
314 template<class BasicTurbulenceModel>
316 {
318  {
319  Cmu_.readIfPresent(this->coeffDict());
320  C1_.readIfPresent(this->coeffDict());
321  C2_.readIfPresent(this->coeffDict());
322  C3_.readIfPresent(this->coeffDict());
323  Cp_.readIfPresent(this->coeffDict());
324  sigmak_.readIfPresent(this->coeffDict());
325  sigmaEps_.readIfPresent(this->coeffDict());
326 
327  return true;
328  }
329 
330  return false;
331 }
332 
333 
334 template<class BasicTurbulenceModel>
336 {
337  this->nut_ = Cmu_*sqr(k_)/epsilon_;
338  this->nut_.correctBoundaryConditions();
339  fv::options::New(this->mesh_).correct(this->nut_);
340 
341  BasicTurbulenceModel::correctNut();
342 }
343 
344 
345 template<class BasicTurbulenceModel>
348 {
349  if (!liquidTurbulencePtr_)
350  {
351  const volVectorField& U = this->U_;
352 
353  const transportModel& gas = this->transport();
354  const twoPhaseSystem& fluid =
355  refCast<const twoPhaseSystem>(gas.fluid());
356  const transportModel& liquid = fluid.otherPhase(gas);
357 
358  liquidTurbulencePtr_ =
360  (
361  U.db().lookupObject<mixtureKEpsilon<BasicTurbulenceModel>>
362  (
364  (
366  liquid.name()
367  )
368  )
369  );
370  }
371 
372  return *liquidTurbulencePtr_;
373 }
374 
375 
376 template<class BasicTurbulenceModel>
378 {
379  const mixtureKEpsilon<BasicTurbulenceModel>& liquidTurbulence =
380  this->liquidTurbulence();
381 
382  const transportModel& gas = this->transport();
383  const twoPhaseSystem& fluid = refCast<const twoPhaseSystem>(gas.fluid());
384  const transportModel& liquid = fluid.otherPhase(gas);
385 
386  const volScalarField& alphag = this->alpha_;
387 
388  volScalarField magUr(mag(liquidTurbulence.U() - this->U()));
389 
391  (
392  (6*this->Cmu_/(4*sqrt(3.0/2.0)))
393  *fluid.Kd()/liquid.rho()
394  *(liquidTurbulence.k_/liquidTurbulence.epsilon_)
395  );
396  volScalarField Ct0((3 + beta)/(1 + beta + 2*gas.rho()/liquid.rho()));
397  volScalarField fAlphad((180 + (-4.71e3 + 4.26e4*alphag)*alphag)*alphag);
398 
399  return sqr(1 + (Ct0 - 1)*exp(-fAlphad));
400 }
401 
402 
403 template<class BasicTurbulenceModel>
405 {
406  const transportModel& gas = this->transport();
407  const twoPhaseSystem& fluid = refCast<const twoPhaseSystem>(gas.fluid());
408  return fluid.otherPhase(gas).rho();
409 }
410 
411 
412 template<class BasicTurbulenceModel>
414 {
415  const transportModel& gas = this->transport();
416  const twoPhaseSystem& fluid = refCast<const twoPhaseSystem>(gas.fluid());
417  const virtualMassModel& virtualMass =
418  fluid.lookupSubModel<virtualMassModel>(gas, fluid.otherPhase(gas));
419  return
420  gas.rho()
421  + virtualMass.Cvm()*fluid.otherPhase(gas).rho();
422 }
423 
424 
425 template<class BasicTurbulenceModel>
427 {
428  const volScalarField& alphag = this->alpha_;
429  const volScalarField& alphal = this->liquidTurbulence().alpha_;
430 
431  return alphal*rholEff() + alphag*rhogEff();
432 }
433 
434 
435 template<class BasicTurbulenceModel>
437 (
438  const volScalarField& fc,
439  const volScalarField& fd
440 ) const
441 {
442  const volScalarField& alphag = this->alpha_;
443  const volScalarField& alphal = this->liquidTurbulence().alpha_;
444 
445  return (alphal*rholEff()*fc + alphag*rhogEff()*fd)/rhom_();
446 }
447 
448 
449 template<class BasicTurbulenceModel>
451 (
452  const volScalarField& fc,
453  const volScalarField& fd
454 ) const
455 {
456  const volScalarField& alphag = this->alpha_;
457  const volScalarField& alphal = this->liquidTurbulence().alpha_;
458 
459  return
460  (alphal*rholEff()*fc + alphag*rhogEff()*Ct2_()*fd)
461  /(alphal*rholEff() + alphag*rhogEff()*Ct2_());
462 }
463 
464 
465 template<class BasicTurbulenceModel>
467 (
468  const surfaceScalarField& fc,
469  const surfaceScalarField& fd
470 ) const
471 {
472  const volScalarField& alphag = this->alpha_;
473  const volScalarField& alphal = this->liquidTurbulence().alpha_;
474 
476  surfaceScalarField alphagf(fvc::interpolate(alphag));
477 
478  surfaceScalarField rholEfff(fvc::interpolate(rholEff()));
479  surfaceScalarField rhogEfff(fvc::interpolate(rhogEff()));
480 
481  return
482  (alphalf*rholEfff*fc + alphagf*rhogEfff*fvc::interpolate(Ct2_())*fd)
483  /(alphalf*rholEfff + alphagf*rhogEfff*fvc::interpolate(Ct2_()));
484 }
485 
486 
487 template<class BasicTurbulenceModel>
489 {
490  const mixtureKEpsilon<BasicTurbulenceModel>& liquidTurbulence =
491  this->liquidTurbulence();
492 
493  const transportModel& gas = this->transport();
494  const twoPhaseSystem& fluid = refCast<const twoPhaseSystem>(gas.fluid());
495  const transportModel& liquid = fluid.otherPhase(gas);
496 
497  const dragModel& drag = fluid.lookupSubModel<dragModel>(gas, liquid);
498 
499  volScalarField magUr(mag(liquidTurbulence.U() - this->U()));
500 
501  // Lahey model
502  tmp<volScalarField> bubbleG
503  (
504  Cp_
505  *liquid*liquid.rho()
506  *(
507  pow3(magUr)
508  + pow(drag.CdRe()*liquid.nu()/gas.d(), 4.0/3.0)
509  *pow(magUr, 5.0/3.0)
510  )
511  *gas
512  /gas.d()
513  );
514 
515  // Simple model
516  // tmp<volScalarField> bubbleG
517  // (
518  // Cp_*liquid*drag.K()*sqr(magUr)
519  // );
520 
521  return bubbleG;
522 }
523 
524 
525 template<class BasicTurbulenceModel>
527 {
528  return fvm::Su(bubbleG()/rhom_(), km_());
529 }
530 
531 
532 template<class BasicTurbulenceModel>
534 {
535  return fvm::Su(C3_*epsilonm_()*bubbleG()/(rhom_()*km_()), epsilonm_());
536 }
537 
538 
539 template<class BasicTurbulenceModel>
541 {
542  const transportModel& gas = this->transport();
543  const twoPhaseSystem& fluid = refCast<const twoPhaseSystem>(gas.fluid());
544 
545  // Only solve the mixture turbulence for the gas-phase
546  if (&gas != &fluid.phase1())
547  {
548  // This is the liquid phase but check the model for the gas-phase
549  // is consistent
550  this->liquidTurbulence();
551 
552  return;
553  }
554 
555  if (!this->turbulence_)
556  {
557  return;
558  }
559 
560  // Initialise the mixture fields if they have not yet been constructed
561  initMixtureFields();
562 
563  // Local references to gas-phase properties
564  tmp<surfaceScalarField> phig = this->phi();
565  const volVectorField& Ug = this->U_;
566  const volScalarField& alphag = this->alpha_;
567  volScalarField& kg = this->k_;
568  volScalarField& epsilong = this->epsilon_;
569  volScalarField& nutg = this->nut_;
570 
571  // Local references to liquid-phase properties
572  mixtureKEpsilon<BasicTurbulenceModel>& liquidTurbulence =
573  this->liquidTurbulence();
574  tmp<surfaceScalarField> phil = liquidTurbulence.phi();
575  const volVectorField& Ul = liquidTurbulence.U_;
576  const volScalarField& alphal = liquidTurbulence.alpha_;
577  volScalarField& kl = liquidTurbulence.k_;
578  volScalarField& epsilonl = liquidTurbulence.epsilon_;
579  volScalarField& nutl = liquidTurbulence.nut_;
580 
581  // Local references to mixture properties
582  volScalarField& rhom = rhom_();
583  volScalarField& km = km_();
584  volScalarField& epsilonm = epsilonm_();
585 
586  fv::options& fvOptions(fv::options::New(this->mesh_));
587 
589 
590  // Update the effective mixture density
591  rhom = this->rhom();
592 
593  // Mixture flux
594  surfaceScalarField phim("phim", mixFlux(phil, phig));
595 
596  // Mixture velocity divergence
597  volScalarField divUm
598  (
599  mixU
600  (
601  fvc::div(fvc::absolute(phil, Ul)),
603  )
604  );
605 
607  {
608  tmp<volTensorField> tgradUl = fvc::grad(Ul);
610  (
611  new volScalarField
612  (
613  this->GName(),
614  nutl*(tgradUl() && dev(twoSymm(tgradUl())))
615  )
616  );
617  tgradUl.clear();
618 
619  // Update k, epsilon and G at the wall
620  kl.boundaryFieldRef().updateCoeffs();
621  epsilonl.boundaryFieldRef().updateCoeffs();
622 
623  Gc.ref().checkOut();
624  }
625 
627  {
628  tmp<volTensorField> tgradUg = fvc::grad(Ug);
630  (
631  new volScalarField
632  (
633  this->GName(),
634  nutg*(tgradUg() && dev(twoSymm(tgradUg())))
635  )
636  );
637  tgradUg.clear();
638 
639  // Update k, epsilon and G at the wall
640  kg.boundaryFieldRef().updateCoeffs();
641  epsilong.boundaryFieldRef().updateCoeffs();
642 
643  Gd.ref().checkOut();
644  }
645 
646  // Mixture turbulence generation
647  volScalarField Gm(mix(Gc, Gd));
648 
649  // Mixture turbulence viscosity
650  volScalarField nutm(mixU(nutl, nutg));
651 
652  // Update the mixture k and epsilon boundary conditions
653  km == mix(kl, kg);
654  bound(km, this->kMin_);
655  epsilonm == mix(epsilonl, epsilong);
656  bound(epsilonm, this->epsilonMin_);
657 
658  // Dissipation equation
659  tmp<fvScalarMatrix> epsEqn
660  (
661  fvm::ddt(epsilonm)
662  + fvm::div(phim, epsilonm)
663  - fvm::Sp(fvc::div(phim), epsilonm)
664  - fvm::laplacian(DepsilonEff(nutm), epsilonm)
665  ==
666  C1_*Gm*epsilonm/km
667  - fvm::SuSp(((2.0/3.0)*C1_)*divUm, epsilonm)
668  - fvm::Sp(C2_*epsilonm/km, epsilonm)
669  + epsilonSource()
670  + fvOptions(epsilonm)
671  );
672 
673  epsEqn.ref().relax();
674  fvOptions.constrain(epsEqn.ref());
675  epsEqn.ref().boundaryManipulate(epsilonm.boundaryFieldRef());
676  solve(epsEqn);
677  fvOptions.correct(epsilonm);
678  bound(epsilonm, this->epsilonMin_);
679 
680 
681  // Turbulent kinetic energy equation
682  tmp<fvScalarMatrix> kmEqn
683  (
684  fvm::ddt(km)
685  + fvm::div(phim, km)
686  - fvm::Sp(fvc::div(phim), km)
687  - fvm::laplacian(DkEff(nutm), km)
688  ==
689  Gm
690  - fvm::SuSp((2.0/3.0)*divUm, km)
691  - fvm::Sp(epsilonm/km, km)
692  + kSource()
693  + fvOptions(km)
694  );
695 
696  kmEqn.ref().relax();
697  fvOptions.constrain(kmEqn.ref());
698  solve(kmEqn);
699  fvOptions.correct(km);
700  bound(km, this->kMin_);
702 
703  volScalarField Cc2(rhom/(alphal*rholEff() + alphag*rhogEff()*Ct2_()));
704  kl = Cc2*km;
706  epsilonl = Cc2*epsilonm;
707  epsilonl.correctBoundaryConditions();
708  liquidTurbulence.correctNut();
709 
710  Ct2_() = Ct2();
711  kg = Ct2_()*kl;
713  epsilong = Ct2_()*epsilonl;
714  epsilong.correctBoundaryConditions();
715  nutg = Ct2_()*(liquidTurbulence.nu()/this->nu())*nutl;
716 }
717 
718 
719 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
720 
721 } // End namespace RASModels
722 } // End namespace Foam
723 
724 // ************************************************************************* //
Foam::RASModels::mixtureKEpsilon::mixFlux
tmp< surfaceScalarField > mixFlux(const surfaceScalarField &fc, const surfaceScalarField &fd) const
Definition: mixtureKEpsilon.C:467
Foam::virtualMassModel
Definition: virtualMassModel.H:54
Foam::RASModels::mixtureKEpsilon::mix
tmp< volScalarField > mix(const volScalarField &fc, const volScalarField &fd) const
Definition: mixtureKEpsilon.C:437
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:104
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:451
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:129
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:62
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:321
Foam::tmp::clear
void clear() const noexcept
Definition: tmpI.H:325
Foam::tmp< volScalarField >
Foam::RASModels::mixtureKEpsilon::initMixtureFields
void initMixtureFields()
Definition: mixtureKEpsilon.C:224
Foam::RASModel
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:51
Foam::fv::options::New
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:104
Foam::RASModels::mixtureKEpsilon::correctInletOutlet
void correctInletOutlet(volScalarField &vsf, const volScalarField &refVsf) const
Definition: mixtureKEpsilon.C:197
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:175
Foam::RASModels::mixtureKEpsilon::bubbleG
tmp< volScalarField > bubbleG() const
Definition: mixtureKEpsilon.C:488
Foam::RASModels::mixtureKEpsilon::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: mixtureKEpsilon.C:540
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::RASModels::mixtureKEpsilon::correctNut
virtual void correctNut()
Definition: mixtureKEpsilon.C:335
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:96
Foam::RASModels::mixtureKEpsilon::epsilonSource
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: mixtureKEpsilon.C:533
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::RASModels::mixtureKEpsilon::rhom
tmp< volScalarField > rhom() const
Definition: mixtureKEpsilon.C:426
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:258
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:122
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:404
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:43
Foam::RASModels::mixtureKEpsilon::Ct2
tmp< volScalarField > Ct2() const
Definition: mixtureKEpsilon.C:377
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:315
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:909
Foam::RASModels::mixtureKEpsilon::rhogEff
tmp< volScalarField > rhogEff() const
Definition: mixtureKEpsilon.C:413
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:752
Foam::List< word >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::IOobject::groupName
static word groupName(StringType name, const word &group)
Create dot-delimited name.group.
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:526
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:55
epsilon
scalar epsilon
Definition: evaluateNearWall.H:7
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:120
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:106