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 #include "fvmLaplacian.H"
28 #include "fvcLaplacian.H"
29 #include "coordinateSystem.H"
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 template<class BasicSolidThermo, class MixtureType>
35 {
36  const scalarField& hCells = this->he_;
37  const auto& pCells = this->p_;
38 
39  scalarField& TCells = this->T_.primitiveFieldRef();
40  scalarField& CpCells = this->Cp_.primitiveFieldRef();
41  scalarField& CvCells = this->Cv_.primitiveFieldRef();
42  scalarField& rhoCells = this->rho_.primitiveFieldRef();
43  scalarField& 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  kappaCells[celli] =
65  transportMixture.kappa(pCells[celli], TCells[celli]);
66  }
67 
68 
69  volScalarField::Boundary& heBf =
70  this->he().boundaryFieldRef();
71 
72  const auto& pBf = this->p_.boundaryField();
73 
74  volScalarField::Boundary& TBf =
75  this->T_.boundaryFieldRef();
76 
77  volScalarField::Boundary& CpBf =
78  this->Cp_.boundaryFieldRef();
79 
80  volScalarField::Boundary& CvBf =
81  this->Cv_.boundaryFieldRef();
82 
83  volScalarField::Boundary& rhoBf =
84  this->rho_.boundaryFieldRef();
85 
86  volScalarField::Boundary& kappaBf =
87  this->kappa_.boundaryFieldRef();
88 
89  forAll(this->T_.boundaryField(), patchi)
90  {
91  fvPatchScalarField& phe = heBf[patchi];
92  const auto& pp = pBf[patchi];
93  fvPatchScalarField& pT = TBf[patchi];
94  fvPatchScalarField& pCp = CpBf[patchi];
95  fvPatchScalarField& pCv = CvBf[patchi];
96  fvPatchScalarField& prho = rhoBf[patchi];
97  fvPatchScalarField& pkappa = kappaBf[patchi];
98 
99  if (pT.fixesValue())
100  {
101  forAll(pT, facei)
102  {
103  const typename MixtureType::thermoMixtureType&
104  thermoMixture = this->patchFaceThermoMixture(patchi, facei);
105 
106  const typename MixtureType::transportMixtureType&
107  transportMixture =
108  this->patchFaceTransportMixture
109  (patchi, facei, thermoMixture);
110 
111  phe[facei] = thermoMixture.HE(pp[facei], pT[facei]);
112 
113  prho[facei] = thermoMixture.rho(pp[facei], pT[facei]);
114  pCp[facei] = thermoMixture.Cp(pp[facei], pT[facei]);
115  pCv[facei] = thermoMixture.Cv(pp[facei], pT[facei]);
116 
117  pkappa[facei] = transportMixture.kappa(pp[facei], pT[facei]);
118  }
119  }
120  else
121  {
122  forAll(pT, facei)
123  {
124  const typename MixtureType::thermoMixtureType&
125  thermoMixture = this->patchFaceThermoMixture(patchi, facei);
126 
127  const typename MixtureType::transportMixtureType&
128  transportMixture =
129  this->patchFaceTransportMixture
130  (patchi, facei, thermoMixture);
131 
132  pT[facei] = thermoMixture.THE(phe[facei], pp[facei] ,pT[facei]);
133 
134  prho[facei] = thermoMixture.rho(pp[facei], pT[facei]);
135  pCp[facei] = thermoMixture.Cp(pp[facei], pT[facei]);
136  pCv[facei] = thermoMixture.Cv(pp[facei], pT[facei]);
137 
138  pkappa[facei] = transportMixture.kappa(pp[facei], pT[facei]);
139  }
140  }
141  }
142 
143  this->kappa_.correctBoundaryConditions();
144 }
145 
146 
147 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
148 
149 template<class BasicSolidThermo, class MixtureType>
152 (
153  const fvMesh& mesh,
154  const word& phaseName
155 )
156 :
158 {
159  calculate();
160 }
161 
162 
163 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
164 
165 template<class BasicSolidThermo, class MixtureType>
167 {}
168 
169 
170 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
171 
172 template<class BasicSolidThermo, class MixtureType>
174 {
175  if (debug)
176  {
177  InfoInFunction << endl;
178  }
179 
180  calculate();
181 
182  if (debug)
183  {
184  Info<< " Finished" << endl;
185  }
186 }
187 
188 
189 template<class BasicSolidThermo, class MixtureType>
192 {
193  const fvMesh& mesh = this->T_.mesh();
194 
195  tmp<volVectorField> tKappa
196  (
198  (
199  "Kappa",
200  mesh,
202  )
203  );
204 
205  const auto& pCells = this->p_;
206  const scalarField& TCells = this->T_;
207 
208  volVectorField& Kappa = tKappa.ref();
209  vectorField& KappaCells = Kappa.primitiveFieldRef();
210 
211  forAll(KappaCells, celli)
212  {
213  Kappa[celli] =
214  this->cellTransportMixture
215  (
216  celli
217  ).Kappa(pCells[celli], TCells[celli]);
218  }
219 
220  volVectorField::Boundary& KappaBf = Kappa.boundaryFieldRef();
221 
222  forAll(KappaBf, patchi)
223  {
224  const auto& pp = this->p_.boundaryField()[patchi];
225  const scalarField& pT = this->T_.boundaryField()[patchi];
226 
227  vectorField& Kappap = KappaBf[patchi];
228 
229  forAll(Kappap, facei)
230  {
231  Kappap[facei] =
232  this->patchFaceTransportMixture
233  (
234  patchi,
235  facei
236  ).Kappa(pp[facei], pT[facei]);
237  }
238  }
239 
240  return tKappa;
241 }
242 
243 
244 template<class BasicSolidThermo, class MixtureType>
247 (
248  const label patchi
249 ) const
250 {
251  const auto& pp = this->p_.boundaryField()[patchi];
252  const scalarField& Tp = this->T_.boundaryField()[patchi];
253 
254  tmp<vectorField> tKappa(new vectorField(Tp.size()));
255  vectorField& Kappap = tKappa.ref();
256 
257  forAll(Tp, facei)
258  {
259  Kappap[facei] =
260  this->patchFaceTransportMixture
261  (
262  patchi,
263  facei
264  ).Kappa(pp[patchi], Tp[facei]);
265  }
266 
267  return tKappa;
268 }
269 
270 
271 template<class BasicSolidThermo, class MixtureType>
274 {
275  const fvMesh& mesh = this->T_.mesh();
276 
277  const coordinateSystem coordinates
278  (
279  coordinateSystem::New(mesh, this->properties())
280  );
281 
282  const tmp<volVectorField> tKappa(Kappa());
283  const volVectorField& Kappa = tKappa();
284 
285  tmp<volSymmTensorField> tKappaLocal
286  (
288  (
289  "KappaLocal",
290  mesh,
292  )
293  );
294  volSymmTensorField& KappaLocal = tKappaLocal.ref();
295 
296  KappaLocal.primitiveFieldRef() =
297  coordinates.R(mesh.C()).transformVector(Kappa);
298 
299  forAll(KappaLocal.boundaryField(), patchi)
300  {
301  KappaLocal.boundaryFieldRef()[patchi] =
302  coordinates.R(mesh.boundary()[patchi].Cf())
303  .transformVector(Kappa.boundaryField()[patchi]);
304  }
305 
306  return tKappaLocal;
307 }
308 
309 
310 template<class BasicSolidThermo, class MixtureType>
313 (
314  const label patchi
315 ) const
316 {
317  const fvMesh& mesh = this->T_.mesh();
318 
319  const coordinateSystem coordinates
320  (
321  coordinateSystem::New(mesh, this->properties())
322  );
323 
324  return
325  coordinates.R(mesh.boundary()[patchi].Cf())
326  .transformVector(Kappa(patchi));
327 }
328 
329 
330 template<class BasicSolidThermo, class MixtureType>
333 {
334  const fvMesh& mesh = this->T_.mesh();
335  mesh.schemes().setFluxRequired(this->T_.name());
336 
337  return
338  - (
339  isotropic()
340  ? fvm::laplacian(this->kappa(), this->T_)().flux()
341  : fvm::laplacian(KappaLocal(), this->T_)().flux()
342  )/mesh.magSf();
343 }
344 
345 
346 template<class BasicSolidThermo, class MixtureType>
349 (
350  volScalarField& e
351 ) const
352 {
353  return
354  - (
355  isotropic()
356  ? fvc::laplacian(this->kappa(), this->T_)
357  + correction
358  (
360  (
361  this->kappa()/this->Cv(),
362  e,
363  "laplacian(alphae,e)"
364  )
365  )
366  : fvc::laplacian(KappaLocal(), this->T_)
367  + correction
368  (
370  (
371  KappaLocal()/this->Cv(),
372  e,
373  "laplacian(alphae,e)"
374  )
375  )
376  );
377 }
378 
379 
380 // ************************************************************************* //
Base class for other coordinate system specifications.
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
scalar Cv(const scalar p, const scalar T) const
Definition: HtoEthermo.H:2
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1672
const Boundary & boundaryField() const
Return const-reference to the boundary field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Calculate the matrix for the laplacian of the field.
void setFluxRequired(const word &name) const
Definition: fvSchemes.C:487
volVectorField vectorField(fieldObject, mesh)
fvMesh & mesh
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
const dimensionSet dimLength
const dimensionSet dimTime
const dimensionSet & dimensions() const
Return dimensions.
virtual void correct()
Update properties.
virtual tmp< volVectorField > Kappa() const
Anisotropic thermal conductivity [W/m/K].
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
label calculate(const fvMesh &mesh, const labelHashSet &patchIDs, const scalar minFaceFraction, GeometricField< scalar, PatchField, GeoMesh > &distance)
Calculate distance data from patches.
fvPatchField< scalar > fvPatchScalarField
A class for handling words, derived from string.
Definition: word.H:59
Calculate the laplacian of the given field.
static const zero Zero
Definition: zero.H:97
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
volScalarField scalarField(fieldObject, mesh)
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
thermo he()
heSolidThermo(const fvMesh &, const word &phaseName)
Construct from mesh and phase name.
const dimensionSet dimEnergy
dimensioned< symmTensor > dimensionedSymmTensor
Dimensioned tensor obtained from generic dimensioned type.
virtual tmp< surfaceScalarField > q() const
Return the heat flux [W/m^2].
Internal & ref()
Return a reference to the dimensioned internal field.
label patchi
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
Energy for a solid mixture.
Definition: heSolidThermo.H:49
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
messageStream Info
const volVectorField & C() const
Return cell centres.
A class for managing temporary objects.
Definition: PtrList.H:53
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition: fvcFlux.C:32
virtual ~heSolidThermo()
Destructor.
const dimensionSet dimTemperature
const coordinateRotation & R() const
Return const reference to co-ordinate rotation.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:800
#define InfoInFunction
Report an information message using Foam::Info.