PaSR.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2021 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "PaSR.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace combustionModels
34 {
35  defineTypeNameAndDebug(PaSR, 0);
36  addToRunTimeSelectionTable(combustionModel, PaSR, dictionary);
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
44 (
45  const word& modelType,
46  const fluidReactionThermo& thermo,
48  const word& combustionProperties
49 )
50 :
51  laminar(modelType, thermo, turb, combustionProperties),
52  Cmix_(this->coeffs().template lookup<scalar>("Cmix")),
53  kappa_
54  (
55  IOobject
56  (
57  thermo.phasePropertyName(typeName + ":kappa"),
58  this->mesh().time().timeName(),
59  this->mesh(),
62  ),
63  this->mesh(),
65  )
66 {}
67 
68 
69 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
70 
72 {}
73 
74 
75 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
76 
78 {
80 
81  tmp<volScalarField> tepsilon(this->turbulence().epsilon());
82  const scalarField& epsilon = tepsilon();
83 
84  tmp<volScalarField> tnuEff(this->turbulence().nuEff());
85  const scalarField& nuEff = tnuEff();
86 
87  tmp<volScalarField> ttc(this->chemistryPtr_->tc());
88  const scalarField& tc = ttc();
89 
90  forAll(epsilon, i)
91  {
92  const scalar tk =
93  Cmix_*sqrt(max(nuEff[i]/(epsilon[i] + small), 0));
94 
95  if (tk > small)
96  {
97  kappa_[i] = tc[i]/(tc[i] + tk);
98  }
99  else
100  {
101  kappa_[i] = 1.0;
102  }
103  }
104 }
105 
106 
109 {
110  return kappa_*laminar::R(Y);
111 }
112 
113 
116 {
117  return volScalarField::New
118  (
119  this->thermo().phasePropertyName(typeName + ":Qdot"),
120  kappa_*laminar::Qdot()
121  );
122 }
123 
124 
126 {
127  if (laminar::read())
128  {
129  this->coeffs().lookup("Cmix") >> Cmix_;
130  return true;
131  }
132  else
133  {
134  return false;
135  }
136 }
137 
138 
139 // ************************************************************************* //
virtual void correct()
Correct combustion rate.
Definition: laminar.C:94
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
fluidReactionThermo & thermo
Definition: createFields.H:28
static word phasePropertyName(const word &name, const word &phaseName)
Name of a property for a given phase.
Definition: basicThermo.H:150
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
Definition: PaSR.C:108
virtual tmp< volScalarField > Qdot() const
Heat release rate [kg/m/s^3].
Definition: PaSR.C:115
virtual bool read()
Update properties from given dictionary.
Definition: PaSR.C:125
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
dimensionedScalar sqrt(const dimensionedScalar &ds)
Base-class for multi-component fluid thermodynamic properties.
const dimensionSet dimless
fvMesh & mesh
virtual tmp< volScalarField > Qdot() const
Heat release rate [kg/m/s^3].
Definition: laminar.C:139
Macros for easy insertion into run-time selection tables.
Laminar combustion model.
Definition: laminar.H:52
A class for handling words, derived from string.
Definition: word.H:59
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(dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::momentumTransportModel > turbulence(compressible::momentumTransportModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
virtual void correct()
Correct combustion rate.
Definition: PaSR.C:77
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
PtrList< volScalarField > & Y
scalar epsilon
virtual ~PaSR()
Destructor.
Definition: PaSR.C:71
PaSR(const word &modelType, const fluidReactionThermo &thermo, const compressibleMomentumTransportModel &turb, const word &combustionProperties)
Construct from components.
Definition: PaSR.C:44
A class for managing temporary objects.
Definition: PtrList.H:53
defineTypeNameAndDebug(diffusion, 0)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
Definition: laminar.C:127
Base class for single-phase compressible turbulence models.
addToRunTimeSelectionTable(combustionModel, diffusion, dictionary)
Namespace for OpenFOAM.
virtual bool read()
Update properties from given dictionary.
Definition: laminar.C:145