incompressibleThreePhaseMixture.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  Copyright (C) 2019-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
31 #include "surfaceFields.H"
32 #include "fvc.H"
33 
34 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
35 
36 void Foam::incompressibleThreePhaseMixture::calcNu()
37 {
38  nuModel1_->correct();
39  nuModel2_->correct();
40  nuModel3_->correct();
41 
42  // Average kinematic viscosity calculated from dynamic viscosity
43  nu_ = mu()/(alpha1_*rho1_ + alpha2_*rho2_ + alpha3_*rho3_);
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
50 (
51  const volVectorField& U,
52  const surfaceScalarField& phi
53 )
54 :
55  IOdictionary
56  (
57  IOobject
58  (
59  "transportProperties",
60  U.time().constant(),
61  U.db(),
62  IOobject::MUST_READ_IF_MODIFIED,
63  IOobject::NO_WRITE
64  )
65  ),
66 
67  phase1Name_(get<wordList>("phases")[0]),
68  phase2Name_(get<wordList>("phases")[1]),
69  phase3Name_(get<wordList>("phases")[2]),
70 
71  alpha1_
72  (
73  IOobject
74  (
75  IOobject::groupName("alpha", phase1Name_),
76  U.time().timeName(),
77  U.mesh(),
78  IOobject::MUST_READ,
79  IOobject::AUTO_WRITE
80  ),
81  U.mesh()
82  ),
83 
84  alpha2_
85  (
86  IOobject
87  (
88  IOobject::groupName("alpha", phase2Name_),
89  U.time().timeName(),
90  U.mesh(),
91  IOobject::MUST_READ,
92  IOobject::AUTO_WRITE
93  ),
94  U.mesh()
95  ),
96 
97  alpha3_
98  (
99  IOobject
100  (
101  IOobject::groupName("alpha", phase3Name_),
102  U.time().timeName(),
103  U.mesh(),
104  IOobject::MUST_READ,
105  IOobject::AUTO_WRITE
106  ),
107  U.mesh()
108  ),
109 
110  U_(U),
111  phi_(phi),
112 
113  nu_
114  (
115  IOobject
116  (
117  "nu",
118  U.time().timeName(),
119  U.db()
120  ),
121  U.mesh(),
122  dimensionedScalar(dimensionSet(0, 2, -1, 0, 0), Zero)
123  ),
124 
125  nuModel1_
126  (
127  viscosityModel::New
128  (
129  "nu1",
130  subDict(phase1Name_),
131  U,
132  phi
133  )
134  ),
135  nuModel2_
136  (
137  viscosityModel::New
138  (
139  "nu2",
140  subDict(phase2Name_),
141  U,
142  phi
143  )
144  ),
145  nuModel3_
146  (
147  viscosityModel::New
148  (
149  "nu3",
150  subDict(phase3Name_),
151  U,
152  phi
153  )
154  ),
155 
156  rho1_("rho", dimDensity, nuModel1_->viscosityProperties()),
157  rho2_("rho", dimDensity, nuModel2_->viscosityProperties()),
158  rho3_("rho", dimDensity, nuModel3_->viscosityProperties())
159 {
160  alpha3_ == 1.0 - alpha1_ - alpha2_;
161  calcNu();
162 }
163 
164 
165 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
166 
169 {
170  return tmp<volScalarField>
171  (
172  new volScalarField
173  (
174  "mu",
175  alpha1_*rho1_*nuModel1_->nu()
176  + alpha2_*rho2_*nuModel2_->nu()
177  + alpha3_*rho3_*nuModel3_->nu()
178  )
179  );
180 }
181 
182 
185 {
186  surfaceScalarField alpha1f(fvc::interpolate(alpha1_));
187  surfaceScalarField alpha2f(fvc::interpolate(alpha2_));
188  surfaceScalarField alpha3f(fvc::interpolate(alpha3_));
189 
190  return tmp<surfaceScalarField>
191  (
193  (
194  "mu",
195  alpha1f*rho1_*fvc::interpolate(nuModel1_->nu())
196  + alpha2f*rho2_*fvc::interpolate(nuModel2_->nu())
197  + alpha3f*rho3_*fvc::interpolate(nuModel3_->nu())
198  )
199  );
200 }
201 
202 
205 {
206  surfaceScalarField alpha1f(fvc::interpolate(alpha1_));
207  surfaceScalarField alpha2f(fvc::interpolate(alpha2_));
208  surfaceScalarField alpha3f(fvc::interpolate(alpha3_));
209 
210  return tmp<surfaceScalarField>
211  (
213  (
214  "nu",
215  (
216  alpha1f*rho1_*fvc::interpolate(nuModel1_->nu())
217  + alpha2f*rho2_*fvc::interpolate(nuModel2_->nu())
218  + alpha3f*rho3_*fvc::interpolate(nuModel3_->nu())
219  )/(alpha1f*rho1_ + alpha2f*rho2_ + alpha3f*rho3_)
220  )
221  );
222 }
223 
224 
226 {
227  if (transportModel::read())
228  {
229  if
230  (
231  nuModel1_().read(*this)
232  && nuModel2_().read(*this)
233  && nuModel3_().read(*this)
234  )
235  {
236  nuModel1_->viscosityProperties().readEntry("rho", rho1_);
237  nuModel2_->viscosityProperties().readEntry("rho", rho2_);
238  nuModel3_->viscosityProperties().readEntry("rho", rho3_);
239 
240  return true;
241  }
242  }
243 
244  return false;
245 }
246 
247 
248 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
Foam::surfaceFields.
incompressibleThreePhaseMixture(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
tmp< volScalarField > mu() const
Return the dynamic laminar viscosity.
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
Definition: volFieldsFwd.H:82
Macros for easy insertion into run-time selection tables.
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:127
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
word timeName
Definition: getTimeIndex.H:3
dynamicFvMesh & mesh
tmp< surfaceScalarField > muf() const
Return the face-interpolated dynamic laminar viscosity.
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
List< word > wordList
List of word.
Definition: fileName.H:58
U
Definition: pEqn.H:72
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual bool read()=0
Read transportProperties dictionary.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
tmp< surfaceScalarField > nuf() const
Return the face-interpolated dynamic laminar viscosity.
bool read()
Read base transportProperties dictionary.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133