incompressibleInterPhaseTransportModel.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) 2021-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 
33 (
34  const volVectorField& U,
35  const surfaceScalarField& phi,
36  const surfaceScalarField& alphaPhi1,
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(),
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  alphaPhi2_ =
75  (
77  (
78  IOobject::groupName("alphaPhi", alpha2.group()),
79  (phi_ - alphaPhi1_)
80  )
81  );
82 
83  turbulence1_ =
84  (
86  (
87  alpha1,
88  U,
89  alphaPhi1_,
90  phi,
91  mixture.nuModel1()
92  )
93  );
94 
95  turbulence2_ =
96  (
98  (
99  alpha2,
100  U,
101  alphaPhi2_(),
102  phi,
103  mixture.nuModel2()
104  )
105  );
106  }
107  else
108  {
110  (
111  U,
112  phi,
113  mixture
114  );
115 
116  turbulence_->validate();
117  }
118 }
119 
120 
121 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
122 
125 (
126  const volScalarField& rho,
128 ) const
129 {
130  if (twoPhaseTransport_)
131  {
132  return
133  mixture_.rho1()*turbulence1_->divDevSigma(U)
134  + mixture_.rho2()*turbulence2_->divDevSigma(U);
135  }
136  else
137  {
138  return turbulence_->divDevTau(rho, U);
139  }
140 }
141 
142 
144 {
145  if (twoPhaseTransport_)
146  {
147  alphaPhi2_.ref() = (phi_ - alphaPhi1_);
148  }
149 }
150 
151 
153 {
154  if (twoPhaseTransport_)
155  {
156  turbulence1_->predict();
157  turbulence2_->predict();
158  }
159  else
160  {
161  turbulence_->predict();
162  }
163 }
164 
165 
167 {
168  if (twoPhaseTransport_)
169  {
170  turbulence1_->correct();
171  turbulence2_->correct();
172  }
173  else
174  {
175  turbulence_->correct();
176  }
177 }
178 
179 
180 // ************************************************************************* //
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
static word group(const word &name)
Return group (extension part of name)
Definition: IOobject.C:134
static word groupName(Name name, const word &group)
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:860
void correct()
Correct the phase or mixture transport models.
incompressibleInterPhaseTransportModel(const volVectorField &U, const surfaceScalarField &phi, const surfaceScalarField &alphaPhi1, const incompressibleTwoPhaseVoFMixture &mixture)
Construct from components.
void predict()
Predict the phase or mixture transport models.
tmp< fvVectorMatrix > divDevTau(const volScalarField &rho, volVectorField &U) const
Return the effective momentum stress divergence.
static autoPtr< incompressibleMomentumTransportModel > New(const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
Return a reference to the selected turbulence model.
Class to represent a mixture of two constant density phases.
static autoPtr< phaseIncompressibleMomentumTransportModel > New(const alphaField &alpha, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
Return a reference to the selected turbulence model.
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
SurfaceField< scalar > surfaceScalarField