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-2025 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 
46 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
47 
48 template<class BasicMomentumTransportModel>
50 (
51  const volScalarField& Cc2
52 )
53 {
54  epsilonm_() = max
55  (
56  epsilonm_(),
57  Cmu_*sqr(k_)/(this->nutMaxCoeff_*Cc2*this->nu())
58  );
59 }
60 
61 
62 template<class BasicMomentumTransportModel>
64 {
65  this->nut_ = Cmu_*sqr(k_)/epsilon_;
66  this->nut_.correctBoundaryConditions();
67  fvConstraints::New(this->mesh_).constrain(this->nut_);
68 }
69 
70 
71 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
72 
73 template<class BasicMomentumTransportModel>
75 (
76  const alphaField& alpha,
77  const rhoField& rho,
78  const volVectorField& U,
79  const surfaceScalarField& alphaRhoPhi,
80  const surfaceScalarField& phi,
81  const viscosity& viscosity,
82  const word& type
83 )
84 :
85  eddyViscosity<RASModel<BasicMomentumTransportModel>>
86  (
87  type,
88  alpha,
89  rho,
90  U,
91  alphaRhoPhi,
92  phi,
93  viscosity
94  ),
95 
96  gasTurbulencePtr_(nullptr),
97 
98  Cmu_("Cmu", this->coeffDict(), 0.09),
99  C1_("C1", this->coeffDict(), 1.44),
100  C2_("C2", this->coeffDict(), 1.92),
101  C3_("C3", this->coeffDict(), C2_.value()),
102  Cp_("Cp", this->coeffDict(), 0.25),
103  alphap_("alphap", this->coeffDict(), 1),
104  sigmak_("sigmak", this->coeffDict(), 1.0),
105  sigmaEps_("sigmaEps", this->coeffDict(), 1.3),
106 
107  k_
108  (
109  IOobject
110  (
111  this->groupName("k"),
112  this->runTime_.name(),
113  this->mesh_,
114  IOobject::MUST_READ,
115  IOobject::AUTO_WRITE
116  ),
117  this->mesh_
118  ),
119  epsilon_
120  (
121  IOobject
122  (
123  this->groupName("epsilon"),
124  this->runTime_.name(),
125  this->mesh_,
126  IOobject::MUST_READ,
127  IOobject::AUTO_WRITE
128  ),
129  this->mesh_
130  )
131 {
132  bound(k_, this->kMin_);
133 
134  const phaseModel& phase = refCast<const phaseModel>(this->properties());
135 
136  // Construct mixture properties only for liquid phase (phase 1)
137  if (phase.index() == 1)
138  {
139  km_.set
140  (
141  new volScalarField
142  (
143  IOobject
144  (
145  "km",
146  this->runTime_.name(),
147  this->mesh_,
150  ),
151  this->mesh_
152  )
153  );
154 
155  epsilonm_.set
156  (
157  new volScalarField
158  (
159  IOobject
160  (
161  "epsilonm",
162  this->runTime_.name(),
163  this->mesh_,
166  ),
167  this->mesh_
168  )
169  );
170 
171  Ct2_.set
172  (
173  new volScalarField
174  (
175  IOobject
176  (
177  "Ct2",
178  this->runTime_.name(),
179  this->mesh_,
182  ),
183  this->mesh_,
185  )
186  );
187 
188  rhom_.set
189  (
190  new volScalarField
191  (
192  IOobject
193  (
194  "rhom",
195  this->runTime_.name(),
196  this->mesh_,
199  ),
200  this->mesh_,
202  )
203  );
204  }
205 }
206 
207 
208 template<class BasicMomentumTransportModel>
210 {
211  static bool initialised = false;
212 
213  if (!initialised)
214  {
215  Ct2_() = Ct2();
216  rhom_() = rhom();
217  initialised = true;
218  }
219 }
220 
221 
222 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
223 
224 template<class BasicMomentumTransportModel>
226 {
228  {
229  Cmu_.readIfPresent(this->coeffDict());
230  C1_.readIfPresent(this->coeffDict());
231  C2_.readIfPresent(this->coeffDict());
232  C3_.readIfPresent(this->coeffDict());
233  Cp_.readIfPresent(this->coeffDict());
234  sigmak_.readIfPresent(this->coeffDict());
235  sigmaEps_.readIfPresent(this->coeffDict());
236 
237  return true;
238  }
239  else
240  {
241  return false;
242  }
243 }
244 
245 
246 template<class BasicMomentumTransportModel>
249 {
250  if (!gasTurbulencePtr_)
251  {
252  const volVectorField& U = this->U_;
253 
254  const phaseModel& liquid =
255  refCast<const phaseModel>(this->properties());
256  const phaseSystem& fluid = liquid.fluid();
257  const phaseModel& gas = fluid.otherPhase(liquid);
258 
259  gasTurbulencePtr_ =
261  (
262  U.db().lookupObject
263  <
265  >
266  (
268  (
269  momentumTransportModel::typeName,
270  gas.name()
271  )
272  )
273  );
274  }
275 
276  return *gasTurbulencePtr_;
277 }
278 
279 
280 template<class BasicMomentumTransportModel>
282 {
283  const mixtureKEpsilon<BasicMomentumTransportModel>& gasTurbulence =
284  this->gasTurbulence();
285 
286  const phaseModel& liquid = refCast<const phaseModel>(this->properties());
287  const phaseSystem& fluid = liquid.fluid();
288  const phaseModel& gas = fluid.otherPhase(liquid);
289 
293 
294  const volScalarField& alphag = gasTurbulence.alpha();
295  const volScalarField magUr(mag(this->U() - gasTurbulence.U()));
296 
297  volScalarField beta
298  (
299  (6*this->Cmu_/(4*sqrt(3.0/2.0)))
300  *drag.K()/liquid.rho()
301  *(k_/epsilon_)
302  );
303  volScalarField Ct0((3 + beta)/(1 + beta + 2*gas.rho()/liquid.rho()));
304  volScalarField fAlphad((180 + (-4.71e3 + 4.26e4*alphag)*alphag)*alphag);
305 
306  return sqr(1 + (Ct0 - 1)*exp(-fAlphad));
307 }
308 
309 
310 template<class BasicMomentumTransportModel>
313 {
314  const phaseModel& liquid = refCast<const phaseModel>(this->properties());
315  return liquid.rho();
316 }
317 
318 
319 template<class BasicMomentumTransportModel>
322 {
323  const phaseModel& liquid = refCast<const phaseModel>(this->properties());
324  const phaseSystem& fluid = liquid.fluid();
325  const phaseModel& gas = fluid.otherPhase(liquid);
326 
331 
332  return gas.rho() + virtualMass.Cvm()*liquid.rho();
333 }
334 
335 
336 template<class BasicMomentumTransportModel>
338 {
339  const volScalarField& alphal = this->alpha_;
340  const volScalarField& alphag = this->gasTurbulence().alpha_;
341 
342  return alphal*rholEff() + alphag*rhogEff();
343 }
344 
345 
346 template<class BasicMomentumTransportModel>
348 (
349  const volScalarField& fc,
350  const volScalarField& fd
351 ) const
352 {
353  const volScalarField& alphal = this->alpha_;
354  const volScalarField& alphag = this->gasTurbulence().alpha_;
355 
356  return (alphal*rholEff()*fc + alphag*rhogEff()*fd)/rhom_();
357 }
358 
359 
360 template<class BasicMomentumTransportModel>
362 (
363  const volScalarField& fc,
364  const volScalarField& fd
365 ) const
366 {
367  const volScalarField& alphal = this->alpha_;
368  const volScalarField& alphag = this->gasTurbulence().alpha_;
369 
370  return
371  (alphal*rholEff()*fc + alphag*rhogEff()*Ct2_()*fd)
372  /(alphal*rholEff() + alphag*rhogEff()*Ct2_());
373 }
374 
375 
376 template<class BasicMomentumTransportModel>
378 (
379  const surfaceScalarField& fc,
380  const surfaceScalarField& fd
381 ) const
382 {
383  const volScalarField& alphal = this->alpha_;
384  const volScalarField& alphag = this->gasTurbulence().alpha_;
385 
386  surfaceScalarField alphalf(fvc::interpolate(alphal));
387  surfaceScalarField alphagf(fvc::interpolate(alphag));
388 
389  surfaceScalarField rholEfff(fvc::interpolate(rholEff()));
390  surfaceScalarField rhogEfff(fvc::interpolate(rhogEff()));
391 
392  return
393  (alphalf*rholEfff*fc + alphagf*rhogEfff*fvc::interpolate(Ct2_())*fd)
394  /(alphalf*rholEfff + alphagf*rhogEfff*fvc::interpolate(Ct2_()));
395 }
396 
397 
398 template<class BasicMomentumTransportModel>
401 {
402  const mixtureKEpsilon<BasicMomentumTransportModel>& gasTurbulence =
403  this->gasTurbulence();
404 
405  const phaseModel& liquid = refCast<const phaseModel>(this->properties());
406  const phaseSystem& fluid = liquid.fluid();
407  const phaseModel& gas = fluid.otherPhase(liquid);
408 
412 
413  volScalarField magUr(mag(this->U() - gasTurbulence.U()));
414 
415  // Lahey model
416  tmp<volScalarField> bubbleG
417  (
418  Cp_
419  *pos(alphap_ - gas)*liquid*liquid.rho()
420  *(
421  pow3(magUr)
422  + pow(drag.CdRe()*liquid.fluidThermo().nu()/gas.d(), 4.0/3.0)
423  *pow(magUr, 5.0/3.0)
424  )
425  *gas
426  /gas.d()
427  );
428 
429  // Simple model
430  // tmp<volScalarField> bubbleG
431  // (
432  // Cp_*liquid*drag.K()*sqr(magUr)
433  // );
434 
435  return bubbleG;
436 }
437 
438 
439 template<class BasicMomentumTransportModel>
442 {
443  return fvm::Su(bubbleG()/rhom_(), km_());
444 }
445 
446 
447 template<class BasicMomentumTransportModel>
450 {
451  return fvm::Su(C3_*epsilonm_()*bubbleG()/(rhom_()*km_()), epsilonm_());
452 }
453 
454 
455 template<class BasicMomentumTransportModel>
457 {
458  const phaseModel& phase = refCast<const phaseModel>(this->properties());
459 
460  // Only solve the mixture turbulence for the liquid phase (phase 1)
461  if (phase.index() == 0)
462  {
463  // This is the liquid phase but check the model for the gas-phase
464  // is consistent
465  this->gasTurbulence();
466 
467  return;
468  }
469  else
470  {
471  initMixtureFields();
472  }
473 
474  if (!this->turbulence_)
475  {
476  return;
477  }
478 
479  // Local references to liquid-phase properties
480  tmp<surfaceScalarField> phil = this->phi();
481  const volVectorField& Ul = this->U_;
482  const volScalarField& alphal = this->alpha_;
483  volScalarField& kl = this->k_;
484  volScalarField& epsilonl = this->epsilon_;
485  volScalarField& nutl = this->nut_;
486 
487  // Local references to gas-phase properties
489  this->gasTurbulence();
490  tmp<surfaceScalarField> phig = gasTurbulence.phi();
491  const volVectorField& Ug = gasTurbulence.U_;
492  const volScalarField& alphag = gasTurbulence.alpha_;
493  volScalarField& kg = gasTurbulence.k_;
494  volScalarField& epsilong = gasTurbulence.epsilon_;
495  volScalarField& nutg = gasTurbulence.nut_;
496 
497  // Local references to mixture properties
498  volScalarField& rhom = rhom_();
499  volScalarField& km = km_();
500  volScalarField& epsilonm = epsilonm_();
501 
502  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
504  (
505  Foam::fvConstraints::New(this->mesh_)
506  );
507 
509 
510  // Update the effective mixture density
511  rhom = this->rhom();
512 
513  // Mixture flux
514  const surfaceScalarField phim("phim", mixFlux(phil, phig));
515 
516  // Mixture velocity divergence
517  const volScalarField divUm
518  (
519  mixU
520  (
521  fvc::div(fvc::absolute(phil, Ul)),
522  fvc::div(fvc::absolute(phig, Ug))
523  )
524  );
525 
527  {
528  tmp<volTensorField> tgradUl = fvc::grad(Ul);
530  (
531  new volScalarField
532  (
533  this->GName(),
534  nutl*(tgradUl() && dev(twoSymm(tgradUl())))
535  )
536  );
537  tgradUl.clear();
538 
539  // Update k, epsilon and G at the wall
541  epsilonl.boundaryFieldRef().updateCoeffs();
542 
543  Gc.ref().checkOut();
544  }
545 
547  {
548  tmp<volTensorField> tgradUg = fvc::grad(Ug);
550  (
551  new volScalarField
552  (
553  this->GName(),
554  nutg*(tgradUg() && dev(twoSymm(tgradUg())))
555  )
556  );
557  tgradUg.clear();
558 
559  // Update k, epsilon and G at the wall
561  epsilong.boundaryFieldRef().updateCoeffs();
562 
563  Gd.ref().checkOut();
564  }
565 
566  // Mixture turbulence generation
567  const volScalarField Gm(mix(Gc, Gd));
568 
569  // Mixture turbulence viscosity
570  const volScalarField nutm(mixU(nutl, nutg));
571 
572  // Update the mixture k and epsilon boundary conditions
573  km == mix(kl, kg);
574  epsilonm == mix(epsilonl, epsilong);
575 
576  // Dissipation equation
577  tmp<fvScalarMatrix> epsEqn
578  (
579  fvm::ddt(epsilonm)
580  + fvm::div(phim, epsilonm)
581  + fvm::SuSp(-fvc::div(phim), epsilonm)
582  - fvm::laplacian(DepsilonEff(nutm), epsilonm)
583  ==
584  C1_*Gm*epsilonm/km
585  - fvm::SuSp(((2.0/3.0)*C1_)*divUm, epsilonm)
586  - fvm::Sp(C2_*epsilonm/km, epsilonm)
587  + epsilonSource()
588  + fvModels.source(rhom, epsilonm)/rhom
589  );
590 
591  epsEqn.ref().relax();
592  fvConstraints.constrain(epsEqn.ref());
593  epsEqn.ref().boundaryManipulate(epsilonm.boundaryFieldRef());
594  solve(epsEqn);
595  fvConstraints.constrain(epsilonm);
596 
597  const volScalarField Cc2(rhom/(alphal*rholEff() + alphag*rhogEff()*Ct2_()));
598  boundEpsilonm(Cc2);
599 
600  // Turbulent kinetic energy equation
601  tmp<fvScalarMatrix> kmEqn
602  (
603  fvm::ddt(km)
604  + fvm::div(phim, km)
605  + fvm::SuSp(-fvc::div(phim), km)
606  - fvm::laplacian(DkEff(nutm), km)
607  ==
608  Gm
609  - fvm::SuSp((2.0/3.0)*divUm, km)
610  - fvm::Sp(epsilonm/km, km)
611  + kSource()
612  + fvModels.source(rhom, km)/rhom
613  );
614 
615  kmEqn.ref().relax();
616  fvConstraints.constrain(kmEqn.ref());
617  solve(kmEqn);
619  bound(km, this->kMin_);
621  boundEpsilonm(Cc2);
622 
623  kl = Cc2*km;
625  epsilonl = Cc2*epsilonm;
626  epsilonl.correctBoundaryConditions();
627  correctNut();
628 
629  Ct2_() = Ct2();
630  kg = Ct2_()*kl;
632  epsilong = Ct2_()*epsilonl;
633  epsilong.correctBoundaryConditions();
634  nutg = Ct2_()*(this->nu()/gasTurbulence.nu())*nutl;
635 }
636 
637 
638 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
639 
640 } // End namespace RASModels
641 } // End namespace Foam
642 
643 // ************************************************************************* //
Bound the given scalar field where it is below the specified minimum.
void updateCoeffs()
Update the boundary condition coefficients.
Generic GeometricField class.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
void correctBoundaryConditions()
Correct 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 void correctNut()
Correct the eddy-viscosity nut.
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
void boundEpsilonm(const volScalarField &Cc2)
Bound epsilonm.
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:103
Class to represent a interface between phases where one phase is considered dispersed within the othe...
Base class for Lagrangian drag models.
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
volScalarField nut_
Definition: eddyViscosity.H:60
Finite volume constraints.
Definition: fvConstraints.H:67
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
Finite volume models.
Definition: fvModels.H:65
tmp< fvMatrix< Type > > source(const VolField< Type > &field) const
Return source for an equation.
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:181
label index() const
Return the index of the phase.
Definition: phaseModel.C:157
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:145
virtual const volScalarField & rho() const =0
Return the density field.
Class to represent a system of phases.
Definition: phaseSystem.H:74
const phaseModel & otherPhase(const phaseModel &phase) const
Return the phase not given as an argument in a two-phase system.
Definition: phaseSystemI.H:174
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:253
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
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.
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)
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 > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
Namespace for OpenFOAM.
dimensionedScalar pos(const dimensionedScalar &ds)
dimensionedScalar exp(const dimensionedScalar &ds)
bool read(const char *, int32_t &)
Definition: int32IO.C:85
void dev(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void twoSymm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
const dimensionSet dimless
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimDensity
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
bool bound(volScalarField &, const dimensionedScalar &min)
Bound the given scalar field where it is below the specified min value.
Definition: bound.C:31
void pow3(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
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.