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