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 -------------------------------------------------------------------------------
11 License
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 
33 template<class ReactionThermo>
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  (
46  IOobject
47  (
48  thermo.phasePropertyName(typeName + ":kappa"),
49  this->mesh().time().timeName(),
50  this->mesh(),
51  IOobject::NO_READ,
52  IOobject::AUTO_WRITE
53  ),
54  this->mesh(),
56  )
57 {}
58 
59 
60 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
61 
62 template<class ReactionThermo>
64 {}
65 
66 
67 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
68 
69 template<class ReactionThermo>
71 {
72  if (this->active())
73  {
75 
76  tmp<volScalarField> tepsilon(this->turbulence().epsilon());
77  const scalarField& epsilon = tepsilon();
78 
79  tmp<volScalarField> tmuEff(this->turbulence().muEff());
80  const scalarField& muEff = tmuEff();
81 
82  tmp<volScalarField> ttc(this->tc());
83  const scalarField& tc = ttc();
84 
85  tmp<volScalarField> trho(this->rho());
86  const scalarField& rho = trho();
87 
88  forAll(epsilon, i)
89  {
90  const scalar tk =
91  Cmix_*sqrt(max(muEff[i]/rho[i]/(epsilon[i] + SMALL), 0));
92 
93  if (tk > SMALL)
94  {
95  kappa_[i] = tc[i]/(tc[i] + tk);
96  }
97  else
98  {
99  kappa_[i] = 1.0;
100  }
101  }
102 
103  // Evaluate bcs
104  kappa_.correctBoundaryConditions();
105  }
106 }
107 
108 
109 template<class ReactionThermo>
112 {
113  return kappa_*laminar<ReactionThermo>::R(Y);
114 }
115 
116 
117 template<class ReactionThermo>
120 {
121  return tmp<volScalarField>
122  (
123  new volScalarField
124  (
125  this->thermo().phasePropertyName(typeName + ":Qdot"),
127  )
128  );
129 }
130 
131 
132 template<class ReactionThermo>
134 {
136  {
137  this->coeffs().readEntry("Cmix", Cmix_);
138  return true;
139  }
140 
141  return false;
142 }
143 
144 
145 // ************************************************************************* //
virtual void correct()
Correct combustion rate.
Definition: laminar.C:78
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
Definition: PaSR.C:104
virtual tmp< volScalarField > Qdot() const
Heat release rate [kg/m/s3].
Definition: PaSR.C:112
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
virtual bool read()
Update properties from given dictionary.
Definition: PaSR.C:126
compressible::turbulenceModel & turb
dimensionedScalar sqrt(const dimensionedScalar &ds)
tmp< volScalarField > trho
const dimensionSet dimless
Dimensionless.
Partially stirred reactor turbulent combustion model.
Definition: PaSR.H:56
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
psiReactionThermo & thermo
Definition: createFields.H:28
word timeName
Definition: getTimeIndex.H:3
Laminar combustion model.
Definition: laminar.H:52
dynamicFvMesh & mesh
A class for handling words, derived from Foam::string.
Definition: word.H:63
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
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
virtual void correct()
Correct combustion rate.
Definition: PaSR.C:63
PtrList< volScalarField > & Y
scalar epsilon
virtual ~PaSR()
Destructor.
Definition: PaSR.C:56
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
Definition: laminar.C:117
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127