eddyDissipationDiffusionModel.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) 2016-2019 OpenCFD Ltd
9 -------------------------------------------------------------------------------
10 License
11 
12  OpenFOAM is free software: you can redistribute it and/or modify it
13  under the terms of the GNU General Public License as published by
14  the Free Software Foundation, either version 3 of the License, or
15  (at your option) any later version.
16 
17  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
18  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
19  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
20  for more details.
21 
22  You should have received a copy of the GNU General Public License
23  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
24 
25 \*---------------------------------------------------------------------------*/
26 
28 
29 namespace Foam
30 {
31 namespace combustionModels
32 {
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
36 template<class ReactionThermo, class ThermoType>
39 (
40  const word& modelType,
41  ReactionThermo& thermo,
43  const word& combustionProperties
44 )
45 :
46  eddyDissipationModelBase<ReactionThermo, ThermoType>
47  (
48  modelType,
49  thermo,
50  turb,
51  combustionProperties
52  ),
53  Cd_(this->coeffs().getScalar("Cd"))
54 {}
55 
56 
57 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
58 
59 template<class ReactionThermo, class ThermoType>
62 {}
63 
64 
65 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
66 
67 
68 template<class ReactionThermo, class ThermoType>
71 {
72  return (max(this->rtTurb(), this->rtDiff()));
73 }
74 
75 
76 template<class ReactionThermo, class ThermoType>
79 {
80  auto tdelta = volScalarField::New
81  (
82  "tdelta",
84  this->mesh(),
87  );
88  auto& delta = tdelta.ref();
89 
90  delta.ref() = cbrt(this->mesh().V());
91  delta.correctBoundaryConditions();
92 
93  // NOTE: Assume Prt = 1
94  return Cd_*this->turbulence().nuEff()/sqr(delta);
95 }
96 
97 
98 template<class ReactionThermo, class ThermoType>
100 {
102  {
103  this->coeffs().readEntry("Cd", Cd_);
104  return true;
105  }
106 
107  return false;
108 }
109 
110 
111 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
112 
113 } // End namespace combustionModels
114 } // End namespace Foam
115 
116 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:221
scalar delta
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
compressible::turbulenceModel & turb
dynamicFvMesh & mesh
A class for handling words, derived from Foam::string.
Definition: word.H:63
dimensionedScalar cbrt(const dimensionedScalar &ds)
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
Abstract base class for turbulence models (RAS, LES and laminar).
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
Standard Eddy Dissipation Model based on the assumption that the reaction rates are controlled by the...
tmp< volScalarField > rtDiff() const
Return the reciprocal of the diffusion time scale.
Eddy dissipation model based on the principle of mixed is burnt.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
virtual tmp< volScalarField > timeScale()
Calculate time scale.
Do not request registration (bool: false)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127