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-2024 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 "fvMatrix.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
32 (
33  const volScalarField& rho,
34  const volVectorField& U,
35  const surfaceScalarField& phi,
36  const surfaceScalarField& rhoPhi,
37  const surfaceScalarField& alphaPhi1,
38  const surfaceScalarField& alphaPhi2,
39  const surfaceScalarField& alphaRhoPhi1,
40  const surfaceScalarField& alphaRhoPhi2,
42 )
43 :
44  twoPhaseTransport_(false),
45  mixture_(mixture),
46  phi_(phi),
47  alphaPhi1_(alphaPhi1),
48  alphaPhi2_(alphaPhi2),
49  alphaRhoPhi1_(alphaRhoPhi1),
50  alphaRhoPhi2_(alphaRhoPhi2)
51 {
52  {
53  IOdictionary momentumTransport
54  (
55  IOobject
56  (
57  momentumTransportModel::typeName,
58  U.time().constant(),
59  U.db(),
62  )
63  );
64 
65  word simulationType
66  (
67  momentumTransport.lookup("simulationType")
68  );
69 
70  if (simulationType == "twoPhaseTransport")
71  {
72  twoPhaseTransport_ = true;
73  }
74  }
75 
76  if (twoPhaseTransport_)
77  {
78  momentumTransport1_ =
79  (
81  (
82  mixture_.alpha1(),
83  mixture_.thermo1().rho(),
84  U,
85  alphaRhoPhi1_,
86  phi_,
87  mixture_.thermo1()
88  )
89  );
90 
91  momentumTransport2_ =
92  (
94  (
95  mixture_.alpha2(),
96  mixture_.thermo2().rho(),
97  U,
98  alphaRhoPhi2_,
99  phi_,
100  mixture_.thermo2()
101  )
102  );
103  }
104  else
105  {
106  mixtureMomentumTransport_ = compressible::momentumTransportModel::New
107  (
108  rho,
109  U,
110  rhoPhi,
111  mixture_
112  );
113 
114  mixtureMomentumTransport_->validate();
115  }
116 }
117 
118 
119 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
120 
123 (
125 ) const
126 {
127  if (twoPhaseTransport_)
128  {
129  return
130  momentumTransport1_->divDevTau(U)
131  + momentumTransport2_->divDevTau(U);
132  }
133  else
134  {
135  return mixtureMomentumTransport_->divDevTau(U);
136  }
137 }
138 
139 
141 {
142  if (twoPhaseTransport_)
143  {
144  momentumTransport1_->predict();
145  momentumTransport2_->predict();
146  }
147  else
148  {
149  mixtureMomentumTransport_->predict();
150  }
151 }
152 
153 
155 {
156  if (twoPhaseTransport_)
157  {
158  momentumTransport1_->correct();
159  momentumTransport2_->correct();
160  }
161  else
162  {
163  mixtureMomentumTransport_->correct();
164  }
165 }
166 
167 
168 // ************************************************************************* //
Generic GeometricField class.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
compressibleInterPhaseTransportModel(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const surfaceScalarField &rhoPhi, const surfaceScalarField &alphaPhi1, const surfaceScalarField &alphaPhi2, const surfaceScalarField &alphaRhoPhi1, const surfaceScalarField &alphaRhoPhi2, const compressibleTwoPhaseVoFMixture &mixture)
Construct from components.
void correct()
Correct the phase or mixture transport models.
void predict()
Predict the phase or mixture transport models.
tmp< fvVectorMatrix > divDevTau(volVectorField &U) const
Return the effective momentum stress divergence.
static autoPtr< compressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
Return a reference to the selected turbulence model.
Class to represent a mixture of two rhoFluidThermo-based phases.
const rhoFluidThermo & thermo1() const
Return the thermo for phase 1.
const rhoFluidThermo & thermo2() const
Return the thermo for phase 2.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
static autoPtr< phaseCompressibleMomentumTransportModel > New(const alphaField &alpha, const volScalarField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
Return a reference to the selected turbulence model.
virtual tmp< volScalarField > rho() const =0
Density [kg/m^3].
A class for managing temporary objects.
Definition: tmp.H:55
const volScalarField & alpha1() const
Return the phase-fraction of phase 1.
const volScalarField & alpha2() const
Return the phase-fraction of phase 2.
A class for handling words, derived from string.
Definition: word.H:62
U
Definition: pEqn.H:72