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 | Copyright (C) 2011-2016 OpenFOAM Foundation
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 "heRhoThermo.H"
27 
28 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
30 template<class BasicPsiThermo, class MixtureType>
32 {
33  const scalarField& hCells = this->he();
34  const scalarField& pCells = this->p_;
35 
36  scalarField& TCells = this->T_.primitiveFieldRef();
37  scalarField& psiCells = this->psi_.primitiveFieldRef();
38  scalarField& rhoCells = this->rho_.primitiveFieldRef();
39  scalarField& muCells = this->mu_.primitiveFieldRef();
40  scalarField& alphaCells = this->alpha_.primitiveFieldRef();
41 
42  forAll(TCells, celli)
43  {
44  const typename MixtureType::thermoType& mixture_ =
45  this->cellMixture(celli);
46 
47  TCells[celli] = mixture_.THE
48  (
49  hCells[celli],
50  pCells[celli],
51  TCells[celli]
52  );
53 
54  psiCells[celli] = mixture_.psi(pCells[celli], TCells[celli]);
55  rhoCells[celli] = mixture_.rho(pCells[celli], TCells[celli]);
56 
57  muCells[celli] = mixture_.mu(pCells[celli], TCells[celli]);
58  alphaCells[celli] = mixture_.alphah(pCells[celli], TCells[celli]);
59  }
60 
61  volScalarField::Boundary& pBf =
62  this->p_.boundaryFieldRef();
63 
64  volScalarField::Boundary& TBf =
65  this->T_.boundaryFieldRef();
66 
67  volScalarField::Boundary& psiBf =
68  this->psi_.boundaryFieldRef();
69 
70  volScalarField::Boundary& rhoBf =
71  this->rho_.boundaryFieldRef();
72 
73  volScalarField::Boundary& heBf =
74  this->he().boundaryFieldRef();
75 
76  volScalarField::Boundary& muBf =
77  this->mu_.boundaryFieldRef();
78 
79  volScalarField::Boundary& alphaBf =
80  this->alpha_.boundaryFieldRef();
81 
82  forAll(this->T_.boundaryField(), patchi)
83  {
84  fvPatchScalarField& pp = pBf[patchi];
85  fvPatchScalarField& pT = TBf[patchi];
86  fvPatchScalarField& ppsi = psiBf[patchi];
87  fvPatchScalarField& prho = rhoBf[patchi];
88  fvPatchScalarField& phe = heBf[patchi];
89  fvPatchScalarField& pmu = muBf[patchi];
90  fvPatchScalarField& palpha = alphaBf[patchi];
91 
92  if (pT.fixesValue())
93  {
94  forAll(pT, facei)
95  {
96  const typename MixtureType::thermoType& mixture_ =
97  this->patchFaceMixture(patchi, facei);
98 
99  phe[facei] = mixture_.HE(pp[facei], pT[facei]);
100 
101  ppsi[facei] = mixture_.psi(pp[facei], pT[facei]);
102  prho[facei] = mixture_.rho(pp[facei], pT[facei]);
103  pmu[facei] = mixture_.mu(pp[facei], pT[facei]);
104  palpha[facei] = mixture_.alphah(pp[facei], pT[facei]);
105  }
106  }
107  else
108  {
109  forAll(pT, facei)
110  {
111  const typename MixtureType::thermoType& mixture_ =
112  this->patchFaceMixture(patchi, facei);
113 
114  pT[facei] = mixture_.THE(phe[facei], pp[facei], pT[facei]);
115 
116  ppsi[facei] = mixture_.psi(pp[facei], pT[facei]);
117  prho[facei] = mixture_.rho(pp[facei], pT[facei]);
118  pmu[facei] = mixture_.mu(pp[facei], pT[facei]);
119  palpha[facei] = mixture_.alphah(pp[facei], pT[facei]);
120  }
121  }
122  }
123 }
124 
125 
126 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
127 
128 template<class BasicPsiThermo, class MixtureType>
130 (
131  const fvMesh& mesh,
132  const word& phaseName
133 )
134 :
136 {
137  calculate();
138 }
139 
140 
141 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
142 
143 template<class BasicPsiThermo, class MixtureType>
145 {}
146 
147 
148 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
149 
150 template<class BasicPsiThermo, class MixtureType>
152 {
153  if (debug)
154  {
155  InfoInFunction << endl;
156  }
157 
158  calculate();
159 
160  if (debug)
161  {
162  Info<< " Finished" << endl;
163  }
164 }
165 
166 
167 // ************************************************************************* //
volScalarField & he
Definition: YEEqn.H:56
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
virtual ~heRhoThermo()
Destructor.
Definition: heRhoThermo.C:144
virtual void correct()
Update properties.
Definition: heRhoThermo.C:151
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
dynamicFvMesh & mesh
fvPatchField< scalar > fvPatchScalarField
A class for handling words, derived from string.
Definition: word.H:59
volScalarField scalarField(fieldObject, mesh)
label patchi
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Energy for a mixture based on density.
Definition: heRhoThermo.H:50
messageStream Info
#define InfoInFunction
Report an information message using Foam::Info.