All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
incompressibleThreePhaseMixture.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 
28 #include "surfaceFields.H"
29 #include "fvc.H"
30 
31 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
32 
33 void Foam::incompressibleThreePhaseMixture::calcNu()
34 {
35  nuModel1_->correct();
36  nuModel2_->correct();
37  nuModel3_->correct();
38 
39  // Average kinematic viscosity calculated from dynamic viscosity
40  nu_ = mu()/(alpha1_*rho1_ + alpha2_*rho2_ + alpha3_*rho3_);
41 }
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
47 (
48  const volVectorField& U,
49  const surfaceScalarField& phi
50 )
51 :
53  (
54  IOobject
55  (
56  "phaseProperties",
57  U.time().constant(),
58  U.db(),
61  )
62  ),
63 
64  phase1Name_(wordList(lookup("phases"))[0]),
65  phase2Name_(wordList(lookup("phases"))[1]),
66  phase3Name_(wordList(lookup("phases"))[2]),
67 
68  alpha1_
69  (
70  IOobject
71  (
72  IOobject::groupName("alpha", phase1Name_),
73  U.time().timeName(),
74  U.mesh(),
77  ),
78  U.mesh()
79  ),
80 
81  alpha2_
82  (
83  IOobject
84  (
85  IOobject::groupName("alpha", phase2Name_),
86  U.time().timeName(),
87  U.mesh(),
90  ),
91  U.mesh()
92  ),
93 
94  alpha3_
95  (
96  IOobject
97  (
98  IOobject::groupName("alpha", phase3Name_),
99  U.time().timeName(),
100  U.mesh(),
103  ),
104  U.mesh()
105  ),
106 
107  U_(U),
108  phi_(phi),
109 
110  nu_
111  (
112  IOobject
113  (
114  "nu",
115  U.time().timeName(),
116  U.db()
117  ),
118  U.mesh(),
119  dimensionedScalar(dimensionSet(0, 2, -1, 0, 0), 0),
120  calculatedFvPatchScalarField::typeName
121  ),
122 
123  nuModel1_(viscosityModel::New(U.mesh(), phase1Name_)),
124  nuModel2_(viscosityModel::New(U.mesh(), phase2Name_)),
125  nuModel3_(viscosityModel::New(U.mesh(), phase3Name_)),
126 
127  rho1_("rho", dimDensity, nuModel1_()),
128  rho2_("rho", dimDensity, nuModel2_()),
129  rho3_("rho", dimDensity, nuModel3_())
130 {
131  alpha3_ == 1.0 - alpha1_ - alpha2_;
132  calcNu();
133 }
134 
135 
136 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
137 
140 {
141  return volScalarField::New
142  (
143  "mu",
144  alpha1_*rho1_*nuModel1_->nu()
145  + alpha2_*rho2_*nuModel2_->nu()
146  + alpha3_*rho3_*nuModel3_->nu()
147  );
148 }
149 
150 
153 {
154  surfaceScalarField alpha1f(fvc::interpolate(alpha1_));
155  surfaceScalarField alpha2f(fvc::interpolate(alpha2_));
156  surfaceScalarField alpha3f(fvc::interpolate(alpha3_));
157 
159  (
160  "mu",
161  alpha1f*rho1_*fvc::interpolate(nuModel1_->nu())
162  + alpha2f*rho2_*fvc::interpolate(nuModel2_->nu())
163  + alpha3f*rho3_*fvc::interpolate(nuModel3_->nu())
164  );
165 }
166 
167 
170 {
171  surfaceScalarField alpha1f(fvc::interpolate(alpha1_));
172  surfaceScalarField alpha2f(fvc::interpolate(alpha2_));
173  surfaceScalarField alpha3f(fvc::interpolate(alpha3_));
174 
176  (
177  "nu",
178  (
179  alpha1f*rho1_*fvc::interpolate(nuModel1_->nu())
180  + alpha2f*rho2_*fvc::interpolate(nuModel2_->nu())
181  + alpha3f*rho3_*fvc::interpolate(nuModel3_->nu())
182  )/(alpha1f*rho1_ + alpha2f*rho2_ + alpha3f*rho3_)
183  );
184 }
185 
186 
188 {
189  if (regIOobject::read())
190  {
191  nuModel1_->lookup("rho") >> rho1_;
192  nuModel2_->lookup("rho") >> rho2_;
193  nuModel3_->lookup("rho") >> rho3_;
194 
195  return true;
196  }
197  else
198  {
199  return false;
200  }
201 }
202 
203 
204 // ************************************************************************* //
Foam::surfaceFields.
const surfaceScalarField & phi() const
Return the flux.
incompressibleThreePhaseMixture(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
tmp< volScalarField > mu() const
Return the dynamic laminar viscosity.
virtual bool read()
Read object.
IOobject(const word &name, const fileName &instance, const objectRegistry &registry, readOption r=NO_READ, writeOption w=NO_WRITE, bool registerObject=true)
Construct from name, instance, registry, io options.
Definition: IOobject.C:167
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:59
fvMesh & mesh
Macros for easy insertion into run-time selection tables.
static word groupName(Name name, const word &group)
const dimensionSet dimDensity
word timeName
Definition: getTimeIndex.H:3
tmp< surfaceScalarField > muf() const
Return the face-interpolated dynamic laminar viscosity.
IOdictionary(const IOobject &io, const word &wantedType)
Construct given an IOobject, supply wanted typeName.
Definition: IOdictionary.C:47
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:96
const volVectorField & U() const
Return the velocity.
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.
List< word > wordList
A List of words.
Definition: fileName.H:54
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const Time & time() const
Return time.
Definition: IOobject.C:318
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:312
tmp< surfaceScalarField > nuf() const
Return the face-interpolated dynamic laminar viscosity.
bool read()
Read base phaseProperties dictionary.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864