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