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-------------------------------------------------------------------------------
10License
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
29namespace Foam
30{
31namespace combustionModels
32{
33
34// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35
36template<class ReactionThermo, class ThermoType>
37eddyDissipationDiffusionModel<ReactionThermo, ThermoType>::
38eddyDissipationDiffusionModel
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
59template<class ReactionThermo, class ThermoType>
62{}
63
64
65// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
66
67
68template<class ReactionThermo, class ThermoType>
72 return (max(this->rtTurb(), this->rtDiff()));
73}
74
75
76template<class ReactionThermo, class ThermoType>
81 (
83 (
85 (
86 "tdelta",
87 this->mesh().time().timeName(),
88 this->mesh(),
91 ),
92 this->mesh(),
94 zeroGradientFvPatchScalarField::typeName
95 )
96 );
97
98 volScalarField& delta = tdelta.ref();
99 delta.ref() = cbrt(this->mesh().V());
100 delta.correctBoundaryConditions();
101
102 // NOTE: Assume Prt = 1
103 return Cd_*this->turbulence().nuEff()/sqr(delta);
104}
105
106
107template<class ReactionThermo, class ThermoType>
111 {
112 this->coeffs().readEntry("Cd", Cd_);
113 return true;
114 }
115
116 return false;
117}
118
119
120// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
121
122} // End namespace combustionModels
123} // End namespace Foam
124
125// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
scalar delta
compressible::turbulenceModel & turb
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
Eddy dissipation model based on the principle of mixed is burnt.
tmp< volScalarField > rtDiff() const
Return the reciprocal of the diffusion time scale.
virtual tmp< volScalarField > timeScale()
Calculate time scale.
Standard Eddy Dissipation Model based on the assumption that the reaction rates are controlled by the...
Abstract base class for turbulence models (RAS, LES and laminar).
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
A class for handling words, derived from Foam::string.
Definition: word.H:68
dynamicFvMesh & mesh
compressible::turbulenceModel & turbulence
word timeName
Definition: getTimeIndex.H:3
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
dimensionedScalar cbrt(const dimensionedScalar &ds)