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-2020 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& alphaPhi10,
37  const twoPhaseMixtureThermo& mixture
38 )
39 :
40  twoPhaseTransport_(false),
41  mixture_(mixture),
42  phi_(phi),
43  alphaPhi10_(alphaPhi10)
44 {
45  {
46  IOdictionary momentumTransport
47  (
48  IOobject
49  (
50  momentumTransportModel::typeName,
51  U.time().constant(),
52  U.db(),
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)*alphaPhi10_
83  )
84  );
85 
86  alphaRhoPhi2_ =
87  (
89  (
90  IOobject::groupName("alphaRhoPhi", alpha2.group()),
91  fvc::interpolate(rho2)*(phi_ - alphaPhi10_)
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  turbulence1_->mut()
166  )
167  + mixture_.alpha2()*mixture_.thermo2().alphaEff
168  (
169  turbulence2_->mut()
170  );
171  }
172  else
173  {
174  return mixture_.alphaEff(turbulence_->mut());
175  }
176 }
177 
178 
181 (
182  volVectorField& U
183 ) const
184 {
185  if (twoPhaseTransport_)
186  {
187  return
188  turbulence1_->divDevTau(U)
189  + turbulence2_->divDevTau(U);
190  }
191  else
192  {
193  return turbulence_->divDevTau(U);
194  }
195 }
196 
197 
199 {
200  if (twoPhaseTransport_)
201  {
202  const volScalarField& rho1 = mixture_.thermo1().rho();
203  const volScalarField& rho2 = mixture_.thermo2().rho();
204 
205  alphaRhoPhi1_.ref() = fvc::interpolate(rho1)*alphaPhi10_;
206  alphaRhoPhi2_.ref() = fvc::interpolate(rho2)*(phi_ - alphaPhi10_);
207  }
208 }
209 
210 
212 {
213  if (twoPhaseTransport_)
214  {
215  turbulence1_->correct();
216  turbulence2_->correct();
217  }
218  else
219  {
220  turbulence_->correct();
221  }
222 }
223 
224 
225 // ************************************************************************* //
static word group(const word &name)
Return group (extension part of name)
Definition: IOobject.C:144
const volScalarField & alpha2() const
Return the phase-fraction of phase 2.
static autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const transportModel &transportModel)
Return a reference to the selected turbulence model.
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:177
static autoPtr< PhaseCompressibleMomentumTransportModel > New(const alphaField &alpha, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transportModel)
Return a reference to the selected turbulence model.
void correctPhasePhi()
Correct the phase mass-fluxes.
void correct()
Correct the phase or mixture transport models.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
tmp< fvVectorMatrix > divDevTau(volVectorField &U) const
Return the effective momentum stress divergence.
tmp< volScalarField > alphaEff() const
Return the effective temperature transport coefficient.
const volScalarField & alpha1() const
Return the phase-fraction of phase 1.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
static word groupName(Name name, const word &group)
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.
Internal & ref()
Return a reference to the dimensioned internal field.
compressibleInterPhaseTransportModel(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const surfaceScalarField &rhoPhi, const surfaceScalarField &alphaPhi10, const twoPhaseMixtureThermo &mixture)
Construct from components.
IOdictionary(const IOobject &)
Construct given an IOobject.
Definition: IOdictionary.C:30
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField