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-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 
33 (
34  const volVectorField& U,
35  const surfaceScalarField& phi,
36  const surfaceScalarField& alphaPhi1,
37  const surfaceScalarField& alphaPhi2,
39 )
40 :
41  twoPhaseTransport_(false),
42  mixture_(mixture),
43  phi_(phi),
44  alphaPhi1_(alphaPhi1),
45  alphaPhi2_(alphaPhi2)
46 {
47  {
48  IOdictionary momentumTransport
49  (
50  IOobject
51  (
52  momentumTransportModel::typeName,
53  U.time().constant(),
54  U.db(),
57  )
58  );
59 
60  word simulationType
61  (
62  momentumTransport.lookup("simulationType")
63  );
64 
65  if (simulationType == "twoPhaseTransport")
66  {
67  twoPhaseTransport_ = true;
68  }
69  }
70 
71  if (twoPhaseTransport_)
72  {
73  const volScalarField& alpha1(mixture_.alpha1());
74  const volScalarField& alpha2(mixture_.alpha2());
75 
76  turbulence1_ =
77  (
79  (
80  alpha1,
81  U,
82  alphaPhi1_,
83  phi_,
84  mixture_.nuModel1()
85  )
86  );
87 
88  turbulence2_ =
89  (
91  (
92  alpha2,
93  U,
94  alphaPhi2_,
95  phi_,
96  mixture_.nuModel2()
97  )
98  );
99  }
100  else
101  {
103  (
104  U,
105  phi_,
106  mixture_
107  );
108 
109  turbulence_->validate();
110  }
111 }
112 
113 
114 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
115 
118 (
119  const volScalarField& rho,
121 ) const
122 {
123  if (twoPhaseTransport_)
124  {
125  return
126  mixture_.rho1()*turbulence1_->divDevSigma(U)
127  + mixture_.rho2()*turbulence2_->divDevSigma(U);
128  }
129  else
130  {
131  return turbulence_->divDevTau(rho, U);
132  }
133 }
134 
135 
137 {
138  if (twoPhaseTransport_)
139  {
140  turbulence1_->predict();
141  turbulence2_->predict();
142  }
143  else
144  {
145  turbulence_->predict();
146  }
147 }
148 
149 
151 {
152  if (twoPhaseTransport_)
153  {
154  turbulence1_->correct();
155  turbulence2_->correct();
156  }
157  else
158  {
159  turbulence_->correct();
160  }
161 }
162 
163 
164 // ************************************************************************* //
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
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
void correct()
Correct the phase or mixture transport models.
incompressibleInterPhaseTransportModel(const volVectorField &U, const surfaceScalarField &phi, const surfaceScalarField &alphaPhi1, const surfaceScalarField &alphaPhi2, 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.
const viscosityModel & nuModel1() const
Return const-access to phase1 viscosityModel.
const viscosityModel & nuModel2() const
Return const-access to phase2 viscosityModel.
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