incompressibleTwoPhaseMixture.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 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 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  defineTypeNameAndDebug(incompressibleTwoPhaseMixture, 0);
36 }
37 
38 
39 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
40 
41 //- Calculate and return the laminar viscosity
43 {
44  nuModel1_->correct();
45  nuModel2_->correct();
46 
47  const volScalarField limitedAlpha1
48  (
49  "limitedAlpha1",
50  min(max(alpha1_, scalar(0)), scalar(1))
51  );
52 
53  // Average kinematic viscosity calculated from dynamic viscosity
54  nu_ = mu()/(limitedAlpha1*rho1_ + (scalar(1) - limitedAlpha1)*rho2_);
55 }
56 
57 
58 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
59 
61 (
62  const volVectorField& U,
63  const surfaceScalarField& phi
64 )
65 :
67  (
68  IOobject
69  (
70  "transportProperties",
71  U.time().constant(),
72  U.db(),
75  )
76  ),
77  twoPhaseMixture(U.mesh(), *this),
78 
79  nuModel1_
80  (
82  (
83  "nu1",
84  subDict(phase1Name_),
85  U,
86  phi
87  )
88  ),
89  nuModel2_
90  (
92  (
93  "nu2",
94  subDict(phase2Name_),
95  U,
96  phi
97  )
98  ),
99 
100  rho1_("rho", dimDensity, nuModel1_->viscosityProperties()),
101  rho2_("rho", dimDensity, nuModel2_->viscosityProperties()),
102 
103  U_(U),
104  phi_(phi),
105 
106  nu_
107  (
108  IOobject
109  (
110  "nu",
111  U_.time().timeName(),
112  U_.db()
113  ),
114  U_.mesh(),
116  calculatedFvPatchScalarField::typeName
117  )
118 {
119  calcNu();
120 }
121 
122 
123 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
124 
127 {
128  const volScalarField limitedAlpha1
129  (
130  min(max(alpha1_, scalar(0)), scalar(1))
131  );
132 
133  return tmp<volScalarField>
134  (
135  new volScalarField
136  (
137  "mu",
138  limitedAlpha1*rho1_*nuModel1_->nu()
139  + (scalar(1) - limitedAlpha1)*rho2_*nuModel2_->nu()
140  )
141  );
142 }
143 
144 
147 {
148  const surfaceScalarField alpha1f
149  (
150  min(max(fvc::interpolate(alpha1_), scalar(0)), scalar(1))
151  );
152 
154  (
156  (
157  "muf",
158  alpha1f*rho1_*fvc::interpolate(nuModel1_->nu())
159  + (scalar(1) - alpha1f)*rho2_*fvc::interpolate(nuModel2_->nu())
160  )
161  );
162 }
163 
164 
167 {
168  const surfaceScalarField alpha1f
169  (
170  min(max(fvc::interpolate(alpha1_), scalar(0)), scalar(1))
171  );
172 
174  (
176  (
177  "nuf",
178  (
179  alpha1f*rho1_*fvc::interpolate(nuModel1_->nu())
180  + (scalar(1) - alpha1f)*rho2_*fvc::interpolate(nuModel2_->nu())
181  )/(alpha1f*rho1_ + (scalar(1) - alpha1f)*rho2_)
182  )
183  );
184 }
185 
186 
188 {
189  if (regIOobject::read())
190  {
191  if
192  (
193  nuModel1_().read
194  (
195  subDict(phase1Name_ == "1" ? "phase1": phase1Name_)
196  )
197  && nuModel2_().read
198  (
199  subDict(phase2Name_ == "2" ? "phase2": phase2Name_)
200  )
201  )
202  {
203  nuModel1_->viscosityProperties().lookup("rho") >> rho1_;
204  nuModel2_->viscosityProperties().lookup("rho") >> rho2_;
205 
206  return true;
207  }
208  else
209  {
210  return false;
211  }
212  }
213  else
214  {
215  return false;
216  }
217 }
218 
219 
220 // ************************************************************************* //
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::incompressibleTwoPhaseMixture::mu
tmp< volScalarField > mu() const
Return the dynamic laminar viscosity.
Definition: incompressibleTwoPhaseMixture.C:126
Foam::incompressibleTwoPhaseMixture::nu_
volScalarField nu_
Definition: incompressibleTwoPhaseMixture.H:71
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::dimDensity
const dimensionSet dimDensity
Foam::incompressibleTwoPhaseMixture::read
virtual bool read()
Read base transportProperties dictionary.
Definition: incompressibleTwoPhaseMixture.C:187
Foam::incompressibleTwoPhaseMixture::nuModel2_
autoPtr< viscosityModel > nuModel2_
Definition: incompressibleTwoPhaseMixture.H:63
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::read
bool read(const char *, int32_t &)
Definition: int32IO.C:87
Foam::incompressibleTwoPhaseMixture::muf
tmp< surfaceScalarField > muf() const
Return the face-interpolated dynamic laminar viscosity.
Definition: incompressibleTwoPhaseMixture.C:146
Foam::fvc::interpolate
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Definition: surfaceInterpolate.C:110
Foam::regIOobject::read
virtual bool read()
Read object.
Definition: regIOobjectRead.C:171
incompressibleTwoPhaseMixture.H
Foam::incompressibleTwoPhaseMixture::nuModel1_
autoPtr< viscosityModel > nuModel1_
Definition: incompressibleTwoPhaseMixture.H:62
surfaceFields.H
Foam::surfaceFields.
Foam::incompressibleTwoPhaseMixture::rho1_
dimensionedScalar rho1_
Definition: incompressibleTwoPhaseMixture.H:65
Foam::twoPhaseMixture::alpha1_
volScalarField alpha1_
Definition: twoPhaseMixture.H:57
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
U
U
Definition: pEqn.H:46
Foam::IOobject::MUST_READ_IF_MODIFIED
@ MUST_READ_IF_MODIFIED
Definition: IOobject.H:109
Foam::twoPhaseMixture
A two-phase mixture model.
Definition: twoPhaseMixture.H:48
Foam::dimViscosity
const dimensionSet dimViscosity
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::viscosityModel::New
static autoPtr< viscosityModel > New(const word &name, const dictionary &viscosityProperties, const volVectorField &U, const surfaceScalarField &phi)
Return a reference to the selected viscosity model.
Definition: viscosityModelNew.C:33
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::incompressibleTwoPhaseMixture::rho2_
dimensionedScalar rho2_
Definition: incompressibleTwoPhaseMixture.H:66
Foam::incompressibleTwoPhaseMixture::nuf
tmp< surfaceScalarField > nuf() const
Return the face-interpolated kinematic laminar viscosity.
Definition: incompressibleTwoPhaseMixture.C:166
Foam::incompressibleTwoPhaseMixture::calcNu
void calcNu()
Calculate and return the laminar viscosity.
Definition: incompressibleTwoPhaseMixture.C:42
fvc.H
Foam::incompressibleTwoPhaseMixture::incompressibleTwoPhaseMixture
incompressibleTwoPhaseMixture(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
Definition: incompressibleTwoPhaseMixture.C:61
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)