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-------------------------------------------------------------------------------
10License
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
35namespace 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",
51 min(max(alpha1_, scalar(0)), scalar(1))
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,
65)
66:
68 (
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 (
110 (
111 "nu",
112 U_.time().timeName(),
113 U_.db()
114 ),
115 U_.mesh(),
117 calculatedFvPatchScalarField::typeName
118 )
119{
120 calcNu();
121}
122
123
124// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
125
128{
129 const volScalarField limitedAlpha1
130 (
131 min(max(alpha1_, scalar(0)), scalar(1))
132 );
133
135 (
136 "mu",
137 limitedAlpha1*rho1_*nuModel1_->nu()
138 + (scalar(1) - limitedAlpha1)*rho2_*nuModel2_->nu()
139 );
140}
141
142
145{
146
147 return mu()().boundaryField()[patchI];
148}
149
150
153{
154 const surfaceScalarField alpha1f
155 (
156 min(max(fvc::interpolate(alpha1_), scalar(0)), scalar(1))
157 );
158
160 (
161 "muf",
162 alpha1f*rho1_*fvc::interpolate(nuModel1_->nu())
163 + (scalar(1) - alpha1f)*rho2_*fvc::interpolate(nuModel2_->nu())
164 );
165}
166
167
170{
171 const surfaceScalarField alpha1f
172 (
173 min(max(fvc::interpolate(alpha1_), scalar(0)), scalar(1))
174 );
175
177 (
178 "nuf",
179 (
180 alpha1f*rho1_*fvc::interpolate(nuModel1_->nu())
181 + (scalar(1) - alpha1f)*rho2_*fvc::interpolate(nuModel2_->nu())
182 )/(alpha1f*rho1_ + (scalar(1) - alpha1f)*rho2_)
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().readEntry("rho", rho1_);
204 nuModel2_->viscosityProperties().readEntry("rho", rho2_);
205
206 return true;
207 }
208 }
209
210 return false;
211}
212
213
214// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
surfaceScalarField & phi
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
A two-phase incompressible transportModel.
tmp< surfaceScalarField > muf() const
Return the face-interpolated dynamic laminar viscosity.
tmp< surfaceScalarField > nuf() const
Return the face-interpolated kinematic laminar viscosity.
tmp< volScalarField > mu() const
Return the dynamic laminar viscosity.
void calcNu()
Calculate and return the laminar viscosity.
virtual bool read()
Read base transportProperties dictionary.
constant condensation/saturation model.
virtual bool read()
Read object.
A class for managing temporary objects.
Definition: tmp.H:65
A two-phase mixture model.
volScalarField alpha1_
An abstract base class for incompressible viscosityModels.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
U
Definition: pEqn.H:72
const volScalarField & mu
dynamicFvMesh & mesh
word timeName
Definition: getTimeIndex.H:3
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.
Namespace for OpenFOAM.
const dimensionSet dimViscosity
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
const dimensionSet dimDensity
Foam::surfaceFields.