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-2012 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_.internalField();
35 
36  const scalarField& hCells = this->he_.internalField();
37  const scalarField& pCells = this->p_.internalField();
38  scalarField& rhoCells = this->rho_.internalField();
39  scalarField& alphaCells = this->alpha_.internalField();
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  forAll(this->T_.boundaryField(), patchi)
65  {
66  fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
67  fvPatchScalarField& pT = this->T_.boundaryField()[patchi];
68  fvPatchScalarField& prho = this->rho_.boundaryField()[patchi];
69  fvPatchScalarField& palpha = this->alpha_.boundaryField()[patchi];
70 
71  fvPatchScalarField& ph = this->he_.boundaryField()[patchi];
72 
73  if (pT.fixesValue())
74  {
75  forAll(pT, facei)
76  {
77  const typename MixtureType::thermoType& mixture_ =
78  this->patchFaceMixture(patchi, facei);
79 
80  const typename MixtureType::thermoType& volMixture_ =
81  this->patchFaceVolMixture
82  (
83  pp[facei],
84  pT[facei],
85  patchi,
86  facei
87  );
88 
89 
90  ph[facei] = mixture_.HE(pp[facei], pT[facei]);
91  prho[facei] = volMixture_.rho(pp[facei], pT[facei]);
92 
93  palpha[facei] =
94  volMixture_.kappa(pp[facei], pT[facei])
95  / mixture_.Cpv(pp[facei], pT[facei]);
96  }
97  }
98  else
99  {
100  forAll(pT, facei)
101  {
102  const typename MixtureType::thermoType& mixture_ =
103  this->patchFaceMixture(patchi, facei);
104 
105  const typename MixtureType::thermoType& volMixture_ =
106  this->patchFaceVolMixture
107  (
108  pp[facei],
109  pT[facei],
110  patchi,
111  facei
112  );
113 
114  pT[facei] = mixture_.THE(ph[facei], pp[facei] ,pT[facei]);
115  prho[facei] = volMixture_.rho(pp[facei], pT[facei]);
116 
117  palpha[facei] =
118  volMixture_.kappa(pp[facei], pT[facei])
119  / mixture_.Cpv(pp[facei], pT[facei]);
120  }
121  }
122  }
123 }
124 
125 
126 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
127 
128 template<class BasicSolidThermo, class MixtureType>
131 (
132  const fvMesh& mesh,
133  const word& phaseName
134 )
135 :
137 {
138  calculate();
139 }
140 
141 
142 template<class BasicSolidThermo, class MixtureType>
145 (
146  const fvMesh& mesh,
147  const dictionary& dict,
148  const word& phaseName
149 )
150 :
152 {
153  calculate();
154 }
155 
156 
157 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
158 
159 template<class BasicSolidThermo, class MixtureType>
161 {}
162 
163 
164 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
165 
166 template<class BasicSolidThermo, class MixtureType>
168 {
169  if (debug)
170  {
171  Info<< "entering heSolidThermo<MixtureType>::correct()" << endl;
172  }
173 
174  calculate();
175 
176  if (debug)
177  {
178  Info<< "exiting heSolidThermo<MixtureType>::correct()" << endl;
179  }
180 }
181 
182 
183 template<class BasicSolidThermo, class MixtureType>
186 {
187  const fvMesh& mesh = this->T_.mesh();
188 
189  tmp<volVectorField> tKappa
190  (
191  new volVectorField
192  (
193  IOobject
194  (
195  "Kappa",
196  mesh.time().timeName(),
197  mesh,
198  IOobject::NO_READ,
199  IOobject::NO_WRITE
200  ),
201  mesh,
203  )
204  );
205 
206  volVectorField& Kappa = tKappa();
207  vectorField& KappaCells = Kappa.internalField();
208  const scalarField& TCells = this->T_.internalField();
209  const scalarField& pCells = this->p_.internalField();
210 
211  forAll(KappaCells, celli)
212  {
213  Kappa[celli] =
214  this->cellVolMixture
215  (
216  pCells[celli],
217  TCells[celli],
218  celli
219  ).Kappa(pCells[celli], TCells[celli]);
220  }
221 
222  forAll(Kappa.boundaryField(), patchi)
223  {
224  vectorField& Kappap = Kappa.boundaryField()[patchi];
225  const scalarField& pT = this->T_.boundaryField()[patchi];
226  const scalarField& pp = this->p_.boundaryField()[patchi];
227 
228  forAll(Kappap, facei)
229  {
230  Kappap[facei] =
231  this->patchFaceVolMixture
232  (
233  pp[facei],
234  pT[facei],
235  patchi,
236  facei
237  ).Kappa(pp[facei], pT[facei]);
238  }
239  }
240 
241  return tKappa;
242 }
243 
244 
245 template<class BasicSolidThermo, class MixtureType>
248 (
249  const label patchi
250 ) const
251 {
252  const scalarField& pp = this->p_.boundaryField()[patchi];
253  const scalarField& Tp = this->T_.boundaryField()[patchi];
254  tmp<vectorField> tKappa(new vectorField(pp.size()));
255 
256  vectorField& Kappap = tKappa();
257 
258  forAll(Tp, facei)
259  {
260  Kappap[facei] =
261  this->patchFaceVolMixture
262  (
263  pp[facei],
264  Tp[facei],
265  patchi,
266  facei
267  ).Kappa(pp[facei], Tp[facei]);
268  }
269 
270  return tKappa;
271 }
272 
273 
274 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
volFields.H
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::heThermo< BasicSolidThermo, MixtureType >
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::dimEnergy
const dimensionSet dimEnergy
Foam::heSolidThermo::correct
virtual void correct()
Update properties.
Definition: heSolidThermo.C:167
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::fvPatchField::fixesValue
virtual bool fixesValue() const
Return true if this patch field fixes a value.
Definition: fvPatchField.H:336
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
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:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::Info
messageStream Info
Foam::heSolidThermo::~heSolidThermo
virtual ~heSolidThermo()
Destructor.
Definition: heSolidThermo.C:160
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
Foam::heSolidThermo::Kappa
virtual tmp< volVectorField > Kappa() const
Anisotropic thermal conductivity [W/m/K].
Definition: heSolidThermo.C:185
vectorField
volVectorField vectorField(fieldObject, mesh)
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::heSolidThermo::heSolidThermo
heSolidThermo(const heSolidThermo< BasicSolidThermo, MixtureType > &)
Construct as copy (not implemented)
heSolidThermo.H
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::dimTemperature
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::heSolidThermo::calculate
void calculate()
Calculate the thermo variables.
Definition: heSolidThermo.C:32