heRhoThermo.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) 2015-2022 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 "heRhoThermo.H"
30
31// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32
33template<class BasicPsiThermo, class MixtureType>
35(
36 const volScalarField& p,
37 volScalarField& T,
38 volScalarField& he,
39 volScalarField& psi,
40 volScalarField& rho,
41 volScalarField& mu,
42 volScalarField& alpha,
43 const bool doOldTimes
44)
45{
46 // Note: update oldTimes before current time so that if T.oldTime() is
47 // created from T, it starts from the unconverted T
48 if (doOldTimes && (p.nOldTimes() || T.nOldTimes()))
49 {
50 calculate
51 (
52 p.oldTime(),
53 T.oldTime(),
54 he.oldTime(),
55 psi.oldTime(),
56 rho.oldTime(),
57 mu.oldTime(),
58 alpha.oldTime(),
59 true
60 );
61 }
62
63 const scalarField& hCells = he.primitiveField();
64 const scalarField& pCells = p.primitiveField();
65
66 scalarField& TCells = T.primitiveFieldRef();
67 scalarField& psiCells = psi.primitiveFieldRef();
68 scalarField& rhoCells = rho.primitiveFieldRef();
69 scalarField& muCells = mu.primitiveFieldRef();
70 scalarField& alphaCells = alpha.primitiveFieldRef();
71
72 forAll(TCells, celli)
73 {
74 const typename MixtureType::thermoType& mixture_ =
75 this->cellMixture(celli);
76
77 if (this->updateT())
78 {
79 TCells[celli] = mixture_.THE
80 (
81 hCells[celli],
82 pCells[celli],
83 TCells[celli]
84 );
85 }
86
87 psiCells[celli] = mixture_.psi(pCells[celli], TCells[celli]);
88 rhoCells[celli] = mixture_.rho(pCells[celli], TCells[celli]);
89
90 muCells[celli] = mixture_.mu(pCells[celli], TCells[celli]);
91 alphaCells[celli] = mixture_.alphah(pCells[celli], TCells[celli]);
92 }
93
94 const volScalarField::Boundary& pBf = p.boundaryField();
95 volScalarField::Boundary& TBf = T.boundaryFieldRef();
96 volScalarField::Boundary& psiBf = psi.boundaryFieldRef();
97 volScalarField::Boundary& rhoBf = rho.boundaryFieldRef();
98 volScalarField::Boundary& heBf = he.boundaryFieldRef();
99 volScalarField::Boundary& muBf = mu.boundaryFieldRef();
100 volScalarField::Boundary& alphaBf = alpha.boundaryFieldRef();
101
102 forAll(pBf, patchi)
103 {
104 const fvPatchScalarField& pp = pBf[patchi];
105 fvPatchScalarField& pT = TBf[patchi];
106 fvPatchScalarField& ppsi = psiBf[patchi];
107 fvPatchScalarField& prho = rhoBf[patchi];
108 fvPatchScalarField& phe = heBf[patchi];
109 fvPatchScalarField& pmu = muBf[patchi];
110 fvPatchScalarField& palpha = alphaBf[patchi];
111
112 if (pT.fixesValue())
113 {
114 forAll(pT, facei)
115 {
116 const typename MixtureType::thermoType& mixture_ =
117 this->patchFaceMixture(patchi, facei);
118
119 phe[facei] = mixture_.HE(pp[facei], pT[facei]);
120
121 ppsi[facei] = mixture_.psi(pp[facei], pT[facei]);
122 prho[facei] = mixture_.rho(pp[facei], pT[facei]);
123 pmu[facei] = mixture_.mu(pp[facei], pT[facei]);
124 palpha[facei] = mixture_.alphah(pp[facei], pT[facei]);
125 }
126 }
127 else
128 {
129 forAll(pT, facei)
130 {
131 const typename MixtureType::thermoType& mixture_ =
132 this->patchFaceMixture(patchi, facei);
133
134 if (this->updateT())
135 {
136 pT[facei] = mixture_.THE(phe[facei], pp[facei], pT[facei]);
137 }
138
139 ppsi[facei] = mixture_.psi(pp[facei], pT[facei]);
140 prho[facei] = mixture_.rho(pp[facei], pT[facei]);
141 pmu[facei] = mixture_.mu(pp[facei], pT[facei]);
142 palpha[facei] = mixture_.alphah(pp[facei], pT[facei]);
143 }
144 }
145 }
146}
147
148// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
149
150template<class BasicPsiThermo, class MixtureType>
152(
153 const fvMesh& mesh,
154 const word& phaseName
155)
156:
157 heThermo<BasicPsiThermo, MixtureType>(mesh, phaseName)
158{
159 calculate
160 (
161 this->p_,
162 this->T_,
163 this->he_,
164 this->psi_,
165 this->rho_,
166 this->mu_,
167 this->alpha_,
168 true // Create old time fields
169 );
170}
171
172
173template<class BasicPsiThermo, class MixtureType>
175(
176 const fvMesh& mesh,
177 const word& phaseName,
178 const word& dictName
179)
180:
181 heThermo<BasicPsiThermo, MixtureType>(mesh, phaseName, dictName)
182{
183 calculate
184 (
185 this->p_,
186 this->T_,
187 this->he_,
188 this->psi_,
189 this->rho_,
190 this->mu_,
191 this->alpha_,
192 true // Create old time fields
193 );
194}
195
196
197// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
198
199template<class BasicPsiThermo, class MixtureType>
201{}
202
203
204// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
205
206template<class BasicPsiThermo, class MixtureType>
208{
210
211 calculate
212 (
213 this->p_,
214 this->T_,
215 this->he_,
216 this->psi_,
217 this->rho_,
218 this->mu_,
219 this->alpha_,
220 false // No need to update old times
221 );
222
223 DebugInFunction << "Finished" << endl;
224}
225
226// ************************************************************************* //
volScalarField & he
Definition: YEEqn.H:52
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Energy for a mixture based on density.
Definition: heRhoThermo.H:56
virtual ~heRhoThermo()
Destructor.
Definition: heRhoThermo.C:200
virtual void correct()
Update properties.
Definition: heRhoThermo.C:207
Enthalpy/Internal energy for a mixture.
Definition: heThermo.H:56
volScalarField he_
Energy field.
Definition: heThermo.H:62
ThermoType thermoType
The type of thermodynamics this mixture is instantiated for.
Definition: pureMixture.H:66
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
const volScalarField & psi
const volScalarField & T
const volScalarField & mu
dynamicFvMesh & mesh
const word dictName("faMeshDefinition")
#define DebugInFunction
Report an information message using Foam::Info.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
fvPatchField< scalar > fvPatchScalarField
volScalarField & alpha
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333