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-2020 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& T,
158  const labelList& cells
159 ) const
160 {
161  return
162  scalarField(alpha1(), cells)*thermo1_->he(T, cells)
163  + scalarField(alpha2(), cells)*thermo2_->he(T, cells);
164 }
165 
166 
168 (
169  const scalarField& T,
170  const label patchi
171 ) const
172 {
173  return
174  alpha1().boundaryField()[patchi]*thermo1_->he(T, patchi)
175  + alpha2().boundaryField()[patchi]*thermo2_->he(T, patchi);
176 }
177 
178 
180 {
181  return alpha1()*thermo1_->hs() + alpha2()*thermo2_->hs();
182 }
183 
184 
186 (
187  const volScalarField& p,
188  const volScalarField& T
189 ) const
190 {
191  return alpha1()*thermo1_->hs(p, T) + alpha2()*thermo2_->hs(p, T);
192 }
193 
194 
196 (
197  const scalarField& T,
198  const labelList& cells
199 ) const
200 {
201  return
202  scalarField(alpha1(), cells)*thermo1_->hs(T, cells)
203  + scalarField(alpha2(), cells)*thermo2_->hs(T, cells);
204 }
205 
206 
208 (
209  const scalarField& T,
210  const label patchi
211 ) const
212 {
213  return
214  alpha1().boundaryField()[patchi]*thermo1_->hs(T, patchi)
215  + alpha2().boundaryField()[patchi]*thermo2_->hs(T, patchi);
216 }
217 
218 
220 {
221  return alpha1()*thermo1_->ha() + alpha2()*thermo2_->ha();
222 }
223 
224 
226 (
227  const volScalarField& p,
228  const volScalarField& T
229 ) const
230 {
231  return alpha1()*thermo1_->ha(p, T) + alpha2()*thermo2_->ha(p, T);
232 }
233 
234 
236 (
237  const scalarField& T,
238  const labelList& cells
239 ) const
240 {
241  return
242  scalarField(alpha1(), cells)*thermo1_->ha(T, cells)
243  + scalarField(alpha2(), cells)*thermo2_->ha(T, cells);
244 }
245 
246 
248 (
249  const scalarField& T,
250  const label patchi
251 ) const
252 {
253  return
254  alpha1().boundaryField()[patchi]*thermo1_->ha(T, patchi)
255  + alpha2().boundaryField()[patchi]*thermo2_->ha(T, patchi);
256 }
257 
258 
260 {
261  return alpha1()*thermo1_->hc() + alpha2()*thermo2_->hc();
262 }
263 
264 
266 (
267  const scalarField& h,
268  const scalarField& T0,
269  const labelList& cells
270 ) const
271 {
273  return T0;
274 }
275 
276 
278 (
279  const scalarField& h,
280  const scalarField& T0,
281  const label patchi
282 ) const
283 {
285  return T0;
286 }
287 
288 
290 {
291  return alpha1()*thermo1_->Cp() + alpha2()*thermo2_->Cp();
292 }
293 
294 
296 (
297  const scalarField& T,
298  const label patchi
299 ) const
300 {
301  return
302  alpha1().boundaryField()[patchi]*thermo1_->Cp(T, patchi)
303  + alpha2().boundaryField()[patchi]*thermo2_->Cp(T, patchi);
304 }
305 
306 
308 {
309  return alpha1()*thermo1_->Cv() + alpha2()*thermo2_->Cv();
310 }
311 
312 
314 (
315  const scalarField& T,
316  const label patchi
317 ) const
318 {
319  return
320  alpha1().boundaryField()[patchi]*thermo1_->Cv(T, patchi)
321  + alpha2().boundaryField()[patchi]*thermo2_->Cv(T, patchi);
322 }
323 
324 
326 {
327  return alpha1()*thermo1_->gamma() + alpha2()*thermo2_->gamma();
328 }
329 
330 
332 (
333  const scalarField& T,
334  const label patchi
335 ) const
336 {
337  return
338  alpha1().boundaryField()[patchi]*thermo1_->gamma(T, patchi)
339  + alpha2().boundaryField()[patchi]*thermo2_->gamma(T, patchi);
340 }
341 
342 
344 {
345  return alpha1()*thermo1_->Cpv() + alpha2()*thermo2_->Cpv();
346 }
347 
348 
350 (
351  const scalarField& T,
352  const label patchi
353 ) const
354 {
355  return
356  alpha1().boundaryField()[patchi]*thermo1_->Cpv(T, patchi)
357  + alpha2().boundaryField()[patchi]*thermo2_->Cpv(T, patchi);
358 }
359 
360 
362 {
363  return
364  alpha1()*thermo1_->CpByCpv()
365  + alpha2()*thermo2_->CpByCpv();
366 }
367 
368 
370 (
371  const scalarField& T,
372  const label patchi
373 ) const
374 {
375  return
376  alpha1().boundaryField()[patchi]*thermo1_->CpByCpv(T, patchi)
377  + alpha2().boundaryField()[patchi]*thermo2_->CpByCpv(T, patchi);
378 }
379 
380 
382 {
383  return alpha1()*thermo1_->W() + alpha2()*thermo1_->W();
384 }
385 
386 
388 (
389  const label patchi
390 ) const
391 {
392  return
393  alpha1().boundaryField()[patchi]*thermo1_->W(patchi)
394  + alpha2().boundaryField()[patchi]*thermo1_->W(patchi);
395 }
396 
397 
399 {
400  return mu()/(alpha1()*thermo1_->rho() + alpha2()*thermo2_->rho());
401 }
402 
403 
405 (
406  const label patchi
407 ) const
408 {
409  return
410  mu(patchi)
411  /(
412  alpha1().boundaryField()[patchi]*thermo1_->rho(patchi)
413  + alpha2().boundaryField()[patchi]*thermo2_->rho(patchi)
414  );
415 }
416 
417 
419 {
420  return alpha1()*thermo1_->kappa() + alpha2()*thermo2_->kappa();
421 }
422 
423 
425 (
426  const label patchi
427 ) const
428 {
429  return
430  alpha1().boundaryField()[patchi]*thermo1_->kappa(patchi)
431  + alpha2().boundaryField()[patchi]*thermo2_->kappa(patchi);
432 }
433 
434 
436 {
437  return
438  alpha1()*thermo1_->alphahe()
439  + alpha2()*thermo2_->alphahe();
440 }
441 
442 
444 (
445  const label patchi
446 ) const
447 {
448  return
449  alpha1().boundaryField()[patchi]*thermo1_->alphahe(patchi)
450  + alpha2().boundaryField()[patchi]*thermo2_->alphahe(patchi);
451 }
452 
453 
455 (
456  const volScalarField& alphat
457 ) const
458 {
459  return
460  alpha1()*thermo1_->kappaEff(alphat)
461  + alpha2()*thermo2_->kappaEff(alphat);
462 }
463 
464 
466 (
467  const scalarField& alphat,
468  const label patchi
469 ) const
470 {
471  return
472  alpha1().boundaryField()[patchi]*thermo1_->kappaEff(alphat, patchi)
473  + alpha2().boundaryField()[patchi]*thermo2_->kappaEff(alphat, patchi);
474 }
475 
476 
478 (
479  const volScalarField& alphat
480 ) const
481 {
482  return
483  alpha1()*thermo1_->alphaEff(alphat)
484  + alpha2()*thermo2_->alphaEff(alphat);
485 }
486 
487 
489 (
490  const scalarField& alphat,
491  const label patchi
492 ) const
493 {
494  return
495  alpha1().boundaryField()[patchi]*thermo1_->alphaEff(alphat, patchi)
496  + alpha2().boundaryField()[patchi]*thermo2_->alphaEff(alphat, patchi);
497 }
498 
499 
501 {
502  if (psiThermo::read())
503  {
504  return interfaceProperties::read();
505  }
506  else
507  {
508  return false;
509  }
510 }
511 
512 
513 // ************************************************************************* //
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 > ha() const
Absolute enthalpy [J/kg].
virtual tmp< volScalarField > alphaEff(const volScalarField &alphat) const
Effective thermal diffusivity of mixture [W/m/K].
virtual tmp< volScalarField > hc() const
Enthalpy of formation [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:177
virtual tmp< volScalarField > hs() const
Sensible 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 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:58
alpha2
Definition: alphaEqn.H:115
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
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)
volScalarField & alpha1(mixture.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.
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-8/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.
label patchi
const dimensionedScalar & mu
Atomic mass unit.
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:506
const dimensionedScalar & h
Planck constant.
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
virtual tmp< scalarField > THE(const scalarField &h, const scalarField &T0, const labelList &cells) const
Temperature from enthalpy/internal energy for cell-set.
Namespace for OpenFOAM.
scalar T0
Definition: createFields.H:22