psiuReactionThermo.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 | Copyright (C) 2011-2016 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 "psiuReactionThermo.H"
27 #include "fvMesh.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 
38 /* * * * * * * * * * * * * * * private static data * * * * * * * * * * * * * */
39 
40 defineTypeNameAndDebug(psiuReactionThermo, 0);
41 defineRunTimeSelectionTable(psiuReactionThermo, fvMesh);
42 
43 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
44 
46 {
47  const volScalarField::Boundary& tbf =
48  this->Tu().boundaryField();
49 
50  wordList hbt = tbf.types();
51 
52  forAll(tbf, patchi)
53  {
54  if (isA<fixedValueFvPatchScalarField>(tbf[patchi]))
55  {
56  hbt[patchi] = fixedUnburntEnthalpyFvPatchScalarField::typeName;
57  }
58  else if
59  (
60  isA<zeroGradientFvPatchScalarField>(tbf[patchi])
61  || isA<fixedGradientFvPatchScalarField>(tbf[patchi])
62  )
63  {
64  hbt[patchi] = gradientUnburntEnthalpyFvPatchScalarField::typeName;
65  }
66  else if (isA<mixedFvPatchScalarField>(tbf[patchi]))
67  {
68  hbt[patchi] = mixedUnburntEnthalpyFvPatchScalarField::typeName;
69  }
70  }
71 
72  return hbt;
73 }
74 
76 {
77  volScalarField::Boundary& hbf = heu.boundaryFieldRef();
78 
79  forAll(hbf, patchi)
80  {
81  if
82  (
83  isA<gradientUnburntEnthalpyFvPatchScalarField>(hbf[patchi])
84  )
85  {
86  refCast<gradientUnburntEnthalpyFvPatchScalarField>(hbf[patchi])
87  .gradient() = hbf[patchi].fvPatchField::snGrad();
88  }
89  else if
90  (
91  isA<mixedUnburntEnthalpyFvPatchScalarField>(hbf[patchi])
92  )
93  {
94  refCast<mixedUnburntEnthalpyFvPatchScalarField>(hbf[patchi])
95  .refGrad() = hbf[patchi].fvPatchField::snGrad();
96  }
97  }
98 }
99 
100 
101 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
102 
104 (
105  const fvMesh& mesh,
106  const word& phaseName
107 )
108 :
109  psiReactionThermo(mesh, phaseName)
110 {}
111 
112 
113 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
114 
116 (
117  const fvMesh& mesh,
118  const word& phaseName
119 )
120 {
121  return basicThermo::New<psiuReactionThermo>(mesh, phaseName);
122 }
123 
124 
125 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
126 
128 {}
129 
130 
131 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
132 
133 } // End namespace Foam
134 
135 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
const Boundary & boundaryField() const
Return const-reference to the boundary field.
virtual ~psiuReactionThermo()
Destructor.
virtual const volScalarField & Tu() const =0
Unburnt gas temperature [K].
psiReactionThermo(const fvMesh &, const word &phaseName)
Construct from mesh and phase name.
psiuReactionThermo(const fvMesh &, const word &phaseName)
Construct from mesh and phase name.
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
defineTypeNameAndDebug(combustionModel, 0)
virtual volScalarField & heu()=0
Unburnt gas enthalpy [J/kg].
label patchi
void heuBoundaryCorrection(volScalarField &heu)
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
static autoPtr< psiuReactionThermo > New(const fvMesh &, const word &phaseName=word::null)
Namespace for OpenFOAM.