Loading...
Searching...
No Matches
compressibleInterPhaseTransportModel.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) 2017-2018 OpenFOAM Foundation
9 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
33Foam::compressibleInterPhaseTransportModel::compressibleInterPhaseTransportModel
34(
35 const volScalarField& rho,
36 const volVectorField& U,
37 const surfaceScalarField& phi,
38 const surfaceScalarField& rhoPhi,
39 const surfaceScalarField& alphaPhi10,
40 const twoPhaseMixtureThermo& mixture
41)
42:
43 twoPhaseTransport_(false),
44 mixture_(mixture),
45 phi_(phi),
46 alphaPhi10_(alphaPhi10)
47{
48 {
49 IOdictionary turbulenceProperties
50 (
51 IOobject
52 (
53 turbulenceModel::propertiesName,
54 U.time().constant(),
55 U.db(),
56 IOobject::MUST_READ,
57 IOobject::NO_WRITE
58 )
59 );
60
61 const word simulationType
62 (
63 turbulenceProperties.get<word>("simulationType")
64 );
65
66 if (simulationType == "twoPhaseTransport")
67 {
68 twoPhaseTransport_ = true;
69 }
70 }
71
72 if (twoPhaseTransport_)
73 {
74 const volScalarField& alpha1(mixture_.alpha1());
75 const volScalarField& alpha2(mixture_.alpha2());
76
77 const tmp<volScalarField> trho1 = mixture_.thermo1().rho();
78 const tmp<volScalarField> trho2 = mixture_.thermo2().rho();
79
80 const auto& rho1 = trho1();
81 const auto& rho2 = trho2();
82
83 alphaRhoPhi1_ =
84 (
86 (
87 IOobject::groupName("alphaRhoPhi", alpha1.group()),
88 fvc::interpolate(rho1)*alphaPhi10_
89 )
90 );
91
92 alphaRhoPhi2_ =
93 (
95 (
96 IOobject::groupName("alphaRhoPhi", alpha2.group()),
97 fvc::interpolate(rho2)*(phi_ - alphaPhi10_)
98 )
99 );
100
101 turbulence1_ =
102 (
103 ThermalDiffusivity<PhaseCompressibleTurbulenceModel<fluidThermo>>
104 ::New
105 (
106 alpha1,
107 rho1,
108 U,
109 alphaRhoPhi1_(),
110 phi,
111 mixture.thermo1()
112 )
113 );
114
115 turbulence2_ =
116 (
117 ThermalDiffusivity<PhaseCompressibleTurbulenceModel<fluidThermo>>
118 ::New
119 (
120 alpha2,
121 rho2,
122 U,
123 alphaRhoPhi2_(),
124 phi,
125 mixture.thermo2()
126 )
127 );
128 }
129 else
130 {
131 turbulence_ = compressible::turbulenceModel::New
132 (
133 rho,
134 U,
135 rhoPhi,
136 mixture
137 );
138
139 turbulence_->validate();
140 }
141}
142
143
144// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
145
148{
149 if (twoPhaseTransport_)
150 {
151 return
152 mixture_.alpha1()*mixture_.thermo1().alphaEff
153 (
154 turbulence1_->alphat()
155 )
156 + mixture_.alpha2()*mixture_.thermo2().alphaEff
157 (
158 turbulence2_->alphat()
159 );
160 }
161 else
162 {
163 return mixture_.alphaEff(turbulence_->alphat());
164 }
165}
166
167
168Foam::tmp<Foam::fvVectorMatrix>
170(
172) const
173{
174 if (twoPhaseTransport_)
175 {
176 return
177 turbulence1_->divDevRhoReff(U)
178 + turbulence2_->divDevRhoReff(U);
179 }
180 else
181 {
182 return turbulence_->divDevRhoReff(U);
183 }
184}
185
186
188{
189 if (twoPhaseTransport_)
190 {
191 const tmp<volScalarField> rho1 = mixture_.thermo1().rho();
192 const tmp<volScalarField> rho2 = mixture_.thermo2().rho();
193
194 alphaRhoPhi1_.ref() = fvc::interpolate(rho1)*alphaPhi10_;
195 alphaRhoPhi2_.ref() = fvc::interpolate(rho2)*(phi_ - alphaPhi10_);
196 }
197}
198
199
201{
202 if (twoPhaseTransport_)
203 {
204 turbulence1_->correct();
205 turbulence2_->correct();
206 }
207 else
208 {
209 turbulence_->correct();
210 }
211}
212
213
214// ************************************************************************* //
rhoPhi
Definition rhoEqn.H:10
const volScalarField & alpha1
volScalarField & rho2
const volScalarField & alpha2
volScalarField & rho1
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef().
void correct()
Correct the phase or mixture transport models.
tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the effective momentum stress divergence.
tmp< volScalarField > alphaEff() const
Return the effective temperature transport coefficient.
void correctPhasePhi()
Correct the phase mass-fluxes.
A class for managing temporary objects.
Definition tmp.H:75
U
Definition pEqn.H:72
alphaPhi10
Definition alphaEqn.H:7
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.
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture