All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
hePsiQGDThermo.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) 2019 OpenCFD Ltd.
10  Copyright (C) 2016-2019 ISP RAS (www.ispras.ru) UniCFD Group (www.unicfd.ru)
11 -------------------------------------------------------------------------------
12 
13 License
14  This file is part of QGDsolver, based on OpenFOAM library.
15 
16  OpenFOAM is free software: you can redistribute it and/or modify it
17  under the terms of the GNU General Public License as published by
18  the Free Software Foundation, either version 3 of the License, or
19  (at your option) any later version.
20 
21  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
22  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
23  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
24  for more details.
25 
26  You should have received a copy of the GNU General Public License
27  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
28 
29 Group
30  grpPsiQGDThermo
31 \*---------------------------------------------------------------------------*/
32 
33 #include "hePsiQGDThermo.H"
34 
35 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36 
37 template<class BasicPsiThermo, class MixtureType>
39 {
40  const scalarField& hCells = this->he_;
41  const scalarField& pCells = this->p_;
42 
43  scalarField& TCells = this->T_.primitiveFieldRef();
44  scalarField& psiCells = this->psi_.primitiveFieldRef();
45  scalarField& muCells = this->mu_.primitiveFieldRef();
46  scalarField& alphaCells = this->alpha_.primitiveFieldRef();
47 
48  forAll(TCells, celli)
49  {
50  const typename MixtureType::thermoType& mixture_ =
51  this->cellMixture(celli);
52 
53  TCells[celli] = mixture_.THE
54  (
55  hCells[celli],
56  pCells[celli],
57  TCells[celli]
58  );
59 
60  psiCells[celli] = mixture_.psi(pCells[celli], TCells[celli]);
61 
62  muCells[celli] = mixture_.mu(pCells[celli], TCells[celli]);
63  alphaCells[celli] = mixture_.alphah(pCells[celli], TCells[celli]);
64  }
65 
66  volScalarField::Boundary& pBf =
67  this->p_.boundaryFieldRef();
68 
69  volScalarField::Boundary& TBf =
70  this->T_.boundaryFieldRef();
71 
72  volScalarField::Boundary& psiBf =
73  this->psi_.boundaryFieldRef();
74 
75  volScalarField::Boundary& heBf =
76  this->he().boundaryFieldRef();
77 
78  volScalarField::Boundary& muBf =
79  this->mu_.boundaryFieldRef();
80 
81  volScalarField::Boundary& alphaBf =
82  this->alpha_.boundaryFieldRef();
83 
84  forAll(this->T_.boundaryField(), patchi)
85  {
86  fvPatchScalarField& pp = pBf[patchi];
87  fvPatchScalarField& pT = TBf[patchi];
88  fvPatchScalarField& ppsi = psiBf[patchi];
89  fvPatchScalarField& phe = heBf[patchi];
90  fvPatchScalarField& pmu = muBf[patchi];
91  fvPatchScalarField& palpha = alphaBf[patchi];
92 
93  if (pT.fixesValue())
94  {
95  forAll(pT, facei)
96  {
97  const typename MixtureType::thermoType& mixture_ =
98  this->patchFaceMixture(patchi, facei);
99 
100  phe[facei] = mixture_.HE(pp[facei], pT[facei]);
101 
102  ppsi[facei] = mixture_.psi(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  pmu[facei] = mixture_.mu(pp[facei], pT[facei]);
118  palpha[facei] = mixture_.alphah(pp[facei], pT[facei]);
119  }
120  }
121  }
122 
123  this->gamma_ == (this->Cp() / this->Cv());
124  this->c_ = sqrt(this->gamma_ / this->psi());
125  this->correctQGD(this->mu_, this->alpha_);
126 }
127 
128 
129 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
130 
131 template<class BasicPsiThermo, class MixtureType>
133 (
134  const fvMesh& mesh,
135  const word& phaseName
136 )
137 :
138  heThermo<BasicPsiThermo, MixtureType>(mesh, phaseName)
139 {
140  calculate();
141 
142  // Switch on saving old time
143  this->psi_.oldTime();
144 }
145 
146 
147 template<class BasicPsiThermo, class MixtureType>
149 (
150  const fvMesh& mesh,
151  const word& phaseName,
152  const word& dictionaryName
153 )
154 :
155  heThermo<BasicPsiThermo, MixtureType>(mesh, phaseName, dictionaryName)
156 {
157  calculate();
158 
159  // Switch on saving old time
160  this->psi_.oldTime();
161 }
162 
163 
164 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
165 
166 template<class BasicPsiThermo, class MixtureType>
168 {}
169 
170 
171 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
172 
173 template<class BasicPsiThermo, class MixtureType>
175 {
176  if (debug)
177  {
178  InfoInFunction << endl;
179  }
180 
181  // force the saving of the old-time values
182  this->psi_.oldTime();
183 
184  calculate();
185 
186  if (debug)
187  {
188  Info<< " Finished" << endl;
189  }
190 }
191 
192 template<class BasicPsiThermo, class MixtureType>
193 Foam::tmp<Foam::volScalarField> Foam::hePsiQGDThermo<BasicPsiThermo, MixtureType>::gamma() const
194 {
195  return this->gamma_;
196 }
197 
198 
199 // ************************************************************************* //
virtual ~hePsiQGDThermo()
Destructor.
forAll(Y, i)
Definition: QGDYEqn.H:36
const volScalarField & psi
Definition: createFields.H:20
virtual void correct()
Update properties.
virtual tmp< volScalarField > gamma() const
Energy for a perfect gas mixture with QGD equations.