heSolidThermo.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2017-2018 OpenCFD Ltd
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "heSolidThermo.H"
30 #include "volFields.H"
31 
32 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33 
34 template<class BasicSolidThermo, class MixtureType>
36 {
37  scalarField& TCells = this->T_.primitiveFieldRef();
38 
39  const scalarField& hCells = this->he_;
40  const scalarField& pCells = this->p_;
41  scalarField& rhoCells = this->rho_.primitiveFieldRef();
42  //scalarField& psiCells = this->psi_.primitiveFieldRef();
43  //scalarField& muCells = this->mu_.primitiveFieldRef();
44  scalarField& alphaCells = this->alpha_.primitiveFieldRef();
45 
46  forAll(TCells, celli)
47  {
48  const typename MixtureType::thermoType& mixture_ =
49  this->cellMixture(celli);
50 
51  const typename MixtureType::thermoType& volMixture_ =
52  this->cellVolMixture(pCells[celli], TCells[celli], celli);
53 
54  if (this->updateT())
55  {
56  TCells[celli] = mixture_.THE
57  (
58  hCells[celli],
59  pCells[celli],
60  TCells[celli]
61  );
62  }
63 
64  rhoCells[celli] = volMixture_.rho(pCells[celli], TCells[celli]);
65  //psiCells[celli] = volMixture_.psi(pCells[celli], TCells[celli]);
66  //muCells[celli] = volMixture_.mu(pCells[celli], TCells[celli]);
67 
68  alphaCells[celli] =
69  volMixture_.kappa(pCells[celli], TCells[celli])
70  /
71  mixture_.Cpv(pCells[celli], TCells[celli]);
72  }
73 
74  volScalarField::Boundary& pBf = this->p_.boundaryFieldRef();
75  volScalarField::Boundary& TBf = this->T_.boundaryFieldRef();
76  volScalarField::Boundary& rhoBf = this->rho_.boundaryFieldRef();
77  //volScalarField::Boundary& psiBf = this->psi_.boundaryFieldRef();
78  //volScalarField::Boundary& muBf = this->mu_.boundaryFieldRef();
79  volScalarField::Boundary& heBf = this->he().boundaryFieldRef();
80  volScalarField::Boundary& alphaBf = this->alpha_.boundaryFieldRef();
81 
82  forAll(this->T_.boundaryField(), patchi)
83  {
84  fvPatchScalarField& pp = pBf[patchi];
85  fvPatchScalarField& pT = TBf[patchi];
86  fvPatchScalarField& prho = rhoBf[patchi];
87  //fvPatchScalarField& ppsi = psiBf[patchi];
88  //fvPatchScalarField& pmu = muBf[patchi];
89  fvPatchScalarField& phe = heBf[patchi];
90  fvPatchScalarField& palpha = alphaBf[patchi];
91 
92  if (pT.fixesValue())
93  {
94  forAll(pT, facei)
95  {
96  const typename MixtureType::thermoType& mixture_ =
97  this->patchFaceMixture(patchi, facei);
98 
99  const typename MixtureType::thermoType& volMixture_ =
100  this->patchFaceVolMixture
101  (
102  pp[facei],
103  pT[facei],
104  patchi,
105  facei
106  );
107 
108 
109  phe[facei] = mixture_.HE(pp[facei], pT[facei]);
110  prho[facei] = volMixture_.rho(pp[facei], pT[facei]);
111  //ppsi[facei] = volMixture_.psi(pp[facei], pT[facei]);
112  //pmu[facei] = volMixture_.mu(pp[facei], pT[facei]);
113 
114  palpha[facei] =
115  volMixture_.kappa(pp[facei], pT[facei])
116  / mixture_.Cpv(pp[facei], pT[facei]);
117  }
118  }
119  else
120  {
121  forAll(pT, facei)
122  {
123  const typename MixtureType::thermoType& mixture_ =
124  this->patchFaceMixture(patchi, facei);
125 
126  const typename MixtureType::thermoType& volMixture_ =
127  this->patchFaceVolMixture
128  (
129  pp[facei],
130  pT[facei],
131  patchi,
132  facei
133  );
134 
135  if (this->updateT())
136  {
137  pT[facei] = mixture_.THE(phe[facei], pp[facei] ,pT[facei]);
138  }
139 
140  prho[facei] = volMixture_.rho(pp[facei], pT[facei]);
141  //ppsi[facei] = volMixture_.psi(pp[facei], pT[facei]);
142  //pmu[facei] = volMixture_.mu(pp[facei], pT[facei]);
143 
144  palpha[facei] =
145  volMixture_.kappa(pp[facei], pT[facei])
146  / mixture_.Cpv(pp[facei], pT[facei]);
147  }
148  }
149  }
150 
151  this->alpha_.correctBoundaryConditions();
152 }
153 
154 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
155 
156 template<class BasicSolidThermo, class MixtureType>
159 (
160  const fvMesh& mesh,
161  const word& phaseName
162 )
163 :
165 {
166  calculate();
167  this->mu_ == dimensionedScalar(this->mu_.dimensions(), Zero);
168  this->psi_ == dimensionedScalar(this->psi_.dimensions(), Zero);
169 }
170 
171 
172 template<class BasicSolidThermo, class MixtureType>
175 (
176  const fvMesh& mesh,
177  const dictionary& dict,
178  const word& phaseName
179 )
180 :
182 {
183  calculate();
184  this->mu_ == dimensionedScalar(this->mu_.dimensions(), Zero);
185  this->psi_ == dimensionedScalar(this->psi_.dimensions(), Zero);
186 }
187 
188 
189 template<class BasicSolidThermo, class MixtureType>
192 (
193  const fvMesh& mesh,
194  const word& phaseName,
195  const word& dictName
196 )
197 :
199 {
200  calculate();
201 
202  // TBD. initialise psi, mu (at heThermo level) since these do not
203  // get initialised. Move to heThermo constructor?
204  this->mu_ == dimensionedScalar(this->mu_.dimensions(), Zero);
205  this->psi_ == dimensionedScalar(this->psi_.dimensions(), Zero);
206 }
207 
208 
209 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
210 
211 template<class BasicSolidThermo, class MixtureType>
213 {}
214 
215 
216 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
217 
218 template<class BasicSolidThermo, class MixtureType>
220 {
221  if (debug)
222  {
223  InfoInFunction << endl;
224  }
225 
226  calculate();
227 
228  if (debug)
229  {
230  Info<< " Finished" << endl;
231  }
232 }
233 
234 
235 template<class BasicSolidThermo, class MixtureType>
238 {
239  const fvMesh& mesh = this->T_.mesh();
240 
241  tmp<volVectorField> tKappa
242  (
243  new volVectorField
244  (
245  IOobject
246  (
247  "Kappa",
248  mesh.time().timeName(),
249  mesh,
250  IOobject::NO_READ,
251  IOobject::NO_WRITE
252  ),
253  mesh,
255  )
256  );
257 
258  volVectorField& Kappa = tKappa.ref();
259  vectorField& KappaCells = Kappa.primitiveFieldRef();
260  const scalarField& TCells = this->T_;
261  const scalarField& pCells = this->p_;
262 
263  forAll(KappaCells, celli)
264  {
265  Kappa[celli] =
266  this->cellVolMixture
267  (
268  pCells[celli],
269  TCells[celli],
270  celli
271  ).Kappa(pCells[celli], TCells[celli]);
272  }
273 
274  volVectorField::Boundary& KappaBf = Kappa.boundaryFieldRef();
275 
276  forAll(KappaBf, patchi)
277  {
278  vectorField& Kappap = KappaBf[patchi];
279  const scalarField& pT = this->T_.boundaryField()[patchi];
280  const scalarField& pp = this->p_.boundaryField()[patchi];
281 
282  forAll(Kappap, facei)
283  {
284  Kappap[facei] =
285  this->patchFaceVolMixture
286  (
287  pp[facei],
288  pT[facei],
289  patchi,
290  facei
291  ).Kappa(pp[facei], pT[facei]);
292  }
293  }
294 
295  return tKappa;
296 }
297 
298 
299 template<class BasicSolidThermo, class MixtureType>
302 (
303  const label patchi
304 ) const
305 {
306  const scalarField& pp = this->p_.boundaryField()[patchi];
307  const scalarField& Tp = this->T_.boundaryField()[patchi];
308  tmp<vectorField> tKappa(new vectorField(pp.size()));
309 
310  vectorField& Kappap = tKappa.ref();
311 
312  forAll(Tp, facei)
313  {
314  Kappap[facei] =
315  this->patchFaceVolMixture
316  (
317  pp[facei],
318  Tp[facei],
319  patchi,
320  facei
321  ).Kappa(pp[facei], Tp[facei]);
322  }
323 
324  return tKappa;
325 }
326 
327 
328 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::fvPatchScalarField
fvPatchField< scalar > fvPatchScalarField
Definition: fvPatchFieldsFwd.H:40
volFields.H
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
InfoInFunction
#define InfoInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:316
Foam::heThermo< BasicSolidThermo, MixtureType >
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::dimEnergy
const dimensionSet dimEnergy
Foam::heSolidThermo::correct
virtual void correct()
Update properties.
Definition: heSolidThermo.C:219
dictName
const word dictName("blockMeshDict")
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:258
Foam::label
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:62
Foam::Field< vector >
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::heSolidThermo::~heSolidThermo
virtual ~heSolidThermo()
Destructor.
Definition: heSolidThermo.C:212
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
Foam::heSolidThermo::Kappa
virtual tmp< volVectorField > Kappa() const
Anisotropic thermal conductivity [W/m/K].
Definition: heSolidThermo.C:237
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::heSolidThermo
Energy for a solid mixture.
Definition: heSolidThermo.H:52
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
heSolidThermo.H
he
volScalarField & he
Definition: YEEqn.H:52
Foam::GeometricField::ref
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
Definition: GeometricField.C:718
Foam::dimTemperature
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:55
Foam::dimensioned::dimensions
const dimensionSet & dimensions() const
Return const reference to dimensions.
Definition: dimensionedType.C:390
Foam::GeometricField< vector, fvPatchField, volMesh >