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-2025 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 
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  nu_
65  (
66  IOobject
67  (
68  "nu",
69  mesh.time().name(),
70  mesh
71  ),
72  mesh,
74  calculatedFvPatchScalarField::typeName
75  )
76 {
77  correct();
78 }
79 
80 
81 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
82 
84 {
85  bool incompressible = true;
86 
87  forAll(phases_, phasei)
88  {
89  incompressible =
90  incompressible && phases_[phasei].thermo().incompressible();
91  }
92 
93  return incompressible;
94 }
95 
96 
99 {
100  tmp<volScalarField> tpsiByRho
101  (
102  phases_[0]*phases_[0].thermo().psi()/phases_[0].thermo().rho()
103  );
104 
105  for (label phasei=1; phasei<phases_.size(); phasei++)
106  {
107  tpsiByRho.ref() +=
108  phases_[phasei]*phases_[phasei].thermo().psi()
109  /phases_[phasei].thermo().rho();
110  }
111 
112  return tpsiByRho;
113 }
114 
115 
117 (
118  const volScalarField& nut
119 ) const
120 {
121  tmp<volScalarField> talphaEff
122  (
123  phases_[0]
124  *(
125  phases_[0].thermo().kappa()
126  + phases_[0].thermo().rho()*phases_[0].thermo().Cp()*nut
127  )/phases_[0].thermo().Cv()
128  );
129 
130  for (label phasei=1; phasei<phases_.size(); phasei++)
131  {
132  talphaEff.ref() +=
133  phases_[phasei]
134  *(
135  phases_[phasei].thermo().kappa()
136  + phases_[phasei].thermo().rho()*phases_[phasei].thermo().Cp()*nut
137  )/phases_[phasei].thermo().Cv();
138  }
139 
140  return talphaEff;
141 }
142 
143 
146 {
147  tmp<volScalarField> trCv(phases_[0]/phases_[0].thermo().Cv());
148 
149  for (label phasei=1; phasei<phases_.size(); phasei++)
150  {
151  trCv.ref() += phases_[phasei]/phases_[phasei].thermo().Cv();
152  }
153 
154  return trCv;
155 }
156 
157 
159 {
160  forAll(phases_, phasei)
161  {
162  phases_[phasei].correct(p(), T());
163  }
164 }
165 
166 
168 {
169  rho_ = phases_[0]*phases_[0].thermo().rho();
170  volScalarField mu(phases_[0]*phases_[0].thermo().mu());
171 
172  for (label phasei=1; phasei<phases_.size(); phasei++)
173  {
174  rho_ += phases_[phasei]*phases_[phasei].thermo().rho();
175  mu += phases_[phasei]*phases_[phasei].thermo().mu();
176  }
177 
178  // Update the mixture kinematic viscosity
179  nu_ = mu/rho_;
180 
181  calcAlphas();
182 }
183 
184 
186 (
187  const volScalarField& dp
188 )
189 {
190  forAll(phases_, phasei)
191  {
192  phases_[phasei].thermo().rho() += phases_[phasei].thermo().psi()*dp;
193  }
194 }
195 
196 
197 // ************************************************************************* //
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:433
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 > 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:96
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:197
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
const scalar nut
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 dimKinematicViscosity
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimDensity
defineTypeNameAndDebug(combustionModel, 0)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
volScalarField & p
fluidMulticomponentThermo & thermo
Definition: createFields.H:31