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-2020 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 "heRhoThermo.H"
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 template<class BasicPsiThermo, class MixtureType>
35 (
36  const volScalarField& p,
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  const typename MixtureType::thermoType& volMixture_ =
78  this->cellVolMixture(pCells[celli], TCells[celli], celli);
79 
80  if (this->updateT())
81  {
82  TCells[celli] = mixture_.THE
83  (
84  hCells[celli],
85  pCells[celli],
86  TCells[celli]
87  );
88  }
89 
90  psiCells[celli] = mixture_.psi(pCells[celli], TCells[celli]);
91  rhoCells[celli] = volMixture_.rho(pCells[celli], TCells[celli]);
92 
93  muCells[celli] = mixture_.mu(pCells[celli], TCells[celli]);
94  alphaCells[celli] = mixture_.alphah(pCells[celli], TCells[celli]);
95  }
96 
97  const volScalarField::Boundary& pBf = p.boundaryField();
98  volScalarField::Boundary& TBf = T.boundaryFieldRef();
99  volScalarField::Boundary& psiBf = psi.boundaryFieldRef();
100  volScalarField::Boundary& rhoBf = rho.boundaryFieldRef();
101  volScalarField::Boundary& heBf = he.boundaryFieldRef();
102  volScalarField::Boundary& muBf = mu.boundaryFieldRef();
103  volScalarField::Boundary& alphaBf = alpha.boundaryFieldRef();
104 
105  forAll(pBf, patchi)
106  {
107  const fvPatchScalarField& pp = pBf[patchi];
108  fvPatchScalarField& pT = TBf[patchi];
109  fvPatchScalarField& ppsi = psiBf[patchi];
110  fvPatchScalarField& prho = rhoBf[patchi];
111  fvPatchScalarField& phe = heBf[patchi];
112  fvPatchScalarField& pmu = muBf[patchi];
113  fvPatchScalarField& palpha = alphaBf[patchi];
114 
115  if (pT.fixesValue())
116  {
117  forAll(pT, facei)
118  {
119  const typename MixtureType::thermoType& mixture_ =
120  this->patchFaceMixture(patchi, facei);
121 
122  const typename MixtureType::thermoType& volMixture_ =
123  this->patchFaceVolMixture
124  (
125  pp[facei],
126  pT[facei],
127  patchi,
128  facei
129  );
130 
131  phe[facei] = mixture_.HE(pp[facei], pT[facei]);
132 
133  ppsi[facei] = mixture_.psi(pp[facei], pT[facei]);
134  prho[facei] = volMixture_.rho(pp[facei], pT[facei]);
135  pmu[facei] = mixture_.mu(pp[facei], pT[facei]);
136  palpha[facei] = mixture_.alphah(pp[facei], pT[facei]);
137  }
138  }
139  else
140  {
141  forAll(pT, facei)
142  {
143  const typename MixtureType::thermoType& mixture_ =
144  this->patchFaceMixture(patchi, facei);
145 
146  const typename MixtureType::thermoType& volMixture_ =
147  this->patchFaceVolMixture
148  (
149  pp[facei],
150  pT[facei],
151  patchi,
152  facei
153  );
154 
155  if (this->updateT())
156  {
157  pT[facei] = mixture_.THE(phe[facei], pp[facei], pT[facei]);
158  }
159 
160  ppsi[facei] = mixture_.psi(pp[facei], pT[facei]);
161  prho[facei] = volMixture_.rho(pp[facei], pT[facei]);
162  pmu[facei] = mixture_.mu(pp[facei], pT[facei]);
163  palpha[facei] = mixture_.alphah(pp[facei], pT[facei]);
164  }
165  }
166  }
167 }
168 
169 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
170 
171 template<class BasicPsiThermo, class MixtureType>
173 (
174  const fvMesh& mesh,
175  const word& phaseName
176 )
177 :
179 {
180  calculate
181  (
182  this->p_,
183  this->T_,
184  this->he_,
185  this->psi_,
186  this->rho_,
187  this->mu_,
188  this->alpha_,
189  true // Create old time fields
190  );
191 }
192 
193 
194 template<class BasicPsiThermo, class MixtureType>
196 (
197  const fvMesh& mesh,
198  const word& phaseName,
199  const word& dictName
200 )
201 :
203 {
204  calculate
205  (
206  this->p_,
207  this->T_,
208  this->he_,
209  this->psi_,
210  this->rho_,
211  this->mu_,
212  this->alpha_,
213  true // Create old time fields
214  );
215 }
216 
217 
218 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
219 
220 template<class BasicPsiThermo, class MixtureType>
222 {}
223 
224 
225 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
226 
227 template<class BasicPsiThermo, class MixtureType>
229 {
231 
232  calculate
233  (
234  this->p_,
235  this->T_,
236  this->he_,
237  this->psi_,
238  this->rho_,
239  this->mu_,
240  this->alpha_,
241  false // No need to update old times
242  );
243 
244  DebugInFunction << "Finished" << endl;
245 }
246 
247 // ************************************************************************* //
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
dictName
const word dictName("faMeshDefinition")
Foam::heRhoThermo
Energy for a mixture based on density.
Definition: heRhoThermo.H:53
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
rho
rho
Definition: readInitialConditions.H:88
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
heRhoThermo.H
Foam::heRhoThermo::~heRhoThermo
virtual ~heRhoThermo()
Destructor.
Definition: heRhoThermo.C:221
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
he
volScalarField & he
Definition: YEEqn.H:52
psi
const volScalarField & psi
Definition: createFieldRefs.H:1
Foam::heRhoThermo::correct
virtual void correct()
Update properties.
Definition: heRhoThermo.C:228