heThermo.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-2015 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 "heThermo.H"
29 
30 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
31 
32 template<class BasicThermo, class MixtureType>
35 {
36  volScalarField::GeometricBoundaryField& hbf = h.boundaryField();
37 
38  forAll(hbf, patchi)
39  {
40  if (isA<gradientEnergyFvPatchScalarField>(hbf[patchi]))
41  {
42  refCast<gradientEnergyFvPatchScalarField>(hbf[patchi]).gradient()
43  = hbf[patchi].fvPatchField::snGrad();
44  }
45  else if (isA<mixedEnergyFvPatchScalarField>(hbf[patchi]))
46  {
47  refCast<mixedEnergyFvPatchScalarField>(hbf[patchi]).refGrad()
48  = hbf[patchi].fvPatchField::snGrad();
49  }
50  }
51 }
52 
53 
54 template<class BasicThermo, class MixtureType>
56 {
57  scalarField& heCells = he_.internalField();
58  const scalarField& pCells = this->p_.internalField();
59  const scalarField& TCells = this->T_.internalField();
60 
61  forAll(heCells, celli)
62  {
63  heCells[celli] =
64  this->cellMixture(celli).HE(pCells[celli], TCells[celli]);
65  }
66 
67  forAll(he_.boundaryField(), patchi)
68  {
69  he_.boundaryField()[patchi] == he
70  (
71  this->p_.boundaryField()[patchi],
72  this->T_.boundaryField()[patchi],
73  patchi
74  );
75  }
76 
77  this->heBoundaryCorrection(he_);
78 }
79 
80 
81 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
82 
83 
84 template<class BasicThermo, class MixtureType>
86 (
87  const fvMesh& mesh,
88  const word& phaseName
89 )
90 :
91  BasicThermo(mesh, phaseName),
92  MixtureType(*this, mesh, phaseName),
93 
94  he_
95  (
96  IOobject
97  (
98  BasicThermo::phasePropertyName
99  (
100  MixtureType::thermoType::heName()
101  ),
102  mesh.time().timeName(),
103  mesh,
104  IOobject::NO_READ,
105  IOobject::NO_WRITE
106  ),
107  mesh,
109  this->heBoundaryTypes(),
110  this->heBoundaryBaseTypes()
111  )
112 {
113  init();
114 }
115 
116 
117 template<class BasicThermo, class MixtureType>
119 (
120  const fvMesh& mesh,
121  const dictionary& dict,
122  const word& phaseName
123 )
124 :
125  BasicThermo(mesh, dict, phaseName),
126  MixtureType(*this, mesh, phaseName),
127 
128  he_
129  (
130  IOobject
131  (
132  BasicThermo::phasePropertyName
133  (
134  MixtureType::thermoType::heName()
135  ),
136  mesh.time().timeName(),
137  mesh,
138  IOobject::NO_READ,
139  IOobject::NO_WRITE
140  ),
141  mesh,
143  this->heBoundaryTypes(),
144  this->heBoundaryBaseTypes()
145  )
146 {
147  init();
148 }
149 
150 
151 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
152 
153 template<class BasicThermo, class MixtureType>
155 {}
156 
157 
158 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
159 
160 template<class BasicThermo, class MixtureType>
162 (
163  const volScalarField& p,
164  const volScalarField& T
165 ) const
166 {
167  const fvMesh& mesh = this->T_.mesh();
168 
170  (
171  new volScalarField
172  (
173  IOobject
174  (
175  "he",
176  mesh.time().timeName(),
177  mesh,
178  IOobject::NO_READ,
179  IOobject::NO_WRITE,
180  false
181  ),
182  mesh,
183  he_.dimensions()
184  )
185  );
186 
187  volScalarField& he = the();
188  scalarField& heCells = he.internalField();
189  const scalarField& pCells = p.internalField();
190  const scalarField& TCells = T.internalField();
191 
192  forAll(heCells, celli)
193  {
194  heCells[celli] =
195  this->cellMixture(celli).HE(pCells[celli], TCells[celli]);
196  }
197 
198  forAll(he.boundaryField(), patchi)
199  {
200  scalarField& hep = he.boundaryField()[patchi];
201  const scalarField& pp = p.boundaryField()[patchi];
202  const scalarField& Tp = T.boundaryField()[patchi];
203 
204  forAll(hep, facei)
205  {
206  hep[facei] =
207  this->patchFaceMixture(patchi, facei).HE(pp[facei], Tp[facei]);
208  }
209  }
210 
211  return the;
212 }
213 
214 
215 template<class BasicThermo, class MixtureType>
217 (
218  const scalarField& p,
219  const scalarField& T,
220  const labelList& cells
221 ) const
222 {
223  tmp<scalarField> the(new scalarField(T.size()));
224  scalarField& he = the();
225 
226  forAll(T, celli)
227  {
228  he[celli] = this->cellMixture(cells[celli]).HE(p[celli], T[celli]);
229  }
230 
231  return the;
232 }
233 
234 
235 template<class BasicThermo, class MixtureType>
237 (
238  const scalarField& p,
239  const scalarField& T,
240  const label patchi
241 ) const
242 {
243  tmp<scalarField> the(new scalarField(T.size()));
244  scalarField& he = the();
245 
246  forAll(T, facei)
247  {
248  he[facei] =
249  this->patchFaceMixture(patchi, facei).HE(p[facei], T[facei]);
250  }
251 
252  return the;
253 }
254 
255 
256 template<class BasicThermo, class MixtureType>
259 {
260  const fvMesh& mesh = this->T_.mesh();
261 
263  (
264  new volScalarField
265  (
266  IOobject
267  (
268  "hc",
269  mesh.time().timeName(),
270  mesh,
271  IOobject::NO_READ,
272  IOobject::NO_WRITE,
273  false
274  ),
275  mesh,
276  he_.dimensions()
277  )
278  );
279 
280  volScalarField& hcf = thc();
281  scalarField& hcCells = hcf.internalField();
282 
283  forAll(hcCells, celli)
284  {
285  hcCells[celli] = this->cellMixture(celli).Hc();
286  }
287 
288  forAll(hcf.boundaryField(), patchi)
289  {
290  scalarField& hcp = hcf.boundaryField()[patchi];
291 
292  forAll(hcp, facei)
293  {
294  hcp[facei] = this->patchFaceMixture(patchi, facei).Hc();
295  }
296  }
297 
298  return thc;
299 }
300 
301 
302 template<class BasicThermo, class MixtureType>
304 (
305  const scalarField& p,
306  const scalarField& T,
307  const label patchi
308 ) const
309 {
310  tmp<scalarField> tCp(new scalarField(T.size()));
311  scalarField& cp = tCp();
312 
313  forAll(T, facei)
314  {
315  cp[facei] =
316  this->patchFaceMixture(patchi, facei).Cp(p[facei], T[facei]);
317  }
318 
319  return tCp;
320 }
321 
322 
323 template<class BasicThermo, class MixtureType>
326 {
327  const fvMesh& mesh = this->T_.mesh();
328 
330  (
331  new volScalarField
332  (
333  IOobject
334  (
335  "Cp",
336  mesh.time().timeName(),
337  mesh,
338  IOobject::NO_READ,
339  IOobject::NO_WRITE,
340  false
341  ),
342  mesh,
344  )
345  );
346 
347  volScalarField& cp = tCp();
348 
349  forAll(this->T_, celli)
350  {
351  cp[celli] =
352  this->cellMixture(celli).Cp(this->p_[celli], this->T_[celli]);
353  }
354 
355  forAll(this->T_.boundaryField(), patchi)
356  {
357  const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
358  const fvPatchScalarField& pT = this->T_.boundaryField()[patchi];
359  fvPatchScalarField& pCp = cp.boundaryField()[patchi];
360 
361  forAll(pT, facei)
362  {
363  pCp[facei] =
364  this->patchFaceMixture(patchi, facei).Cp(pp[facei], pT[facei]);
365  }
366  }
367 
368  return tCp;
369 }
370 
371 
372 template<class BasicThermo, class MixtureType>
375 (
376  const scalarField& p,
377  const scalarField& T,
378  const label patchi
379 ) const
380 {
381  tmp<scalarField> tCv(new scalarField(T.size()));
382  scalarField& cv = tCv();
383 
384  forAll(T, facei)
385  {
386  cv[facei] =
387  this->patchFaceMixture(patchi, facei).Cv(p[facei], T[facei]);
388  }
389 
390  return tCv;
391 }
392 
393 
394 template<class BasicThermo, class MixtureType>
397 {
398  const fvMesh& mesh = this->T_.mesh();
399 
401  (
402  new volScalarField
403  (
404  IOobject
405  (
406  "Cv",
407  mesh.time().timeName(),
408  mesh,
409  IOobject::NO_READ,
410  IOobject::NO_WRITE,
411  false
412  ),
413  mesh,
415  )
416  );
417 
418  volScalarField& cv = tCv();
419 
420  forAll(this->T_, celli)
421  {
422  cv[celli] =
423  this->cellMixture(celli).Cv(this->p_[celli], this->T_[celli]);
424  }
425 
426  forAll(this->T_.boundaryField(), patchi)
427  {
428  cv.boundaryField()[patchi] = Cv
429  (
430  this->p_.boundaryField()[patchi],
431  this->T_.boundaryField()[patchi],
432  patchi
433  );
434  }
435 
436  return tCv;
437 }
438 
439 
440 template<class BasicThermo, class MixtureType>
442 (
443  const scalarField& p,
444  const scalarField& T,
445  const label patchi
446 ) const
447 {
448  tmp<scalarField> tgamma(new scalarField(T.size()));
449  scalarField& cpv = tgamma();
450 
451  forAll(T, facei)
452  {
453  cpv[facei] =
454  this->patchFaceMixture(patchi, facei).gamma(p[facei], T[facei]);
455  }
456 
457  return tgamma;
458 }
459 
460 
461 template<class BasicThermo, class MixtureType>
464 {
465  const fvMesh& mesh = this->T_.mesh();
466 
467  tmp<volScalarField> tgamma
468  (
469  new volScalarField
470  (
471  IOobject
472  (
473  "gamma",
474  mesh.time().timeName(),
475  mesh,
476  IOobject::NO_READ,
477  IOobject::NO_WRITE,
478  false
479  ),
480  mesh,
481  dimless
482  )
483  );
484 
485  volScalarField& cpv = tgamma();
486 
487  forAll(this->T_, celli)
488  {
489  cpv[celli] =
490  this->cellMixture(celli).gamma(this->p_[celli], this->T_[celli]);
491  }
492 
493  forAll(this->T_.boundaryField(), patchi)
494  {
495  const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
496  const fvPatchScalarField& pT = this->T_.boundaryField()[patchi];
497  fvPatchScalarField& pgamma = cpv.boundaryField()[patchi];
498 
499  forAll(pT, facei)
500  {
501  pgamma[facei] = this->patchFaceMixture(patchi, facei).gamma
502  (
503  pp[facei],
504  pT[facei]
505  );
506  }
507  }
508 
509  return tgamma;
510 }
511 
512 
513 template<class BasicThermo, class MixtureType>
515 (
516  const scalarField& p,
517  const scalarField& T,
518  const label patchi
519 ) const
520 {
521  tmp<scalarField> tCpv(new scalarField(T.size()));
522  scalarField& cpv = tCpv();
523 
524  forAll(T, facei)
525  {
526  cpv[facei] =
527  this->patchFaceMixture(patchi, facei).Cpv(p[facei], T[facei]);
528  }
529 
530  return tCpv;
531 }
532 
533 
534 template<class BasicThermo, class MixtureType>
537 {
538  const fvMesh& mesh = this->T_.mesh();
539 
541  (
542  new volScalarField
543  (
544  IOobject
545  (
546  "Cpv",
547  mesh.time().timeName(),
548  mesh,
549  IOobject::NO_READ,
550  IOobject::NO_WRITE,
551  false
552  ),
553  mesh,
555  )
556  );
557 
558  volScalarField& cpv = tCpv();
559 
560  forAll(this->T_, celli)
561  {
562  cpv[celli] =
563  this->cellMixture(celli).Cpv(this->p_[celli], this->T_[celli]);
564  }
565 
566  forAll(this->T_.boundaryField(), patchi)
567  {
568  const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
569  const fvPatchScalarField& pT = this->T_.boundaryField()[patchi];
570  fvPatchScalarField& pCpv = cpv.boundaryField()[patchi];
571 
572  forAll(pT, facei)
573  {
574  pCpv[facei] =
575  this->patchFaceMixture(patchi, facei).Cpv(pp[facei], pT[facei]);
576  }
577  }
578 
579  return tCpv;
580 }
581 
582 
583 template<class BasicThermo, class MixtureType>
585 (
586  const scalarField& p,
587  const scalarField& T,
588  const label patchi
589 ) const
590 {
591  tmp<scalarField> tCpByCpv(new scalarField(T.size()));
592  scalarField& cpByCpv = tCpByCpv();
593 
594  forAll(T, facei)
595  {
596  cpByCpv[facei] =
597  this->patchFaceMixture(patchi, facei).cpBycpv(p[facei], T[facei]);
598  }
599 
600  return tCpByCpv;
601 }
602 
603 
604 template<class BasicThermo, class MixtureType>
607 {
608  const fvMesh& mesh = this->T_.mesh();
609 
610  tmp<volScalarField> tCpByCpv
611  (
612  new volScalarField
613  (
614  IOobject
615  (
616  "CpByCpv",
617  mesh.time().timeName(),
618  mesh,
619  IOobject::NO_READ,
620  IOobject::NO_WRITE,
621  false
622  ),
623  mesh,
624  dimless
625  )
626  );
627 
628  volScalarField& cpByCpv = tCpByCpv();
629 
630  forAll(this->T_, celli)
631  {
632  cpByCpv[celli] = this->cellMixture(celli).cpBycpv
633  (
634  this->p_[celli],
635  this->T_[celli]
636  );
637  }
638 
639  forAll(this->T_.boundaryField(), patchi)
640  {
641  const fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
642  const fvPatchScalarField& pT = this->T_.boundaryField()[patchi];
643  fvPatchScalarField& pCpByCpv = cpByCpv.boundaryField()[patchi];
644 
645  forAll(pT, facei)
646  {
647  pCpByCpv[facei] = this->patchFaceMixture(patchi, facei).cpBycpv
648  (
649  pp[facei],
650  pT[facei]
651  );
652  }
653  }
654 
655  return tCpByCpv;
656 }
657 
658 
659 template<class BasicThermo, class MixtureType>
661 (
662  const scalarField& h,
663  const scalarField& p,
664  const scalarField& T0,
665  const labelList& cells
666 ) const
667 {
668  tmp<scalarField> tT(new scalarField(h.size()));
669  scalarField& T = tT();
670 
671  forAll(h, celli)
672  {
673  T[celli] =
674  this->cellMixture(cells[celli]).THE(h[celli], p[celli], T0[celli]);
675  }
676 
677  return tT;
678 }
679 
680 
681 template<class BasicThermo, class MixtureType>
683 (
684  const scalarField& h,
685  const scalarField& p,
686  const scalarField& T0,
687  const label patchi
688 ) const
689 {
690 
691  tmp<scalarField> tT(new scalarField(h.size()));
692  scalarField& T = tT();
693  forAll(h, facei)
694  {
695  T[facei] = this->patchFaceMixture
696  (
697  patchi,
698  facei
699  ).THE(h[facei], p[facei], T0[facei]);
700  }
701 
702  return tT;
703 }
704 
705 
706 template<class BasicThermo, class MixtureType>
709 {
710  tmp<Foam::volScalarField> kappa(Cp()*this->alpha_);
711  kappa().rename("kappa");
712  return kappa;
713 }
714 
715 
716 template<class BasicThermo, class MixtureType>
718 (
719  const label patchi
720 ) const
721 {
722  return
723  Cp
724  (
725  this->p_.boundaryField()[patchi],
726  this->T_.boundaryField()[patchi],
727  patchi
728  )*this->alpha_.boundaryField()[patchi];
729 }
730 
731 
732 template<class BasicThermo, class MixtureType>
735 (
736  const volScalarField& alphat
737 ) const
738 {
739  tmp<Foam::volScalarField> kappaEff(Cp()*(this->alpha_ + alphat));
740  kappaEff().rename("kappaEff");
741  return kappaEff;
742 }
743 
744 
745 template<class BasicThermo, class MixtureType>
748 (
749  const scalarField& alphat,
750  const label patchi
751 ) const
752 {
753  return
754  Cp
755  (
756  this->p_.boundaryField()[patchi],
757  this->T_.boundaryField()[patchi],
758  patchi
759  )
760  *(
761  this->alpha_.boundaryField()[patchi]
762  + alphat
763  );
764 }
765 
766 
767 template<class BasicThermo, class MixtureType>
770 (
771  const volScalarField& alphat
772 ) const
773 {
774  tmp<Foam::volScalarField> alphaEff(this->CpByCpv()*(this->alpha_ + alphat));
775  alphaEff().rename("alphaEff");
776  return alphaEff;
777 }
778 
779 
780 template<class BasicThermo, class MixtureType>
783 (
784  const scalarField& alphat,
785  const label patchi
786 ) const
787 {
788  return
789  this->CpByCpv
790  (
791  this->p_.boundaryField()[patchi],
792  this->T_.boundaryField()[patchi],
793  patchi
794  )
795  *(
796  this->alpha_.boundaryField()[patchi]
797  + alphat
798  );
799 }
800 
801 
802 template<class BasicThermo, class MixtureType>
804 {
805  if (BasicThermo::read())
806  {
807  MixtureType::read(*this);
808  return true;
809  }
810  else
811  {
812  return false;
813  }
814 }
815 
816 
817 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::heThermo::gamma
virtual tmp< volScalarField > gamma() const
Gamma = Cp/Cv [].
Definition: heThermo.C:463
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::heThermo::he
virtual volScalarField & he()
Enthalpy/Internal energy [J/kg].
Definition: heThermo.H:141
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::dimEnergy
const dimensionSet dimEnergy
Foam::heThermo::hc
virtual tmp< volScalarField > hc() const
Chemical enthalpy [J/kg].
Definition: heThermo.C:258
Foam::read
bool read(const char *, int32_t &)
Definition: int32IO.C:87
gradientEnergyFvPatchScalarField.H
alphat
alphat
Definition: TEqn.H:2
Foam::heThermo::kappa
virtual tmp< volScalarField > kappa() const
Thermal diffusivity for temperature of mixture [J/m/s/K].
Definition: heThermo.C:708
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::heThermo::CpByCpv
virtual tmp< volScalarField > CpByCpv() const
Heat capacity ratio [].
Definition: heThermo.C:606
Foam::constant::electromagnetic::kappa
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Foam::heThermo::Cp
virtual tmp< volScalarField > Cp() const
Heat capacity at constant pressure [J/kg/K].
Definition: heThermo.C:325
Foam::heThermo::~heThermo
virtual ~heThermo()
Destructor.
Definition: heThermo.C:154
Foam::heThermo::Cpv
virtual tmp< volScalarField > Cpv() const
Heat capacity at constant pressure/volume [J/kg/K].
Definition: heThermo.C:536
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::constant::universal::h
const dimensionedScalar h
Planck constant.
Definition: createFields.H:6
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
cp
const volScalarField & cp
Definition: setRegionSolidFields.H:8
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::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
alphaEff
const volScalarField & alphaEff
Definition: setAlphaEff.H:93
Foam::heThermo::read
virtual bool read()
Read thermophysical properties dictionary.
Definition: heThermo.C:803
Foam::heThermo::Cv
virtual tmp< volScalarField > Cv() const
Heat capacity at constant volume [J/kg/K].
Definition: heThermo.C:396
mixedEnergyFvPatchScalarField.H
Foam::heThermo::init
void init()
Initialize heThermo.
Definition: heThermo.C:55
T
const volScalarField & T
Definition: createFields.H:25
he
volScalarField & he
Definition: YEEqn.H:56
Foam::heThermo::heThermo
heThermo(const heThermo< BasicThermo, MixtureType > &)
Construct as copy (not implemented)
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
scalarField
volScalarField scalarField(fieldObject, mesh)
Foam::heThermo::alphaEff
virtual tmp< volScalarField > alphaEff(const volScalarField &alphat) const
Effective thermal diffusivity of mixture [kg/m/s].
Definition: heThermo.C:770
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::dimTemperature
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::heThermo::THE
virtual tmp< scalarField > THE(const scalarField &he, const scalarField &p, const scalarField &T0, const labelList &cells) const
Temperature from enthalpy/internal energy for cell-set.
Definition: heThermo.C:661
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
heThermo.H
Foam::heThermo::kappaEff
virtual tmp< volScalarField > kappaEff(const volScalarField &) const
Effective thermal diffusivity for temperature.
Definition: heThermo.C:735
Foam::GeometricField::GeometricBoundaryField
Definition: GeometricField.H:105
Foam::heThermo::heBoundaryCorrection
void heBoundaryCorrection(volScalarField &he)
Correct the enthalpy/internal energy field boundaries.
Definition: heThermo.C:34