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-2023 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& alphaRhoPhi1,
39  const surfaceScalarField& alphaRhoPhi2,
41 )
42 :
43  twoPhaseTransport_(false),
44  mixture_(mixture),
45  phi_(phi),
46  alphaPhi1_(alphaPhi1),
47  alphaRhoPhi1_(alphaRhoPhi1),
48  alphaRhoPhi2_(alphaRhoPhi2)
49 {
50  {
51  IOdictionary momentumTransport
52  (
53  IOobject
54  (
55  momentumTransportModel::typeName,
56  U.time().constant(),
57  U.db(),
60  )
61  );
62 
63  word simulationType
64  (
65  momentumTransport.lookup("simulationType")
66  );
67 
68  if (simulationType == "twoPhaseTransport")
69  {
70  twoPhaseTransport_ = true;
71  }
72  }
73 
74  if (twoPhaseTransport_)
75  {
76  momentumTransport1_ =
77  (
79  (
80  mixture_.alpha1(),
81  mixture_.thermo1().rho(),
82  U,
83  alphaRhoPhi1_,
84  phi,
85  mixture.thermo1()
86  )
87  );
88 
89  momentumTransport2_ =
90  (
92  (
93  mixture_.alpha2(),
94  mixture_.thermo2().rho(),
95  U,
96  alphaRhoPhi2_,
97  phi,
98  mixture.thermo2()
99  )
100  );
101  }
102  else
103  {
104  mixtureMomentumTransport_ = compressible::momentumTransportModel::New
105  (
106  rho,
107  U,
108  rhoPhi,
109  mixture
110  );
111 
112  mixtureMomentumTransport_->validate();
113  }
114 }
115 
116 
117 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
118 
121 (
123 ) const
124 {
125  if (twoPhaseTransport_)
126  {
127  return
128  momentumTransport1_->divDevTau(U)
129  + momentumTransport2_->divDevTau(U);
130  }
131  else
132  {
133  return mixtureMomentumTransport_->divDevTau(U);
134  }
135 }
136 
137 
139 {
140  if (twoPhaseTransport_)
141  {
142  momentumTransport1_->predict();
143  momentumTransport2_->predict();
144  }
145  else
146  {
147  mixtureMomentumTransport_->predict();
148  }
149 }
150 
151 
153 {
154  if (twoPhaseTransport_)
155  {
156  momentumTransport1_->correct();
157  momentumTransport2_->correct();
158  }
159  else
160  {
161  mixtureMomentumTransport_->correct();
162  }
163 }
164 
165 
166 // ************************************************************************* //
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 &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 rhoThermo-based phases.
const rhoThermo & thermo1() const
Return the thermo for phase 1.
const rhoThermo & 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:860
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