Loading...
Searching...
No Matches
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
31#include "fvc.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
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,
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 )
118{
120}
121
122
123// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
124
127{
128 const volScalarField limitedAlpha1
129 (
130 clamp(alpha1_, zero_one{})
131 );
132
134 (
135 "mu",
137 limitedAlpha1*rho1_*nuModel1_->nu()
138 + (scalar(1) - limitedAlpha1)*rho2_*nuModel2_->nu()
139 );
140}
141
142
144Foam::incompressibleTwoPhaseMixture::mu(const label patchI) const
146
147 return mu()().boundaryField()[patchI];
148}
149
150
153{
154 const surfaceScalarField alpha1f
155 (
156 clamp(fvc::interpolate(alpha1_), zero_one{})
157 );
158
160 (
161 "muf",
164 + (scalar(1) - alpha1f)*rho2_*fvc::interpolate(nuModel2_->nu())
165 );
166}
167
168
171{
172 const surfaceScalarField alpha1f
173 (
174 clamp(fvc::interpolate(alpha1_), zero_one{})
175 );
176
178 (
179 "nuf",
181 (
183 + (scalar(1) - alpha1f)*rho2_*fvc::interpolate(nuModel2_->nu())
184 )/(alpha1f*rho1_ + (scalar(1) - alpha1f)*rho2_)
185 );
186}
187
188
190{
191 if (regIOobject::read())
192 {
193 if
194 (
195 nuModel1_().read
196 (
197 subDict(phase1Name_ == "1" ? "phase1": phase1Name_)
198 )
199 && nuModel2_().read
200 (
201 subDict(phase2Name_ == "2" ? "phase2": phase2Name_)
202 )
203 )
204 {
205 nuModel1_->viscosityProperties().readEntry("rho", rho1_);
206 nuModel2_->viscosityProperties().readEntry("rho", rho2_);
207
208 return true;
209 }
210 }
211
212 return false;
213}
214
215
216// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
IOdictionary(const IOobject &io, const dictionary *fallback=nullptr)
Construct given an IOobject and optional fallback dictionary content.
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition IOobject.C:450
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition dictionary.C:441
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.
const volVectorField & U() const
Return const-access to the mixture velocity.
tmp< volScalarField > mu() const
Return the dynamic laminar viscosity.
void calcNu()
Calculate and return the laminar viscosity.
incompressibleTwoPhaseMixture(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
virtual bool read()
Read base transportProperties dictionary.
virtual bool read()
Read object.
A class for managing temporary objects.
Definition tmp.H:75
A two-phase mixture model.
twoPhaseMixture(const fvMesh &mesh, const dictionary &dict)
Construct from components.
An abstract base class for incompressible viscosityModels.
Represents 0/1 range or concept. Used for tagged dispatch or clamping.
Definition pTraits.H:51
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
dynamicFvMesh & mesh
word timeName
Definition getTimeIndex.H:3
Different types of constants.
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
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.
GeometricField< vector, fvPatchField, volMesh > volVectorField
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition int32.H:127
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimDensity
Foam::surfaceFields.