All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
heRhoQGDThermo.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  grpRhoQGDThermo
31 
32 \*---------------------------------------------------------------------------*/
33 
34 #include "heRhoQGDThermo.H"
35 
36 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
37 
38 template<class BasicPsiThermo, class MixtureType>
40 {
41  const scalarField& hCells = this->he();
42  const scalarField& pCells = this->p_;
43 
44  scalarField& TCells = this->T_.primitiveFieldRef();
45  scalarField& psiCells = this->psi_.primitiveFieldRef();
46  scalarField& rhoCells = this->rho_.primitiveFieldRef();
47  scalarField& muCells = this->mu_.primitiveFieldRef();
48  scalarField& alphaCells = this->alpha_.primitiveFieldRef();
49 
50  forAll(TCells, celli)
51  {
52  const typename MixtureType::thermoType& mixture_ =
53  this->cellMixture(celli);
54 
55  TCells[celli] = mixture_.THE
56  (
57  hCells[celli],
58  pCells[celli],
59  TCells[celli]
60  );
61 
62  psiCells[celli] = mixture_.psi(pCells[celli], TCells[celli]);
63  rhoCells[celli] = mixture_.rho(pCells[celli], TCells[celli]);
64 
65  muCells[celli] = mixture_.mu(pCells[celli], TCells[celli]);
66  alphaCells[celli] = mixture_.alphah(pCells[celli], TCells[celli]);
67  }
68 
69  volScalarField::Boundary& pBf =
70  this->p_.boundaryFieldRef();
71 
72  volScalarField::Boundary& TBf =
73  this->T_.boundaryFieldRef();
74 
75  volScalarField::Boundary& psiBf =
76  this->psi_.boundaryFieldRef();
77 
78  volScalarField::Boundary& rhoBf =
79  this->rho_.boundaryFieldRef();
80 
81  volScalarField::Boundary& heBf =
82  this->he().boundaryFieldRef();
83 
84  volScalarField::Boundary& muBf =
85  this->mu_.boundaryFieldRef();
86 
87  volScalarField::Boundary& alphaBf =
88  this->alpha_.boundaryFieldRef();
89 
90 
91  forAll(this->T_.boundaryField(), patchi)
92  {
93  fvPatchScalarField& pp = pBf[patchi];
94  fvPatchScalarField& pT = TBf[patchi];
95  fvPatchScalarField& ppsi = psiBf[patchi];
96  fvPatchScalarField& prho = rhoBf[patchi];
97  fvPatchScalarField& phe = heBf[patchi];
98  fvPatchScalarField& pmu = muBf[patchi];
99  fvPatchScalarField& palpha = alphaBf[patchi];
100 
101  if (pT.fixesValue())
102  {
103  forAll(pT, facei)
104  {
105  const typename MixtureType::thermoType& mixture_ =
106  this->patchFaceMixture(patchi, facei);
107 
108  phe[facei] = mixture_.HE(pp[facei], pT[facei]);
109 
110  ppsi[facei] = mixture_.psi(pp[facei], pT[facei]);
111  prho[facei] = mixture_.rho(pp[facei], pT[facei]);
112  pmu[facei] = mixture_.mu(pp[facei], pT[facei]);
113  palpha[facei] = mixture_.alphah(pp[facei], pT[facei]);
114  }
115  }
116  else
117  {
118  forAll(pT, facei)
119  {
120  const typename MixtureType::thermoType& mixture_ =
121  this->patchFaceMixture(patchi, facei);
122 
123  pT[facei] = mixture_.THE(phe[facei], pp[facei], pT[facei]);
124 
125  ppsi[facei] = mixture_.psi(pp[facei], pT[facei]);
126  prho[facei] = mixture_.rho(pp[facei], pT[facei]);
127  pmu[facei] = mixture_.mu(pp[facei], pT[facei]);
128  palpha[facei] = mixture_.alphah(pp[facei], pT[facei]);
129  }
130  }
131  }
132 
133  if (!this->isochoric())
134  {
135  this->gamma_ == (this->Cp() / this->Cv());
136  this->c_ = sqrt(this->gamma_ / this->psi());
137  }
138  this->correctQGD(this->mu_, this->alpha_);
139 }
140 
141 
142 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
143 
144 template<class BasicPsiThermo, class MixtureType>
146 (
147  const fvMesh& mesh,
148  const word& phaseName
149 )
150 :
151  heThermo<BasicPsiThermo, MixtureType>(mesh, phaseName)
152 {
153  calculate();
154  // Switch on saving old time
155  this->psi_.oldTime();
156 }
157 
158 template<class BasicPsiThermo, class MixtureType>
160 (
161  const fvMesh& mesh,
162  const word& phaseName,
163  const word& dictionaryName
164 )
165 :
166  heThermo<BasicPsiThermo, MixtureType>(mesh, phaseName, dictionaryName)
167 {
168  calculate();
169  // Switch on saving old time
170  this->psi_.oldTime();
171 }
172 
173 
174 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
175 
176 template<class BasicPsiThermo, class MixtureType>
178 {}
179 
180 
181 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
182 
183 template<class BasicPsiThermo, class MixtureType>
185 {
186  if (debug)
187  {
188  InfoInFunction << endl;
189  }
190 
191  // force the saving of the old-time values
192  this->psi_.oldTime();
193 
194  calculate();
195 
196  if (debug)
197  {
198  Info<< " Finished" << endl;
199  }
200 }
201 
202 //template<class BasicPsiThermo, class MixtureType>
203 //Foam::tmp<Foam::volScalarField> Foam::heRhoQGDThermo<BasicPsiThermo, MixtureType>::gamma() const
204 //{
205 // return this->gamma_;
206 //}
207 
208 
209 // ************************************************************************* //
virtual ~heRhoQGDThermo()
Destructor.
forAll(Y, i)
Definition: QGDYEqn.H:36
const volScalarField & psi
Definition: createFields.H:20
virtual void correct()
Update properties.