Loading...
Searching...
No Matches
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 Copyright (C) 2023 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
29#include "waxSolventViscosity.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
37{
38namespace regionModels
39{
41{
42
43// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44
46
48(
52);
53
54
55// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
56
57void waxSolventViscosity::correctMu()
58{
60
61 const auto& Wwax =
63 (
64 IOobject::scopedName(waxSolventEvaporation::typeName, "Wwax")
65 );
66
67 const auto& Wsolvent =
69 (
70 IOobject::scopedName(waxSolventEvaporation::typeName, "Wsolvent")
71 );
72
73 const auto& Ysolvent0 =
75 (
76 IOobject::scopedName(waxSolventEvaporation::typeName, "Ysolvent0")
77 );
78
79 const auto& Ysolvent =
81 (
82 IOobject::scopedName(waxSolventEvaporation::typeName, "Ysolvent")
83 );
84
85 const volScalarField Xsolvent
86 (
87 Ysolvent*Wsolvent/((1 - Ysolvent)*Wwax + Ysolvent*Wsolvent)
88 );
89
90 const dimensionedScalar Xsolvent0
91 (
92 Ysolvent0*Wsolvent/((1 - Ysolvent0)*Wwax + Ysolvent0*Wsolvent)
93 );
94
95 mu_ = pow(muWax_/muSolvent_, (1 - Xsolvent)/(1 - Xsolvent0))*muSolvent_;
96 mu_.correctBoundaryConditions();
97}
98
99
100// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
101
102waxSolventViscosity::waxSolventViscosity
103(
105 const dictionary& dict,
107)
108:
110 muWax_
111 (
113 (
114 IOobject::scopedName(typeName, "muWax"),
115 film.regionMesh().time().timeName(),
116 film.regionMesh().thisDb(),
117 IOobject::NO_READ,
118 IOobject::AUTO_WRITE,
119 IOobject::REGISTER
120 ),
121 film.regionMesh(),
123 fvPatchFieldBase::zeroGradientType()
124 ),
125 muWaxModel_
126 (
128 (
129 film,
130 coeffDict_.subDict("muWax"),
131 muWax_
132 )
133 ),
134 muSolvent_
135 (
137 (
138 IOobject::scopedName(typeName, "muSolvent"),
139 film.regionMesh().time().timeName(),
140 film.regionMesh().thisDb(),
141 IOobject::NO_READ,
142 IOobject::AUTO_WRITE,
143 IOobject::REGISTER
144 ),
145 film.regionMesh(),
147 fvPatchFieldBase::zeroGradientType()
148 ),
149 muSolventModel_
150 (
152 (
153 film,
154 coeffDict_.subDict("muSolvent"),
155 muSolvent_
157 )
158{}
159
160
161// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
162
164(
165 const volScalarField& p,
166 const volScalarField& T
167)
168{
170 muSolventModel_->correct(p, T);
171
172 correctMu();
173}
174
175
176// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
177
178} // End namespace surfaceFilmModels
179} // End namespace regionModels
180} // End namespace Foam
181
182// ************************************************************************* //
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:191
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Template invariant parts for fvPatchField.
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
const fvMesh & regionMesh() const
Return the region mesh database.
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.
static autoPtr< filmViscosityModel > New(surfaceFilmRegionModel &film, const dictionary &dict, volScalarField &mu)
Return a reference to the selected phase change model.
virtual void correct(const volScalarField &p, const volScalarField &T)=0
Correct.
Kinematic form of single-cell layer surface film model.
virtual void correct(const volScalarField &p, const volScalarField &T)
Correct.
autoPtr< filmViscosityModel > muSolventModel_
Solvent viscosity model.
autoPtr< filmViscosityModel > muWaxModel_
Wax viscosity model.
const dictionary coeffDict_
Coefficients dictionary.
const dictionary & dict() const
Return const access to the cloud dictionary.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
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.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
UniformDimensionedField< scalar > uniformDimensionedScalarField
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
const dimensionSet dimDynamicViscosity
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict