Loading...
Searching...
No Matches
PaSR.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-2017 OpenFOAM Foundation
9 Copyright (C) 2019-2023 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 "PaSR.H"
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
33template<class ReactionThermo>
34Foam::combustionModels::PaSR<ReactionThermo>::PaSR
35(
36 const word& modelType,
37 ReactionThermo& thermo,
39 const word& combustionProperties
40)
41:
42 laminar<ReactionThermo>(modelType, thermo, turb, combustionProperties),
43 Cmix_(this->coeffs().getScalar("Cmix")),
44 kappa_
45 (
47 (
48 thermo.phaseScopedName(typeName, "kappa"),
49 this->mesh().time().timeName(),
50 this->mesh(),
54 ),
55 this->mesh(),
57 )
58{}
59
60
61// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
62
63template<class ReactionThermo>
65{}
66
67
68// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
69
70template<class ReactionThermo>
72{
73 if (this->active())
74 {
76
77 tmp<volScalarField> tepsilon(this->turbulence().epsilon());
78 const auto& epsilon = tepsilon();
79
80 tmp<volScalarField> tmuEff(this->turbulence().muEff());
81 const auto& muEff = tmuEff();
82
83 tmp<volScalarField> ttc(this->tc());
84 const auto& tc = ttc();
85
87 const auto& rho = trho();
88
89 forAll(epsilon, i)
90 {
91 const scalar tk =
92 Cmix_*sqrt(max(muEff[i]/rho[i]/(epsilon[i] + SMALL), 0));
93
94 if (tk > SMALL)
95 {
96 kappa_[i] = tc[i]/(tc[i] + tk);
97 }
98 else
99 {
100 kappa_[i] = 1.0;
101 }
102 }
103
104 // Evaluate bcs
105 kappa_.correctBoundaryConditions();
106 }
107}
108
109
110template<class ReactionThermo>
114 return kappa_*laminar<ReactionThermo>::R(Y);
115}
116
117
118template<class ReactionThermo>
121{
123 (
124 this->thermo().phaseScopedName(typeName, "Qdot"),
127 );
128}
129
130
131template<class ReactionThermo>
133{
135 {
136 this->coeffs().readEntry("Cmix", Cmix_);
137 return true;
138 }
139
140 return false;
141}
142
143
144// ************************************************************************* //
compressible::turbulenceModel & turb
virtual ReactionThermo & thermo()
Return access to the thermo package.
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).
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
const dictionary & coeffs() const
Return const dictionary of the model.
const volScalarField & rho() const
Return const access to rho.
Switch active() const noexcept
Is combustion active?
const fvMesh & mesh() const
Return const access to the mesh database.
virtual void correct()
Correct combustion rate.
Definition PaSR.C:64
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
Definition PaSR.C:105
virtual ~PaSR()
Destructor.
Definition PaSR.C:57
virtual tmp< volScalarField > Qdot() const
Heat release rate [kg/m/s3].
Definition PaSR.C:113
virtual bool read()
Update properties from given dictionary.
Definition PaSR.C:125
tmp< volScalarField > tc() const
Return the chemical time scale.
Definition laminar.C:71
virtual void correct()
Correct combustion rate.
Definition laminar.C:78
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
Definition laminar.C:117
virtual tmp< volScalarField > Qdot() const
Heat release rate [kg/m/s3].
Definition laminar.C:136
virtual bool read()
Update properties from given dictionary.
Definition laminar.C:156
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).
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
PtrList< volScalarField > & Y
scalar epsilon
word timeName
Definition getTimeIndex.H:3
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
const dimensionSet dimless
Dimensionless.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
dimensionedScalar sqrt(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< volScalarField > trho
psiReactionThermo & thermo
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299