heRhoThermo.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 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 BasicRhoThermo, class MixtureType>
32 {
33  const scalarField& hCells = this->he();
34  const scalarField& pCells = this->p_;
35 
36  scalarField& TCells = this->T_.primitiveFieldRef();
37  scalarField& CpCells = this->Cp_.primitiveFieldRef();
38  scalarField& CvCells = this->Cv_.primitiveFieldRef();
39  scalarField& psiCells = this->psi_.primitiveFieldRef();
40  scalarField& rhoCells = this->rho_.primitiveFieldRef();
41  scalarField& muCells = this->mu_.primitiveFieldRef();
42  scalarField& alphaCells = this->alpha_.primitiveFieldRef();
43 
44  forAll(TCells, celli)
45  {
46  const typename MixtureType::thermoMixtureType& thermoMixture =
47  this->cellThermoMixture(celli);
48 
49  const typename MixtureType::transportMixtureType& transportMixture =
50  this->cellTransportMixture(celli, thermoMixture);
51 
52  TCells[celli] = thermoMixture.THE
53  (
54  hCells[celli],
55  pCells[celli],
56  TCells[celli]
57  );
58 
59  CpCells[celli] = thermoMixture.Cp(pCells[celli], TCells[celli]);
60  CvCells[celli] = thermoMixture.Cv(pCells[celli], TCells[celli]);
61  psiCells[celli] = thermoMixture.psi(pCells[celli], TCells[celli]);
62  rhoCells[celli] = thermoMixture.rho(pCells[celli], TCells[celli]);
63 
64  muCells[celli] = transportMixture.mu(pCells[celli], TCells[celli]);
65  alphaCells[celli] =
66  transportMixture.kappa(pCells[celli], TCells[celli])
67  /thermoMixture.Cp(pCells[celli], TCells[celli]);
68  }
69 
70  volScalarField::Boundary& pBf =
71  this->p_.boundaryFieldRef();
72 
73  volScalarField::Boundary& TBf =
74  this->T_.boundaryFieldRef();
75 
76  volScalarField::Boundary& CpBf =
77  this->Cp_.boundaryFieldRef();
78 
79  volScalarField::Boundary& CvBf =
80  this->Cv_.boundaryFieldRef();
81 
82  volScalarField::Boundary& psiBf =
83  this->psi_.boundaryFieldRef();
84 
85  volScalarField::Boundary& rhoBf =
86  this->rho_.boundaryFieldRef();
87 
88  volScalarField::Boundary& heBf =
89  this->he().boundaryFieldRef();
90 
91  volScalarField::Boundary& muBf =
92  this->mu_.boundaryFieldRef();
93 
94  volScalarField::Boundary& alphaBf =
95  this->alpha_.boundaryFieldRef();
96 
97  forAll(this->T_.boundaryField(), patchi)
98  {
99  fvPatchScalarField& pp = pBf[patchi];
100  fvPatchScalarField& pT = TBf[patchi];
101  fvPatchScalarField& pCp = CpBf[patchi];
102  fvPatchScalarField& pCv = CvBf[patchi];
103  fvPatchScalarField& ppsi = psiBf[patchi];
104  fvPatchScalarField& prho = rhoBf[patchi];
105  fvPatchScalarField& phe = heBf[patchi];
106  fvPatchScalarField& pmu = muBf[patchi];
107  fvPatchScalarField& palpha = alphaBf[patchi];
108 
109  if (pT.fixesValue())
110  {
111  forAll(pT, facei)
112  {
113  const typename MixtureType::thermoMixtureType& thermoMixture =
114  this->patchFaceThermoMixture(patchi, facei);
115 
116  const typename MixtureType::transportMixtureType&
117  transportMixture =
118  this->patchFaceTransportMixture
119  (patchi, facei, thermoMixture);
120 
121  phe[facei] = thermoMixture.HE(pp[facei], pT[facei]);
122 
123  pCp[facei] = thermoMixture.Cp(pp[facei], pT[facei]);
124  pCv[facei] = thermoMixture.Cv(pp[facei], pT[facei]);
125  ppsi[facei] = thermoMixture.psi(pp[facei], pT[facei]);
126  prho[facei] = thermoMixture.rho(pp[facei], pT[facei]);
127 
128  pmu[facei] = transportMixture.mu(pp[facei], pT[facei]);
129  palpha[facei] =
130  transportMixture.kappa(pp[facei], pT[facei])
131  /thermoMixture.Cp(pp[facei], pT[facei]);
132  }
133  }
134  else
135  {
136  forAll(pT, facei)
137  {
138  const typename MixtureType::thermoMixtureType& thermoMixture =
139  this->patchFaceThermoMixture(patchi, facei);
140 
141  const typename MixtureType::transportMixtureType&
142  transportMixture =
143  this->patchFaceTransportMixture
144  (patchi, facei, thermoMixture);
145 
146  pT[facei] = thermoMixture.THE(phe[facei], pp[facei], pT[facei]);
147 
148  pCp[facei] = thermoMixture.Cp(pp[facei], pT[facei]);
149  pCv[facei] = thermoMixture.Cv(pp[facei], pT[facei]);
150  ppsi[facei] = thermoMixture.psi(pp[facei], pT[facei]);
151  prho[facei] = thermoMixture.rho(pp[facei], pT[facei]);
152 
153  pmu[facei] = transportMixture.mu(pp[facei], pT[facei]);
154  palpha[facei] =
155  transportMixture.kappa(pp[facei], pT[facei])
156  /thermoMixture.Cp(pp[facei], pT[facei]);
157  }
158  }
159  }
160 }
161 
162 
163 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
164 
165 template<class BasicRhoThermo, class MixtureType>
167 (
168  const fvMesh& mesh,
169  const word& phaseName
170 )
171 :
173 {
174  calculate();
175 }
176 
177 
178 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
179 
180 template<class BasicRhoThermo, class MixtureType>
182 {}
183 
184 
185 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
186 
187 template<class BasicRhoThermo, class MixtureType>
189 {
190  if (debug)
191  {
192  InfoInFunction << endl;
193  }
194 
195  calculate();
196 
197  if (debug)
198  {
199  Info<< " Finished" << endl;
200  }
201 }
202 
203 
204 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
heRhoThermo(const fvMesh &, const word &phaseName)
Construct from mesh and phase name.
Definition: heRhoThermo.C:167
virtual ~heRhoThermo()
Destructor.
Definition: heRhoThermo.C:181
virtual void correct()
Update properties.
Definition: heRhoThermo.C:188
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
dynamicFvMesh & mesh
fvPatchField< scalar > fvPatchScalarField
A class for handling words, derived from string.
Definition: word.H:59
volScalarField scalarField(fieldObject, mesh)
thermo he()
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.