thixotropicViscosity.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) 2013-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 "thixotropicViscosity.H"
29 #include "kinematicSingleLayer.H"
31 
32 #include "fvmDdt.H"
33 #include "fvmDiv.H"
34 #include "fvcDiv.H"
35 #include "fvmSup.H"
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41 namespace regionModels
42 {
43 namespace surfaceFilmModels
44 {
45 
46 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
47 
48 defineTypeNameAndDebug(thixotropicViscosity, 0);
49 
51 (
52  filmViscosityModel,
53  thixotropicViscosity,
54  dictionary
55 );
56 
57 
58 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
59 
60 thixotropicViscosity::thixotropicViscosity
61 (
63  const dictionary& dict,
65 )
66 :
67  filmViscosityModel(typeName, film, dict, mu),
68  a_("a", dimless/dimTime, coeffDict_),
69  b_("b", dimless, coeffDict_),
70  d_("d", dimless, coeffDict_),
71  c_("c", pow(dimTime, d_.value() - scalar(1)), coeffDict_),
72  mu0_("mu0", dimPressure*dimTime, coeffDict_),
73  muInf_("muInf", mu0_.dimensions(), coeffDict_),
74  K_(1 - sqrt(muInf_/mu0_)),
75  lambda_
76  (
77  IOobject
78  (
79  typeName + ":lambda",
80  film.regionMesh().time().timeName(),
81  film.regionMesh(),
84  ),
85  film.regionMesh()
86  )
87 {
88  lambda_.min(1);
89  lambda_.max(0);
90 
91  // Initialise viscosity to inf value because it cannot be evaluated yet
92  mu_ = muInf_;
93  mu_.correctBoundaryConditions();
94 }
95 
96 
97 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
98 
100 {}
101 
102 
103 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
104 
106 (
107  const volScalarField& p,
108  const volScalarField& T
109 )
110 {
111  const kinematicSingleLayer& film = filmType<kinematicSingleLayer>();
112 
113  const volVectorField& U = film.U();
114  const volVectorField& Uw = film.Uw();
115  const volScalarField& delta = film.delta();
116  const volScalarField& deltaRho = film.deltaRho();
117  const surfaceScalarField& phi = film.phi();
118  const volScalarField& alpha = film.alpha();
119  const Time& runTime = this->film().regionMesh().time();
120 
121  // Shear rate
122  const volScalarField gDot
123  (
124  "gDot",
125  alpha*mag(U - Uw)/(delta + film.deltaSmall())
126  );
127 
128  if (debug && runTime.writeTime())
129  {
130  gDot.write();
131  }
132 
133  const dimensionedScalar deltaRho0
134  (
135  "deltaRho0",
136  deltaRho.dimensions(),
137  ROOTVSMALL
138  );
139 
140  const surfaceScalarField phiU(phi/fvc::interpolate(deltaRho + deltaRho0));
141 
142  const dimensionedScalar c0("c0", dimless/dimTime, ROOTVSMALL);
143  const volScalarField coeff("coeff", -c_*pow(gDot, d_) + c0);
144 
145  fvScalarMatrix lambdaEqn
146  (
147  fvm::ddt(lambda_)
148  + fvm::div(phiU, lambda_)
149  - fvm::Sp(fvc::div(phiU), lambda_)
150  ==
151  a_*pow((1 - lambda_), b_)
152  + fvm::SuSp(coeff, lambda_)
153 
154  // Include the effect of the impinging droplets added with lambda = 0
155  - fvm::Sp
156  (
157  max
158  (
159  -film.rhoSp(),
160  dimensionedScalar(film.rhoSp().dimensions(), Zero)
161  )/(deltaRho + deltaRho0),
162  lambda_
163  )
164  );
165 
166  lambdaEqn.relax();
167  lambdaEqn.solve();
168 
169  lambda_.min(1);
170  lambda_.max(0);
171 
172  mu_ = muInf_/(sqr(1 - K_*lambda_) + ROOTVSMALL);
173  mu_.correctBoundaryConditions();
174 }
175 
176 
177 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
178 
179 } // End namespace surfaceFilmModels
180 } // End namespace regionModels
181 } // End namespace Foam
182 
183 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::dimPressure
const dimensionSet dimPressure
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:194
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::constant::physicoChemical::mu
const dimensionedScalar mu
Atomic mass unit.
Definition: createFieldRefs.H:4
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::regionModels::surfaceFilmModels::thixotropicViscosity::~thixotropicViscosity
virtual ~thixotropicViscosity()
Destructor.
Definition: thixotropicViscosity.C:99
thixotropicViscosity.H
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::U
virtual const volVectorField & U() const
Return the film velocity [m/s].
Definition: kinematicSingleLayer.C:952
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
fvcDiv.H
Calculate the divergence of the given field.
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer
Kinematic form of single-cell layer surface film model.
Definition: kinematicSingleLayer.H:67
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::rhoSp
volScalarField & rhoSp()
Mass [kg/m2/s].
Definition: kinematicSingleLayerI.H:128
Foam::regionModels::surfaceFilmModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(surfaceFilmRegionModel, kinematicSingleLayer, mesh)
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::deltaSmall
const dimensionedScalar & deltaSmall() const
Return small delta.
Definition: kinematicSingleLayerI.H:68
fvmDiv.H
Calculate the matrix for the divergence of the given field and flux.
Foam::regionModels::surfaceFilmModels::filmViscosityModel
Base class for surface film viscosity models.
Definition: filmViscosityModel.H:57
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::fvm::SuSp
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::regionModels::surfaceFilmModels::surfaceFilmRegionModel
Base class for surface film models.
Definition: surfaceFilmRegionModel.H:55
Foam::regionModels::surfaceFilmModels::thixotropicViscosity::correct
virtual void correct(const volScalarField &p, const volScalarField &T)
Correct.
Definition: thixotropicViscosity.C:106
Foam::fvMatrix::relax
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1183
Foam::regionModels::regionModel::regionMesh
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:64
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::phi
virtual const surfaceScalarField & phi() const
Return the film flux [kg.m/s].
Definition: kinematicSingleLayer.C:976
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
kinematicSingleLayer.H
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned< scalar >
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::fvc::interpolate
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.
fvmSup.H
Calculate the matrix for implicit and explicit sources.
Foam::regionModels::surfaceFilmModels::defineTypeNameAndDebug
defineTypeNameAndDebug(kinematicSingleLayer, 0)
U
U
Definition: pEqn.H:72
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::Uw
virtual const volVectorField & Uw() const
Return the film wall velocity [m/s].
Definition: kinematicSingleLayer.C:964
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::delta
const volScalarField & delta() const
Return const access to the film thickness [m].
Definition: kinematicSingleLayerI.H:86
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::alpha
const volScalarField & alpha() const
Return the film coverage, 1 = covered, 0 = uncovered [].
Definition: kinematicSingleLayerI.H:92
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:280
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::deltaRho
virtual const volScalarField & deltaRho() const
Return the film thickness*density (helper field) [kg/m3].
Definition: kinematicSingleLayer.C:970
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
fvmDdt.H
Calculate the matrix for the first temporal derivative.
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
Foam::IOobject::MUST_READ
Definition: IOobject.H:185