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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2017 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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 "psiuReactionThermo.H"
30#include "fvMesh.H"
35
36// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37
38namespace Foam
39{
40
41/* * * * * * * * * * * * * * * private static data * * * * * * * * * * * * * */
42
46
47// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
48
50{
51 const volScalarField::Boundary& tbf =
52 this->Tu().boundaryField();
53
54 wordList hbt = tbf.types();
55
56 forAll(tbf, patchi)
57 {
58 if (isA<fixedValueFvPatchScalarField>(tbf[patchi]))
59 {
61 }
62 else if
63 (
64 isA<zeroGradientFvPatchScalarField>(tbf[patchi])
65 || isA<fixedGradientFvPatchScalarField>(tbf[patchi])
66 )
67 {
69 }
70 else if (isA<mixedFvPatchScalarField>(tbf[patchi]))
71 {
73 }
74 }
75
76 return hbt;
77}
78
80{
82
83 forAll(hbf, patchi)
84 {
85 if
86 (
87 isA<gradientUnburntEnthalpyFvPatchScalarField>(hbf[patchi])
88 )
89 {
90 refCast<gradientUnburntEnthalpyFvPatchScalarField>(hbf[patchi])
91 .gradient() = hbf[patchi].fvPatchField::snGrad();
92 }
93 else if
94 (
95 isA<mixedUnburntEnthalpyFvPatchScalarField>(hbf[patchi])
96 )
97 {
98 refCast<mixedUnburntEnthalpyFvPatchScalarField>(hbf[patchi])
99 .refGrad() = hbf[patchi].fvPatchField::snGrad();
100 }
101 }
102}
103
104
105// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
106
108(
109 const fvMesh& mesh,
110 const word& phaseName
111)
112:
113 psiReactionThermo(mesh, phaseName)
114{}
115
116
118(
119 const fvMesh& mesh,
120 const word& phaseName,
121 const word& dictName
122)
123:
124 psiReactionThermo(mesh, phaseName, dictName)
125{}
126
127
128// * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
129
131(
132 const fvMesh& mesh,
133 const word& phaseName
134)
135{
136 return basicThermo::New<psiuReactionThermo>(mesh, phaseName);
137}
138
139
141(
142 const fvMesh& mesh,
143 const word& phaseName,
144 const word& dictName
145)
146{
147 return basicThermo::New<psiuReactionThermo>(mesh, phaseName, dictName);
148}
149
150
151// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
152
154{}
155
156
157// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
158
159} // End namespace Foam
160
161// ************************************************************************* //
wordList types() const
Return a list of the patch types.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Foam::psiReactionThermo.
Foam::psiuReactionThermo.
void heuBoundaryCorrection(volScalarField &heu)
virtual ~psiuReactionThermo()
Destructor.
virtual volScalarField & heu()=0
Unburnt gas enthalpy [J/kg].
virtual const volScalarField & Tu() const =0
Unburnt gas temperature [K].
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
dynamicFvMesh & mesh
const word dictName("faMeshDefinition")
Namespace for OpenFOAM.
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
static const char *const typeName
The type name used in ensight case files.