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-2021 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"
27 
28 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32  defineTypeNameAndDebug(twoPhaseMixtureThermo, 0);
33 }
34 
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
39 (
40  const volVectorField& U,
41  const surfaceScalarField& phi
42 )
43 :
44  rhoThermo::composite(U.mesh(), word::null),
45  twoPhaseMixture(U.mesh(), *this),
46  interfaceProperties(alpha1(), alpha2(), U, *this),
47  thermo1_(nullptr),
48  thermo2_(nullptr),
49  Alpha1_
50  (
51  IOobject
52  (
53  IOobject::groupName("Alpha", phase1Name()),
54  U.mesh().time().timeName(),
55  U.mesh()
56  ),
57  alpha1(),
58  calculatedFvPatchScalarField::typeName
59  ),
60  Alpha2_
61  (
62  IOobject
63  (
64  IOobject::groupName("Alpha", phase2Name()),
65  U.mesh().time().timeName(),
66  U.mesh()
67  ),
68  alpha2(),
69  calculatedFvPatchScalarField::typeName
70  )
71 {
72  {
74  (
75  IOobject
76  (
77  IOobject::groupName("T", phase1Name()),
78  U.mesh().time().timeName(),
79  U.mesh()
80  ),
81  T_,
82  calculatedFvPatchScalarField::typeName
83  );
84  T1.write();
85  }
86 
87  {
89  (
90  IOobject
91  (
92  IOobject::groupName("T", phase2Name()),
93  U.mesh().time().timeName(),
94  U.mesh()
95  ),
96  T_,
97  calculatedFvPatchScalarField::typeName
98  );
99  T2.write();
100  }
101 
102  // Note: we're writing files to be read in immediately afterwards.
103  // Avoid any thread-writing problems.
104  fileHandler().flush();
105 
106  thermo1_ = rhoThermo::New(U.mesh(), phase1Name());
107  thermo2_ = rhoThermo::New(U.mesh(), phase2Name());
108 
109  // thermo1_->validate(phase1Name(), "e");
110  // thermo2_->validate(phase2Name(), "e");
111 
112  correct();
113 }
114 
115 
116 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
117 
119 {}
120 
121 
122 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
123 
125 {
126  thermo1_->T() = T_;
127  thermo1_->he() = thermo1_->he(p_, T_);
128  thermo1_->correct();
129 
130  thermo2_->T() = T_;
131  thermo2_->he() = thermo2_->he(p_, T_);
132  thermo2_->correct();
133 }
134 
135 
137 {
138  {
139  const volScalarField alphaRho1(alpha1()*thermo1_->rho());
140  const volScalarField alphaRho2(alpha2()*thermo2_->rho());
141 
142  rho_ = alphaRho1 + alphaRho2;
143  Alpha1_ = alphaRho1/rho_;
144  Alpha2_ = alphaRho2/rho_;
145  }
146 
147  psi_ = alpha1()*thermo1_->psi() + alpha2()*thermo2_->psi();
148  mu_ = alpha1()*thermo1_->mu() + alpha2()*thermo2_->mu();
149  alpha_ = alpha1()*thermo1_->alpha() + alpha2()*thermo2_->alpha();
150 
152 }
153 
154 
156 {
157  return thermo1_->thermoName() + ',' + thermo2_->thermoName();
158 }
159 
160 
162 {
163  return thermo1_->incompressible() && thermo2_->incompressible();
164 }
165 
166 
168 {
169  return thermo1_->isochoric() && thermo2_->isochoric();
170 }
171 
172 
174 (
175  const volScalarField& p,
176  const volScalarField& T
177 ) const
178 {
179  return Alpha1_*thermo1_->he(p, T) + Alpha2_*thermo2_->he(p, T);
180 }
181 
182 
184 (
185  const scalarField& T,
186  const labelList& cells
187 ) const
188 {
189  return
190  scalarField(Alpha1_, cells)*thermo1_->he(T, cells)
191  + scalarField(Alpha2_, cells)*thermo2_->he(T, cells);
192 }
193 
194 
196 (
197  const scalarField& T,
198  const label patchi
199 ) const
200 {
201  return
202  Alpha1_.boundaryField()[patchi]*thermo1_->he(T, patchi)
203  + Alpha2_.boundaryField()[patchi]*thermo2_->he(T, patchi);
204 }
205 
206 
208 {
209  return Alpha1_*thermo1_->hs() + Alpha2_*thermo2_->hs();
210 }
211 
212 
214 (
215  const volScalarField& p,
216  const volScalarField& T
217 ) const
218 {
219  return Alpha1_*thermo1_->hs(p, T) + Alpha2_*thermo2_->hs(p, T);
220 }
221 
222 
224 (
225  const scalarField& T,
226  const labelList& cells
227 ) const
228 {
229  return
230  scalarField(Alpha1_, cells)*thermo1_->hs(T, cells)
231  + scalarField(Alpha2_, cells)*thermo2_->hs(T, cells);
232 }
233 
234 
236 (
237  const scalarField& T,
238  const label patchi
239 ) const
240 {
241  return
242  Alpha1_.boundaryField()[patchi]*thermo1_->hs(T, patchi)
243  + Alpha2_.boundaryField()[patchi]*thermo2_->hs(T, patchi);
244 }
245 
246 
248 {
249  return Alpha1_*thermo1_->ha() + Alpha2_*thermo2_->ha();
250 }
251 
252 
254 (
255  const volScalarField& p,
256  const volScalarField& T
257 ) const
258 {
259  return Alpha1_*thermo1_->ha(p, T) + Alpha2_*thermo2_->ha(p, T);
260 }
261 
262 
264 (
265  const scalarField& T,
266  const labelList& cells
267 ) const
268 {
269  return
270  scalarField(Alpha1_, cells)*thermo1_->ha(T, cells)
271  + scalarField(Alpha2_, cells)*thermo2_->ha(T, cells);
272 }
273 
274 
276 (
277  const scalarField& T,
278  const label patchi
279 ) const
280 {
281  return
282  Alpha1_.boundaryField()[patchi]*thermo1_->ha(T, patchi)
283  + Alpha2_.boundaryField()[patchi]*thermo2_->ha(T, patchi);
284 }
285 
286 
288 {
289  return Alpha1_*thermo1_->hc() + Alpha2_*thermo2_->hc();
290 }
291 
292 
294 (
295  const volScalarField& h,
296  const volScalarField& p,
297  const volScalarField& T0
298 ) const
299 {
301  return T0;
302 }
303 
304 
306 (
307  const scalarField& h,
308  const scalarField& T0,
309  const labelList& cells
310 ) const
311 {
313  return T0;
314 }
315 
316 
318 (
319  const scalarField& h,
320  const scalarField& T0,
321  const label patchi
322 ) const
323 {
325  return T0;
326 }
327 
328 
330 {
331  return Alpha1_*thermo1_->Cp() + Alpha2_*thermo2_->Cp();
332 }
333 
334 
336 (
337  const scalarField& T,
338  const label patchi
339 ) const
340 {
341  return
342  Alpha1_.boundaryField()[patchi]*thermo1_->Cp(T, patchi)
343  + Alpha2_.boundaryField()[patchi]*thermo2_->Cp(T, patchi);
344 }
345 
346 
348 {
349  return Alpha1_*thermo1_->Cv() + Alpha2_*thermo2_->Cv();
350 }
351 
352 
354 (
355  const scalarField& T,
356  const label patchi
357 ) const
358 {
359  return
360  Alpha1_.boundaryField()[patchi]*thermo1_->Cv(T, patchi)
361  + Alpha2_.boundaryField()[patchi]*thermo2_->Cv(T, patchi);
362 }
363 
364 
366 {
367  return Cp()/Cv();
368 }
369 
370 
372 (
373  const scalarField& T,
374  const label patchi
375 ) const
376 {
377  return Cp(T, patchi)/Cv(T, patchi);
378 }
379 
380 
382 {
383  return Alpha1_*thermo1_->Cpv() + Alpha2_*thermo2_->Cpv();
384 }
385 
386 
388 (
389  const scalarField& T,
390  const label patchi
391 ) const
392 {
393  return
394  Alpha1_.boundaryField()[patchi]*thermo1_->Cpv(T, patchi)
395  + Alpha2_.boundaryField()[patchi]*thermo2_->Cpv(T, patchi);
396 }
397 
398 
400 {
401  return 1/(Alpha1_/thermo1_->W() + Alpha2_/thermo2_->W());
402 }
403 
404 
406 (
407  const label patchi
408 ) const
409 {
410  return
411  1/(
412  Alpha1_.boundaryField()[patchi]/thermo1_->W(patchi)
413  + Alpha2_.boundaryField()[patchi]/thermo2_->W(patchi)
414  );
415 }
416 
417 
419 {
420  return mu()/rho_;
421 }
422 
423 
425 (
426  const label patchi
427 ) const
428 {
429  return mu(patchi)/rho_.boundaryField()[patchi];
430 }
431 
432 
434 {
435  return alpha1()*thermo1_->kappa() + alpha2()*thermo2_->kappa();
436 }
437 
438 
440 (
441  const label patchi
442 ) const
443 {
444  return
445  alpha1().boundaryField()[patchi]*thermo1_->kappa(patchi)
446  + alpha2().boundaryField()[patchi]*thermo2_->kappa(patchi);
447 }
448 
449 
451 {
452  return alpha1()*thermo1_->alphahe() + alpha2()*thermo2_->alphahe();
453 }
454 
455 
457 (
458  const label patchi
459 ) const
460 {
461  return
462  alpha1().boundaryField()[patchi]*thermo1_->alphahe(patchi)
463  + alpha2().boundaryField()[patchi]*thermo2_->alphahe(patchi);
464 }
465 
466 
468 (
469  const volScalarField& alphat
470 ) const
471 {
472  return
473  alpha1()*thermo1_->kappaEff(alphat)
474  + alpha2()*thermo2_->kappaEff(alphat);
475 }
476 
477 
479 (
480  const scalarField& alphat,
481  const label patchi
482 ) const
483 {
484  return
485  alpha1().boundaryField()[patchi]*thermo1_->kappaEff(alphat, patchi)
486  + alpha2().boundaryField()[patchi]*thermo2_->kappaEff(alphat, patchi);
487 }
488 
489 
491 (
492  const volScalarField& alphat
493 ) const
494 {
495  return
496  alpha1()*thermo1_->alphaEff(alphat)
497  + alpha2()*thermo2_->alphaEff(alphat);
498 }
499 
500 
502 (
503  const scalarField& alphat,
504  const label patchi
505 ) const
506 {
507  return
508  alpha1().boundaryField()[patchi]*thermo1_->alphaEff(alphat, patchi)
509  + alpha2().boundaryField()[patchi]*thermo2_->alphaEff(alphat, patchi);
510 }
511 
512 
514 {
516  {
517  return interfaceProperties::read();
518  }
519  else
520  {
521  return false;
522  }
523 }
524 
525 
526 // ************************************************************************* //
scalar Cv(const scalar p, const scalar T) const
Definition: HtoEthermo.H:2
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].
const Boundary & boundaryField() const
Return const-reference to the boundary field.
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:167
volScalarField & alpha1(mixture.alpha1())
virtual tmp< volScalarField > hs() const
Sensible enthalpy [J/kg].
static autoPtr< rhoThermo > New(const fvMesh &, const word &phaseName=word::null)
Standard selection based on fvMesh.
Definition: rhoThermo.C:143
bool read()
Read transportProperties dictionary.
virtual tmp< volScalarField > alphahe() const
Thermal diffusivity for energy of mixture [kg/m/s].
virtual tmp< volScalarField > THE(const volScalarField &h, const volScalarField &p, const volScalarField &T0) const
Temperature from enthalpy/internal energy.
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
const dimensionedScalar h
Planck constant.
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.
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-9/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
static word groupName(Name name, const word &group)
virtual tmp< volScalarField > W() const
Molecular weight [kg/kmol].
virtual ~twoPhaseMixtureThermo()
Destructor.
word timeName
Definition: getTimeIndex.H:3
List< label > labelList
A List of labels.
Definition: labelList.H:56
const fileOperation & fileHandler()
Get current file handler.
const dimensionedScalar mu
Atomic mass unit.
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
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].
scalar Cp(const scalar p, const scalar T) const
Definition: EtoHthermo.H:2
virtual tmp< volScalarField > nu() const
Kinematic viscosity of mixture [m^2/s].
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
virtual bool read()
Read thermophysical properties dictionary.
Definition: basicThermo.C:476
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:370
Namespace for OpenFOAM.
scalar T0
Definition: createFields.H:22