Loading...
Searching...
No Matches
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{
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
65// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
66
67
68template<class ReactionThermo, class ThermoType>
71{
72 return (max(this->rtTurb(), this->rtDiff()));
73}
74
75
76template<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
98template<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// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
scalar delta
compressible::turbulenceModel & turb
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())
@ NO_REGISTER
Do not request registration (bool: false).
virtual ReactionThermo & thermo()
Return access to the thermo package.
const dictionary & coeffs() const
Return const dictionary of the model.
tmp< volScalarField > rtDiff() const
Return the reciprocal of the diffusion time scale.
virtual tmp< volScalarField > timeScale()
Calculate time scale.
tmp< volScalarField > rtTurb() const
Return the reciprocal of the turbulent mixing time scale.
Abstract base class for turbulence models (RAS, LES and laminar).
scalar getScalar(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Same as get<scalar>(const word&, keyType::option).
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
dynamicFvMesh & mesh
Namespace for OpenFOAM.
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)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensionedScalar cbrt(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.