mixtureKEpsilon.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2013-2021 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 "fvModels.H"
28 #include "bound.H"
29 #include "phaseSystem.H"
30 #include "dragModel.H"
31 #include "virtualMassModel.H"
34 #include "fvmSup.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace RASModels
41 {
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
45 template<class BasicMomentumTransportModel>
47 (
48  const alphaField& alpha,
49  const rhoField& rho,
50  const volVectorField& U,
51  const surfaceScalarField& alphaRhoPhi,
52  const surfaceScalarField& phi,
53  const transportModel& transport,
54  const word& type
55 )
56 :
58  (
59  type,
60  alpha,
61  rho,
62  U,
63  alphaRhoPhi,
64  phi,
65  transport
66  ),
67 
68  liquidTurbulencePtr_(nullptr),
69 
70  Cmu_
71  (
73  (
74  "Cmu",
75  this->coeffDict_,
76  0.09
77  )
78  ),
79  C1_
80  (
82  (
83  "C1",
84  this->coeffDict_,
85  1.44
86  )
87  ),
88  C2_
89  (
91  (
92  "C2",
93  this->coeffDict_,
94  1.92
95  )
96  ),
97  C3_
98  (
100  (
101  "C3",
102  this->coeffDict_,
103  C2_.value()
104  )
105  ),
106  Cp_
107  (
109  (
110  "Cp",
111  this->coeffDict_,
112  0.25
113  )
114  ),
115  alphap_
116  (
118  (
119  "alphap",
120  this->coeffDict_,
121  1
122  )
123  ),
124  sigmak_
125  (
127  (
128  "sigmak",
129  this->coeffDict_,
130  1.0
131  )
132  ),
133  sigmaEps_
134  (
136  (
137  "sigmaEps",
138  this->coeffDict_,
139  1.3
140  )
141  ),
142 
143  k_
144  (
145  IOobject
146  (
147  IOobject::groupName("k", alphaRhoPhi.group()),
148  this->runTime_.timeName(),
149  this->mesh_,
152  ),
153  this->mesh_
154  ),
155  epsilon_
156  (
157  IOobject
158  (
159  IOobject::groupName("epsilon", alphaRhoPhi.group()),
160  this->runTime_.timeName(),
161  this->mesh_,
164  ),
165  this->mesh_
166  )
167 {
168  bound(k_, this->kMin_);
169  bound(epsilon_, this->epsilonMin_);
170 
171  if (type == typeName)
172  {
173  this->printCoeffs(type);
174  }
175 }
176 
177 
178 template<class BasicMomentumTransportModel>
180 (
181  const volScalarField& epsilon
182 ) const
183 {
184  const volScalarField::Boundary& ebf = epsilon.boundaryField();
185 
186  wordList ebt = ebf.types();
187 
188  forAll(ebf, patchi)
189  {
190  if (isA<fixedValueFvPatchScalarField>(ebf[patchi]))
191  {
192  ebt[patchi] = fixedValueFvPatchScalarField::typeName;
193  }
194  }
195 
196  return ebt;
197 }
198 
199 
200 template<class BasicMomentumTransportModel>
202 (
203  volScalarField& vsf,
204  const volScalarField& refVsf
205 ) const
206 {
207  volScalarField::Boundary& bf = vsf.boundaryFieldRef();
208  const volScalarField::Boundary& refBf =
209  refVsf.boundaryField();
210 
211  forAll(bf, patchi)
212  {
213  if
214  (
215  isA<inletOutletFvPatchScalarField>(bf[patchi])
216  && isA<inletOutletFvPatchScalarField>(refBf[patchi])
217  )
218  {
219  refCast<inletOutletFvPatchScalarField>
220  (bf[patchi]).refValue() =
221  refCast<const inletOutletFvPatchScalarField>
222  (refBf[patchi]).refValue();
223  }
224  }
225 }
226 
227 
228 template<class BasicMomentumTransportModel>
230 {
231  if (rhom_.valid()) return;
232 
233  // Local references to gas-phase properties
234  const volScalarField& kg = this->k_;
235  const volScalarField& epsilong = this->epsilon_;
236 
237  // Local references to liquid-phase properties
239  this->liquidTurbulence();
240  const volScalarField& kl = turbc.k_;
241  const volScalarField& epsilonl = turbc.epsilon_;
242 
243  word startTimeName
244  (
245  this->runTime_.timeName(this->runTime_.startTime().value())
246  );
247 
248  Ct2_.set
249  (
250  new volScalarField
251  (
252  IOobject
253  (
254  "Ct2",
255  startTimeName,
256  this->mesh_,
259  ),
260  Ct2()
261  )
262  );
263 
264  rhom_.set
265  (
266  new volScalarField
267  (
268  IOobject
269  (
270  "rhom",
271  startTimeName,
272  this->mesh_,
275  ),
276  rhom()
277  )
278  );
279 
280  km_.set
281  (
282  new volScalarField
283  (
284  IOobject
285  (
286  "km",
287  startTimeName,
288  this->mesh_,
291  ),
292  mix(kl, kg),
293  kl.boundaryField().types()
294  )
295  );
296  correctInletOutlet(km_(), kl);
297 
298  epsilonm_.set
299  (
300  new volScalarField
301  (
302  IOobject
303  (
304  "epsilonm",
305  startTimeName,
306  this->mesh_,
309  ),
310  mix(epsilonl, epsilong),
311  epsilonBoundaryTypes(epsilonl)
312  )
313  );
314  correctInletOutlet(epsilonm_(), epsilonl);
315 }
316 
317 
318 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
319 
320 template<class BasicMomentumTransportModel>
322 {
324  {
325  Cmu_.readIfPresent(this->coeffDict());
326  C1_.readIfPresent(this->coeffDict());
327  C2_.readIfPresent(this->coeffDict());
328  C3_.readIfPresent(this->coeffDict());
329  Cp_.readIfPresent(this->coeffDict());
330  sigmak_.readIfPresent(this->coeffDict());
331  sigmaEps_.readIfPresent(this->coeffDict());
332 
333  return true;
334  }
335  else
336  {
337  return false;
338  }
339 }
340 
341 
342 template<class BasicMomentumTransportModel>
344 {
345  this->nut_ = Cmu_*sqr(k_)/epsilon_;
346  this->nut_.correctBoundaryConditions();
347  fvConstraints::New(this->mesh_).constrain(this->nut_);
348 }
349 
350 
351 template<class BasicMomentumTransportModel>
354 {
355  if (!liquidTurbulencePtr_)
356  {
357  const volVectorField& U = this->U_;
358 
359  const phaseModel& gas = refCast<const phaseModel>(this->transport());
360  const phaseSystem& fluid = gas.fluid();
361  const phaseModel& liquid = fluid.otherPhase(gas);
362 
363  liquidTurbulencePtr_ =
365  (
366  U.db().lookupObject
367  <
369  >
370  (
372  (
373  momentumTransportModel::typeName,
374  liquid.name()
375  )
376  )
377  );
378  }
379 
380  return *liquidTurbulencePtr_;
381 }
382 
383 
384 template<class BasicMomentumTransportModel>
386 {
387  const mixtureKEpsilon<BasicMomentumTransportModel>& liquidTurbulence =
388  this->liquidTurbulence();
389 
390  const phaseModel& gas = refCast<const phaseModel>(this->transport());
391  const phaseSystem& fluid = gas.fluid();
392  const phaseModel& liquid = fluid.otherPhase(gas);
393 
394  const dragModel& drag = fluid.lookupSubModel<dragModel>(gas, liquid);
395 
396  const volScalarField& alphag = this->alpha_;
397 
398  volScalarField magUr(mag(liquidTurbulence.U() - this->U()));
399 
400  volScalarField beta
401  (
402  (6*this->Cmu_/(4*sqrt(3.0/2.0)))
403  *drag.K()/liquid.rho()
404  *(liquidTurbulence.k_/liquidTurbulence.epsilon_)
405  );
406  volScalarField Ct0((3 + beta)/(1 + beta + 2*gas.rho()/liquid.rho()));
407  volScalarField fAlphad((180 + (-4.71e3 + 4.26e4*alphag)*alphag)*alphag);
408 
409  return sqr(1 + (Ct0 - 1)*exp(-fAlphad));
410 }
411 
412 
413 template<class BasicMomentumTransportModel>
416 {
417  const phaseModel& gas = refCast<const phaseModel>(this->transport());
418  const phaseSystem& fluid = gas.fluid();
419  return fluid.otherPhase(gas).rho();
420 }
421 
422 
423 template<class BasicMomentumTransportModel>
426 {
427  const phaseModel& gas = refCast<const phaseModel>(this->transport());
428  const phaseSystem& fluid = gas.fluid();
429  const virtualMassModel& virtualMass =
430  fluid.lookupSubModel<virtualMassModel>(gas, fluid.otherPhase(gas));
431  return gas.rho() + virtualMass.Cvm()*fluid.otherPhase(gas).rho();
432 }
433 
434 
435 template<class BasicMomentumTransportModel>
437 {
438  const volScalarField& alphag = this->alpha_;
439  const volScalarField& alphal = this->liquidTurbulence().alpha_;
440 
441  return alphal*rholEff() + alphag*rhogEff();
442 }
443 
444 
445 template<class BasicMomentumTransportModel>
447 (
448  const volScalarField& fc,
449  const volScalarField& fd
450 ) const
451 {
452  const volScalarField& alphag = this->alpha_;
453  const volScalarField& alphal = this->liquidTurbulence().alpha_;
454 
455  return (alphal*rholEff()*fc + alphag*rhogEff()*fd)/rhom_();
456 }
457 
458 
459 template<class BasicMomentumTransportModel>
461 (
462  const volScalarField& fc,
463  const volScalarField& fd
464 ) const
465 {
466  const volScalarField& alphag = this->alpha_;
467  const volScalarField& alphal = this->liquidTurbulence().alpha_;
468 
469  return
470  (alphal*rholEff()*fc + alphag*rhogEff()*Ct2_()*fd)
471  /(alphal*rholEff() + alphag*rhogEff()*Ct2_());
472 }
473 
474 
475 template<class BasicMomentumTransportModel>
477 (
478  const surfaceScalarField& fc,
479  const surfaceScalarField& fd
480 ) const
481 {
482  const volScalarField& alphag = this->alpha_;
483  const volScalarField& alphal = this->liquidTurbulence().alpha_;
484 
485  surfaceScalarField alphalf(fvc::interpolate(alphal));
486  surfaceScalarField alphagf(fvc::interpolate(alphag));
487 
488  surfaceScalarField rholEfff(fvc::interpolate(rholEff()));
489  surfaceScalarField rhogEfff(fvc::interpolate(rhogEff()));
490 
491  return
492  (alphalf*rholEfff*fc + alphagf*rhogEfff*fvc::interpolate(Ct2_())*fd)
493  /(alphalf*rholEfff + alphagf*rhogEfff*fvc::interpolate(Ct2_()));
494 }
495 
496 
497 template<class BasicMomentumTransportModel>
500 {
501  const mixtureKEpsilon<BasicMomentumTransportModel>& liquidTurbulence =
502  this->liquidTurbulence();
503 
504  const phaseModel& gas = refCast<const phaseModel>(this->transport());
505  const phaseSystem& fluid = gas.fluid();
506  const phaseModel& liquid = fluid.otherPhase(gas);
507 
508  const dragModel& drag = fluid.lookupSubModel<dragModel>(gas, liquid);
509 
510  volScalarField magUr(mag(liquidTurbulence.U() - this->U()));
511 
512  // Lahey model
513  tmp<volScalarField> bubbleG
514  (
515  Cp_
516  *pos(alphap_ - gas)*liquid*liquid.rho()
517  *(
518  pow3(magUr)
519  + pow(drag.CdRe()*liquid.thermo().nu()/gas.d(), 4.0/3.0)
520  *pow(magUr, 5.0/3.0)
521  )
522  *gas
523  /gas.d()
524  );
525 
526  // Simple model
527  // tmp<volScalarField> bubbleG
528  // (
529  // Cp_*liquid*drag.K()*sqr(magUr)
530  // );
531 
532  return bubbleG;
533 }
534 
535 
536 template<class BasicMomentumTransportModel>
539 {
540  return fvm::Su(bubbleG()/rhom_(), km_());
541 }
542 
543 
544 template<class BasicMomentumTransportModel>
547 {
548  return fvm::Su(C3_*epsilonm_()*bubbleG()/(rhom_()*km_()), epsilonm_());
549 }
550 
551 
552 template<class BasicMomentumTransportModel>
554 {
555  const phaseModel& gas = refCast<const phaseModel>(this->transport());
556  const phaseSystem& fluid = gas.fluid();
557 
558  // Only solve the mixture turbulence for the gas-phase
559  if (&gas != &fluid.phases()[0])
560  {
561  // This is the liquid phase but check the model for the gas-phase
562  // is consistent
563  this->liquidTurbulence();
564 
565  return;
566  }
567 
568  if (!this->turbulence_)
569  {
570  return;
571  }
572 
573  // Initialise the mixture fields if they have not yet been constructed
574  initMixtureFields();
575 
576  // Local references to gas-phase properties
577  tmp<surfaceScalarField> phig = this->phi();
578  const volVectorField& Ug = this->U_;
579  const volScalarField& alphag = this->alpha_;
580  volScalarField& kg = this->k_;
581  volScalarField& epsilong = this->epsilon_;
582  volScalarField& nutg = this->nut_;
583 
584  // Local references to liquid-phase properties
586  this->liquidTurbulence();
587  tmp<surfaceScalarField> phil = liquidTurbulence.phi();
588  const volVectorField& Ul = liquidTurbulence.U_;
589  const volScalarField& alphal = liquidTurbulence.alpha_;
590  volScalarField& kl = liquidTurbulence.k_;
591  volScalarField& epsilonl = liquidTurbulence.epsilon_;
592  volScalarField& nutl = liquidTurbulence.nut_;
593 
594  // Local references to mixture properties
595  volScalarField& rhom = rhom_();
596  volScalarField& km = km_();
597  volScalarField& epsilonm = epsilonm_();
598 
599  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
601  (
602  Foam::fvConstraints::New(this->mesh_)
603  );
604 
606 
607  // Update the effective mixture density
608  rhom = this->rhom();
609 
610  // Mixture flux
611  const surfaceScalarField phim("phim", mixFlux(phil, phig));
612 
613  // Mixture velocity divergence
614  const volScalarField divUm
615  (
616  mixU
617  (
618  fvc::div(fvc::absolute(phil, Ul)),
619  fvc::div(fvc::absolute(phig, Ug))
620  )
621  );
622 
624  {
625  tmp<volTensorField> tgradUl = fvc::grad(Ul);
627  (
628  new volScalarField
629  (
630  this->GName(),
631  nutl*(tgradUl() && dev(twoSymm(tgradUl())))
632  )
633  );
634  tgradUl.clear();
635 
636  // Update k, epsilon and G at the wall
637  kl.boundaryFieldRef().updateCoeffs();
638  epsilonl.boundaryFieldRef().updateCoeffs();
639 
640  Gc.ref().checkOut();
641  }
642 
644  {
645  tmp<volTensorField> tgradUg = fvc::grad(Ug);
647  (
648  new volScalarField
649  (
650  this->GName(),
651  nutg*(tgradUg() && dev(twoSymm(tgradUg())))
652  )
653  );
654  tgradUg.clear();
655 
656  // Update k, epsilon and G at the wall
657  kg.boundaryFieldRef().updateCoeffs();
658  epsilong.boundaryFieldRef().updateCoeffs();
659 
660  Gd.ref().checkOut();
661  }
662 
663  // Mixture turbulence generation
664  const volScalarField Gm(mix(Gc, Gd));
665 
666  // Mixture turbulence viscosity
667  const volScalarField nutm(mixU(nutl, nutg));
668 
669  // Update the mixture k and epsilon boundary conditions
670  km == mix(kl, kg);
671  bound(km, this->kMin_);
672  epsilonm == mix(epsilonl, epsilong);
673  bound(epsilonm, this->epsilonMin_);
674 
675  // Dissipation equation
676  tmp<fvScalarMatrix> epsEqn
677  (
678  fvm::ddt(epsilonm)
679  + fvm::div(phim, epsilonm)
680  - fvm::Sp(fvc::div(phim), epsilonm)
681  - fvm::laplacian(DepsilonEff(nutm), epsilonm)
682  ==
683  C1_*Gm*epsilonm/km
684  - fvm::SuSp(((2.0/3.0)*C1_)*divUm, epsilonm)
685  - fvm::Sp(C2_*epsilonm/km, epsilonm)
686  + epsilonSource()
687  + fvModels.source(epsilonm)
688  );
689 
690  epsEqn.ref().relax();
691  fvConstraints.constrain(epsEqn.ref());
692  epsEqn.ref().boundaryManipulate(epsilonm.boundaryFieldRef());
693  solve(epsEqn);
694  fvConstraints.constrain(epsilonm);
695  bound(epsilonm, this->epsilonMin_);
696 
697 
698  // Turbulent kinetic energy equation
699  tmp<fvScalarMatrix> kmEqn
700  (
701  fvm::ddt(km)
702  + fvm::div(phim, km)
703  - fvm::Sp(fvc::div(phim), km)
704  - fvm::laplacian(DkEff(nutm), km)
705  ==
706  Gm
707  - fvm::SuSp((2.0/3.0)*divUm, km)
708  - fvm::Sp(epsilonm/km, km)
709  + kSource()
710  + fvModels.source(km)
711  );
712 
713  kmEqn.ref().relax();
714  fvConstraints.constrain(kmEqn.ref());
715  solve(kmEqn);
716  fvConstraints.constrain(km);
717  bound(km, this->kMin_);
719 
720  const volScalarField Cc2(rhom/(alphal*rholEff() + alphag*rhogEff()*Ct2_()));
721  kl = Cc2*km;
723  epsilonl = Cc2*epsilonm;
724  epsilonl.correctBoundaryConditions();
725  liquidTurbulence.correctNut();
726 
727  Ct2_() = Ct2();
728  kg = Ct2_()*kl;
730  epsilong = Ct2_()*epsilonl;
731  epsilong.correctBoundaryConditions();
732  nutg = Ct2_()*(liquidTurbulence.nu()/this->nu())*nutl;
733 }
734 
735 
736 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
737 
738 } // End namespace RASModels
739 } // End namespace Foam
740 
741 // ************************************************************************* //
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Definition: liquid.H:50
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:434
tmp< surfaceScalarField > mixFlux(const surfaceScalarField &fc, const surfaceScalarField &fd) const
tmp< volScalarField > mixU(const volScalarField &fc, const volScalarField &fd) const
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:237
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
tmp< fvMatrix< Type > > Su(const DimensionedField< Type, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
tmp< volScalarField > mix(const volScalarField &fc, const volScalarField &fd) const
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
dimensionedScalar sqrt(const dimensionedScalar &ds)
virtual tmp< volScalarField > rho() const =0
Return the density field.
wordList epsilonBoundaryTypes(const volScalarField &epsilon) const
const word & name() const
Definition: phaseModel.H:110
Generic dimensioned Type class.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:49
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
BasicMomentumTransportModel::alphaField alphaField
Definition: RASModel.H:88
Info<< "Reading strained laminar flame speed field Su\"<< endl;volScalarField Su(IOobject("Su", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field betav\"<< endl;volScalarField betav(IOobject("betav", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field Lobs\"<< endl;volScalarField Lobs(IOobject("Lobs", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field CT\"<< endl;volSymmTensorField CT(IOobject("CT", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field Nv\"<< endl;volScalarField Nv(IOobject("Nv", mesh.facesInstance(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE), mesh);Info<< "Reading field nsv\"<< 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:192
tmp< volScalarField > rhom() const
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
tmp< volScalarField > bubbleG() const
tmp< volScalarField > rhogEff() const
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:49
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
virtual tmp< volScalarField > K() const
Return the drag coefficient K.
dimensionedScalar pos(const dimensionedScalar &ds)
bool read(const char *, int32_t &)
Definition: int32IO.C:85
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedScalar exp(const dimensionedScalar &ds)
BasicMomentumTransportModel::transportModel transportModel
Definition: RASModel.H:90
Foam::fvConstraints & fvConstraints
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:70
alphal
Definition: alphavPsi.H:12
A class for handling words, derived from string.
Definition: word.H:59
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-9/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:68
const rhoThermo & thermo() const
Return const-access to phase rhoThermo.
Definition: phaseModel.H:121
static word groupName(Name name, const word &group)
Finite volume models.
Definition: fvModels.H:60
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 phaseModel & otherPhase(const phaseModel &phase) const
Return the phase not given as an argument in a two-phase system.
Definition: phaseSystemI.H:112
virtual tmp< fvScalarMatrix > epsilonSource() const
phaseSystem & fluid
Definition: createFields.H:11
phi
Definition: correctPhi.H:3
virtual tmp< fvScalarMatrix > kSource() const
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.
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:96
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.
surfaceScalarField phig(-rhorAUf *ghf *fvc::snGrad(rho) *mesh.magSf())
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
Foam::fvModels & fvModels
tmp< volScalarField > Ct2() const
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
BasicMomentumTransportModel::rhoField rhoField
Definition: RASModel.H:89
label patchi
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:188
U
Definition: pEqn.H:72
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
Finite volume constraints.
Definition: fvConstraints.H:57
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
const phaseSystem & fluid() const
Return the system to which this phase belongs.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
virtual tmp< volScalarField > nu() const
Kinematic viscosity of mixture [m^2/s].
Definition: fluidThermo.C:84
void correctBoundaryConditions()
Correct boundary field.
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
dimensioned< scalar > mag(const dimensioned< Type > &)
mixtureKEpsilon(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
Construct from components.
A class for managing temporary objects.
Definition: PtrList.H:53
virtual tmp< volScalarField > CdRe() const =0
Drag coefficient.
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:312
void correctInletOutlet(volScalarField &vsf, const volScalarField &refVsf) const
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
tmp< volScalarField > rholEff() const
Calculate the matrix for implicit and explicit sources.
Namespace for OpenFOAM.