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