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-2016 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 "fvOptions.H"
28 #include "bound.H"
29 #include "twoPhaseSystem.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>
45 mixtureKEpsilon<BasicTurbulenceModel>::mixtureKEpsilon
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::Boundary& 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 {
199  volScalarField::Boundary& bf = vsf.boundaryFieldRef();
200  const volScalarField::Boundary& refBf =
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  fv::options::New(this->mesh_).correct(this->nut_);
339 
340  BasicTurbulenceModel::correctNut();
341 }
342 
343 
344 template<class BasicTurbulenceModel>
347 {
348  if (!liquidTurbulencePtr_)
349  {
350  const volVectorField& U = this->U_;
351 
352  const transportModel& gas = this->transport();
353  const twoPhaseSystem& fluid =
354  refCast<const twoPhaseSystem>(gas.fluid());
355  const transportModel& liquid = fluid.otherPhase(gas);
356 
357  liquidTurbulencePtr_ =
359  (
361  (
363  (
365  liquid.name()
366  )
367  )
368  );
369  }
370 
371  return *liquidTurbulencePtr_;
372 }
373 
374 
375 template<class BasicTurbulenceModel>
377 {
378  const mixtureKEpsilon<BasicTurbulenceModel>& liquidTurbulence =
379  this->liquidTurbulence();
380 
381  const transportModel& gas = this->transport();
382  const twoPhaseSystem& fluid = refCast<const twoPhaseSystem>(gas.fluid());
383  const transportModel& liquid = fluid.otherPhase(gas);
384 
385  const volScalarField& alphag = this->alpha_;
386 
387  volScalarField magUr(mag(liquidTurbulence.U() - this->U()));
388 
390  (
391  (6*this->Cmu_/(4*sqrt(3.0/2.0)))
392  *fluid.drag(gas).K()/liquid.rho()
393  *(liquidTurbulence.k_/liquidTurbulence.epsilon_)
394  );
395  volScalarField Ct0((3 + beta)/(1 + beta + 2*gas.rho()/liquid.rho()));
396  volScalarField fAlphad((180 + (-4.71e3 + 4.26e4*alphag)*alphag)*alphag);
397 
398  return sqr(1 + (Ct0 - 1)*exp(-fAlphad));
399 }
400 
401 
402 template<class BasicTurbulenceModel>
404 {
405  const transportModel& gas = this->transport();
406  const twoPhaseSystem& fluid = refCast<const twoPhaseSystem>(gas.fluid());
407  return fluid.otherPhase(gas).rho();
408 }
409 
410 
411 template<class BasicTurbulenceModel>
413 {
414  const transportModel& gas = this->transport();
415  const twoPhaseSystem& fluid = refCast<const twoPhaseSystem>(gas.fluid());
416  return
417  gas.rho()
418  + fluid.virtualMass(gas).Cvm()*fluid.otherPhase(gas).rho();
419 }
420 
421 
422 template<class BasicTurbulenceModel>
424 {
425  const volScalarField& alphag = this->alpha_;
426  const volScalarField& alphal = this->liquidTurbulence().alpha_;
427 
428  return alphal*rholEff() + alphag*rhogEff();
429 }
430 
431 
432 template<class BasicTurbulenceModel>
434 (
435  const volScalarField& fc,
436  const volScalarField& fd
437 ) const
438 {
439  const volScalarField& alphag = this->alpha_;
440  const volScalarField& alphal = this->liquidTurbulence().alpha_;
441 
442  return (alphal*rholEff()*fc + alphag*rhogEff()*fd)/rhom_();
443 }
444 
445 
446 template<class BasicTurbulenceModel>
448 (
449  const volScalarField& fc,
450  const volScalarField& fd
451 ) const
452 {
453  const volScalarField& alphag = this->alpha_;
454  const volScalarField& alphal = this->liquidTurbulence().alpha_;
455 
456  return
457  (alphal*rholEff()*fc + alphag*rhogEff()*Ct2_()*fd)
458  /(alphal*rholEff() + alphag*rhogEff()*Ct2_());
459 }
460 
461 
462 template<class BasicTurbulenceModel>
464 (
465  const surfaceScalarField& fc,
466  const surfaceScalarField& fd
467 ) const
468 {
469  const volScalarField& alphag = this->alpha_;
470  const volScalarField& alphal = this->liquidTurbulence().alpha_;
471 
472  surfaceScalarField alphalf(fvc::interpolate(alphal));
473  surfaceScalarField alphagf(fvc::interpolate(alphag));
474 
475  surfaceScalarField rholEfff(fvc::interpolate(rholEff()));
476  surfaceScalarField rhogEfff(fvc::interpolate(rhogEff()));
477 
478  return
479  (alphalf*rholEfff*fc + alphagf*rhogEfff*fvc::interpolate(Ct2_())*fd)
480  /(alphalf*rholEfff + alphagf*rhogEfff*fvc::interpolate(Ct2_()));
481 }
482 
483 
484 template<class BasicTurbulenceModel>
486 {
487  const mixtureKEpsilon<BasicTurbulenceModel>& liquidTurbulence =
488  this->liquidTurbulence();
489 
490  const transportModel& gas = this->transport();
491  const twoPhaseSystem& fluid = refCast<const twoPhaseSystem>(gas.fluid());
492  const transportModel& liquid = fluid.otherPhase(gas);
493 
494  volScalarField magUr(mag(liquidTurbulence.U() - this->U()));
495 
496  // Lahey model
497  tmp<volScalarField> bubbleG
498  (
499  Cp_
500  *liquid*liquid.rho()
501  *(
502  pow3(magUr)
503  + pow(fluid.drag(gas).CdRe()*liquid.nu()/gas.d(), 4.0/3.0)
504  *pow(magUr, 5.0/3.0)
505  )
506  *gas
507  /gas.d()
508  );
509 
510  // Simple model
511  // tmp<volScalarField> bubbleG
512  // (
513  // Cp_*liquid*fluid.drag(gas).K()*sqr(magUr)
514  // );
515 
516  return bubbleG;
517 }
518 
519 
520 template<class BasicTurbulenceModel>
522 {
523  return fvm::Su(bubbleG()/rhom_(), km_());
524 }
525 
526 
527 template<class BasicTurbulenceModel>
529 {
530  return fvm::Su(C3_*epsilonm_()*bubbleG()/(rhom_()*km_()), epsilonm_());
531 }
532 
533 
534 template<class BasicTurbulenceModel>
536 {
537  const transportModel& gas = this->transport();
538  const twoPhaseSystem& fluid = refCast<const twoPhaseSystem>(gas.fluid());
539 
540  // Only solve the mixture turbulence for the gas-phase
541  if (&gas != &fluid.phase1())
542  {
543  // This is the liquid phase but check the model for the gas-phase
544  // is consistent
545  this->liquidTurbulence();
546 
547  return;
548  }
549 
550  if (!this->turbulence_)
551  {
552  return;
553  }
554 
555  // Initialise the mixture fields if they have not yet been constructed
556  initMixtureFields();
557 
558  // Local references to gas-phase properties
559  tmp<surfaceScalarField> phig = this->phi();
560  const volVectorField& Ug = this->U_;
561  const volScalarField& alphag = this->alpha_;
562  volScalarField& kg = this->k_;
563  volScalarField& epsilong = this->epsilon_;
564  volScalarField& nutg = this->nut_;
565 
566  // Local references to liquid-phase properties
567  mixtureKEpsilon<BasicTurbulenceModel>& liquidTurbulence =
568  this->liquidTurbulence();
569  tmp<surfaceScalarField> phil = liquidTurbulence.phi();
570  const volVectorField& Ul = liquidTurbulence.U_;
571  const volScalarField& alphal = liquidTurbulence.alpha_;
572  volScalarField& kl = liquidTurbulence.k_;
573  volScalarField& epsilonl = liquidTurbulence.epsilon_;
574  volScalarField& nutl = liquidTurbulence.nut_;
575 
576  // Local references to mixture properties
577  volScalarField& rhom = rhom_();
578  volScalarField& km = km_();
579  volScalarField& epsilonm = epsilonm_();
580 
581  fv::options& fvOptions(fv::options::New(this->mesh_));
582 
584 
585  // Update the effective mixture density
586  rhom = this->rhom();
587 
588  // Mixture flux
589  surfaceScalarField phim("phim", mixFlux(phil, phig));
590 
591  // Mixture velocity divergence
592  volScalarField divUm
593  (
594  mixU
595  (
596  fvc::div(fvc::absolute(phil, Ul)),
597  fvc::div(fvc::absolute(phig, Ug))
598  )
599  );
600 
602  {
603  tmp<volTensorField> tgradUl = fvc::grad(Ul);
605  (
606  new volScalarField
607  (
608  this->GName(),
609  nutl*(tgradUl() && dev(twoSymm(tgradUl())))
610  )
611  );
612  tgradUl.clear();
613 
614  // Update k, epsilon and G at the wall
615  kl.boundaryFieldRef().updateCoeffs();
616  epsilonl.boundaryFieldRef().updateCoeffs();
617 
618  Gc.ref().checkOut();
619  }
620 
622  {
623  tmp<volTensorField> tgradUg = fvc::grad(Ug);
625  (
626  new volScalarField
627  (
628  this->GName(),
629  nutg*(tgradUg() && dev(twoSymm(tgradUg())))
630  )
631  );
632  tgradUg.clear();
633 
634  // Update k, epsilon and G at the wall
635  kg.boundaryFieldRef().updateCoeffs();
636  epsilong.boundaryFieldRef().updateCoeffs();
637 
638  Gd.ref().checkOut();
639  }
640 
641  // Mixture turbulence generation
642  volScalarField Gm(mix(Gc, Gd));
643 
644  // Mixture turbulence viscosity
645  volScalarField nutm(mixU(nutl, nutg));
646 
647  // Update the mixture k and epsilon boundary conditions
648  km == mix(kl, kg);
649  bound(km, this->kMin_);
650  epsilonm == mix(epsilonl, epsilong);
651  bound(epsilonm, this->epsilonMin_);
652 
653  // Dissipation equation
654  tmp<fvScalarMatrix> epsEqn
655  (
656  fvm::ddt(epsilonm)
657  + fvm::div(phim, epsilonm)
658  - fvm::Sp(fvc::div(phim), epsilonm)
659  - fvm::laplacian(DepsilonEff(nutm), epsilonm)
660  ==
661  C1_*Gm*epsilonm/km
662  - fvm::SuSp(((2.0/3.0)*C1_)*divUm, epsilonm)
663  - fvm::Sp(C2_*epsilonm/km, epsilonm)
664  + epsilonSource()
665  + fvOptions(epsilonm)
666  );
667 
668  epsEqn.ref().relax();
669  fvOptions.constrain(epsEqn.ref());
670  epsEqn.ref().boundaryManipulate(epsilonm.boundaryFieldRef());
671  solve(epsEqn);
672  fvOptions.correct(epsilonm);
673  bound(epsilonm, this->epsilonMin_);
674 
675 
676  // Turbulent kinetic energy equation
677  tmp<fvScalarMatrix> kmEqn
678  (
679  fvm::ddt(km)
680  + fvm::div(phim, km)
681  - fvm::Sp(fvc::div(phim), km)
682  - fvm::laplacian(DkEff(nutm), km)
683  ==
684  Gm
685  - fvm::SuSp((2.0/3.0)*divUm, km)
686  - fvm::Sp(epsilonm/km, km)
687  + kSource()
688  + fvOptions(km)
689  );
690 
691  kmEqn.ref().relax();
692  fvOptions.constrain(kmEqn.ref());
693  solve(kmEqn);
694  fvOptions.correct(km);
695  bound(km, this->kMin_);
697 
698  volScalarField Cc2(rhom/(alphal*rholEff() + alphag*rhogEff()*Ct2_()));
699  kl = Cc2*km;
701  epsilonl = Cc2*epsilonm;
702  epsilonl.correctBoundaryConditions();
703  liquidTurbulence.correctNut();
704 
705  Ct2_() = Ct2();
706  kg = Ct2_()*kl;
708  epsilong = Ct2_()*epsilonl;
709  epsilong.correctBoundaryConditions();
710  nutg = Ct2_()*(liquidTurbulence.nu()/this->nu())*nutl;
711 }
712 
713 
714 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
715 
716 } // End namespace RASModels
717 } // End namespace Foam
718 
719 // ************************************************************************* //
surfaceScalarField & phi
BasicTurbulenceModel::alphaField alphaField
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
U
Definition: pEqn.H:83
const virtualMassModel & virtualMass(const phaseModel &phase) const
Return the virtual mass model for the given phase.
BasicTurbulenceModel::rhoField rhoField
tmp< volScalarField > rhogEff() const
multiphaseSystem & fluid
Definition: createFields.H:10
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:221
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
BasicTurbulenceModel::transportModel transportModel
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
tmp< fvMatrix< Type > > Su(const DimensionedField< Type, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
dimensionedScalar sqrt(const dimensionedScalar &ds)
tmp< volScalarField > mix(const volScalarField &fc, const volScalarField &fd) const
tmp< volScalarField > mixU(const volScalarField &fc, const volScalarField &fd) const
Generic dimensioned Type class.
virtual tmp< fvScalarMatrix > epsilonSource() const
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
Finite-volume options.
Definition: fvOptions.H:52
tmp< volScalarField > rholEff() const
virtual tmp< fvScalarMatrix > kSource() const
tmp< volScalarField > Ct2() const
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
fv::options & fvOptions
Class which solves the volume fraction equations for two phases.
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:49
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:230
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedScalar exp(const dimensionedScalar &ds)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
tmp< surfaceScalarField > mixFlux(const surfaceScalarField &fc, const surfaceScalarField &fd) const
tmp< fvMatrix< Type > > SuSp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
static const word propertiesName
Default name of the turbulence properties dictionary.
alphal
Definition: alphavPsi.H:12
A class for handling words, derived from string.
Definition: word.H:59
SolverPerformance< Type > solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
static word groupName(Name name, const word &group)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
virtual bool read()
Re-read model coefficients if they have changed.
const dimensionedScalar & rho() const
Definition: phaseModel.H:165
surfaceScalarField phig("phig",-rhorAUf *ghf *fvc::snGrad(rho)*mesh.magSf())
const phaseModel & otherPhase(const phaseModel &phase) const
Constant access the phase not given as an argument.
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
Bound the given scalar field if it has gone unbounded.
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil)-psi)*pSat, rhoMin);#1"/home/ubuntu/OpenFOAM-4.1/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H"1{alphav=max(min((rho-rholSat)/(rhovSat-rholSat), scalar(1)), scalar(0));alphal=1.0-alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:72
tmp< volScalarField > bubbleG() const
virtual tmp< volScalarField > Cvm() const =0
Return the virtual mass coefficient.
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.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
wordList epsilonBoundaryTypes(const volScalarField &epsilon) const
dimensionedScalar pow3(const dimensionedScalar &ds)
label patchi
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:188
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
Mixture k-epsilon turbulence model for two-phase gas-liquid systems.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
void correctBoundaryConditions()
Correct boundary field.
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
dimensioned< scalar > mag(const dimensioned< Type > &)
const dragModel & drag(const phaseModel &phase) const
Return the drag model for the given phase.
const phaseModel & phase1() const
Constant access phase model 1.
A class for managing temporary objects.
Definition: PtrList.H:54
virtual tmp< volScalarField > CdRe() const =0
Drag coefficient.
static options & New(const fvMesh &mesh)
Construct fvOptions and register to datbase if not present.
Definition: fvOptions.C:103
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual tmp< volScalarField > K(const volScalarField &Ur) const =0
The dragfunction K used in the momentum eq.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
volScalarField & nu
tmp< volScalarField > rhom() const
void correctInletOutlet(volScalarField &vsf, const volScalarField &refVsf) const
Calculate the matrix for implicit and explicit sources.
Namespace for OpenFOAM.