kineticTheoryModel.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2018 OpenFOAM Foundation
9  Copyright (C) 2019-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "kineticTheoryModel.H"
30 #include "mathematicalConstants.H"
31 #include "twoPhaseSystem.H"
32 #include "fvOptions.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
36 Foam::RASModels::kineticTheoryModel::kineticTheoryModel
37 (
38  const volScalarField& alpha,
39  const volScalarField& rho,
40  const volVectorField& U,
41  const surfaceScalarField& alphaRhoPhi,
42  const surfaceScalarField& phi,
43  const transportModel& phase,
44  const word& propertiesName,
45  const word& type
46 )
47 :
48  eddyViscosity
49  <
51  >
52  (
53  type,
54  alpha,
55  rho,
56  U,
57  alphaRhoPhi,
58  phi,
59  phase,
60  propertiesName
61  ),
62 
63  phase_(phase),
64 
65  viscosityModel_
66  (
67  kineticTheoryModels::viscosityModel::New
68  (
69  coeffDict_
70  )
71  ),
72  conductivityModel_
73  (
74  kineticTheoryModels::conductivityModel::New
75  (
76  coeffDict_
77  )
78  ),
79  radialModel_
80  (
81  kineticTheoryModels::radialModel::New
82  (
83  coeffDict_
84  )
85  ),
86  granularPressureModel_
87  (
88  kineticTheoryModels::granularPressureModel::New
89  (
90  coeffDict_
91  )
92  ),
93  frictionalStressModel_
94  (
95  kineticTheoryModels::frictionalStressModel::New
96  (
97  coeffDict_
98  )
99  ),
100 
101  equilibrium_(coeffDict_.get<bool>("equilibrium")),
102  e_("e", dimless, coeffDict_),
103  alphaMax_("alphaMax", dimless, coeffDict_),
104  alphaMinFriction_("alphaMinFriction", dimless, coeffDict_),
105  residualAlpha_("residualAlpha", dimless, coeffDict_),
106  maxNut_("maxNut", dimViscosity, 1000, coeffDict_),
107 
108  Theta_
109  (
110  IOobject
111  (
112  IOobject::groupName("Theta", phase.name()),
113  U.time().timeName(),
114  U.mesh(),
115  IOobject::MUST_READ,
116  IOobject::AUTO_WRITE
117  ),
118  U.mesh()
119  ),
120 
121  lambda_
122  (
123  IOobject
124  (
125  IOobject::groupName("lambda", phase.name()),
126  U.time().timeName(),
127  U.mesh(),
128  IOobject::NO_READ,
129  IOobject::NO_WRITE
130  ),
131  U.mesh(),
133  ),
134 
135  gs0_
136  (
137  IOobject
138  (
139  IOobject::groupName("gs0", phase.name()),
140  U.time().timeName(),
141  U.mesh(),
142  IOobject::NO_READ,
143  IOobject::NO_WRITE
144  ),
145  U.mesh(),
147  ),
148 
149  kappa_
150  (
151  IOobject
152  (
153  IOobject::groupName("kappa", phase.name()),
154  U.time().timeName(),
155  U.mesh(),
156  IOobject::NO_READ,
157  IOobject::NO_WRITE
158  ),
159  U.mesh(),
161  ),
162 
163  nuFric_
164  (
165  IOobject
166  (
167  IOobject::groupName("nuFric", phase.name()),
168  U.time().timeName(),
169  U.mesh(),
170  IOobject::NO_READ,
171  IOobject::AUTO_WRITE
172  ),
173  U.mesh(),
175  )
176 {
177  if (type == typeName)
178  {
179  printCoeffs(type);
180  }
181 }
182 
183 
184 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
185 
187 {}
188 
189 
190 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
191 
193 {
194  if
195  (
197  <
199  >::read()
200  )
201  {
202  coeffDict().readEntry("equilibrium", equilibrium_);
203  e_.readIfPresent(coeffDict());
204  alphaMax_.readIfPresent(coeffDict());
205  alphaMinFriction_.readIfPresent(coeffDict());
206 
207  viscosityModel_->read();
208  conductivityModel_->read();
209  radialModel_->read();
210  granularPressureModel_->read();
211  frictionalStressModel_->read();
212 
213  return true;
214  }
215 
216  return false;
217 }
218 
219 
222 {
224  return nut_;
225 }
226 
227 
230 {
232  return nut_;
233 }
234 
235 
238 {
240  return nullptr;
241 }
242 
243 
246 {
248  (
250  (
251  IOobject
252  (
253  IOobject::groupName("R", U_.group()),
254  runTime_.timeName(),
255  mesh_,
258  ),
259  - (nut_)*dev(twoSymm(fvc::grad(U_)))
260  - (lambda_*fvc::div(phi_))*symmTensor::I
261  )
262  );
263 }
264 
265 
268 {
269  const volScalarField& rho = phase_.rho();
270 
271  tmp<volScalarField> tpPrime
272  (
273  Theta_
274  *granularPressureModel_->granularPressureCoeffPrime
275  (
276  alpha_,
277  radialModel_->g0(alpha_, alphaMinFriction_, alphaMax_),
278  radialModel_->g0prime(alpha_, alphaMinFriction_, alphaMax_),
279  rho,
280  e_
281  )
282  + frictionalStressModel_->frictionalPressurePrime
283  (
284  phase_,
285  alphaMinFriction_,
286  alphaMax_
287  )
288  );
289 
290  volScalarField::Boundary& bpPrime =
291  tpPrime.ref().boundaryFieldRef();
292 
293  forAll(bpPrime, patchi)
294  {
295  if (!bpPrime[patchi].coupled())
296  {
297  bpPrime[patchi] == 0;
298  }
299  }
300 
301  return tpPrime;
302 }
303 
304 
307 {
308  return fvc::interpolate(pPrime());
309 }
310 
311 
314 {
316  (
318  (
319  IOobject
320  (
321  IOobject::groupName("devRhoReff", U_.group()),
322  runTime_.timeName(),
323  mesh_,
326  ),
327  - (rho_*nut_)
328  *dev(twoSymm(fvc::grad(U_)))
329  - ((rho_*lambda_)*fvc::div(phi_))*symmTensor::I
330  )
331  );
332 }
333 
334 
337 (
339 ) const
340 {
341  return
342  (
343  - fvm::laplacian(rho_*nut_, U)
344  - fvc::div
345  (
346  (rho_*nut_)*dev2(T(fvc::grad(U)))
347  + ((rho_*lambda_)*fvc::div(phi_))
349  )
350  );
351 }
352 
353 
355 {
356  // Local references
357  volScalarField alpha(max(alpha_, scalar(0)));
358  const volScalarField& rho = phase_.rho();
359  const surfaceScalarField& alphaRhoPhi = alphaRhoPhi_;
360  const volVectorField& U = U_;
361  const volVectorField& Uc_ =
362  refCast<const twoPhaseSystem>(phase_.fluid()).otherPhase(phase_).U();
363 
364  const scalar sqrtPi = sqrt(constant::mathematical::pi);
365  dimensionedScalar ThetaSmall("ThetaSmall", Theta_.dimensions(), 1e-6);
366  dimensionedScalar ThetaSmallSqrt(sqrt(ThetaSmall));
367 
368  tmp<volScalarField> tda(phase_.d());
369  const volScalarField& da = tda();
370 
371  tmp<volTensorField> tgradU(fvc::grad(U_));
372  const volTensorField& gradU(tgradU());
373  volSymmTensorField D(symm(gradU));
374 
375  // Calculating the radial distribution function
376  gs0_ = radialModel_->g0(alpha, alphaMinFriction_, alphaMax_);
377 
378  if (!equilibrium_)
379  {
380  // Particle viscosity (Table 3.2, p.47)
381  nut_ = viscosityModel_->nu(alpha, Theta_, gs0_, rho, da, e_);
382 
383  volScalarField ThetaSqrt("sqrtTheta", sqrt(Theta_));
384 
385  // Bulk viscosity p. 45 (Lun et al. 1984).
386  lambda_ = (4.0/3.0)*sqr(alpha)*da*gs0_*(1 + e_)*ThetaSqrt/sqrtPi;
387 
388  // Stress tensor, Definitions, Table 3.1, p. 43
390  (
391  rho*(2*nut_*D + (lambda_ - (2.0/3.0)*nut_)*tr(D)*I)
392  );
393 
394  // Dissipation (Eq. 3.24, p.50)
395  volScalarField gammaCoeff
396  (
397  "gammaCoeff",
398  12*(1 - sqr(e_))
399  *max(sqr(alpha), residualAlpha_)
400  *rho*gs0_*(1.0/da)*ThetaSqrt/sqrtPi
401  );
402 
403  // Drag
405  (
406  refCast<const twoPhaseSystem>(phase_.fluid()).Kd()
407  );
408 
409  // Eq. 3.25, p. 50 Js = J1 - J2
410  volScalarField J1("J1", 3*beta);
411  volScalarField J2
412  (
413  "J2",
414  0.25*sqr(beta)*da*magSqr(U - Uc_)
415  /(
416  max(alpha, residualAlpha_)*rho
417  *sqrtPi*(ThetaSqrt + ThetaSmallSqrt)
418  )
419  );
420 
421  // particle pressure - coefficient in front of Theta (Eq. 3.22, p. 45)
422  volScalarField PsCoeff
423  (
424  granularPressureModel_->granularPressureCoeff
425  (
426  alpha,
427  gs0_,
428  rho,
429  e_
430  )
431  );
432 
433  // 'thermal' conductivity (Table 3.3, p. 49)
434  kappa_ = conductivityModel_->kappa(alpha, Theta_, gs0_, rho, da, e_);
435 
437 
438  // Construct the granular temperature equation (Eq. 3.20, p. 44)
439  // NB. note that there are two typos in Eq. 3.20:
440  // Ps should be without grad
441  // the laplacian has the wrong sign
442  fvScalarMatrix ThetaEqn
443  (
444  1.5*
445  (
446  fvm::ddt(alpha, rho, Theta_)
447  + fvm::div(alphaRhoPhi, Theta_)
448  - fvc::Sp(fvc::ddt(alpha, rho) + fvc::div(alphaRhoPhi), Theta_)
449  )
450  - fvm::laplacian(kappa_, Theta_, "laplacian(kappa,Theta)")
451  ==
452  - fvm::SuSp((PsCoeff*I) && gradU, Theta_)
453  + (tau && gradU)
454  + fvm::Sp(-gammaCoeff, Theta_)
455  + fvm::Sp(-J1, Theta_)
456  + fvm::Sp(J2/(Theta_ + ThetaSmall), Theta_)
457  + fvOptions(alpha, rho, Theta_)
458  );
459 
460  ThetaEqn.relax();
461  fvOptions.constrain(ThetaEqn);
462  ThetaEqn.solve();
463  fvOptions.correct(Theta_);
464  }
465  else
466  {
467  // Equilibrium => dissipation == production
468  // Eq. 4.14, p.82
469  volScalarField K1("K1", 2*(1 + e_)*rho*gs0_);
471  (
472  "K3",
473  0.5*da*rho*
474  (
475  (sqrtPi/(3*(3.0 - e_)))
476  *(1 + 0.4*(1 + e_)*(3*e_ - 1)*alpha*gs0_)
477  +1.6*alpha*gs0_*(1 + e_)/sqrtPi
478  )
479  );
480 
482  (
483  "K2",
484  4*da*rho*(1 + e_)*alpha*gs0_/(3*sqrtPi) - 2*K3/3.0
485  );
486 
487  volScalarField K4("K4", 12*(1 - sqr(e_))*rho*gs0_/(da*sqrtPi));
488 
489  volScalarField trD
490  (
491  "trD",
492  alpha/(alpha + residualAlpha_)
493  *fvc::div(phi_)
494  );
495  volScalarField tr2D("tr2D", sqr(trD));
496  volScalarField trD2("trD2", tr(D & D));
497 
498  volScalarField t1("t1", K1*alpha + rho);
499  volScalarField l1("l1", -t1*trD);
500  volScalarField l2("l2", sqr(t1)*tr2D);
501  volScalarField l3
502  (
503  "l3",
504  4.0
505  *K4
506  *alpha
507  *(2*K3*trD2 + K2*tr2D)
508  );
509 
510  Theta_ = sqr
511  (
512  (l1 + sqrt(l2 + l3))
513  /(2*max(alpha, residualAlpha_)*K4)
514  );
515 
516  kappa_ = conductivityModel_->kappa(alpha, Theta_, gs0_, rho, da, e_);
517  }
518 
519  Theta_.max(0);
520  Theta_.min(100);
521 
522  {
523  // particle viscosity (Table 3.2, p.47)
524  nut_ = viscosityModel_->nu(alpha, Theta_, gs0_, rho, da, e_);
525 
526  volScalarField ThetaSqrt("sqrtTheta", sqrt(Theta_));
527 
528  // Bulk viscosity p. 45 (Lun et al. 1984).
529  lambda_ = (4.0/3.0)*sqr(alpha)*da*gs0_*(1 + e_)*ThetaSqrt/sqrtPi;
530 
531  // Frictional pressure
532  volScalarField pf
533  (
534  frictionalStressModel_->frictionalPressure
535  (
536  phase_,
537  alphaMinFriction_,
538  alphaMax_
539  )
540  );
541 
542  nuFric_ = frictionalStressModel_->nu
543  (
544  phase_,
545  alphaMinFriction_,
546  alphaMax_,
547  pf/rho,
548  D
549  );
550 
551  // Limit viscosity and add frictional viscosity
552  nut_.min(maxNut_);
553  nuFric_ = min(nuFric_, maxNut_ - nut_);
554  nut_ += nuFric_;
555  }
556 
557  if (debug)
558  {
559  Info<< typeName << ':' << nl
560  << " max(Theta) = " << max(Theta_).value() << nl
561  << " max(nut) = " << max(nut_).value() << endl;
562  }
563 }
564 
565 
566 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::RASModels::kineticTheoryModel::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kineticTheoryModel.C:192
Foam::symm
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:84
mathematicalConstants.H
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
fvOptions.H
K1
#define K1
Definition: SHA1.C:144
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::RASModels::kineticTheoryModel::k
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
Definition: kineticTheoryModel.C:221
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::RASModel
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:51
Foam::fv::options::New
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:104
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
K4
#define K4
Definition: SHA1.C:147
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::phaseCompressibleTurbulenceModel
ThermalDiffusivity< PhaseCompressibleTurbulenceModel< phaseModel > > phaseCompressibleTurbulenceModel
Typedef for phaseCompressibleTurbulenceModel.
Definition: phaseCompressibleTurbulenceModel.H:47
Foam::RASModels::kineticTheoryModel::pPrimef
virtual tmp< surfaceScalarField > pPrimef() const
Return the face-phase-pressure'.
Definition: kineticTheoryModel.C:306
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Kd
const volScalarField Kd(fluid.Kd())
Foam::dev2
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:117
K3
#define K3
Definition: SHA1.C:146
rho
rho
Definition: readInitialConditions.H:88
Foam::GeometricField::min
void min(const dimensioned< Type > &dt)
Use the minimum of the field and specified value.
Definition: GeometricField.C:1132
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::EddyDiffusivity< phaseCompressibleTurbulenceModel >
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
kineticTheoryModel.H
NotImplemented
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:445
Foam::fvm::SuSp
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:228
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:23
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
Foam::SymmTensor::I
static const SymmTensor I
Definition: SymmTensor.H:78
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::fvMatrix::relax
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:559
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
beta
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
Foam::RASModels::kineticTheoryModel::correct
virtual void correct()
Solve the kinetic theory equations and correct the viscosity.
Definition: kineticTheoryModel.C:354
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
K2
#define K2
Definition: SHA1.C:145
Foam::fv::options
Finite-volume options.
Definition: fvOptions.H:55
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::dimViscosity
const dimensionSet dimViscosity
Foam::dimDynamicViscosity
const dimensionSet dimDynamicViscosity
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:43
Foam::fvc::interpolate
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.
Foam::RASModels::kineticTheoryModel::epsilon
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
Definition: kineticTheoryModel.C:229
Foam::RASModels::kineticTheoryModel::pPrime
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure'.
Definition: kineticTheoryModel.C:267
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:62
U
U
Definition: pEqn.H:72
Foam::RASModels::kineticTheoryModel::devRhoReff
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
Definition: kineticTheoryModel.C:313
Foam::RASModels::kineticTheoryModel::~kineticTheoryModel
virtual ~kineticTheoryModel()
Destructor.
Definition: kineticTheoryModel.C:186
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::compressible::RASModel
RASModel< EddyDiffusivity< turbulenceModel > > RASModel
Definition: turbulentFluidThermoModel.H:66
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
bool
bool
Definition: EEqn.H:20
D
const dimensionedScalar & D
Definition: solveBulkSurfactant.H:4
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:76
Foam::RASModel< EddyDiffusivity< phaseCompressibleTurbulenceModel > >::printCoeffs
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: RASModel.C:34
Foam::fvMatrix::solve
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Definition: fvMatrixSolve.C:301
Foam::eddyViscosity
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:55
Foam::roots::type
type
Types of root.
Definition: Roots.H:54
Foam::RASModels::kineticTheoryModel::omega
virtual tmp< volScalarField > omega() const
Return the specific dissipation rate.
Definition: kineticTheoryModel.C:237
Foam::fvc::ddt
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:47
Foam::tr
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:51
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
Foam::fvc::Sp
tmp< GeometricField< Type, fvPatchField, volMesh > > Sp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:69
Foam::RASModels::kineticTheoryModel::R
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
Definition: kineticTheoryModel.C:245
Foam::RASModels::kineticTheoryModel::divDevRhoReff
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
Definition: kineticTheoryModel.C:337
Foam::GeometricField< symmTensor, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
Foam::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:95
Foam::I
static const Identity< scalar > I
Definition: Identity.H:95
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:106