Loading...
Searching...
No Matches
diffusion.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) 2012-2017 OpenFOAM Foundation
9 Copyright (C) 2019-2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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 "diffusion.H"
30#include "fvcGrad.H"
31
32namespace Foam
34namespace combustionModels
35{
36
37// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38
39template<class ReactionThermo, class ThermoType>
40diffusion<ReactionThermo, ThermoType>::diffusion
41(
42 const word& modelType,
43 ReactionThermo& thermo,
45 const word& combustionProperties
46)
47:
48 singleStepCombustion<ReactionThermo, ThermoType>
49 (
50 modelType,
51 thermo,
52 turb,
53 combustionProperties
54 ),
55 C_(this->coeffs().getScalar("C")),
56 oxidantName_(this->coeffs().template getOrDefault<word>("oxidant", "O2"))
57{}
58
59
60// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
61
62template<class ReactionThermo, class ThermoType>
64{}
65
66
67// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
68
69template<class ReactionThermo, class ThermoType>
71{
72 this->wFuel_ == Zero;
73
74 if (this->active())
75 {
76 this->singleMixturePtr_->fresCorrect();
77
78 const label fuelI = this->singleMixturePtr_->fuelIndex();
79
80 const volScalarField& YFuel =
81 this->thermo().composition().Y()[fuelI];
82
83 if (this->thermo().composition().contains(oxidantName_))
84 {
85 const volScalarField& YO2 =
86 this->thermo().composition().Y(oxidantName_);
87
88 this->wFuel_ ==
89 C_*this->turbulence().muEff()
90 *mag(fvc::grad(YFuel) & fvc::grad(YO2))
91 *pos0(YFuel)*pos0(YO2);
92 }
93 }
94}
95
96
97template<class ReactionThermo, class ThermoType>
99{
101 {
102 this->coeffs().readEntry("C", C_);
103 this->coeffs().readIfPresent("oxidant", oxidantName_);
104 return true;
105 }
106
107 return false;
108}
109
110
111// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
112
113} // End namespace combustionModels
114} // End namespace Foam
115
116// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
compressible::turbulenceModel & turb
virtual ReactionThermo & thermo()
Return access to the thermo package.
const dictionary & coeffs() const
Return const dictionary of the model.
Switch active() const noexcept
Is combustion active?
virtual void correct()
Correct combustion rate.
Definition diffusion.C:63
virtual ~diffusion()
Destructor.
Definition diffusion.C:56
virtual bool read()
Update properties.
Definition diffusion.C:91
singleStepReactingMixture< ThermoType > * singleMixturePtr_
Pointer to singleStepReactingMixture mixture.
volScalarField wFuel_
Fuel consumption rate.
virtual bool read()
Update properties from given dictionary.
Abstract base class for turbulence models (RAS, LES and laminar).
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
scalar getScalar(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Same as get<scalar>(const word&, keyType::option).
A class for handling words, derived from Foam::string.
Definition word.H:66
basicSpecieMixture & composition
Calculate the gradient of the given field.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
Namespace for OpenFOAM.
dimensionedScalar pos0(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
psiReactionThermo & thermo