hePsiThermo.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 "hePsiThermo.H"
27 
28 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
30 template<class BasicPsiThermo, class MixtureType>
32 {
33  const scalarField& hCells = this->he_.internalField();
34  const scalarField& pCells = this->p_.internalField();
35 
36  scalarField& TCells = this->T_.internalField();
37  scalarField& psiCells = this->psi_.internalField();
38  scalarField& muCells = this->mu_.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  TCells[celli] = mixture_.THE
47  (
48  hCells[celli],
49  pCells[celli],
50  TCells[celli]
51  );
52 
53  psiCells[celli] = mixture_.psi(pCells[celli], TCells[celli]);
54 
55  muCells[celli] = mixture_.mu(pCells[celli], TCells[celli]);
56  alphaCells[celli] = mixture_.alphah(pCells[celli], TCells[celli]);
57  }
58 
59  forAll(this->T_.boundaryField(), patchi)
60  {
61  fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
62  fvPatchScalarField& pT = this->T_.boundaryField()[patchi];
63  fvPatchScalarField& ppsi = this->psi_.boundaryField()[patchi];
64 
65  fvPatchScalarField& ph = this->he_.boundaryField()[patchi];
66 
67  fvPatchScalarField& pmu = this->mu_.boundaryField()[patchi];
68  fvPatchScalarField& palpha = this->alpha_.boundaryField()[patchi];
69 
70  if (pT.fixesValue())
71  {
72  forAll(pT, facei)
73  {
74  const typename MixtureType::thermoType& mixture_ =
75  this->patchFaceMixture(patchi, facei);
76 
77  ph[facei] = mixture_.HE(pp[facei], pT[facei]);
78 
79  ppsi[facei] = mixture_.psi(pp[facei], pT[facei]);
80  pmu[facei] = mixture_.mu(pp[facei], pT[facei]);
81  palpha[facei] = mixture_.alphah(pp[facei], pT[facei]);
82  }
83  }
84  else
85  {
86  forAll(pT, facei)
87  {
88  const typename MixtureType::thermoType& mixture_ =
89  this->patchFaceMixture(patchi, facei);
90 
91  pT[facei] = mixture_.THE(ph[facei], pp[facei], pT[facei]);
92 
93  ppsi[facei] = mixture_.psi(pp[facei], pT[facei]);
94  pmu[facei] = mixture_.mu(pp[facei], pT[facei]);
95  palpha[facei] = mixture_.alphah(pp[facei], pT[facei]);
96  }
97  }
98  }
99 }
100 
101 
102 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
103 
104 template<class BasicPsiThermo, class MixtureType>
106 (
107  const fvMesh& mesh,
108  const word& phaseName
109 )
110 :
112 {
113  calculate();
114 
115  // Switch on saving old time
116  this->psi_.oldTime();
117 }
118 
119 
120 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
121 
122 template<class BasicPsiThermo, class MixtureType>
124 {}
125 
126 
127 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
128 
129 template<class BasicPsiThermo, class MixtureType>
131 {
132  if (debug)
133  {
134  Info<< "entering hePsiThermo<BasicPsiThermo, MixtureType>::correct()"
135  << endl;
136  }
137 
138  // force the saving of the old-time values
139  this->psi_.oldTime();
140 
141  calculate();
142 
143  if (debug)
144  {
145  Info<< "exiting hePsiThermo<BasicPsiThermo, MixtureType>::correct()"
146  << endl;
147  }
148 }
149 
150 
151 // ************************************************************************* //
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::heThermo< BasicPsiThermo, MixtureType >
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::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::hePsiThermo::hePsiThermo
hePsiThermo(const hePsiThermo< BasicPsiThermo, MixtureType > &)
Construct as copy (not implemented)
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::Info
messageStream Info
Foam::hePsiThermo::calculate
void calculate()
Calculate the thermo variables.
Definition: hePsiThermo.C:31
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::hePsiThermo::~hePsiThermo
virtual ~hePsiThermo()
Destructor.
Definition: hePsiThermo.C:123
patchi
label patchi
Definition: getPatchFieldScalar.H:1
hePsiThermo.H
Foam::hePsiThermo::correct
virtual void correct()
Update properties.
Definition: hePsiThermo.C:130