heRhoThermo.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 "heRhoThermo.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& rhoCells = this->rho_.internalField();
39  scalarField& muCells = this->mu_.internalField();
40  scalarField& alphaCells = this->alpha_.internalField();
41 
42  forAll(TCells, celli)
43  {
44  const typename MixtureType::thermoType& mixture_ =
45  this->cellMixture(celli);
46 
47  TCells[celli] = mixture_.THE
48  (
49  hCells[celli],
50  pCells[celli],
51  TCells[celli]
52  );
53 
54  psiCells[celli] = mixture_.psi(pCells[celli], TCells[celli]);
55  rhoCells[celli] = mixture_.rho(pCells[celli], TCells[celli]);
56 
57  muCells[celli] = mixture_.mu(pCells[celli], TCells[celli]);
58  alphaCells[celli] = mixture_.alphah(pCells[celli], TCells[celli]);
59  }
60 
61  forAll(this->T_.boundaryField(), patchi)
62  {
63  fvPatchScalarField& pp = this->p_.boundaryField()[patchi];
64  fvPatchScalarField& pT = this->T_.boundaryField()[patchi];
65  fvPatchScalarField& ppsi = this->psi_.boundaryField()[patchi];
66  fvPatchScalarField& prho = this->rho_.boundaryField()[patchi];
67 
68  fvPatchScalarField& ph = this->he().boundaryField()[patchi];
69 
70  fvPatchScalarField& pmu = this->mu_.boundaryField()[patchi];
71  fvPatchScalarField& palpha = this->alpha_.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  ph[facei] = mixture_.HE(pp[facei], pT[facei]);
81 
82  ppsi[facei] = mixture_.psi(pp[facei], pT[facei]);
83  prho[facei] = mixture_.rho(pp[facei], pT[facei]);
84  pmu[facei] = mixture_.mu(pp[facei], pT[facei]);
85  palpha[facei] = mixture_.alphah(pp[facei], pT[facei]);
86  }
87  }
88  else
89  {
90  forAll(pT, facei)
91  {
92  const typename MixtureType::thermoType& mixture_ =
93  this->patchFaceMixture(patchi, facei);
94 
95  pT[facei] = mixture_.THE(ph[facei], pp[facei], pT[facei]);
96 
97  ppsi[facei] = mixture_.psi(pp[facei], pT[facei]);
98  prho[facei] = mixture_.rho(pp[facei], pT[facei]);
99  pmu[facei] = mixture_.mu(pp[facei], pT[facei]);
100  palpha[facei] = mixture_.alphah(pp[facei], pT[facei]);
101  }
102  }
103  }
104 }
105 
106 
107 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
108 
109 template<class BasicPsiThermo, class MixtureType>
111 (
112  const fvMesh& mesh,
113  const word& phaseName
114 )
115 :
117 {
118  calculate();
119 }
120 
121 
122 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
123 
124 template<class BasicPsiThermo, class MixtureType>
126 {}
127 
128 
129 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
130 
131 template<class BasicPsiThermo, class MixtureType>
133 {
134  if (debug)
135  {
136  Info<< "entering heRhoThermo<MixtureType>::correct()" << endl;
137  }
138 
139  calculate();
140 
141  if (debug)
142  {
143  Info<< "exiting heRhoThermo<MixtureType>::correct()" << endl;
144  }
145 }
146 
147 
148 // ************************************************************************* //
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
heRhoThermo.H
Foam::heRhoThermo::~heRhoThermo
virtual ~heRhoThermo()
Destructor.
Definition: heRhoThermo.C:125
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::Info
messageStream Info
Foam::heRhoThermo::heRhoThermo
heRhoThermo(const heRhoThermo< BasicPsiThermo, MixtureType > &)
Construct as copy (not implemented)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
he
volScalarField & he
Definition: YEEqn.H:56
Foam::heRhoThermo::calculate
void calculate()
Calculate the thermo variables.
Definition: heRhoThermo.C:31
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::heRhoThermo::correct
virtual void correct()
Update properties.
Definition: heRhoThermo.C:132