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-2020 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 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
40 
42 {
43  nuModel1_->correct();
44  nuModel2_->correct();
45 
46  const volScalarField limitedAlpha1
47  (
48  "limitedAlpha1",
49  min(max(alpha1_, scalar(0)), scalar(1))
50  );
51 
52  // Average kinematic viscosity calculated from dynamic viscosity
53  nu_ = mu()/(limitedAlpha1*rho1_ + (scalar(1) - limitedAlpha1)*rho2_);
54 }
55 
56 
57 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
58 
60 (
61  const volVectorField& U,
62  const surfaceScalarField& phi
63 )
64 :
66  (
67  IOobject
68  (
69  "transportProperties",
70  U.time().constant(),
71  U.db(),
74  )
75  ),
76  twoPhaseMixture(U.mesh(), *this),
77 
78  nuModel1_
79  (
81  (
82  "nu1",
83  subDict(phase1Name_),
84  U,
85  phi
86  )
87  ),
88  nuModel2_
89  (
91  (
92  "nu2",
93  subDict(phase2Name_),
94  U,
95  phi
96  )
97  ),
98 
99  rho1_("rho", dimDensity, nuModel1_->viscosityProperties()),
100  rho2_("rho", dimDensity, nuModel2_->viscosityProperties()),
101 
102  U_(U),
103  phi_(phi),
104 
105  nu_
106  (
107  IOobject
108  (
109  "nu",
110  U_.time().timeName(),
111  U_.db()
112  ),
113  U_.mesh(),
115  calculatedFvPatchScalarField::typeName
116  )
117 {
118  calcNu();
119 }
120 
121 
122 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
123 
126 {
127  const volScalarField limitedAlpha1
128  (
129  min(max(alpha1_, scalar(0)), scalar(1))
130  );
131 
132  return volScalarField::New
133  (
134  "mu",
135  limitedAlpha1*rho1_*nuModel1_->nu()
136  + (scalar(1) - limitedAlpha1)*rho2_*nuModel2_->nu()
137  );
138 }
139 
140 
143 {
144  const surfaceScalarField alpha1f
145  (
146  min(max(fvc::interpolate(alpha1_), scalar(0)), scalar(1))
147  );
148 
150  (
151  "muf",
152  alpha1f*rho1_*fvc::interpolate(nuModel1_->nu())
153  + (scalar(1) - alpha1f)*rho2_*fvc::interpolate(nuModel2_->nu())
154  );
155 }
156 
157 
160 {
161  const surfaceScalarField alpha1f
162  (
163  min(max(fvc::interpolate(alpha1_), scalar(0)), scalar(1))
164  );
165 
167  (
168  "nuf",
169  (
170  alpha1f*rho1_*fvc::interpolate(nuModel1_->nu())
171  + (scalar(1) - alpha1f)*rho2_*fvc::interpolate(nuModel2_->nu())
172  )/(alpha1f*rho1_ + (scalar(1) - alpha1f)*rho2_)
173  );
174 }
175 
176 
178 {
179  if (regIOobject::read())
180  {
181  if
182  (
183  nuModel1_().read
184  (
185  subDict(phase1Name_ == "1" ? "phase1": phase1Name_)
186  )
187  && nuModel2_().read
188  (
189  subDict(phase2Name_ == "2" ? "phase2": phase2Name_)
190  )
191  )
192  {
193  nuModel1_->viscosityProperties().lookup("rho") >> rho1_;
194  nuModel2_->viscosityProperties().lookup("rho") >> rho2_;
195 
196  return true;
197  }
198  else
199  {
200  return false;
201  }
202  }
203  else
204  {
205  return false;
206  }
207 }
208 
209 
210 // ************************************************************************* //
Foam::surfaceFields.
tmp< surfaceScalarField > nuf() const
Return the face-interpolated kinematic laminar viscosity.
virtual bool read()
Read object.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
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,.
const dimensionSet dimViscosity
tmp< surfaceScalarField > muf() const
Return the face-interpolated dynamic laminar viscosity.
Macros for easy insertion into run-time selection tables.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:52
bool read(const char *, int32_t &)
Definition: int32IO.C:85
const word & constant() const
Return constant name.
Definition: TimePaths.H:124
tmp< volScalarField > mu() const
Return the dynamic laminar viscosity.
const Mesh & mesh() const
Return mesh.
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
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
const dimensionedScalar & mu
Atomic mass unit.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const Time & time() const
Return time.
Definition: IOobject.C:328
static autoPtr< viscosityModel > New(const word &name, const dictionary &viscosityProperties, const volVectorField &U, const surfaceScalarField &phi)
Return a reference to the selected viscosity model.
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: PtrList.H:53
const objectRegistry & db() const
Return the local objectRegistry.
Definition: IOobject.C:322
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Namespace for OpenFOAM.