psiuMulticomponentThermo.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-2023 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 
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
38 }
39 
41 (
42  "heheuPsiThermo"
43 );
44 
45 
46 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
47 
49 {
50  const volScalarField::Boundary& tbf =
51  this->Tu().boundaryField();
52 
53  wordList hbt = tbf.types();
54 
55  forAll(tbf, patchi)
56  {
57  if (isA<fixedValueFvPatchScalarField>(tbf[patchi]))
58  {
59  hbt[patchi] = fixedUnburntEnthalpyFvPatchScalarField::typeName;
60  }
61  else if
62  (
63  isA<zeroGradientFvPatchScalarField>(tbf[patchi])
64  || isA<fixedGradientFvPatchScalarField>(tbf[patchi])
65  )
66  {
67  hbt[patchi] = gradientUnburntEnthalpyFvPatchScalarField::typeName;
68  }
69  else if (isA<mixedFvPatchScalarField>(tbf[patchi]))
70  {
71  hbt[patchi] = mixedUnburntEnthalpyFvPatchScalarField::typeName;
72  }
73  }
74 
75  return hbt;
76 }
77 
79 {
81 
82  forAll(hbf, patchi)
83  {
84  if
85  (
86  isA<gradientUnburntEnthalpyFvPatchScalarField>(hbf[patchi])
87  )
88  {
89  refCast<gradientUnburntEnthalpyFvPatchScalarField>(hbf[patchi])
90  .gradient() = hbf[patchi].fvPatchField::snGrad();
91  }
92  else if
93  (
94  isA<mixedUnburntEnthalpyFvPatchScalarField>(hbf[patchi])
95  )
96  {
97  refCast<mixedUnburntEnthalpyFvPatchScalarField>(hbf[patchi])
98  .refGrad() = hbf[patchi].fvPatchField::snGrad();
99  }
100  }
101 }
102 
103 
104 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
105 
107 (
108  const dictionary& dict,
109  const wordList& specieNames,
110  const fvMesh& mesh,
111  const word& phaseName
112 )
113 :
114  species_(specieNames),
115  Y_(species_.size())
116 {
117  forAll(species_, i)
118  {
119  Y_.set
120  (
121  i,
122  new volScalarField
123  (
124  IOobject
125  (
127  mesh.time().name(),
128  mesh,
131  ),
132  mesh
133  )
134  );
135  }
136 }
137 
138 
139 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
140 
143 (
144  const fvMesh& mesh,
145  const word& phaseName
146 )
147 {
148  return basicThermo::New<psiuMulticomponentThermo>(mesh, phaseName);
149 }
150 
151 
152 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
153 
155 {}
156 
157 
159 {}
160 
161 
162 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
163 
165 {
166  return p()*psiu();
167 }
168 
169 
171 {
172  return p()*psib();
173 }
174 
175 
176 const Foam::speciesTable&
178 {
179  return species_;
180 }
181 
182 
185 {
186  return Y_;
187 }
188 
189 
192 {
193  return Y_;
194 }
195 
196 
197 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Generic GeometricBoundaryField class.
wordList types() const
Return a list of the patch field types.
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
static word groupName(Name name, const word &group)
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
virtual const fvMesh & mesh() const =0
Return const access to the mesh.
virtual const word & phaseName() const =0
Phase name.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
const word & name() const
Return const reference to name.
virtual const volScalarField & p() const =0
Pressure [Pa].
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:418
A wordList with hashed indices for faster lookup by name.
virtual const speciesTable & species() const
The table of species.
implementation(const dictionary &, const wordList &, const fvMesh &, const word &)
Construct from dictionary, specie names, mesh and phase name.
virtual PtrList< volScalarField > & Y()
Access the mass-fraction fields.
PtrList< volScalarField > Y_
Species mass fractions.
Base-class for combustion fluid thermodynamic properties based on compressibility.
tmp< volScalarField > rhob() const
Burnt gas density [kg/m^3].
static const word derivedThermoName
The derived name.
void heuBoundaryCorrection(volScalarField &heu)
...
wordList heuBoundaryTypes()
Return the unburnt enthalpy/internal energy field boundary types.
virtual ~psiuMulticomponentThermo()
Destructor.
virtual tmp< volScalarField > psib() const =0
Burnt gas compressibility [s^2/m^2].
virtual const volScalarField & Tu() const =0
Unburnt gas temperature [K].
virtual tmp< volScalarField > psiu() const =0
Unburnt gas compressibility [s^2/m^2].
tmp< volScalarField > rhou() const
Unburnt gas density [kg/m^3].
static autoPtr< psiuMulticomponentThermo > New(const fvMesh &, const word &phaseName=word::null)
Standard selection based on fvMesh.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
label patchi
Namespace for OpenFOAM.
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
defineTypeNameAndDebug(combustionModel, 0)
dictionary dict