psiuReactionThermo.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-2020 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 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  defineTypeNameAndDebug(psiuReactionThermo, 0);
38  defineRunTimeSelectionTable(psiuReactionThermo, fvMesh);
39 }
40 
41 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
42 
44 {
45  const volScalarField::Boundary& tbf =
46  this->Tu().boundaryField();
47 
48  wordList hbt = tbf.types();
49 
50  forAll(tbf, patchi)
51  {
52  if (isA<fixedValueFvPatchScalarField>(tbf[patchi]))
53  {
54  hbt[patchi] = fixedUnburntEnthalpyFvPatchScalarField::typeName;
55  }
56  else if
57  (
58  isA<zeroGradientFvPatchScalarField>(tbf[patchi])
59  || isA<fixedGradientFvPatchScalarField>(tbf[patchi])
60  )
61  {
62  hbt[patchi] = gradientUnburntEnthalpyFvPatchScalarField::typeName;
63  }
64  else if (isA<mixedFvPatchScalarField>(tbf[patchi]))
65  {
66  hbt[patchi] = mixedUnburntEnthalpyFvPatchScalarField::typeName;
67  }
68  }
69 
70  return hbt;
71 }
72 
74 {
76 
77  forAll(hbf, patchi)
78  {
79  if
80  (
81  isA<gradientUnburntEnthalpyFvPatchScalarField>(hbf[patchi])
82  )
83  {
84  refCast<gradientUnburntEnthalpyFvPatchScalarField>(hbf[patchi])
85  .gradient() = hbf[patchi].fvPatchField::snGrad();
86  }
87  else if
88  (
89  isA<mixedUnburntEnthalpyFvPatchScalarField>(hbf[patchi])
90  )
91  {
92  refCast<mixedUnburntEnthalpyFvPatchScalarField>(hbf[patchi])
93  .refGrad() = hbf[patchi].fvPatchField::snGrad();
94  }
95  }
96 }
97 
98 
99 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
100 
102 (
103  const fvMesh& mesh,
104  const word& phaseName
105 )
106 {}
107 
108 
109 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
110 
112 (
113  const fvMesh& mesh,
114  const word& phaseName
115 )
116 {
117  return basicThermo::New<psiuReactionThermo>(mesh, phaseName);
118 }
119 
120 
121 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
122 
124 {}
125 
126 
128 {}
129 
130 
131 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
wordList types() const
Return a list of the patch field types.
virtual ~psiuReactionThermo()
Destructor.
wordList heuBoundaryTypes()
Return the unburnt enthalpy/internal energy field boundary types.
fvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
implementation(const fvMesh &, const word &phaseName)
Construct from mesh and phase name.
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
defineTypeNameAndDebug(combustionModel, 0)
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:95
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)
Standard selection based on fvMesh.
Namespace for OpenFOAM.