FreeFOAM The Cross-Platform CFD Toolkit
hhuCombustionThermo.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) 1991-2010 OpenCFD Ltd.
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 "hhuCombustionThermo.H"
27 #include <finiteVolume/fvMesh.H>
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 
38 /* * * * * * * * * * * * * * * private static data * * * * * * * * * * * * * */
39 
40 defineTypeNameAndDebug(hhuCombustionThermo, 0);
41 defineRunTimeSelectionTable(hhuCombustionThermo, fvMesh);
42 
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 
46 {
47  const volScalarField::GeometricBoundaryField& tbf = Tu_.boundaryField();
48 
49  wordList hbt = tbf.types();
50 
51  forAll(tbf, patchi)
52  {
53  if (isA<fixedValueFvPatchScalarField>(tbf[patchi]))
54  {
55  hbt[patchi] = fixedUnburntEnthalpyFvPatchScalarField::typeName;
56  }
57  else if
58  (
59  isA<zeroGradientFvPatchScalarField>(tbf[patchi])
60  || isA<fixedGradientFvPatchScalarField>(tbf[patchi])
61  )
62  {
63  hbt[patchi] = gradientUnburntEnthalpyFvPatchScalarField::typeName;
64  }
65  else if (isA<mixedFvPatchScalarField>(tbf[patchi]))
66  {
67  hbt[patchi] = mixedUnburntEnthalpyFvPatchScalarField::typeName;
68  }
69  }
70 
71  return hbt;
72 }
73 
75 {
76  volScalarField::GeometricBoundaryField& hbf = hu.boundaryField();
77 
78  forAll(hbf, patchi)
79  {
80  if
81  (
82  isA<gradientUnburntEnthalpyFvPatchScalarField>(hbf[patchi])
83  )
84  {
85  refCast<gradientUnburntEnthalpyFvPatchScalarField>(hbf[patchi])
86  .gradient() = hbf[patchi].fvPatchField::snGrad();
87  }
88  else if
89  (
90  isA<mixedUnburntEnthalpyFvPatchScalarField>(hbf[patchi])
91  )
92  {
93  refCast<mixedUnburntEnthalpyFvPatchScalarField>(hbf[patchi])
94  .refGrad() = hbf[patchi].fvPatchField::snGrad();
95  }
96  }
97 }
98 
99 
100 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
101 
103 :
104  hCombustionThermo(mesh),
105 
106  Tu_
107  (
108  IOobject
109  (
110  "Tu",
111  mesh.time().timeName(),
112  mesh,
113  IOobject::MUST_READ,
114  IOobject::AUTO_WRITE
115  ),
116  mesh
117  ),
118 
119  hu_
120  (
121  IOobject
122  (
123  "hu",
124  mesh.time().timeName(),
125  mesh,
126  IOobject::NO_READ,
127  IOobject::NO_WRITE
128  ),
129  mesh,
130  dimensionSet(0, 2, -2, 0, 0),
131  huBoundaryTypes()
132  )
133 {}
134 
135 
136 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
137 
139 {}
140 
141 
142 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
143 
144 } // End namespace Foam
145 
146 // ************************ vim: set sw=4 sts=4 et: ************************ //