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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2013-2017 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "twoPhaseMixtureThermo.H"
32 #include "collatedFileOperation.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(twoPhaseMixtureThermo, 0);
39 }
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
45 (
46  const volVectorField& U,
47  const surfaceScalarField& phi
48 )
49 :
50  psiThermo(U.mesh(), word::null),
51  twoPhaseMixture(U.mesh(), *this),
52  interfaceProperties(alpha1(), U, *this),
53  thermo1_(nullptr),
54  thermo2_(nullptr)
55 {
56  {
57  volScalarField T1(IOobject::groupName("T", phase1Name()), T_);
58  T1.write();
59  }
60 
61  {
62  volScalarField T2(IOobject::groupName("T", phase2Name()), T_);
63  T2.write();
64  }
65 
66  // Note: we're writing files to be read in immediately afterwards.
67  // Avoid any thread-writing problems.
68  fileHandler().flush();
69 
70  thermo1_ = rhoThermo::New(U.mesh(), phase1Name());
71  thermo2_ = rhoThermo::New(U.mesh(), phase2Name());
72 
73  // thermo1_->validate(phase1Name(), "e");
74  // thermo2_->validate(phase2Name(), "e");
75 
76  correct();
77 }
78 
79 
80 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
81 
83 {}
84 
85 
86 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
87 
89 {
90  thermo1_->he() = thermo1_->he(p_, T_);
91  thermo1_->correct();
92 
93  thermo2_->he() = thermo2_->he(p_, T_);
94  thermo2_->correct();
95 }
96 
97 
99 {
100  psi_ = alpha1()*thermo1_->psi() + alpha2()*thermo2_->psi();
101  mu_ = alpha1()*thermo1_->mu() + alpha2()*thermo2_->mu();
102  alpha_ = alpha1()*thermo1_->alpha() + alpha2()*thermo2_->alpha();
103 
105 }
106 
107 
109 {
110  return thermo1_->thermoName() + ',' + thermo2_->thermoName();
111 }
112 
113 
115 {
116  return thermo1_->incompressible() && thermo2_->incompressible();
117 }
118 
119 
121 {
122  return thermo1_->isochoric() && thermo2_->isochoric();
123 }
124 
125 
127 (
128  const volScalarField& p,
129  const volScalarField& T
130 ) const
131 {
132  return alpha1()*thermo1_->he(p, T) + alpha2()*thermo2_->he(p, T);
133 }
134 
135 
137 (
138  const scalarField& p,
139  const scalarField& T,
140  const labelList& cells
141 ) const
142 {
143  return
144  scalarField(alpha1(), cells)*thermo1_->he(p, T, cells)
145  + scalarField(alpha2(), cells)*thermo2_->he(p, T, cells);
146 }
147 
148 
150 (
151  const scalarField& p,
152  const scalarField& T,
153  const label patchi
154 ) const
155 {
156  return
157  alpha1().boundaryField()[patchi]*thermo1_->he(p, T, patchi)
158  + alpha2().boundaryField()[patchi]*thermo2_->he(p, T, patchi);
159 }
160 
161 
163 {
164  return alpha1()*thermo1_->hc() + alpha2()*thermo2_->hc();
165 }
166 
167 
169 (
170  const scalarField& h,
171  const scalarField& p,
172  const scalarField& T0,
173  const labelList& cells
174 ) const
175 {
177  return T0;
178 }
179 
180 
182 (
183  const scalarField& h,
184  const scalarField& p,
185  const scalarField& T0,
186  const label patchi
187 ) const
188 {
190  return T0;
191 }
192 
193 
195 {
196  return alpha1()*thermo1_->Cp() + alpha2()*thermo2_->Cp();
197 }
198 
199 
201 (
202  const scalarField& p,
203  const scalarField& T,
204  const label patchi
205 ) const
206 {
207  return
208  alpha1().boundaryField()[patchi]*thermo1_->Cp(p, T, patchi)
209  + alpha2().boundaryField()[patchi]*thermo2_->Cp(p, T, patchi);
210 }
211 
212 
214 {
215  return alpha1()*thermo1_->Cv() + alpha2()*thermo2_->Cv();
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_->Cv(p, T, patchi)
228  + alpha2().boundaryField()[patchi]*thermo2_->Cv(p, T, patchi);
229 }
230 
231 
233 {
234  return alpha1()*thermo1_->gamma() + alpha2()*thermo2_->gamma();
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_->gamma(p, T, patchi)
247  + alpha2().boundaryField()[patchi]*thermo2_->gamma(p, T, patchi);
248 }
249 
250 
252 {
253  return alpha1()*thermo1_->Cpv() + alpha2()*thermo2_->Cpv();
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_->Cpv(p, T, patchi)
266  + alpha2().boundaryField()[patchi]*thermo2_->Cpv(p, T, patchi);
267 }
268 
269 
271 {
272  return
273  alpha1()*thermo1_->CpByCpv()
274  + alpha2()*thermo2_->CpByCpv();
275 }
276 
277 
279 (
280  const scalarField& p,
281  const scalarField& T,
282  const label patchi
283 ) const
284 {
285  return
286  alpha1().boundaryField()[patchi]*thermo1_->CpByCpv(p, T, patchi)
287  + alpha2().boundaryField()[patchi]*thermo2_->CpByCpv(p, T, patchi);
288 }
289 
290 
292 {
293  return alpha1()*thermo1_->W() + alpha2()*thermo1_->W();
294 }
295 
296 
298 {
299  return mu()/(alpha1()*thermo1_->rho() + alpha2()*thermo2_->rho());
300 }
301 
302 
304 (
305  const label patchi
306 ) const
307 {
308  return
309  mu(patchi)
310  /(
311  alpha1().boundaryField()[patchi]*thermo1_->rho(patchi)
312  + alpha2().boundaryField()[patchi]*thermo2_->rho(patchi)
313  );
314 }
315 
316 
318 {
319  return alpha1()*thermo1_->kappa() + alpha2()*thermo2_->kappa();
320 }
321 
322 
324 (
325  const label patchi
326 ) const
327 {
328  return
329  alpha1().boundaryField()[patchi]*thermo1_->kappa(patchi)
330  + alpha2().boundaryField()[patchi]*thermo2_->kappa(patchi);
331 }
332 
333 
335 {
336  return
337  alpha1()*thermo1_->alphahe()
338  + alpha2()*thermo2_->alphahe();
339 }
340 
341 
343 (
344  const label patchi
345 ) const
346 {
347  return
348  alpha1().boundaryField()[patchi]*thermo1_->alphahe(patchi)
349  + alpha2().boundaryField()[patchi]*thermo2_->alphahe(patchi);
350 }
351 
352 
354 (
355  const volScalarField& alphat
356 ) const
357 {
358  return
359  alpha1()*thermo1_->kappaEff(alphat)
360  + alpha2()*thermo2_->kappaEff(alphat);
361 }
362 
363 
365 (
366  const scalarField& alphat,
367  const label patchi
368 ) const
369 {
370  return
371  alpha1().boundaryField()[patchi]*thermo1_->kappaEff(alphat, patchi)
372  + alpha2().boundaryField()[patchi]*thermo2_->kappaEff(alphat, patchi);
373 }
374 
375 
377 (
378  const volScalarField& alphat
379 ) const
380 {
381  return
382  alpha1()*thermo1_->alphaEff(alphat)
383  + alpha2()*thermo2_->alphaEff(alphat);
384 }
385 
386 
388 (
389  const scalarField& alphat,
390  const label patchi
391 ) const
392 {
393  return
394  alpha1().boundaryField()[patchi]*thermo1_->alphaEff(alphat, patchi)
395  + alpha2().boundaryField()[patchi]*thermo2_->alphaEff(alphat, patchi);
396 }
397 
398 
400 {
401  if (psiThermo::read())
402  {
403  return interfaceProperties::read();
404  }
405 
406  return false;
407 }
408 
409 
410 // ************************************************************************* //
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.
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:199
refPtr< fileOperation > fileHandler(std::nullptr_t)
Delete current file handler - forwards to fileOperation::handler()
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.
const volScalarField & alpha2
virtual bool read()
Read base transportProperties dictionary.
virtual tmp< volScalarField > gamma() const
Gamma = Cp/Cv [].
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:76
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
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].
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
dynamicFvMesh & mesh
const cellShapeList & cells
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual tmp< volScalarField > W() const
Molecular weight [kg/kmol].
virtual ~twoPhaseMixtureThermo()
Destructor.
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/chef2/andy/OpenFOAM/release/v2406/OpenFOAM-v2406/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=clamp((rho - rholSat)/(rhovSat - rholSat), zero_one{});alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:63
defineTypeNameAndDebug(combustionModel, 0)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
virtual tmp< volScalarField > Cpv() const
Heat capacity at constant pressure/volume [J/kg/K].
virtual volScalarField & he()
Enthalpy/Internal energy [J/kg].
const dimensionedScalar h
Planck constant.
virtual word thermoName() const
Return the name of the thermo physics.
const dimensionedScalar mu
Atomic mass unit.
U
Definition: pEqn.H:72
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:643
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:696
Namespace for OpenFOAM.
scalar T0
Definition: createFields.H:22
const volScalarField & alpha1