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-2019 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 "volFields.H"
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 template<class BasicSolidThermo, class MixtureType>
33 {
34  scalarField& TCells = this->T_.primitiveFieldRef();
35 
36  const scalarField& hCells = this->he_;
37  const scalarField& pCells = this->p_;
38  scalarField& rhoCells = this->rho_.primitiveFieldRef();
39  scalarField& alphaCells = this->alpha_.primitiveFieldRef();
40 
41  forAll(TCells, celli)
42  {
43  const typename MixtureType::thermoType& mixture_ =
44  this->cellMixture(celli);
45 
46  rhoCells[celli] = mixture_.rho(pCells[celli], TCells[celli]);
47 
48  TCells[celli] = mixture_.THE
49  (
50  hCells[celli],
51  pCells[celli],
52  TCells[celli]
53  );
54 
55  const typename MixtureType::thermoType& volMixture_ =
56  this->cellVolMixture(pCells[celli], TCells[celli], celli);
57 
58  alphaCells[celli] =
59  volMixture_.kappa(pCells[celli], TCells[celli])
60  /mixture_.Cpv(pCells[celli], TCells[celli]);
61  }
62 
63  volScalarField::Boundary& pBf =
64  this->p_.boundaryFieldRef();
65 
66  volScalarField::Boundary& TBf =
67  this->T_.boundaryFieldRef();
68 
69  volScalarField::Boundary& rhoBf =
70  this->rho_.boundaryFieldRef();
71 
72  volScalarField::Boundary& heBf =
73  this->he().boundaryFieldRef();
74 
75  volScalarField::Boundary& alphaBf =
76  this->alpha_.boundaryFieldRef();
77 
78  forAll(this->T_.boundaryField(), patchi)
79  {
80  fvPatchScalarField& pp = pBf[patchi];
81  fvPatchScalarField& pT = TBf[patchi];
82  fvPatchScalarField& prho = rhoBf[patchi];
83  fvPatchScalarField& phe = heBf[patchi];
84  fvPatchScalarField& palpha = alphaBf[patchi];
85 
86  if (pT.fixesValue())
87  {
88  forAll(pT, facei)
89  {
90  const typename MixtureType::thermoType& mixture_ =
91  this->patchFaceMixture(patchi, facei);
92 
93  const typename MixtureType::thermoType& volMixture_ =
94  this->patchFaceVolMixture
95  (
96  pp[facei],
97  pT[facei],
98  patchi,
99  facei
100  );
101 
102 
103  phe[facei] = mixture_.HE(pp[facei], pT[facei]);
104  prho[facei] = volMixture_.rho(pp[facei], pT[facei]);
105 
106  palpha[facei] =
107  volMixture_.kappa(pp[facei], pT[facei])
108  / mixture_.Cpv(pp[facei], pT[facei]);
109  }
110  }
111  else
112  {
113  forAll(pT, facei)
114  {
115  const typename MixtureType::thermoType& mixture_ =
116  this->patchFaceMixture(patchi, facei);
117 
118  const typename MixtureType::thermoType& volMixture_ =
119  this->patchFaceVolMixture
120  (
121  pp[facei],
122  pT[facei],
123  patchi,
124  facei
125  );
126 
127  pT[facei] = mixture_.THE(phe[facei], pp[facei] ,pT[facei]);
128  prho[facei] = volMixture_.rho(pp[facei], pT[facei]);
129 
130  palpha[facei] =
131  volMixture_.kappa(pp[facei], pT[facei])
132  / mixture_.Cpv(pp[facei], pT[facei]);
133  }
134  }
135  }
136 
137  this->alpha_.correctBoundaryConditions();
138 }
139 
140 
141 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
142 
143 template<class BasicSolidThermo, class MixtureType>
146 (
147  const fvMesh& mesh,
148  const word& phaseName
149 )
150 :
152 {
153  calculate();
154 }
155 
156 
157 template<class BasicSolidThermo, class MixtureType>
160 (
161  const fvMesh& mesh,
162  const dictionary& dict,
163  const word& phaseName
164 )
165 :
167 {
168  calculate();
169 }
170 
171 
172 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
173 
174 template<class BasicSolidThermo, class MixtureType>
176 {}
177 
178 
179 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
180 
181 template<class BasicSolidThermo, class MixtureType>
183 {
184  if (debug)
185  {
186  InfoInFunction << endl;
187  }
188 
189  calculate();
190 
191  if (debug)
192  {
193  Info<< " Finished" << endl;
194  }
195 }
196 
197 
198 template<class BasicSolidThermo, class MixtureType>
201 {
202  const fvMesh& mesh = this->T_.mesh();
203 
204  tmp<volVectorField> tKappa
205  (
207  (
208  "Kappa",
209  mesh,
211  )
212  );
213 
214  volVectorField& Kappa = tKappa.ref();
215  vectorField& KappaCells = Kappa.primitiveFieldRef();
216  const scalarField& TCells = this->T_;
217  const scalarField& pCells = this->p_;
218 
219  forAll(KappaCells, celli)
220  {
221  Kappa[celli] =
222  this->cellVolMixture
223  (
224  pCells[celli],
225  TCells[celli],
226  celli
227  ).Kappa(pCells[celli], TCells[celli]);
228  }
229 
230  volVectorField::Boundary& KappaBf = Kappa.boundaryFieldRef();
231 
232  forAll(KappaBf, patchi)
233  {
234  vectorField& Kappap = KappaBf[patchi];
235  const scalarField& pT = this->T_.boundaryField()[patchi];
236  const scalarField& pp = this->p_.boundaryField()[patchi];
237 
238  forAll(Kappap, facei)
239  {
240  Kappap[facei] =
241  this->patchFaceVolMixture
242  (
243  pp[facei],
244  pT[facei],
245  patchi,
246  facei
247  ).Kappa(pp[facei], pT[facei]);
248  }
249  }
250 
251  return tKappa;
252 }
253 
254 
255 template<class BasicSolidThermo, class MixtureType>
258 (
259  const label patchi
260 ) const
261 {
262  const scalarField& pp = this->p_.boundaryField()[patchi];
263  const scalarField& Tp = this->T_.boundaryField()[patchi];
264  tmp<vectorField> tKappa(new vectorField(pp.size()));
265 
266  vectorField& Kappap = tKappa.ref();
267 
268  forAll(Tp, facei)
269  {
270  Kappap[facei] =
271  this->patchFaceVolMixture
272  (
273  pp[facei],
274  Tp[facei],
275  patchi,
276  facei
277  ).Kappa(pp[facei], Tp[facei]);
278  }
279 
280  return tKappa;
281 }
282 
283 
284 // ************************************************************************* //
volScalarField & he
Definition: YEEqn.H:50
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:158
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
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
volVectorField vectorField(fieldObject, mesh)
virtual void correct()
Update properties.
virtual tmp< volVectorField > Kappa() const
Anisotropic thermal conductivity [W/m/K].
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
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
volScalarField scalarField(fieldObject, mesh)
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
const dimensionSet dimEnergy
label patchi
Energy for a solid mixture.
Definition: heSolidThermo.H:49
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
A class for managing temporary objects.
Definition: PtrList.H:53
virtual ~heSolidThermo()
Destructor.
#define InfoInFunction
Report an information message using Foam::Info.