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 | Copyright (C) 2011-2016 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  const typename MixtureType::thermoType& volMixture_ =
47  this->cellVolMixture(pCells[celli], TCells[celli], celli);
48 
49  TCells[celli] = mixture_.THE
50  (
51  hCells[celli],
52  pCells[celli],
53  TCells[celli]
54  );
55 
56  rhoCells[celli] = volMixture_.rho(pCells[celli], TCells[celli]);
57 
58  alphaCells[celli] =
59  volMixture_.kappa(pCells[celli], TCells[celli])
60  /
61  mixture_.Cpv(pCells[celli], TCells[celli]);
62  }
63 
64  volScalarField::Boundary& pBf =
65  this->p_.boundaryFieldRef();
66 
67  volScalarField::Boundary& TBf =
68  this->T_.boundaryFieldRef();
69 
70  volScalarField::Boundary& rhoBf =
71  this->rho_.boundaryFieldRef();
72 
73  volScalarField::Boundary& heBf =
74  this->he().boundaryFieldRef();
75 
76  volScalarField::Boundary& alphaBf =
77  this->alpha_.boundaryFieldRef();
78 
79  forAll(this->T_.boundaryField(), patchi)
80  {
81  fvPatchScalarField& pp = pBf[patchi];
82  fvPatchScalarField& pT = TBf[patchi];
83  fvPatchScalarField& prho = rhoBf[patchi];
84  fvPatchScalarField& phe = heBf[patchi];
85  fvPatchScalarField& palpha = alphaBf[patchi];
86 
87  if (pT.fixesValue())
88  {
89  forAll(pT, facei)
90  {
91  const typename MixtureType::thermoType& mixture_ =
92  this->patchFaceMixture(patchi, facei);
93 
94  const typename MixtureType::thermoType& volMixture_ =
95  this->patchFaceVolMixture
96  (
97  pp[facei],
98  pT[facei],
99  patchi,
100  facei
101  );
102 
103 
104  phe[facei] = mixture_.HE(pp[facei], pT[facei]);
105  prho[facei] = volMixture_.rho(pp[facei], pT[facei]);
106 
107  palpha[facei] =
108  volMixture_.kappa(pp[facei], pT[facei])
109  / mixture_.Cpv(pp[facei], pT[facei]);
110  }
111  }
112  else
113  {
114  forAll(pT, facei)
115  {
116  const typename MixtureType::thermoType& mixture_ =
117  this->patchFaceMixture(patchi, facei);
118 
119  const typename MixtureType::thermoType& volMixture_ =
120  this->patchFaceVolMixture
121  (
122  pp[facei],
123  pT[facei],
124  patchi,
125  facei
126  );
127 
128  pT[facei] = mixture_.THE(phe[facei], pp[facei] ,pT[facei]);
129  prho[facei] = volMixture_.rho(pp[facei], pT[facei]);
130 
131  palpha[facei] =
132  volMixture_.kappa(pp[facei], pT[facei])
133  / mixture_.Cpv(pp[facei], pT[facei]);
134  }
135  }
136  }
137 }
138 
139 
140 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
141 
142 template<class BasicSolidThermo, class MixtureType>
145 (
146  const fvMesh& mesh,
147  const word& phaseName
148 )
149 :
151 {
152  calculate();
153 }
154 
155 
156 template<class BasicSolidThermo, class MixtureType>
159 (
160  const fvMesh& mesh,
161  const dictionary& dict,
162  const word& phaseName
163 )
164 :
166 {
167  calculate();
168 }
169 
170 
171 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
172 
173 template<class BasicSolidThermo, class MixtureType>
175 {}
176 
177 
178 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
179 
180 template<class BasicSolidThermo, class MixtureType>
182 {
183  if (debug)
184  {
185  InfoInFunction << endl;
186  }
187 
188  calculate();
189 
190  if (debug)
191  {
192  Info<< " Finished" << endl;
193  }
194 }
195 
196 
197 template<class BasicSolidThermo, class MixtureType>
200 {
201  const fvMesh& mesh = this->T_.mesh();
202 
203  tmp<volVectorField> tKappa
204  (
205  new volVectorField
206  (
207  IOobject
208  (
209  "Kappa",
210  mesh.time().timeName(),
211  mesh,
212  IOobject::NO_READ,
213  IOobject::NO_WRITE
214  ),
215  mesh,
217  )
218  );
219 
220  volVectorField& Kappa = tKappa.ref();
221  vectorField& KappaCells = Kappa.primitiveFieldRef();
222  const scalarField& TCells = this->T_;
223  const scalarField& pCells = this->p_;
224 
225  forAll(KappaCells, celli)
226  {
227  Kappa[celli] =
228  this->cellVolMixture
229  (
230  pCells[celli],
231  TCells[celli],
232  celli
233  ).Kappa(pCells[celli], TCells[celli]);
234  }
235 
236  volVectorField::Boundary& KappaBf = Kappa.boundaryFieldRef();
237 
238  forAll(KappaBf, patchi)
239  {
240  vectorField& Kappap = KappaBf[patchi];
241  const scalarField& pT = this->T_.boundaryField()[patchi];
242  const scalarField& pp = this->p_.boundaryField()[patchi];
243 
244  forAll(Kappap, facei)
245  {
246  Kappap[facei] =
247  this->patchFaceVolMixture
248  (
249  pp[facei],
250  pT[facei],
251  patchi,
252  facei
253  ).Kappa(pp[facei], pT[facei]);
254  }
255  }
256 
257  return tKappa;
258 }
259 
260 
261 template<class BasicSolidThermo, class MixtureType>
264 (
265  const label patchi
266 ) const
267 {
268  const scalarField& pp = this->p_.boundaryField()[patchi];
269  const scalarField& Tp = this->T_.boundaryField()[patchi];
270  tmp<vectorField> tKappa(new vectorField(pp.size()));
271 
272  vectorField& Kappap = tKappa.ref();
273 
274  forAll(Tp, facei)
275  {
276  Kappap[facei] =
277  this->patchFaceVolMixture
278  (
279  pp[facei],
280  Tp[facei],
281  patchi,
282  facei
283  ).Kappa(pp[facei], Tp[facei]);
284  }
285 
286  return tKappa;
287 }
288 
289 
290 // ************************************************************************* //
volScalarField & he
Definition: YEEqn.H:51
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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:137
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:163
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
volVectorField vectorField(fieldObject, mesh)
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:644
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
virtual void correct()
Update properties.
virtual tmp< volVectorField > Kappa() const
Anisotropic thermal conductivity [W/m/K].
dynamicFvMesh & mesh
fvPatchField< scalar > fvPatchScalarField
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)
const dimensionSet dimEnergy
Internal & ref()
Return a reference to the dimensioned internal field.
label patchi
Energy for a solid mixture.
Definition: heSolidThermo.H:49
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
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
virtual ~heSolidThermo()
Destructor.
#define InfoInFunction
Report an information message using Foam::Info.