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 -------------------------------------------------------------------------------
10 License
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"
29 #include "kinematicSingleLayer.H"
30 #include "waxSolventEvaporation.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace regionModels
38 {
39 namespace surfaceFilmModels
40 {
41 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
44 defineTypeNameAndDebug(waxSolventViscosity, 0);
45 
47 (
48  filmViscosityModel,
49  waxSolventViscosity,
50  dictionary
51 );
52 
53 
54 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
55 
56 void waxSolventViscosity::correctMu()
57 {
58  const kinematicSingleLayer& film = filmType<kinematicSingleLayer>();
59 
61  (
63  (
64  waxSolventEvaporation::typeName + ":Wwax"
65  )
66  );
67 
68  const uniformDimensionedScalarField Wsolvent
69  (
71  (
72  waxSolventEvaporation::typeName + ":Wsolvent"
73  )
74  );
75 
76  const uniformDimensionedScalarField Ysolvent0
77  (
79  (
80  waxSolventEvaporation::typeName + ":Ysolvent0"
81  )
82  );
83 
84  const volScalarField& Ysolvent
85  (
87  (
88  waxSolventEvaporation::typeName + ":Ysolvent"
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 
109 waxSolventViscosity::waxSolventViscosity
110 (
112  const dictionary& dict,
114 )
115 :
116  filmViscosityModel(typeName, film, dict, mu),
117  muWax_
118  (
119  IOobject
120  (
121  typeName + ":muWax",
122  film.regionMesh().time().timeName(),
123  film.regionMesh().thisDb(),
124  IOobject::NO_READ,
125  IOobject::AUTO_WRITE
126  ),
127  film.regionMesh(),
130  ),
131  muWaxModel_
132  (
134  (
135  film,
136  coeffDict_.subDict("muWax"),
137  muWax_
138  )
139  ),
140  muSolvent_
141  (
142  IOobject
143  (
144  typeName + ":muSolvent",
145  film.regionMesh().time().timeName(),
146  film.regionMesh().thisDb(),
147  IOobject::NO_READ,
148  IOobject::AUTO_WRITE
149  ),
150  film.regionMesh(),
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 // ************************************************************************* //
dictionary dict
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 word zeroGradientType
A zeroGradient patch field type.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
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.
Template invariant parts for fvPatchField.
Definition: fvPatchField.H:77
UniformDimensionedField< scalar > uniformDimensionedScalarField
const dimensionSet dimDynamicViscosity
addToRunTimeSelectionTable(surfaceFilmRegionModel, kinematicSingleLayer, mesh)
Macros for easy insertion into run-time selection tables.
const surfaceFilmRegionModel & film() const
Return const access to the film surface film model.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
word timeName
Definition: getTimeIndex.H:3
autoPtr< filmViscosityModel > muWaxModel_
Wax viscosity model.
Base class for surface film viscosity models.
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:26
autoPtr< filmViscosityModel > muSolventModel_
Solvent viscosity model.
volScalarField & mu_
Reference to the viscosity field.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const dimensionedScalar mu
Atomic mass unit.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual void correct(const volScalarField &p, const volScalarField &T)
Correct.
void correctBoundaryConditions()
Correct boundary field.
volScalarField & p
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
defineTypeNameAndDebug(kinematicSingleLayer, 0)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127