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-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 "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& alphaCells = this->alpha_.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  alphaCells[celli] =
65  transportMixture.kappa(pCells[celli], TCells[celli])
66  /thermoMixture.Cv(pCells[celli], TCells[celli]);
67  }
68 
69 
70  volScalarField::Boundary& heBf =
71  this->he().boundaryFieldRef();
72 
73  const auto& pBf = this->p_.boundaryField();
74 
75  volScalarField::Boundary& TBf =
76  this->T_.boundaryFieldRef();
77 
78  volScalarField::Boundary& CpBf =
79  this->Cp_.boundaryFieldRef();
80 
81  volScalarField::Boundary& CvBf =
82  this->Cv_.boundaryFieldRef();
83 
84  volScalarField::Boundary& rhoBf =
85  this->rho_.boundaryFieldRef();
86 
87  volScalarField::Boundary& alphaBf =
88  this->alpha_.boundaryFieldRef();
89 
90  forAll(this->T_.boundaryField(), patchi)
91  {
92  fvPatchScalarField& phe = heBf[patchi];
93  const auto& pp = pBf[patchi];
94  fvPatchScalarField& pT = TBf[patchi];
95  fvPatchScalarField& pCp = CpBf[patchi];
96  fvPatchScalarField& pCv = CvBf[patchi];
97  fvPatchScalarField& prho = rhoBf[patchi];
98  fvPatchScalarField& palpha = alphaBf[patchi];
99 
100  if (pT.fixesValue())
101  {
102  forAll(pT, facei)
103  {
104  const typename MixtureType::thermoMixtureType&
105  thermoMixture = this->patchFaceThermoMixture(patchi, facei);
106 
107  const typename MixtureType::transportMixtureType&
108  transportMixture =
109  this->patchFaceTransportMixture
110  (patchi, facei, thermoMixture);
111 
112  phe[facei] = thermoMixture.HE(pp[facei], pT[facei]);
113 
114  prho[facei] = thermoMixture.rho(pp[facei], pT[facei]);
115  pCp[facei] = thermoMixture.Cp(pp[facei], pT[facei]);
116  pCv[facei] = thermoMixture.Cv(pp[facei], pT[facei]);
117 
118  palpha[facei] =
119  transportMixture.kappa(pp[facei], pT[facei])
120  /thermoMixture.Cv(pp[facei], pT[facei]);
121  }
122  }
123  else
124  {
125  forAll(pT, facei)
126  {
127  const typename MixtureType::thermoMixtureType&
128  thermoMixture = this->patchFaceThermoMixture(patchi, facei);
129 
130  const typename MixtureType::transportMixtureType&
131  transportMixture =
132  this->patchFaceTransportMixture
133  (patchi, facei, thermoMixture);
134 
135  pT[facei] = thermoMixture.THE(phe[facei], pp[facei] ,pT[facei]);
136 
137  prho[facei] = thermoMixture.rho(pp[facei], pT[facei]);
138  pCp[facei] = thermoMixture.Cp(pp[facei], pT[facei]);
139  pCv[facei] = thermoMixture.Cv(pp[facei], pT[facei]);
140 
141  palpha[facei] =
142  transportMixture.kappa(pp[facei], pT[facei])
143  /thermoMixture.Cv(pp[facei], pT[facei]);
144  }
145  }
146  }
147 
148  this->alpha_.correctBoundaryConditions();
149 }
150 
151 
152 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
153 
154 template<class BasicSolidThermo, class MixtureType>
157 (
158  const fvMesh& mesh,
159  const word& phaseName
160 )
161 :
163 {
164  calculate();
165 }
166 
167 
168 template<class BasicSolidThermo, class MixtureType>
171 (
172  const fvMesh& mesh,
173  const dictionary& dict,
174  const word& phaseName
175 )
176 :
178 {
179  calculate();
180 }
181 
182 
183 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
184 
185 template<class BasicSolidThermo, class MixtureType>
187 {}
188 
189 
190 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
191 
192 template<class BasicSolidThermo, class MixtureType>
194 {
195  if (debug)
196  {
197  InfoInFunction << endl;
198  }
199 
200  calculate();
201 
202  if (debug)
203  {
204  Info<< " Finished" << endl;
205  }
206 }
207 
208 
209 template<class BasicSolidThermo, class MixtureType>
212 {
213  const fvMesh& mesh = this->T_.mesh();
214 
215  tmp<volVectorField> tKappa
216  (
218  (
219  "Kappa",
220  mesh,
222  )
223  );
224 
225  const auto& pCells = this->p_;
226  const scalarField& TCells = this->T_;
227 
228  volVectorField& Kappa = tKappa.ref();
229  vectorField& KappaCells = Kappa.primitiveFieldRef();
230 
231  forAll(KappaCells, celli)
232  {
233  Kappa[celli] =
234  this->cellTransportMixture
235  (
236  celli
237  ).Kappa(pCells[celli], TCells[celli]);
238  }
239 
240  volVectorField::Boundary& KappaBf = Kappa.boundaryFieldRef();
241 
242  forAll(KappaBf, patchi)
243  {
244  const auto& pp = this->p_.boundaryField()[patchi];
245  const scalarField& pT = this->T_.boundaryField()[patchi];
246 
247  vectorField& Kappap = KappaBf[patchi];
248 
249  forAll(Kappap, facei)
250  {
251  Kappap[facei] =
252  this->patchFaceTransportMixture
253  (
254  patchi,
255  facei
256  ).Kappa(pp[facei], pT[facei]);
257  }
258  }
259 
260  return tKappa;
261 }
262 
263 
264 template<class BasicSolidThermo, class MixtureType>
267 (
268  const label patchi
269 ) const
270 {
271  const auto& pp = this->p_.boundaryField()[patchi];
272  const scalarField& Tp = this->T_.boundaryField()[patchi];
273 
274  tmp<vectorField> tKappa(new vectorField(Tp.size()));
275  vectorField& Kappap = tKappa.ref();
276 
277  forAll(Tp, facei)
278  {
279  Kappap[facei] =
280  this->patchFaceTransportMixture
281  (
282  patchi,
283  facei
284  ).Kappa(pp[patchi], Tp[facei]);
285  }
286 
287  return tKappa;
288 }
289 
290 
291 template<class BasicSolidThermo, class MixtureType>
294 {
295  const fvMesh& mesh = this->T_.mesh();
296 
297  const coordinateSystem coordinates
298  (
299  coordinateSystem::New(mesh, this->properties())
300  );
301 
302  const tmp<volVectorField> tKappa(Kappa());
303  const volVectorField& Kappa = tKappa();
304 
305  tmp<volSymmTensorField> tKappaLocal
306  (
308  (
309  "KappaLocal",
310  mesh,
312  )
313  );
314  volSymmTensorField& KappaLocal = tKappaLocal.ref();
315 
316  KappaLocal.primitiveFieldRef() =
317  coordinates.R(mesh.C()).transformVector(Kappa);
318 
319  forAll(KappaLocal.boundaryField(), patchi)
320  {
321  KappaLocal.boundaryFieldRef()[patchi] =
322  coordinates.R(mesh.boundary()[patchi].Cf())
323  .transformVector(Kappa.boundaryField()[patchi]);
324  }
325 
326  return tKappaLocal;
327 }
328 
329 
330 template<class BasicSolidThermo, class MixtureType>
333 (
334  const label patchi
335 ) const
336 {
337  const fvMesh& mesh = this->T_.mesh();
338 
339  const coordinateSystem coordinates
340  (
341  coordinateSystem::New(mesh, this->properties())
342  );
343 
344  return
345  coordinates.R(mesh.boundary()[patchi].Cf())
346  .transformVector(Kappa(patchi));
347 }
348 
349 
350 template<class BasicSolidThermo, class MixtureType>
353 {
354  const fvMesh& mesh = this->T_.mesh();
355  mesh.setFluxRequired(this->T_.name());
356 
357  return
358  - (
359  isotropic()
360  ? fvm::laplacian(this->kappa(), this->T_)().flux()
361  : fvm::laplacian(KappaLocal(), this->T_)().flux()
362  )/mesh.magSf();
363 }
364 
365 
366 template<class BasicSolidThermo, class MixtureType>
369 (
370  volScalarField& e
371 ) const
372 {
373  return
374  - (
375  isotropic()
376  ? fvc::laplacian(this->kappa(), this->T_)
377  + correction(fvm::laplacian(this->alpha(), e))
378  : fvc::laplacian(KappaLocal(), this->T_)
379  + correction
380  (
382  (
383  KappaLocal()/this->Cv(),
384  e,
385  "laplacian(" + this->alpha().name() + ",e)"
386  )
387  )
388  );
389 }
390 
391 
392 // ************************************************************************* //
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
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
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
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
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:500
volVectorField vectorField(fieldObject, 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].
dynamicFvMesh & mesh
fvPatchField< scalar > fvPatchScalarField
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
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.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
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:78
messageStream Info
const volVectorField & C() const
Return cell centres as volVectorField.
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:540
#define InfoInFunction
Report an information message using Foam::Info.