Loading...
Searching...
No Matches
incompressibleTwoPhaseInteractingMixture.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) 2014-2015 OpenFOAM Foundation
9 Copyright (C) 2019-2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
39}
40
41
42// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43
46(
47 const volVectorField& U,
48 const surfaceScalarField& phi
49)
50:
51 IOdictionary
52 (
53 IOobject
54 (
55 "transportProperties",
56 U.time().constant(),
57 U.db(),
58 IOobject::MUST_READ_IF_MODIFIED,
59 IOobject::NO_WRITE
60 )
61 ),
62 twoPhaseMixture(U.mesh(), *this),
63
64 muModel_
65 (
66 mixtureViscosityModel::New
67 (
68 "mu",
69 subDict(phase1Name_),
70 U,
71 phi
72 )
73 ),
74
75 nucModel_
76 (
77 viscosityModel::New
78 (
79 "nuc",
80 subDict(phase2Name_),
81 U,
82 phi
83 )
84 ),
85
86 rhod_("rho", dimDensity, muModel_->viscosityProperties()),
87 rhoc_("rho", dimDensity, nucModel_->viscosityProperties()),
88 dd_
89 (
90 "d",
92 muModel_->viscosityProperties().getOrDefault("d", 0.0)
93 ),
94 alphaMax_(muModel_->viscosityProperties().getOrDefault("alphaMax", 1.0)),
95
96 U_(U),
97 phi_(phi),
98
99 mu_
100 (
101 IOobject
102 (
103 "mu",
104 U_.time().timeName(),
105 U_.db()
106 ),
107 U_.mesh(),
108 dimensionedScalar(dimensionSet(1, -1, -1, 0, 0), Zero)
109 )
110{
111 correct();
112}
113
114
115// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
116
118{
119 if (regIOobject::read())
120 {
121 if
122 (
125 )
126 {
127 muModel_->viscosityProperties().readEntry("rho", rhod_);
128 nucModel_->viscosityProperties().readEntry("rho", rhoc_);
129
131 (
132 "d",
133 dimLength,
134 muModel_->viscosityProperties().getOrDefault("d", 0)
135 );
136
137 alphaMax_ =
138 muModel_->viscosityProperties().getOrDefault
139 (
140 "alphaMax",
141 1.0
142 );
143
144 return true;
145 }
146 }
147
148 return false;
149}
150
151
152// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
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 for interacting phases requiring the direct evaluation of t...
incompressibleTwoPhaseInteractingMixture(const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
dimensionedScalar dd_
Optional diameter of the dispersed phase particles.
scalar alphaMax_
Optional maximum dispersed phase-fraction (e.g. packing limit).
virtual bool read()
Read base transportProperties dictionary.
constant condensation/saturation model.
virtual bool read()
Read object.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
thermo correct()
dynamicFvMesh & mesh
word timeName
Definition getTimeIndex.H:3
Namespace for OpenFOAM.
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.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimDensity
Foam::surfaceFields.