incompressibleTwoPhaseMixture.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) 2011-2021 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 "surfaceInterpolate.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34  defineTypeNameAndDebug(incompressibleTwoPhaseMixture, 0);
35 }
36 
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
41 (
42  const volVectorField& U,
43  const surfaceScalarField& phi
44 )
45 :
46  twoPhaseMixture(U.mesh()),
47 
48  nuModel1_(viscosityModel::New(U.mesh(), phase1Name())),
49  nuModel2_(viscosityModel::New(U.mesh(), phase2Name())),
50 
51  rho1_("rho", dimDensity, nuModel1_()),
52  rho2_("rho", dimDensity, nuModel2_()),
53 
54  U_(U),
55 
56  nu_
57  (
58  IOobject
59  (
60  "nu",
61  U_.time().timeName(),
62  U_.db()
63  ),
64  U_.mesh(),
66  calculatedFvPatchScalarField::typeName
67  )
68 {
69  correct();
70 }
71 
72 
73 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
74 
77 {
78  return volScalarField::New
79  (
80  "rho",
81  alpha1()*rho1_ + (scalar(1) - alpha1())*rho2_
82  );
83 }
84 
85 
88 {
89  const volScalarField limitedAlpha1
90  (
91  min(max(alpha1(), scalar(0)), scalar(1))
92  );
93 
94  return volScalarField::New
95  (
96  "mu",
97  limitedAlpha1*rho1_*nuModel1_->nu()
98  + (scalar(1) - limitedAlpha1)*rho2_*nuModel2_->nu()
99  );
100 }
101 
102 
105 {
106  const surfaceScalarField alpha1f
107  (
108  min(max(fvc::interpolate(alpha1()), scalar(0)), scalar(1))
109  );
110 
112  (
113  "muf",
114  alpha1f*rho1_*fvc::interpolate(nuModel1_->nu())
115  + (scalar(1) - alpha1f)*rho2_*fvc::interpolate(nuModel2_->nu())
116  );
117 }
118 
119 
122 {
123  const surfaceScalarField alpha1f
124  (
125  min(max(fvc::interpolate(alpha1()), scalar(0)), scalar(1))
126  );
127 
129  (
130  "nuf",
131  (
132  alpha1f*rho1_*fvc::interpolate(nuModel1_->nu())
133  + (scalar(1) - alpha1f)*rho2_*fvc::interpolate(nuModel2_->nu())
134  )/(alpha1f*rho1_ + (scalar(1) - alpha1f)*rho2_)
135  );
136 }
137 
138 
140 {
141  if (twoPhaseMixture::read())
142  {
143  nuModel1_->lookup("rho") >> rho1_;
144  nuModel2_->lookup("rho") >> rho2_;
145 
146  return true;
147  }
148  else
149  {
150  return false;
151  }
152 }
153 
154 
156 {
157  nuModel1_->correct();
158  nuModel2_->correct();
159 
160  const volScalarField limitedAlpha1
161  (
162  "limitedAlpha1",
163  min(max(alpha1(), scalar(0)), scalar(1))
164  );
165 
166  // Average kinematic viscosity calculated from dynamic viscosity
167  nu_ = mu()/(limitedAlpha1*rho1_ + (scalar(1) - limitedAlpha1)*rho2_);
168 }
169 
170 
171 // ************************************************************************* //
const dimensionSet dimViscosity
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
tmp< surfaceScalarField > nuf() const
Return the face-interpolated kinematic laminar viscosity.
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
incompressibleTwoPhaseMixture(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
volScalarField & alpha1(mixture.alpha1())
tmp< volScalarField > rho() const
Return the mixture density.
tmp< surfaceScalarField > muf() const
Return the face-interpolated dynamic laminar viscosity.
Macros for easy insertion into run-time selection tables.
const dimensionSet dimDensity
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
tmp< volScalarField > mu() const
Return the dynamic laminar viscosity.
static autoPtr< viscosityModel > New(const fvMesh &mesh, const word &group=word::null)
Return a reference to the selected viscosity model.
const dimensionedScalar mu
Atomic mass unit.
const Mesh & mesh() const
Return mesh.
defineTypeNameAndDebug(combustionModel, 0)
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
virtual void correct()
Correct the laminar viscosity.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A two-phase mixture model.
virtual bool read()
Read base phaseProperties dictionary.
virtual bool read()=0
Read base phaseProperties dictionary.
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
Namespace for OpenFOAM.