heSolidThermo.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-2022 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 "heSolidThermo.H"
27 
28 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
30 template<class BasicSolidThermo, class MixtureType>
32 {
33  const bool isotropic = this->isotropic();
34 
35  const scalarField& hCells = this->he_;
36  const auto& pCells = this->p_;
37 
38  scalarField& TCells = this->T_.primitiveFieldRef();
39  scalarField& CpCells = this->Cp_.primitiveFieldRef();
40  scalarField& CvCells = this->Cv_.primitiveFieldRef();
41  scalarField& rhoCells = this->rho_.primitiveFieldRef();
42  scalarField& kappaCells = this->kappa_.primitiveFieldRef();
43  vectorField& KappaCells = this->Kappa_.primitiveFieldRef();
44 
45  forAll(TCells, celli)
46  {
47  const typename MixtureType::thermoMixtureType& thermoMixture =
48  this->cellThermoMixture(celli);
49 
50  const typename MixtureType::transportMixtureType& transportMixture =
51  this->cellTransportMixture(celli, thermoMixture);
52 
53  TCells[celli] = thermoMixture.THE
54  (
55  hCells[celli],
56  pCells[celli],
57  TCells[celli]
58  );
59 
60  CpCells[celli] = thermoMixture.Cp(pCells[celli], TCells[celli]);
61  CvCells[celli] = thermoMixture.Cv(pCells[celli], TCells[celli]);
62  rhoCells[celli] = thermoMixture.rho(pCells[celli], TCells[celli]);
63 
64  if (isotropic)
65  {
66  kappaCells[celli] =
67  transportMixture.kappa(pCells[celli], TCells[celli]);
68  }
69  else
70  {
71  KappaCells[celli] =
72  transportMixture.Kappa(pCells[celli], TCells[celli]);
73  }
74  }
75 
76 
78  this->he().boundaryFieldRef();
79 
80  const auto& pBf = this->p_.boundaryField();
81 
83  this->T_.boundaryFieldRef();
84 
86  this->Cp_.boundaryFieldRef();
87 
89  this->Cv_.boundaryFieldRef();
90 
92  this->rho_.boundaryFieldRef();
93 
94  volScalarField::Boundary& kappaBf =
95  this->kappa_.boundaryFieldRef();
96 
97  volVectorField::Boundary& KappaBf =
98  this->Kappa_.boundaryFieldRef();
99 
100  forAll(this->T_.boundaryField(), patchi)
101  {
102  fvPatchScalarField& phe = heBf[patchi];
103  const auto& pp = pBf[patchi];
104  fvPatchScalarField& pT = TBf[patchi];
105  fvPatchScalarField& pCp = CpBf[patchi];
106  fvPatchScalarField& pCv = CvBf[patchi];
107  fvPatchScalarField& prho = rhoBf[patchi];
108  fvPatchScalarField& pkappa = kappaBf[patchi];
109  fvPatchVectorField& pKappa = KappaBf[patchi];
110 
111  if (pT.fixesValue())
112  {
113  forAll(pT, facei)
114  {
115  const typename MixtureType::thermoMixtureType&
116  thermoMixture = this->patchFaceThermoMixture(patchi, facei);
117 
118  const typename MixtureType::transportMixtureType&
119  transportMixture =
120  this->patchFaceTransportMixture
121  (patchi, facei, thermoMixture);
122 
123  phe[facei] = thermoMixture.HE(pp[facei], pT[facei]);
124 
125  prho[facei] = thermoMixture.rho(pp[facei], pT[facei]);
126  pCp[facei] = thermoMixture.Cp(pp[facei], pT[facei]);
127  pCv[facei] = thermoMixture.Cv(pp[facei], pT[facei]);
128 
129  if (isotropic)
130  {
131  pkappa[facei] =
132  transportMixture.kappa(pp[facei], pT[facei]);
133  }
134  else
135  {
136  pKappa[facei] =
137  transportMixture.Kappa(pp[facei], pT[facei]);
138  }
139  }
140  }
141  else
142  {
143  forAll(pT, facei)
144  {
145  const typename MixtureType::thermoMixtureType&
146  thermoMixture = this->patchFaceThermoMixture(patchi, facei);
147 
148  const typename MixtureType::transportMixtureType&
149  transportMixture =
150  this->patchFaceTransportMixture
151  (patchi, facei, thermoMixture);
152 
153  pT[facei] = thermoMixture.THE(phe[facei], pp[facei] ,pT[facei]);
154 
155  prho[facei] = thermoMixture.rho(pp[facei], pT[facei]);
156  pCp[facei] = thermoMixture.Cp(pp[facei], pT[facei]);
157  pCv[facei] = thermoMixture.Cv(pp[facei], pT[facei]);
158 
159  if (isotropic)
160  {
161  pkappa[facei] =
162  transportMixture.kappa(pp[facei], pT[facei]);
163  }
164  else
165  {
166  pKappa[facei] =
167  transportMixture.Kappa(pp[facei], pT[facei]);
168  }
169  }
170  }
171  }
172 }
173 
174 
175 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
176 
177 template<class BasicSolidThermo, class MixtureType>
180 (
181  const fvMesh& mesh,
182  const word& phaseName
183 )
184 :
185  heThermo<BasicSolidThermo, MixtureType>(mesh, phaseName),
186  Kappa_
187  (
188  IOobject
189  (
190  BasicSolidThermo::phasePropertyName("Kappa", phaseName),
191  mesh.time().name(),
192  mesh,
193  IOobject::NO_READ,
194  IOobject::NO_WRITE
195  ),
196  mesh,
198  )
199 {
200  calculate();
201 }
202 
203 
204 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
205 
206 template<class BasicSolidThermo, class MixtureType>
208 {}
209 
210 
211 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
212 
213 template<class BasicSolidThermo, class MixtureType>
215 {
216  if (debug)
217  {
218  InfoInFunction << endl;
219  }
220 
221  calculate();
222 
223  if (debug)
224  {
225  Info<< " Finished" << endl;
226  }
227 }
228 
229 
230 template<class BasicSolidThermo, class MixtureType>
233 {
234  return Kappa_;
235 }
236 
237 
238 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Generic GeometricBoundaryField class.
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual bool fixesValue() const
Return true if this patch field fixes a value.
Definition: fvPatchField.H:314
heSolidThermo(const fvMesh &, const word &phaseName)
Construct from mesh and phase name.
virtual void correct()
Update properties.
virtual const volVectorField & Kappa() const
Anisotropic thermal conductivity [W/m/K].
virtual ~heSolidThermo()
Destructor.
void calculate()
Calculate the thermo variables.
Definition: heSolidThermo.C:31
Enthalpy/Internal energy for a mixture.
Definition: heThermo.H:55
A class for handling words, derived from string.
Definition: word.H:62
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.
static const zero Zero
Definition: zero.H:97
const dimensionSet dimEnergy
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
const dimensionSet dimLength
const dimensionSet dimTemperature
const dimensionSet dimTime
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
thermo he()