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-2023 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 "dispersedDragModel.H"
36 #include "fvmSup.H"
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42 namespace RASModels
43 {
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 
47 template<class BasicMomentumTransportModel>
49 (
50  const alphaField& alpha,
51  const rhoField& rho,
52  const volVectorField& U,
53  const surfaceScalarField& alphaRhoPhi,
54  const surfaceScalarField& phi,
55  const viscosity& viscosity,
56  const word& type
57 )
58 :
59  eddyViscosity<RASModel<BasicMomentumTransportModel>>
60  (
61  type,
62  alpha,
63  rho,
64  U,
65  alphaRhoPhi,
66  phi,
67  viscosity
68  ),
69 
70  gasTurbulencePtr_(nullptr),
71 
72  Cmu_
73  (
74  dimensioned<scalar>::lookupOrAddToDict
75  (
76  "Cmu",
77  this->coeffDict_,
78  0.09
79  )
80  ),
81  C1_
82  (
83  dimensioned<scalar>::lookupOrAddToDict
84  (
85  "C1",
86  this->coeffDict_,
87  1.44
88  )
89  ),
90  C2_
91  (
92  dimensioned<scalar>::lookupOrAddToDict
93  (
94  "C2",
95  this->coeffDict_,
96  1.92
97  )
98  ),
99  C3_
100  (
101  dimensioned<scalar>::lookupOrAddToDict
102  (
103  "C3",
104  this->coeffDict_,
105  C2_.value()
106  )
107  ),
108  Cp_
109  (
110  dimensioned<scalar>::lookupOrAddToDict
111  (
112  "Cp",
113  this->coeffDict_,
114  0.25
115  )
116  ),
117  alphap_
118  (
119  dimensioned<scalar>::lookupOrAddToDict
120  (
121  "alphap",
122  this->coeffDict_,
123  1
124  )
125  ),
126  sigmak_
127  (
128  dimensioned<scalar>::lookupOrAddToDict
129  (
130  "sigmak",
131  this->coeffDict_,
132  1.0
133  )
134  ),
135  sigmaEps_
136  (
137  dimensioned<scalar>::lookupOrAddToDict
138  (
139  "sigmaEps",
140  this->coeffDict_,
141  1.3
142  )
143  ),
144 
145  k_
146  (
147  IOobject
148  (
149  this->groupName("k"),
150  this->runTime_.name(),
151  this->mesh_,
152  IOobject::MUST_READ,
153  IOobject::AUTO_WRITE
154  ),
155  this->mesh_
156  ),
157  epsilon_
158  (
159  IOobject
160  (
161  this->groupName("epsilon"),
162  this->runTime_.name(),
163  this->mesh_,
164  IOobject::MUST_READ,
165  IOobject::AUTO_WRITE
166  ),
167  this->mesh_
168  )
169 {
170  bound(k_, this->kMin_);
171  bound(epsilon_, this->epsilonMin_);
172 
173  if (type == typeName)
174  {
175  this->printCoeffs(type);
176  }
177 
178  const phaseModel& phase = refCast<const phaseModel>(this->properties());
179 
180  // Construct mixture properties only for liquid phase (phase 1)
181  if (phase.index() == 1)
182  {
183  km_.set
184  (
185  new volScalarField
186  (
187  IOobject
188  (
189  "km",
190  this->runTime_.name(),
191  this->mesh_,
194  ),
195  this->mesh_
196  )
197  );
198 
199  epsilonm_.set
200  (
201  new volScalarField
202  (
203  IOobject
204  (
205  "epsilonm",
206  this->runTime_.name(),
207  this->mesh_,
210  ),
211  this->mesh_
212  )
213  );
214 
215  Ct2_.set
216  (
217  new volScalarField
218  (
219  IOobject
220  (
221  "Ct2",
222  this->runTime_.name(),
223  this->mesh_,
226  ),
227  this->mesh_,
229  )
230  );
231 
232  rhom_.set
233  (
234  new volScalarField
235  (
236  IOobject
237  (
238  "rhom",
239  this->runTime_.name(),
240  this->mesh_,
243  ),
244  this->mesh_,
246  )
247  );
248  }
249 }
250 
251 
252 template<class BasicMomentumTransportModel>
254 {
255  static bool initialised = false;
256 
257  if (!initialised)
258  {
259  Ct2_() = Ct2();
260  rhom_() = rhom();
261  initialised = true;
262  }
263 }
264 
265 
266 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
267 
268 template<class BasicMomentumTransportModel>
270 {
272  {
273  Cmu_.readIfPresent(this->coeffDict());
274  C1_.readIfPresent(this->coeffDict());
275  C2_.readIfPresent(this->coeffDict());
276  C3_.readIfPresent(this->coeffDict());
277  Cp_.readIfPresent(this->coeffDict());
278  sigmak_.readIfPresent(this->coeffDict());
279  sigmaEps_.readIfPresent(this->coeffDict());
280 
281  return true;
282  }
283  else
284  {
285  return false;
286  }
287 }
288 
289 
290 template<class BasicMomentumTransportModel>
292 {
293  this->nut_ = Cmu_*sqr(k_)/epsilon_;
294  this->nut_.correctBoundaryConditions();
295  fvConstraints::New(this->mesh_).constrain(this->nut_);
296 }
297 
298 
299 template<class BasicMomentumTransportModel>
302 {
303  if (!gasTurbulencePtr_)
304  {
305  const volVectorField& U = this->U_;
306 
307  const phaseModel& liquid =
308  refCast<const phaseModel>(this->properties());
309  const phaseSystem& fluid = liquid.fluid();
310  const phaseModel& gas = fluid.phases()[0];
311 
312  gasTurbulencePtr_ =
314  (
315  U.db().lookupObject
316  <
318  >
319  (
321  (
322  momentumTransportModel::typeName,
323  gas.name()
324  )
325  )
326  );
327  }
328 
329  return *gasTurbulencePtr_;
330 }
331 
332 
333 template<class BasicMomentumTransportModel>
335 {
336  const mixtureKEpsilon<BasicMomentumTransportModel>& gasTurbulence =
337  this->gasTurbulence();
338 
339  const phaseModel& liquid = refCast<const phaseModel>(this->properties());
340  const phaseSystem& fluid = liquid.fluid();
341  const phaseModel& gas = fluid.phases()[0];
342 
346 
347  const volScalarField& alphag = gasTurbulence.alpha();
348  const volScalarField magUr(mag(this->U() - gasTurbulence.U()));
349 
350  volScalarField beta
351  (
352  (6*this->Cmu_/(4*sqrt(3.0/2.0)))
353  *drag.K()/liquid.rho()
354  *(k_/epsilon_)
355  );
356  volScalarField Ct0((3 + beta)/(1 + beta + 2*gas.rho()/liquid.rho()));
357  volScalarField fAlphad((180 + (-4.71e3 + 4.26e4*alphag)*alphag)*alphag);
358 
359  return sqr(1 + (Ct0 - 1)*exp(-fAlphad));
360 }
361 
362 
363 template<class BasicMomentumTransportModel>
366 {
367  const phaseModel& liquid = refCast<const phaseModel>(this->properties());
368  return liquid.rho();
369 }
370 
371 
372 template<class BasicMomentumTransportModel>
375 {
376  const phaseModel& liquid = refCast<const phaseModel>(this->properties());
377  const phaseSystem& fluid = liquid.fluid();
378  const phaseModel& gas = fluid.phases()[0];
379 
384 
385  return gas.rho() + virtualMass.Cvm()*liquid.rho();
386 }
387 
388 
389 template<class BasicMomentumTransportModel>
391 {
392  const volScalarField& alphal = this->alpha_;
393  const volScalarField& alphag = this->gasTurbulence().alpha_;
394 
395  return alphal*rholEff() + alphag*rhogEff();
396 }
397 
398 
399 template<class BasicMomentumTransportModel>
401 (
402  const volScalarField& fc,
403  const volScalarField& fd
404 ) const
405 {
406  const volScalarField& alphal = this->alpha_;
407  const volScalarField& alphag = this->gasTurbulence().alpha_;
408 
409  return (alphal*rholEff()*fc + alphag*rhogEff()*fd)/rhom_();
410 }
411 
412 
413 template<class BasicMomentumTransportModel>
415 (
416  const volScalarField& fc,
417  const volScalarField& fd
418 ) const
419 {
420  const volScalarField& alphal = this->alpha_;
421  const volScalarField& alphag = this->gasTurbulence().alpha_;
422 
423  return
424  (alphal*rholEff()*fc + alphag*rhogEff()*Ct2_()*fd)
425  /(alphal*rholEff() + alphag*rhogEff()*Ct2_());
426 }
427 
428 
429 template<class BasicMomentumTransportModel>
431 (
432  const surfaceScalarField& fc,
433  const surfaceScalarField& fd
434 ) const
435 {
436  const volScalarField& alphal = this->alpha_;
437  const volScalarField& alphag = this->gasTurbulence().alpha_;
438 
439  surfaceScalarField alphalf(fvc::interpolate(alphal));
440  surfaceScalarField alphagf(fvc::interpolate(alphag));
441 
442  surfaceScalarField rholEfff(fvc::interpolate(rholEff()));
443  surfaceScalarField rhogEfff(fvc::interpolate(rhogEff()));
444 
445  return
446  (alphalf*rholEfff*fc + alphagf*rhogEfff*fvc::interpolate(Ct2_())*fd)
447  /(alphalf*rholEfff + alphagf*rhogEfff*fvc::interpolate(Ct2_()));
448 }
449 
450 
451 template<class BasicMomentumTransportModel>
454 {
455  const mixtureKEpsilon<BasicMomentumTransportModel>& gasTurbulence =
456  this->gasTurbulence();
457 
458  const phaseModel& liquid = refCast<const phaseModel>(this->properties());
459  const phaseSystem& fluid = liquid.fluid();
460  const phaseModel& gas = fluid.phases()[0];
461 
465 
466  volScalarField magUr(mag(this->U() - gasTurbulence.U()));
467 
468  // Lahey model
469  tmp<volScalarField> bubbleG
470  (
471  Cp_
472  *pos(alphap_ - gas)*liquid*liquid.rho()
473  *(
474  pow3(magUr)
475  + pow(drag.CdRe()*liquid.thermo().nu()/gas.d(), 4.0/3.0)
476  *pow(magUr, 5.0/3.0)
477  )
478  *gas
479  /gas.d()
480  );
481 
482  // Simple model
483  // tmp<volScalarField> bubbleG
484  // (
485  // Cp_*liquid*drag.K()*sqr(magUr)
486  // );
487 
488  return bubbleG;
489 }
490 
491 
492 template<class BasicMomentumTransportModel>
495 {
496  return fvm::Su(bubbleG()/rhom_(), km_());
497 }
498 
499 
500 template<class BasicMomentumTransportModel>
503 {
504  return fvm::Su(C3_*epsilonm_()*bubbleG()/(rhom_()*km_()), epsilonm_());
505 }
506 
507 
508 template<class BasicMomentumTransportModel>
510 {
511  const phaseModel& phase = refCast<const phaseModel>(this->properties());
512 
513  // Only solve the mixture turbulence for the liquid phase (phase 1)
514  if (phase.index() == 0)
515  {
516  // This is the liquid phase but check the model for the gas-phase
517  // is consistent
518  this->gasTurbulence();
519 
520  return;
521  }
522  else
523  {
524  initMixtureFields();
525  }
526 
527  if (!this->turbulence_)
528  {
529  return;
530  }
531 
532  // Local references to liquid-phase properties
533  tmp<surfaceScalarField> phil = this->phi();
534  const volVectorField& Ul = this->U_;
535  const volScalarField& alphal = this->alpha_;
536  volScalarField& kl = this->k_;
537  volScalarField& epsilonl = this->epsilon_;
538  volScalarField& nutl = this->nut_;
539 
540  // Local references to gas-phase properties
542  this->gasTurbulence();
543  tmp<surfaceScalarField> phig = gasTurbulence.phi();
544  const volVectorField& Ug = gasTurbulence.U_;
545  const volScalarField& alphag = gasTurbulence.alpha_;
546  volScalarField& kg = gasTurbulence.k_;
547  volScalarField& epsilong = gasTurbulence.epsilon_;
548  volScalarField& nutg = gasTurbulence.nut_;
549 
550  // Local references to mixture properties
551  volScalarField& rhom = rhom_();
552  volScalarField& km = km_();
553  volScalarField& epsilonm = epsilonm_();
554 
555  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
557  (
558  Foam::fvConstraints::New(this->mesh_)
559  );
560 
562 
563  // Update the effective mixture density
564  rhom = this->rhom();
565 
566  // Mixture flux
567  const surfaceScalarField phim("phim", mixFlux(phil, phig));
568 
569  // Mixture velocity divergence
570  const volScalarField divUm
571  (
572  mixU
573  (
574  fvc::div(fvc::absolute(phil, Ul)),
575  fvc::div(fvc::absolute(phig, Ug))
576  )
577  );
578 
580  {
581  tmp<volTensorField> tgradUl = fvc::grad(Ul);
583  (
584  new volScalarField
585  (
586  this->GName(),
587  nutl*(tgradUl() && dev(twoSymm(tgradUl())))
588  )
589  );
590  tgradUl.clear();
591 
592  // Update k, epsilon and G at the wall
594  epsilonl.boundaryFieldRef().updateCoeffs();
595 
596  Gc.ref().checkOut();
597  }
598 
600  {
601  tmp<volTensorField> tgradUg = fvc::grad(Ug);
603  (
604  new volScalarField
605  (
606  this->GName(),
607  nutg*(tgradUg() && dev(twoSymm(tgradUg())))
608  )
609  );
610  tgradUg.clear();
611 
612  // Update k, epsilon and G at the wall
614  epsilong.boundaryFieldRef().updateCoeffs();
615 
616  Gd.ref().checkOut();
617  }
618 
619  // Mixture turbulence generation
620  const volScalarField Gm(mix(Gc, Gd));
621 
622  // Mixture turbulence viscosity
623  const volScalarField nutm(mixU(nutl, nutg));
624 
625  // Update the mixture k and epsilon boundary conditions
626  km == mix(kl, kg);
627  bound(km, this->kMin_);
628  epsilonm == mix(epsilonl, epsilong);
629  bound(epsilonm, this->epsilonMin_);
630 
631  // Dissipation equation
632  tmp<fvScalarMatrix> epsEqn
633  (
634  fvm::ddt(epsilonm)
635  + fvm::div(phim, epsilonm)
636  + fvm::SuSp(-fvc::div(phim), epsilonm)
637  - fvm::laplacian(DepsilonEff(nutm), epsilonm)
638  ==
639  C1_*Gm*epsilonm/km
640  - fvm::SuSp(((2.0/3.0)*C1_)*divUm, epsilonm)
641  - fvm::Sp(C2_*epsilonm/km, epsilonm)
642  + epsilonSource()
643  + fvModels.source(epsilonm)
644  );
645 
646  epsEqn.ref().relax();
647  fvConstraints.constrain(epsEqn.ref());
648  epsEqn.ref().boundaryManipulate(epsilonm.boundaryFieldRef());
649  solve(epsEqn);
650  fvConstraints.constrain(epsilonm);
651  bound(epsilonm, this->epsilonMin_);
652 
653  // Turbulent kinetic energy equation
654  tmp<fvScalarMatrix> kmEqn
655  (
656  fvm::ddt(km)
657  + fvm::div(phim, km)
658  + fvm::SuSp(-fvc::div(phim), km)
659  - fvm::laplacian(DkEff(nutm), km)
660  ==
661  Gm
662  - fvm::SuSp((2.0/3.0)*divUm, km)
663  - fvm::Sp(epsilonm/km, km)
664  + kSource()
665  + fvModels.source(km)
666  );
667 
668  kmEqn.ref().relax();
669  fvConstraints.constrain(kmEqn.ref());
670  solve(kmEqn);
672  bound(km, this->kMin_);
674 
675  const volScalarField Cc2(rhom/(alphal*rholEff() + alphag*rhogEff()*Ct2_()));
676  kl = Cc2*km;
678  epsilonl = Cc2*epsilonm;
679  epsilonl.correctBoundaryConditions();
680  correctNut();
681 
682  Ct2_() = Ct2();
683  kg = Ct2_()*kl;
685  epsilong = Ct2_()*epsilonl;
686  epsilong.correctBoundaryConditions();
687  nutg = Ct2_()*(this->nu()/gasTurbulence.nu())*nutl;
688 }
689 
690 
691 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
692 
693 } // End namespace RASModels
694 } // End namespace Foam
695 
696 // ************************************************************************* //
Bound the given scalar field where it is below the specified minimum.
void updateCoeffs()
Update the boundary condition coefficients.
Generic GeometricField class.
void correctBoundaryConditions()
Correct boundary field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
static word groupName(Name name, const word &group)
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:56
Mixture k-epsilon turbulence model for two-phase gas-liquid systems.
virtual tmp< fvScalarMatrix > epsilonSource() const
tmp< volScalarField > Ct2() const
tmp< volScalarField > rhogEff() const
tmp< volScalarField > mixU(const volScalarField &fc, const volScalarField &fd) const
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
tmp< volScalarField > rhom() const
tmp< volScalarField > bubbleG() const
autoPtr< volScalarField > rhom_
autoPtr< volScalarField > km_
autoPtr< volScalarField > Ct2_
virtual tmp< fvScalarMatrix > kSource() const
tmp< surfaceScalarField > mixFlux(const surfaceScalarField &fc, const surfaceScalarField &fd) const
tmp< volScalarField > rholEff() const
autoPtr< volScalarField > epsilonm_
tmp< volScalarField > mix(const volScalarField &fc, const volScalarField &fd) const
virtual bool read()
Re-read model coefficients if they have changed.
mixtureKEpsilon(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:100
Generic dimensioned Type class.
Class to represent a interface between phases where one phase is considered dispersed within the othe...
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
volScalarField nut_
Definition: eddyViscosity.H:60
Finite volume constraints.
Definition: fvConstraints.H:62
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
Finite volume models.
Definition: fvModels.H:65
BasicMomentumTransportModel::alphaField alphaField
BasicMomentumTransportModel::rhoField rhoField
Generic thermophysical properties class for a liquid in which the functions and coefficients for each...
Definition: liquid.H:53
scalar rho(scalar p, scalar T) const
Liquid density [kg/m^3].
Definition: liquidI.H:26
tmp< volScalarField > d() const
Return the Sauter-mean diameter.
Definition: phaseModel.C:145
label index() const
Return the index of the phase.
Definition: phaseModel.C:121
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:109
virtual const volScalarField & rho() const =0
Return the density field.
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:73
const phaseModelList & phases() const
Return the phase models.
Definition: phaseSystemI.H:41
const ModelType & lookupInterfacialModel(const phaseInterface &interface) const
Return a sub model for an interface.
A class for managing temporary objects.
Definition: tmp.H:55
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:237
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
virtual tmp< volScalarField > Cvm() const =0
Return the virtual mass coefficient.
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvConstraints & fvConstraints(Foam::fvConstraints::New(mesh))
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
Calculate the matrix for implicit and explicit sources.
Info<< "Reading strained laminar flame speed field Su\n"<< endl;volScalarField Su(IOobject("Su", runTime.name(), 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:194
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:202
tmp< fvMatrix< Type > > Su(const DimensionedField< Type, volMesh > &, const VolField< Type > &)
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedScalar pos(const dimensionedScalar &ds)
dimensionedScalar exp(const dimensionedScalar &ds)
bool read(const char *, int32_t &)
Definition: int32IO.C:85
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
const dimensionSet dimless
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionSet dimDensity
dimensioned< scalar > mag(const dimensioned< Type > &)
bool bound(volScalarField &, const dimensionedScalar &min)
Bound the given scalar field where it is below the specified min value.
Definition: bound.C:31
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.