incompressibleTwoPhaseInteractingMixture.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) 2014-2022 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 "mixtureViscosityModel.H"
28 #include "relativeVelocityModel.H"
29 #include "fvcDiv.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(incompressibleTwoPhaseInteractingMixture, 0);
37 }
38 
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
44 (
46  const surfaceScalarField& phi,
48 )
49 :
50  twoPhaseMixture(U.mesh()),
51 
52  U_(U),
53 
54  muModel_(mixtureViscosityModel::New(*this)),
55  nucModel_(viscosityModel::New(U.mesh(), phase2Name())),
56 
57  rhod_("rho", dimDensity, muModel_()),
58  rhoc_("rho", dimDensity, nucModel_()),
59  dd_
60  (
61  "d",
62  dimLength,
63  muModel_->lookupOrDefault("d", 0.0)
64  ),
65  alphaMax_(lookupOrDefault("alphaMax", 1.0)),
66 
67  MRF_(U.mesh()),
68 
69  mu_
70  (
71  IOobject
72  (
73  "mu",
74  U_.time().timeName(),
75  U_.db()
76  ),
77  U_.mesh(),
78  dimensionedScalar(dimensionSet(1, -1, -1, 0, 0), 0),
79  calculatedFvPatchScalarField::typeName
80  ),
81 
82  UdmModel_(relativeVelocityModel::New(*this, *this, g))
83 {
84  correct();
85 }
86 
87 
88 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
89 
92 {}
93 
94 
95 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
96 
99 {
100  return alpha1();
101 }
102 
103 
106 {
107  return alpha2();
108 }
109 
110 
113 {
114  return muModel_();
115 }
116 
117 
120 {
121  return nucModel_();
122 }
123 
124 
127 {
128  return rhod_;
129 }
130 
131 
134 {
135  return rhoc_;
136 };
137 
138 
141 {
142  return dd_;
143 }
144 
145 
146 Foam::scalar
148 {
149  return alphaMax_;
150 }
151 
152 
155 {
156  return U_;
157 }
158 
159 
160 const Foam::IOMRFZoneList&
162 {
163  return MRF_;
164 }
165 
166 
169 {
170  return mu_;
171 }
172 
173 
176 {
177  return mu_.boundaryField()[patchi];
178 }
179 
180 
183 {
184  return alpha1()*rhod_ + alpha2()*rhoc_;
185 }
186 
187 
190 {
191  return
192  alpha1().boundaryField()[patchi]*rhod_.value()
193  + alpha2().boundaryField()[patchi]*rhoc_.value();
194 }
195 
196 
199 {
200  return mu_/rho();
201 }
202 
203 
206 {
207  return mu_.boundaryField()[patchi]/rho(patchi);
208 }
209 
210 
213 {
214  return UdmModel_->Udm();
215 }
216 
217 
220 {
221  return fvc::div(UdmModel_->tauDm());
222 }
223 
224 
226 {
227  MRF_.correctBoundaryVelocity(U_);
228  mu_ = muModel_->mu(rhoc_*nucModel_->nu(), U_);
229  UdmModel_->correct();
230 }
231 
232 
234 {
235  if (twoPhaseMixture::read())
236  {
237  if (muModel_->read() || nucModel_->read())
238  {
239  muModel_->lookup("rho") >> rhod_;
240  nucModel_->lookup("rho") >> rhoc_;
241 
242  dd_ = dimensionedScalar
243  (
244  "d",
245  dimLength,
246  muModel_->lookupOrDefault("d", 0)
247  );
248 
249  alphaMax_ = muModel_->lookupOrDefault( "alphaMax", 1.0);
250 
251  return true;
252  }
253  else
254  {
255  return false;
256  }
257  }
258  else
259  {
260  return false;
261  }
262 }
263 
264 
265 // ************************************************************************* //
virtual void correct()
Correct the laminar viscosity.
const volScalarField & alpha2() const
Return the phase-fraction of phase 2.
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
virtual bool read()
Read base phaseProperties dictionary.
U
Definition: pEqn.H:72
UniformDimensionedField< vector > uniformDimensionedVectorField
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-10/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:68
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const volScalarField & alphad() const
Return const-access to the dispersed phase-fraction.
const dimensionedScalar & rhoc() const
Return const-access to continuous-phase density.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
const IOMRFZoneList & MRF() const
Return MRF zones.
Generic dimensioned Type class.
virtual ~incompressibleTwoPhaseInteractingMixture()
Destructor.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:59
fvMesh & mesh
const dimensionedScalar & rhod() const
Return const-access to the dispersed-phase density.
An abstract base class for Newtonian viscosity models.
Macros for easy insertion into run-time selection tables.
const volVectorField & U() const
Return const-access to the mixture velocity.
scalar alphaMax() const
Optional maximum phase-fraction (e.g. packing limit)
const dimensionSet dimLength
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
const volScalarField & alpha1() const
Return the phase-fraction of phase 1.
virtual tmp< volScalarField > nu() const
Return the mixture viscosity.
void correctBoundaryVelocity(volVectorField &U) const
Correct the boundary velocity for the rotation of the MRF region.
Definition: MRFZoneList.C:381
incompressibleTwoPhaseInteractingMixture(volVectorField &U, const surfaceScalarField &phi, const uniformDimensionedVectorField &g)
Construct from components.
An abstract base class for incompressible mixtureViscosityModels.
const dimensionSet dimDensity
const Type & value() const
Return const reference to value.
phi
Definition: correctPhi.H:3
word timeName
Definition: getTimeIndex.H:3
Calculate the divergence of the given field.
defineTypeNameAndDebug(combustionModel, 0)
tmp< volVectorField > divTauDm() const
Return the div stress tensor due to the phase transport.
const dimensionedScalar & dd() const
Return the diameter of the dispersed-phase particles.
label patchi
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const volVectorField & Udm() const
Return the diffusion velocity of the dispersed phase.
const mixtureViscosityModel & muModel() const
Return const-access to the mixture viscosityModel.
virtual bool read()=0
Read base phaseProperties dictionary.
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const dimensionedVector & g
const viscosityModel & nucModel() const
Return const-access to the continuous-phase viscosityModel.
List of MRF zones with IO functionality. MRF zones are specified by a list of dictionary entries...
Definition: IOMRFZoneList.H:64
tmp< volScalarField > mu() const
Return the dynamic mixture viscosity.
const volScalarField & alphac() const
Return const-access to the continuous phase-fraction.
virtual tmp< volScalarField > rho() const
Return the mixture density.
Namespace for OpenFOAM.