waxSolventViscosity.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) 2017 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
28#include "waxSolventViscosity.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
36{
37namespace regionModels
38{
39namespace surfaceFilmModels
40{
41
42// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43
45
47(
51);
52
53
54// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
55
56void waxSolventViscosity::correctMu()
57{
58 const kinematicSingleLayer& film = filmType<kinematicSingleLayer>();
59
61 (
63 (
65 )
66 );
67
69 (
71 (
73 )
74 );
75
76 const uniformDimensionedScalarField Ysolvent0
77 (
79 (
81 )
82 );
83
84 const volScalarField& Ysolvent
85 (
87 (
89 )
90 );
91
92 const volScalarField Xsolvent
93 (
94 Ysolvent*Wsolvent/((1 - Ysolvent)*Wwax + Ysolvent*Wsolvent)
95 );
96
97 const dimensionedScalar Xsolvent0
98 (
99 Ysolvent0*Wsolvent/((1 - Ysolvent0)*Wwax + Ysolvent0*Wsolvent)
100 );
101
102 mu_ = pow(muWax_/muSolvent_, (1 - Xsolvent)/(1 - Xsolvent0))*muSolvent_;
104}
105
106
107// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
108
110(
112 const dictionary& dict,
114)
115:
116 filmViscosityModel(typeName, film, dict, mu),
117 muWax_
118 (
120 (
121 typeName + ":muWax",
122 film.regionMesh().time().timeName(),
123 film.regionMesh(),
124 IOobject::NO_READ,
125 IOobject::AUTO_WRITE
126 ),
127 film.regionMesh(),
129 zeroGradientFvPatchScalarField::typeName
130 ),
131 muWaxModel_
132 (
134 (
135 film,
136 coeffDict_.subDict("muWax"),
137 muWax_
138 )
139 ),
140 muSolvent_
141 (
143 (
144 typeName + ":muSolvent",
145 film.regionMesh().time().timeName(),
146 film.regionMesh(),
147 IOobject::NO_READ,
148 IOobject::AUTO_WRITE
149 ),
150 film.regionMesh(),
152 zeroGradientFvPatchScalarField::typeName
153 ),
154 muSolventModel_
155 (
157 (
158 film,
159 coeffDict_.subDict("muSolvent"),
160 muSolvent_
161 )
162 )
163{}
164
165
166// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
167
169(
170 const volScalarField& p,
171 const volScalarField& T
172)
173{
174 muWaxModel_->correct(p, T);
175 muSolventModel_->correct(p, T);
176
177 correctMu();
178}
179
180
181// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
182
183} // End namespace surfaceFilmModels
184} // End namespace regionModels
185} // End namespace Foam
186
187// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
void correctBoundaryConditions()
Correct boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const Type & lookupObject(const word &name, const bool recursive=false) const
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:64
const surfaceFilmRegionModel & film() const
Return const access to the film surface film model.
Base class for surface film viscosity models.
volScalarField & mu_
Reference to the viscosity field.
Kinematic form of single-cell layer surface film model.
autoPtr< filmViscosityModel > muSolventModel_
Solvent viscosity model.
autoPtr< filmViscosityModel > muWaxModel_
Wax viscosity model.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
const volScalarField & T
const volScalarField & mu
word timeName
Definition: getTimeIndex.H:3
Namespace for OpenFOAM.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const dimensionSet dimDynamicViscosity
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.
dictionary dict
static const char *const typeName
The type name used in ensight case files.