twoPhaseMixtureThermo.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) 2013-2017 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 "twoPhaseMixtureThermo.H"
29 
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  defineTypeNameAndDebug(twoPhaseMixtureThermo, 0);
36 }
37 
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
42 (
43  const volVectorField& U,
44  const surfaceScalarField& phi
45 )
46 :
47  psiThermo(U.mesh(), word::null),
48  twoPhaseMixture(U.mesh(), *this),
49  interfaceProperties(alpha1(), U, *this),
50  thermo1_(nullptr),
51  thermo2_(nullptr)
52 {
53  {
54  volScalarField T1(IOobject::groupName("T", phase1Name()), T_);
55  T1.write();
56  }
57 
58  {
59  volScalarField T2(IOobject::groupName("T", phase2Name()), T_);
60  T2.write();
61  }
62 
63  thermo1_ = rhoThermo::New(U.mesh(), phase1Name());
64  thermo2_ = rhoThermo::New(U.mesh(), phase2Name());
65 
66  // thermo1_->validate(phase1Name(), "e");
67  // thermo2_->validate(phase2Name(), "e");
68 
69  correct();
70 }
71 
72 
73 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
74 
76 {}
77 
78 
79 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
80 
82 {
83  thermo1_->he() = thermo1_->he(p_, T_);
84  thermo1_->correct();
85 
86  thermo2_->he() = thermo2_->he(p_, T_);
87  thermo2_->correct();
88 }
89 
90 
92 {
93  psi_ = alpha1()*thermo1_->psi() + alpha2()*thermo2_->psi();
94  mu_ = alpha1()*thermo1_->mu() + alpha2()*thermo2_->mu();
95  alpha_ = alpha1()*thermo1_->alpha() + alpha2()*thermo2_->alpha();
96 
98 }
99 
100 
102 {
103  return thermo1_->incompressible() && thermo2_->incompressible();
104 }
105 
106 
108 {
109  return thermo1_->isochoric() && thermo2_->isochoric();
110 }
111 
112 
114 (
115  const volScalarField& p,
116  const volScalarField& T
117 ) const
118 {
119  return alpha1()*thermo1_->he(p, T) + alpha2()*thermo2_->he(p, T);
120 }
121 
122 
124 (
125  const scalarField& p,
126  const scalarField& T,
127  const labelList& cells
128 ) const
129 {
130  return
131  scalarField(alpha1(), cells)*thermo1_->he(p, T, cells)
132  + scalarField(alpha2(), cells)*thermo2_->he(p, T, cells);
133 }
134 
135 
137 (
138  const scalarField& p,
139  const scalarField& T,
140  const label patchi
141 ) const
142 {
143  return
144  alpha1().boundaryField()[patchi]*thermo1_->he(p, T, patchi)
145  + alpha2().boundaryField()[patchi]*thermo2_->he(p, T, patchi);
146 }
147 
148 
150 {
151  return alpha1()*thermo1_->hc() + alpha2()*thermo2_->hc();
152 }
153 
154 
156 (
157  const scalarField& h,
158  const scalarField& p,
159  const scalarField& T0,
160  const labelList& cells
161 ) const
162 {
164  return T0;
165 }
166 
167 
169 (
170  const scalarField& h,
171  const scalarField& p,
172  const scalarField& T0,
173  const label patchi
174 ) const
175 {
177  return T0;
178 }
179 
180 
182 {
183  return alpha1()*thermo1_->Cp() + alpha2()*thermo2_->Cp();
184 }
185 
186 
188 (
189  const scalarField& p,
190  const scalarField& T,
191  const label patchi
192 ) const
193 {
194  return
195  alpha1().boundaryField()[patchi]*thermo1_->Cp(p, T, patchi)
196  + alpha2().boundaryField()[patchi]*thermo2_->Cp(p, T, patchi);
197 }
198 
199 
201 {
202  return alpha1()*thermo1_->Cv() + alpha2()*thermo2_->Cv();
203 }
204 
205 
207 (
208  const scalarField& p,
209  const scalarField& T,
210  const label patchi
211 ) const
212 {
213  return
214  alpha1().boundaryField()[patchi]*thermo1_->Cv(p, T, patchi)
215  + alpha2().boundaryField()[patchi]*thermo2_->Cv(p, T, patchi);
216 }
217 
218 
220 {
221  return alpha1()*thermo1_->gamma() + alpha2()*thermo2_->gamma();
222 }
223 
224 
226 (
227  const scalarField& p,
228  const scalarField& T,
229  const label patchi
230 ) const
231 {
232  return
233  alpha1().boundaryField()[patchi]*thermo1_->gamma(p, T, patchi)
234  + alpha2().boundaryField()[patchi]*thermo2_->gamma(p, T, patchi);
235 }
236 
237 
239 {
240  return alpha1()*thermo1_->Cpv() + alpha2()*thermo2_->Cpv();
241 }
242 
243 
245 (
246  const scalarField& p,
247  const scalarField& T,
248  const label patchi
249 ) const
250 {
251  return
252  alpha1().boundaryField()[patchi]*thermo1_->Cpv(p, T, patchi)
253  + alpha2().boundaryField()[patchi]*thermo2_->Cpv(p, T, patchi);
254 }
255 
256 
258 {
259  return
260  alpha1()*thermo1_->CpByCpv()
261  + alpha2()*thermo2_->CpByCpv();
262 }
263 
264 
266 (
267  const scalarField& p,
268  const scalarField& T,
269  const label patchi
270 ) const
271 {
272  return
273  alpha1().boundaryField()[patchi]*thermo1_->CpByCpv(p, T, patchi)
274  + alpha2().boundaryField()[patchi]*thermo2_->CpByCpv(p, T, patchi);
275 }
276 
277 
279 {
280  return mu()/(alpha1()*thermo1_->rho() + alpha2()*thermo2_->rho());
281 }
282 
283 
285 (
286  const label patchi
287 ) const
288 {
289  return
290  mu(patchi)
291  /(
292  alpha1().boundaryField()[patchi]*thermo1_->rho(patchi)
293  + alpha2().boundaryField()[patchi]*thermo2_->rho(patchi)
294  );
295 }
296 
297 
299 {
300  return alpha1()*thermo1_->kappa() + alpha2()*thermo2_->kappa();
301 }
302 
303 
305 (
306  const label patchi
307 ) const
308 {
309  return
310  alpha1().boundaryField()[patchi]*thermo1_->kappa(patchi)
311  + alpha2().boundaryField()[patchi]*thermo2_->kappa(patchi);
312 }
313 
314 
316 (
317  const volScalarField& alphat
318 ) const
319 {
320  return
321  alpha1()*thermo1_->kappaEff(alphat)
322  + alpha2()*thermo2_->kappaEff(alphat);
323 }
324 
325 
327 (
328  const scalarField& alphat,
329  const label patchi
330 ) const
331 {
332  return
333  alpha1().boundaryField()[patchi]*thermo1_->kappaEff(alphat, patchi)
334  + alpha2().boundaryField()[patchi]*thermo2_->kappaEff(alphat, patchi)
335  ;
336 }
337 
338 
340 (
341  const volScalarField& alphat
342 ) const
343 {
344  return
345  alpha1()*thermo1_->alphaEff(alphat)
346  + alpha2()*thermo2_->alphaEff(alphat);
347 }
348 
349 
351 (
352  const scalarField& alphat,
353  const label patchi
354 ) const
355 {
356  return
357  alpha1().boundaryField()[patchi]*thermo1_->alphaEff(alphat, patchi)
358  + alpha2().boundaryField()[patchi]*thermo2_->alphaEff(alphat, patchi)
359  ;
360 }
361 
362 
364 {
365  if (psiThermo::read())
366  {
367  return interfaceProperties::read();
368  }
369  else
370  {
371  return false;
372  }
373 }
374 
375 
376 // ************************************************************************* //
virtual tmp< scalarField > THE(const scalarField &h, const scalarField &p, const scalarField &T0, const labelList &cells) const
Temperature from enthalpy/internal energy for cell-set.
virtual tmp< volScalarField > Cp() const
Heat capacity at constant pressure [J/kg/K].
twoPhaseMixtureThermo(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
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
virtual tmp< volScalarField > alphaEff(const volScalarField &alphat) const
Effective thermal diffusivity of mixture [J/m/s/K].
virtual tmp< volScalarField > hc() const
Chemical enthalpy [J/kg].
static autoPtr< rhoThermo > New(const fvMesh &, const word &phaseName=word::null)
Selector.
Definition: rhoThermo.C:141
bool read()
Read transportProperties dictionary.
virtual void correctThermo()
Correct the thermodynamics of each phase.
virtual bool read()
Read base transportProperties dictionary.
virtual tmp< volScalarField > gamma() const
Gamma = Cp/Cv [].
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
virtual bool isochoric() const
Return true if the equation of state is isochoric.
virtual tmp< volScalarField > Cv() const
Heat capacity at constant volume [J/kg/K].
dynamicFvMesh & mesh
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
virtual ~twoPhaseMixtureThermo()
Destructor.
List< label > labelList
A List of labels.
Definition: labelList.H:56
alpha2
Definition: alphaEqn.H:112
virtual tmp< volScalarField > CpByCpv() const
Heat capacity ratio [].
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/ubuntu/OpenFOAM-5.0/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:72
defineTypeNameAndDebug(combustionModel, 0)
virtual tmp< volScalarField > Cpv() const
Heat capacity at constant pressure/volume [J/kg/K].
volScalarField & alpha1
virtual volScalarField & he()
Enthalpy/Internal energy [J/kg].
const dimensionedScalar mu
Atomic mass unit.
label patchi
volScalarField & h
Planck constant.
virtual void correct()
Update mixture properties.
virtual bool incompressible() const
Return true if the equation of state is incompressible.
virtual tmp< volScalarField > kappa() const
Thermal diffusivity for temperature of mixture [J/m/s/K].
virtual tmp< volScalarField > kappaEff(const volScalarField &alphat) const
Effective thermal diffusivity of mixture [J/m/s/K].
virtual tmp< volScalarField > nu() const
Kinematic viscosity of mixture [m^2/s].
virtual bool read()
Read thermophysical properties dictionary.
Definition: basicThermo.C:513
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:366
Namespace for OpenFOAM.