compressibleInterPhaseTransportModel.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) 2017-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 
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
31 (
32  const volScalarField& rho,
33  const volVectorField& U,
34  const surfaceScalarField& phi,
35  const surfaceScalarField& rhoPhi,
36  const surfaceScalarField& alphaPhi1,
37  const compressibleTwoPhaseMixture& mixture
38 )
39 :
40  twoPhaseTransport_(false),
41  mixture_(mixture),
42  phi_(phi),
43  alphaPhi1_(alphaPhi1)
44 {
45  {
46  IOdictionary momentumTransport
47  (
48  IOobject
49  (
50  momentumTransportModel::typeName,
51  U.time().constant(),
52  U.db(),
53  IOobject::MUST_READ,
54  IOobject::NO_WRITE
55  )
56  );
57 
58  word simulationType
59  (
60  momentumTransport.lookup("simulationType")
61  );
62 
63  if (simulationType == "twoPhaseTransport")
64  {
65  twoPhaseTransport_ = true;
66  }
67  }
68 
69  if (twoPhaseTransport_)
70  {
71  const volScalarField& alpha1(mixture_.alpha1());
72  const volScalarField& alpha2(mixture_.alpha2());
73 
74  const volScalarField& rho1 = mixture_.thermo1().rho();
75  const volScalarField& rho2 = mixture_.thermo2().rho();
76 
77  alphaRhoPhi1_ =
78  (
80  (
81  IOobject::groupName("alphaRhoPhi", alpha1.group()),
82  fvc::interpolate(rho1)*alphaPhi1_
83  )
84  );
85 
86  alphaRhoPhi2_ =
87  (
89  (
90  IOobject::groupName("alphaRhoPhi", alpha2.group()),
91  fvc::interpolate(rho2)*(phi_ - alphaPhi1_)
92  )
93  );
94 
95  turbulence1_ =
96  (
98  (
99  alpha1,
100  rho1,
101  U,
102  alphaRhoPhi1_(),
103  phi,
104  mixture.thermo1()
105  )
106  );
107 
108  turbulence2_ =
109  (
111  (
112  alpha2,
113  rho2,
114  U,
115  alphaRhoPhi2_(),
116  phi,
117  mixture.thermo2()
118  )
119  );
120  }
121  else
122  {
124  (
125  rho,
126  U,
127  rhoPhi,
128  mixture
129  );
130 
131  turbulence_->validate();
132  }
133 }
134 
135 
136 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
137 
140 {
141  /* ***HGW
142  if (twoPhaseTransport_)
143  {
144  return
145  mixture_.alpha1()*mixture_.thermo1().alphaEff
146  (
147  turbulence1_->alphat()
148  )
149  + mixture_.alpha2()*mixture_.thermo2().alphaEff
150  (
151  turbulence2_->alphat()
152  );
153  }
154  else
155  {
156  return mixture_.alphaEff(turbulence_->alphat());
157  }
158  */
159 
160  if (twoPhaseTransport_)
161  {
162  return
163  mixture_.alpha1()*mixture_.thermo1().alphaEff
164  (
165  mixture_.thermo1().rho()*turbulence1_->nut()
166  )
167  + mixture_.alpha2()*mixture_.thermo2().alphaEff
168  (
169  mixture_.thermo2().rho()*turbulence2_->nut()
170  );
171  }
172  else
173  {
174  const volScalarField alphat(mixture_.rho()*turbulence_->nut());
175 
176  return
177  mixture_.alpha1()*mixture_.thermo1().alphaEff(alphat)
178  + mixture_.alpha2()*mixture_.thermo2().alphaEff(alphat);
179  }
180 }
181 
182 
185 (
186  volVectorField& U
187 ) const
188 {
189  if (twoPhaseTransport_)
190  {
191  return
192  turbulence1_->divDevTau(U)
193  + turbulence2_->divDevTau(U);
194  }
195  else
196  {
197  return turbulence_->divDevTau(U);
198  }
199 }
200 
201 
203 {
204  if (twoPhaseTransport_)
205  {
206  const volScalarField& rho1 = mixture_.thermo1().rho();
207  const volScalarField& rho2 = mixture_.thermo2().rho();
208 
209  alphaRhoPhi1_.ref() = fvc::interpolate(rho1)*alphaPhi1_;
210  alphaRhoPhi2_.ref() = fvc::interpolate(rho2)*(phi_ - alphaPhi1_);
211  }
212 }
213 
214 
216 {
217  if (twoPhaseTransport_)
218  {
219  turbulence1_->correct();
220  turbulence2_->correct();
221  }
222  else
223  {
224  turbulence_->correct();
225  }
226 }
227 
228 
229 // ************************************************************************* //
const volScalarField & alpha2() const
Return the phase-fraction of phase 2.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
virtual tmp< volScalarField > alphaEff(const volScalarField &alphat) const =0
Effective turbulent thermal diffusivity of energy.
virtual tmp< volScalarField > rho() const =0
Density [kg/m^3].
volScalarField & alpha1(mixture.alpha1())
compressibleInterPhaseTransportModel(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const surfaceScalarField &rhoPhi, const surfaceScalarField &alphaPhi1, const compressibleTwoPhaseMixture &mixture)
Construct from components.
void correctPhasePhi()
Correct the phase mass-fluxes.
void correct()
Correct the phase or mixture transport models.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:59
tmp< fvVectorMatrix > divDevTau(volVectorField &U) const
Return the effective momentum stress divergence.
tmp< volScalarField > alphaEff() const
Return the effective temperature transport coefficient.
alpha2
Definition: alphaEqn.H:115
const volScalarField & alpha1() const
Return the phase-fraction of phase 1.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
const volScalarField & rho() const
Return mixture density [kg/m^3].
const rhoThermo & thermo2() const
Return the thermo for phase 2.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const rhoThermo & thermo1() const
Return the thermo for phase 1.