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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
30 #include "surfaceFields.H"
31 #include "fvc.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
38 }
39 
40 
41 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
42 
44 {
45  nuModel1_->correct();
46  nuModel2_->correct();
47 
48  const volScalarField limitedAlpha1
49  (
50  "limitedAlpha1",
52  );
53 
54  // Average kinematic viscosity calculated from dynamic viscosity
55  nu_ = mu()/(limitedAlpha1*rho1_ + (scalar(1) - limitedAlpha1)*rho2_);
56 }
57 
58 
59 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
60 
62 (
63  const volVectorField& U,
64  const surfaceScalarField& phi
65 )
66 :
68  (
69  IOobject
70  (
71  "transportProperties",
72  U.time().constant(),
73  U.db(),
74  IOobject::MUST_READ_IF_MODIFIED,
75  IOobject::NO_WRITE
76  )
77  ),
78  twoPhaseMixture(U.mesh(), *this),
79 
80  nuModel1_
81  (
83  (
84  "nu1",
85  subDict(phase1Name_),
86  U,
87  phi
88  )
89  ),
90  nuModel2_
91  (
93  (
94  "nu2",
95  subDict(phase2Name_),
96  U,
97  phi
98  )
99  ),
100 
101  rho1_("rho", dimDensity, nuModel1_->viscosityProperties()),
102  rho2_("rho", dimDensity, nuModel2_->viscosityProperties()),
103 
104  U_(U),
105  phi_(phi),
106 
107  nu_
108  (
109  IOobject
110  (
111  "nu",
112  U_.time().timeName(),
113  U_.db()
114  ),
115  U_.mesh(),
117  )
118 {
120 }
121 
122 
123 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
124 
127 {
128  const volScalarField limitedAlpha1
129  (
130  clamp(alpha1_, zero_one{})
131  );
132 
134  (
135  "mu",
136  limitedAlpha1*rho1_*nuModel1_->nu()
137  + (scalar(1) - limitedAlpha1)*rho2_*nuModel2_->nu()
138  );
139 }
140 
141 
143 Foam::incompressibleTwoPhaseMixture::mu(const label patchI) const
144 {
145 
146  return mu()().boundaryField()[patchI];
147 }
148 
149 
152 {
153  const surfaceScalarField alpha1f
154  (
155  clamp(fvc::interpolate(alpha1_), zero_one{})
156  );
157 
159  (
160  "muf",
161  alpha1f*rho1_*fvc::interpolate(nuModel1_->nu())
162  + (scalar(1) - alpha1f)*rho2_*fvc::interpolate(nuModel2_->nu())
163  );
164 }
165 
166 
169 {
170  const surfaceScalarField alpha1f
171  (
172  clamp(fvc::interpolate(alpha1_), zero_one{})
173  );
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 
187 {
188  if (regIOobject::read())
189  {
190  if
191  (
192  nuModel1_().read
193  (
194  subDict(phase1Name_ == "1" ? "phase1": phase1Name_)
195  )
196  && nuModel2_().read
197  (
198  subDict(phase2Name_ == "2" ? "phase2": phase2Name_)
199  )
200  )
201  {
202  nuModel1_->viscosityProperties().readEntry("rho", rho1_);
203  nuModel2_->viscosityProperties().readEntry("rho", rho2_);
204 
205  return true;
206  }
207  }
208 
209  return false;
210 }
211 
212 
213 // ************************************************************************* //
Foam::surfaceFields.
tmp< surfaceScalarField > nuf() const
Return the face-interpolated kinematic laminar viscosity.
A two-phase incompressible transportModel.
virtual bool read()
Read object.
incompressibleTwoPhaseMixture(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
const dimensionSet dimViscosity
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
An abstract base class for incompressible viscosityModels.
tmp< surfaceScalarField > muf() const
Return the face-interpolated dynamic laminar viscosity.
Macros for easy insertion into run-time selection tables.
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:127
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
word timeName
Definition: getTimeIndex.H:3
dynamicFvMesh & mesh
volScalarField alpha1_
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
tmp< volScalarField > mu() const
Return the dynamic laminar viscosity.
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.
const dimensionSet dimDensity
Represents 0/1 range or concept. Used for tagged dispatch or clamping.
Definition: pTraits.H:65
const dimensionedScalar mu
Atomic mass unit.
U
Definition: pEqn.H:72
A two-phase mixture model.
virtual bool read()
Read base transportProperties dictionary.
void calcNu()
Calculate and return the laminar viscosity.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Namespace for OpenFOAM.
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
Definition: dimensionSet.C:271
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127