compressibleMultiphaseVoFMixture.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) 2023 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 
27 #include "surfaceInterpolate.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
34 }
35 
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
40 (
41  const fvMesh& mesh
42 )
43 :
45 
46  multiphaseVoFMixture(mesh, compressibleVoFphase::iNew(mesh, T())),
47 
48  phases_(multiphaseVoFMixture::phases().convert<compressibleVoFphase>()),
49 
50  rho_
51  (
52  IOobject
53  (
54  "rho",
55  mesh.time().name(),
56  mesh,
57  IOobject::NO_READ,
58  IOobject::AUTO_WRITE
59  ),
60  mesh,
61  dimensionedScalar("rho", dimDensity, 0)
62  )
63 {
64  correct();
65 }
66 
67 
68 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
69 
71 {
72  bool incompressible = true;
73 
74  forAll(phases_, phasei)
75  {
76  incompressible =
77  incompressible && phases_[phasei].thermo().incompressible();
78  }
79 
80  return incompressible;
81 }
82 
83 
86 {
87  volScalarField mu(phases_[0].Alpha()*phases_[0].thermo().mu());
88 
89  for (label phasei=1; phasei<phases_.size(); phasei++)
90  {
91  mu += phases_[phasei].Alpha()*phases_[phasei].thermo().mu();
92  }
93 
94  return mu/rho_;
95 }
96 
97 
99 (
100  const label patchi
101 ) const
102 {
104  (
105  phases_[0].Alpha().boundaryField()[patchi]
106  *phases_[0].thermo().mu(patchi)
107  );
108 
109  for (label phasei=1; phasei<phases_.size(); phasei++)
110  {
111  mu +=
112  phases_[phasei].Alpha().boundaryField()[patchi]
113  *phases_[phasei].thermo().mu(patchi);
114  }
115 
116  return mu/rho_.boundaryField()[patchi];
117 }
118 
119 
122 {
123  tmp<volScalarField> tpsiByRho
124  (
125  phases_[0]*phases_[0].thermo().psi()/phases_[0].thermo().rho()
126  );
127 
128  for (label phasei=1; phasei<phases_.size(); phasei++)
129  {
130  tpsiByRho.ref() +=
131  phases_[phasei]*phases_[phasei].thermo().psi()
132  /phases_[phasei].thermo().rho();
133  }
134 
135  return tpsiByRho;
136 }
137 
138 
140 (
141  const volScalarField& nut
142 ) const
143 {
144  tmp<volScalarField> talphaEff
145  (
146  phases_[0]
147  *(
148  phases_[0].thermo().kappa()
149  + phases_[0].thermo().rho()*phases_[0].thermo().Cp()*nut
150  )/phases_[0].thermo().Cv()
151  );
152 
153  for (label phasei=1; phasei<phases_.size(); phasei++)
154  {
155  talphaEff.ref() +=
156  phases_[phasei]
157  *(
158  phases_[phasei].thermo().kappa()
159  + phases_[phasei].thermo().rho()*phases_[phasei].thermo().Cp()*nut
160  )/phases_[phasei].thermo().Cv();
161  }
162 
163  return talphaEff;
164 }
165 
166 
169 {
170  tmp<volScalarField> trCv(phases_[0]/phases_[0].thermo().Cv());
171 
172  for (label phasei=1; phasei<phases_.size(); phasei++)
173  {
174  trCv.ref() += phases_[phasei]/phases_[phasei].thermo().Cv();
175  }
176 
177  return trCv;
178 }
179 
180 
182 {
183  forAll(phases_, phasei)
184  {
185  phases_[phasei].correct(p(), T());
186  }
187 }
188 
189 
191 {
192  rho_ = phases_[0]*phases_[0].thermo().rho();
193 
194  for (label phasei=1; phasei<phases_.size(); phasei++)
195  {
196  rho_ += phases_[phasei]*phases_[phasei].thermo().rho();
197  }
198 
199  forAll(phases_, phasei)
200  {
201  phases_[phasei].Alpha() =
202  phases_[phasei]*phases_[phasei].thermo().rho()/rho_;
203  }
204 
205  calcAlphas();
206 }
207 
208 
210 (
211  const volScalarField& dp
212 )
213 {
214  forAll(phases_, phasei)
215  {
216  phases_[phasei].thermo().rho() += phases_[phasei].thermo().psi()*dp;
217  }
218 }
219 
220 
221 // ************************************************************************* //
scalar Cp(const scalar p, const scalar T) const
Definition: EtoHthermo.H:2
scalar Cv(const scalar p, const scalar T) const
Definition: HtoEthermo.H:2
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
virtual const volScalarField & Cv() const =0
Heat capacity at constant volume [J/kg/K].
Compressible multiphase mixture for interface-capturing simulations.
tmp< volScalarField > rCv() const
Return the phase-averaged reciprocal Cv.
virtual void correctThermo()
Correct the thermodynamics of each phase.
virtual tmp< volScalarField > nu() const
Return the kinematic laminar viscosity.
virtual tmp< volScalarField > alphaEff(const volScalarField &nut) const
Return the effective temperature transport coefficient.
virtual tmp< volScalarField > psiByRho() const
Return the mixture compressibility/density.
compressibleMultiphaseVoFMixture(const fvMesh &mesh)
Construct from fvMesh.
void correctRho(const volScalarField &dp)
Update densities for given pressure change.
bool incompressible() const
Return true if all phases are incompressible.
Single compressible phase derived from the VoFphase.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Multiphase VoF mixture with support for interface properties.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
const scalar nut
label patchi
const volScalarField & psi
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
const dimensionedScalar mu
Atomic mass unit.
Namespace for OpenFOAM.
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
const dimensionSet dimDensity
defineTypeNameAndDebug(combustionModel, 0)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
volScalarField & p
fluidMulticomponentThermo & thermo
Definition: createFields.H:31