twoPhaseMixtureThermo.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2013-2018 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 #include "collatedFileOperation.H"
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  {
55  (
56  IOobject
57  (
58  IOobject::groupName("T", phase1Name()),
59  U.mesh().time().timeName(),
60  U.mesh()
61  ),
62  T_,
63  calculatedFvPatchScalarField::typeName
64  );
65  T1.write();
66  }
67 
68  {
70  (
71  IOobject
72  (
73  IOobject::groupName("T", phase2Name()),
74  U.mesh().time().timeName(),
75  U.mesh()
76  ),
77  T_,
78  calculatedFvPatchScalarField::typeName
79  );
80  T2.write();
81  }
82 
83  // Note: we're writing files to be read in immediately afterwards.
84  // Avoid any thread-writing problems.
85  fileHandler().flush();
86 
87  thermo1_ = rhoThermo::New(U.mesh(), phase1Name());
88  thermo2_ = rhoThermo::New(U.mesh(), phase2Name());
89 
90  // thermo1_->validate(phase1Name(), "e");
91  // thermo2_->validate(phase2Name(), "e");
92 
93  correct();
94 }
95 
96 
97 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
98 
100 {}
101 
102 
103 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
104 
106 {
107  thermo1_->T() = T_;
108  thermo1_->he() = thermo1_->he(p_, T_);
109  thermo1_->correct();
110 
111  thermo2_->T() = T_;
112  thermo2_->he() = thermo2_->he(p_, T_);
113  thermo2_->correct();
114 }
115 
116 
118 {
119  psi_ = alpha1()*thermo1_->psi() + alpha2()*thermo2_->psi();
120  mu_ = alpha1()*thermo1_->mu() + alpha2()*thermo2_->mu();
121  alpha_ = alpha1()*thermo1_->alpha() + alpha2()*thermo2_->alpha();
122 
124 }
125 
126 
128 {
129  return thermo1_->thermoName() + ',' + thermo2_->thermoName();
130 }
131 
132 
134 {
135  return thermo1_->incompressible() && thermo2_->incompressible();
136 }
137 
138 
140 {
141  return thermo1_->isochoric() && thermo2_->isochoric();
142 }
143 
144 
146 (
147  const volScalarField& p,
148  const volScalarField& T
149 ) const
150 {
151  return alpha1()*thermo1_->he(p, T) + alpha2()*thermo2_->he(p, T);
152 }
153 
154 
156 (
157  const scalarField& p,
158  const scalarField& T,
159  const labelList& cells
160 ) const
161 {
162  return
163  scalarField(alpha1(), cells)*thermo1_->he(p, T, cells)
164  + scalarField(alpha2(), cells)*thermo2_->he(p, T, cells);
165 }
166 
167 
169 (
170  const scalarField& p,
171  const scalarField& T,
172  const label patchi
173 ) const
174 {
175  return
176  alpha1().boundaryField()[patchi]*thermo1_->he(p, T, patchi)
177  + alpha2().boundaryField()[patchi]*thermo2_->he(p, T, patchi);
178 }
179 
180 
182 {
183  return alpha1()*thermo1_->hc() + alpha2()*thermo2_->hc();
184 }
185 
186 
188 (
189  const scalarField& h,
190  const scalarField& p,
191  const scalarField& T0,
192  const labelList& cells
193 ) const
194 {
196  return T0;
197 }
198 
199 
201 (
202  const scalarField& h,
203  const scalarField& p,
204  const scalarField& T0,
205  const label patchi
206 ) const
207 {
209  return T0;
210 }
211 
212 
214 {
215  return alpha1()*thermo1_->Cp() + alpha2()*thermo2_->Cp();
216 }
217 
218 
220 (
221  const scalarField& p,
222  const scalarField& T,
223  const label patchi
224 ) const
225 {
226  return
227  alpha1().boundaryField()[patchi]*thermo1_->Cp(p, T, patchi)
228  + alpha2().boundaryField()[patchi]*thermo2_->Cp(p, T, patchi);
229 }
230 
231 
233 {
234  return alpha1()*thermo1_->Cv() + alpha2()*thermo2_->Cv();
235 }
236 
237 
239 (
240  const scalarField& p,
241  const scalarField& T,
242  const label patchi
243 ) const
244 {
245  return
246  alpha1().boundaryField()[patchi]*thermo1_->Cv(p, T, patchi)
247  + alpha2().boundaryField()[patchi]*thermo2_->Cv(p, T, patchi);
248 }
249 
250 
252 {
253  return alpha1()*thermo1_->gamma() + alpha2()*thermo2_->gamma();
254 }
255 
256 
258 (
259  const scalarField& p,
260  const scalarField& T,
261  const label patchi
262 ) const
263 {
264  return
265  alpha1().boundaryField()[patchi]*thermo1_->gamma(p, T, patchi)
266  + alpha2().boundaryField()[patchi]*thermo2_->gamma(p, T, patchi);
267 }
268 
269 
271 {
272  return alpha1()*thermo1_->Cpv() + alpha2()*thermo2_->Cpv();
273 }
274 
275 
277 (
278  const scalarField& p,
279  const scalarField& T,
280  const label patchi
281 ) const
282 {
283  return
284  alpha1().boundaryField()[patchi]*thermo1_->Cpv(p, T, patchi)
285  + alpha2().boundaryField()[patchi]*thermo2_->Cpv(p, T, patchi);
286 }
287 
288 
290 {
291  return
292  alpha1()*thermo1_->CpByCpv()
293  + alpha2()*thermo2_->CpByCpv();
294 }
295 
296 
298 (
299  const scalarField& p,
300  const scalarField& T,
301  const label patchi
302 ) const
303 {
304  return
305  alpha1().boundaryField()[patchi]*thermo1_->CpByCpv(p, T, patchi)
306  + alpha2().boundaryField()[patchi]*thermo2_->CpByCpv(p, T, patchi);
307 }
308 
309 
311 {
312  return alpha1()*thermo1_->W() + alpha2()*thermo1_->W();
313 }
314 
315 
317 (
318  const label patchi
319 ) const
320 {
321  return
322  alpha1().boundaryField()[patchi]*thermo1_->W(patchi)
323  + alpha2().boundaryField()[patchi]*thermo1_->W(patchi);
324 }
325 
326 
328 {
329  return mu()/(alpha1()*thermo1_->rho() + alpha2()*thermo2_->rho());
330 }
331 
332 
334 (
335  const label patchi
336 ) const
337 {
338  return
339  mu(patchi)
340  /(
341  alpha1().boundaryField()[patchi]*thermo1_->rho(patchi)
342  + alpha2().boundaryField()[patchi]*thermo2_->rho(patchi)
343  );
344 }
345 
346 
348 {
349  return alpha1()*thermo1_->kappa() + alpha2()*thermo2_->kappa();
350 }
351 
352 
354 (
355  const label patchi
356 ) const
357 {
358  return
359  alpha1().boundaryField()[patchi]*thermo1_->kappa(patchi)
360  + alpha2().boundaryField()[patchi]*thermo2_->kappa(patchi);
361 }
362 
363 
365 {
366  return
367  alpha1()*thermo1_->alphahe()
368  + alpha2()*thermo2_->alphahe();
369 }
370 
371 
373 (
374  const label patchi
375 ) const
376 {
377  return
378  alpha1().boundaryField()[patchi]*thermo1_->alphahe(patchi)
379  + alpha2().boundaryField()[patchi]*thermo2_->alphahe(patchi);
380 }
381 
382 
384 (
385  const volScalarField& alphat
386 ) const
387 {
388  return
389  alpha1()*thermo1_->kappaEff(alphat)
390  + alpha2()*thermo2_->kappaEff(alphat);
391 }
392 
393 
395 (
396  const scalarField& alphat,
397  const label patchi
398 ) const
399 {
400  return
401  alpha1().boundaryField()[patchi]*thermo1_->kappaEff(alphat, patchi)
402  + alpha2().boundaryField()[patchi]*thermo2_->kappaEff(alphat, patchi);
403 }
404 
405 
407 (
408  const volScalarField& alphat
409 ) const
410 {
411  return
412  alpha1()*thermo1_->alphaEff(alphat)
413  + alpha2()*thermo2_->alphaEff(alphat);
414 }
415 
416 
418 (
419  const scalarField& alphat,
420  const label patchi
421 ) const
422 {
423  return
424  alpha1().boundaryField()[patchi]*thermo1_->alphaEff(alphat, patchi)
425  + alpha2().boundaryField()[patchi]*thermo2_->alphaEff(alphat, patchi);
426 }
427 
428 
430 {
431  if (psiThermo::read())
432  {
433  return interfaceProperties::read();
434  }
435  else
436  {
437  return false;
438  }
439 }
440 
441 
442 // ************************************************************************* //
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 [W/m/K].
virtual tmp< volScalarField > hc() const
Chemical enthalpy [J/kg].
IOobject(const word &name, const fileName &instance, const objectRegistry &registry, readOption r=NO_READ, writeOption w=NO_WRITE, bool registerObject=true)
Construct from name, instance, registry, io options.
Definition: IOobject.C:209
static autoPtr< rhoThermo > New(const fvMesh &, const word &phaseName=word::null)
Selector.
Definition: rhoThermo.C:141
bool read()
Read transportProperties dictionary.
virtual tmp< volScalarField > alphahe() const
Thermal diffusivity for energy of mixture [kg/m/s].
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
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
const volScalarField & alpha1
virtual tmp< volScalarField > W() const
Molecular weight [kg/kmol].
virtual ~twoPhaseMixtureThermo()
Destructor.
List< label > labelList
A List of labels.
Definition: labelList.H:56
const fileOperation & fileHandler()
Get current file handler.
alpha2
Definition: alphaEqn.H:115
virtual tmp< volScalarField > CpByCpv() const
Heat capacity ratio [].
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-7/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:68
defineTypeNameAndDebug(combustionModel, 0)
virtual tmp< volScalarField > Cpv() const
Heat capacity at constant pressure/volume [J/kg/K].
virtual volScalarField & he()
Enthalpy/Internal energy [J/kg].
virtual void flush() const
Forcibly wait until all output done. Flush any cached data.
virtual word thermoName() const
Return the name of the thermo physics.
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 [W/m/K].
virtual tmp< volScalarField > kappaEff(const volScalarField &alphat) const
Effective thermal diffusivity of mixture [W/m/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.