RhoFluidThermo.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-2023 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 "RhoFluidThermo.H"
27 
28 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
30 template<class BaseThermo>
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& kappaCells = this->kappa_.primitiveFieldRef();
43 
44  auto Yslicer = this->Yslicer();
45 
46  forAll(TCells, celli)
47  {
48  auto composition = this->cellComposition(Yslicer, celli);
49 
50  const typename BaseThermo::mixtureType::thermoMixtureType&
51  thermoMixture = this->thermoMixture(composition);
52 
53  const typename BaseThermo::mixtureType::transportMixtureType&
54  transportMixture =
55  this->transportMixture(composition, thermoMixture);
56 
57  TCells[celli] = thermoMixture.The
58  (
59  hCells[celli],
60  pCells[celli],
61  TCells[celli]
62  );
63 
64  CpCells[celli] = thermoMixture.Cp(pCells[celli], TCells[celli]);
65  CvCells[celli] = thermoMixture.Cv(pCells[celli], TCells[celli]);
66  psiCells[celli] = thermoMixture.psi(pCells[celli], TCells[celli]);
67  rhoCells[celli] = thermoMixture.rho(pCells[celli], TCells[celli]);
68 
69  muCells[celli] = transportMixture.mu(pCells[celli], TCells[celli]);
70  kappaCells[celli] =
71  transportMixture.kappa(pCells[celli], TCells[celli]);
72  }
73 
74  volScalarField::Boundary& pBf =
75  this->p_.boundaryFieldRef();
76 
77  volScalarField::Boundary& TBf =
78  this->T_.boundaryFieldRef();
79 
80  volScalarField::Boundary& CpBf =
81  this->Cp_.boundaryFieldRef();
82 
83  volScalarField::Boundary& CvBf =
84  this->Cv_.boundaryFieldRef();
85 
86  volScalarField::Boundary& psiBf =
87  this->psi_.boundaryFieldRef();
88 
89  volScalarField::Boundary& rhoBf =
90  this->rho_.boundaryFieldRef();
91 
92  volScalarField::Boundary& heBf =
93  this->he().boundaryFieldRef();
94 
95  volScalarField::Boundary& muBf =
96  this->mu_.boundaryFieldRef();
97 
98  volScalarField::Boundary& kappaBf =
99  this->kappa_.boundaryFieldRef();
100 
101  forAll(this->T_.boundaryField(), patchi)
102  {
103  fvPatchScalarField& pp = pBf[patchi];
104  fvPatchScalarField& pT = TBf[patchi];
105  fvPatchScalarField& pCp = CpBf[patchi];
106  fvPatchScalarField& pCv = CvBf[patchi];
107  fvPatchScalarField& ppsi = psiBf[patchi];
108  fvPatchScalarField& prho = rhoBf[patchi];
109  fvPatchScalarField& phe = heBf[patchi];
110  fvPatchScalarField& pmu = muBf[patchi];
111  fvPatchScalarField& pkappa = kappaBf[patchi];
112 
113  if (pT.fixesValue())
114  {
115  forAll(pT, facei)
116  {
117  auto composition =
118  this->patchFaceComposition(Yslicer, patchi, facei);
119 
120  const typename BaseThermo::mixtureType::thermoMixtureType&
121  thermoMixture = this->thermoMixture(composition);
122 
123  const typename BaseThermo::mixtureType::transportMixtureType&
124  transportMixture =
125  this->transportMixture(composition, thermoMixture);
126 
127  phe[facei] = thermoMixture.he(pp[facei], pT[facei]);
128 
129  pCp[facei] = thermoMixture.Cp(pp[facei], pT[facei]);
130  pCv[facei] = thermoMixture.Cv(pp[facei], pT[facei]);
131  ppsi[facei] = thermoMixture.psi(pp[facei], pT[facei]);
132  prho[facei] = thermoMixture.rho(pp[facei], pT[facei]);
133 
134  pmu[facei] = transportMixture.mu(pp[facei], pT[facei]);
135  pkappa[facei] = transportMixture.kappa(pp[facei], pT[facei]);
136  }
137  }
138  else
139  {
140  forAll(pT, facei)
141  {
142  auto composition =
143  this->patchFaceComposition(Yslicer, patchi, facei);
144 
145  const typename BaseThermo::mixtureType::thermoMixtureType&
146  thermoMixture = this->thermoMixture(composition);
147 
148  const typename BaseThermo::mixtureType::transportMixtureType&
149  transportMixture =
150  this->transportMixture(composition, thermoMixture);
151 
152  pT[facei] = thermoMixture.The(phe[facei], pp[facei], pT[facei]);
153 
154  pCp[facei] = thermoMixture.Cp(pp[facei], pT[facei]);
155  pCv[facei] = thermoMixture.Cv(pp[facei], pT[facei]);
156  ppsi[facei] = thermoMixture.psi(pp[facei], pT[facei]);
157  prho[facei] = thermoMixture.rho(pp[facei], pT[facei]);
158 
159  pmu[facei] = transportMixture.mu(pp[facei], pT[facei]);
160  pkappa[facei] = transportMixture.kappa(pp[facei], pT[facei]);
161  }
162  }
163  }
164 }
165 
166 
167 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
168 
169 template<class BaseThermo>
171 (
172  const fvMesh& mesh,
173  const word& phaseName
174 )
175 :
176  BaseThermo(mesh, phaseName)
177 {
178  calculate();
179 }
180 
181 
182 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
183 
184 template<class BaseThermo>
186 {}
187 
188 
189 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
190 
191 template<class BaseThermo>
193 {
194  if (BaseThermo::debug)
195  {
196  InfoInFunction << endl;
197  }
198 
199  calculate();
200 
201  if (BaseThermo::debug)
202  {
203  Info<< " Finished" << endl;
204  }
205 }
206 
207 
208 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Thermo implementation based on density.
virtual void correct()
Update properties.
RhoFluidThermo(const fvMesh &, const word &phaseName)
Construct from mesh and phase name.
virtual ~RhoFluidThermo()
Destructor.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
A class for handling words, derived from string.
Definition: word.H:62
volScalarField scalarField(fieldObject, mesh)
label patchi
#define InfoInFunction
Report an information message using Foam::Info.
label calculate(const fvMesh &mesh, const labelHashSet &patchIDs, const scalar minFaceFraction, GeometricField< scalar, PatchField, GeoMesh > &distance)
Calculate distance data from patches.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
messageStream Info
fvPatchField< scalar > fvPatchScalarField
thermo he()