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,
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().thisDb(),
82  IOobject::MUST_READ,
83  IOobject::AUTO_WRITE
84  ),
85  film.regionMesh()
86  )
87 {
89 
90  // Initialise viscosity to inf value because it cannot be evaluated yet
93 }
94 
95 
96 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
97 
99 {}
100 
101 
102 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
103 
105 (
106  const volScalarField& p,
107  const volScalarField& T
108 )
109 {
110  const kinematicSingleLayer& film = filmType<kinematicSingleLayer>();
111 
112  const volVectorField& U = film.U();
113  const volVectorField& Uw = film.Uw();
114  const volScalarField& delta = film.delta();
115  const volScalarField& deltaRho = film.deltaRho();
116  const surfaceScalarField& phi = film.phi();
117  const volScalarField& alpha = film.alpha();
118  const Time& runTime = this->film().regionMesh().time();
119 
120  // Shear rate
121  const volScalarField gDot
122  (
123  "gDot",
124  alpha*mag(U - Uw)/(delta + film.deltaSmall())
125  );
126 
127  if (debug && runTime.writeTime())
128  {
129  gDot.write();
130  }
131 
132  const dimensionedScalar deltaRho0
133  (
134  "deltaRho0",
135  deltaRho.dimensions(),
136  ROOTVSMALL
137  );
138 
139  const surfaceScalarField phiU(phi/fvc::interpolate(deltaRho + deltaRho0));
140 
141  const dimensionedScalar c0("c0", dimless/dimTime, ROOTVSMALL);
142  const volScalarField coeff("coeff", -c_*pow(gDot, d_) + c0);
143 
144  fvScalarMatrix lambdaEqn
145  (
147  + fvm::div(phiU, lambda_)
148  - fvm::Sp(fvc::div(phiU), lambda_)
149  ==
150  a_*pow((1 - lambda_), b_)
151  + fvm::SuSp(coeff, lambda_)
152 
153  // Include the effect of the impinging droplets added with lambda = 0
154  - fvm::Sp
155  (
156  max
157  (
158  -film.rhoSp(),
159  dimensionedScalar(film.rhoSp().dimensions(), Zero)
160  )/(deltaRho + deltaRho0),
161  lambda_
162  )
163  );
164 
165  lambdaEqn.relax();
166  lambdaEqn.solve();
167 
169 
170  mu_ = muInf_/(sqr(1 - K_*lambda_) + ROOTVSMALL);
172 }
173 
174 
175 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
176 
177 } // End namespace surfaceFilmModels
178 } // End namespace regionModels
179 } // End namespace Foam
180 
181 // ************************************************************************* //
scalar delta
void clamp_range(const dimensioned< MinMax< Type >> &range)
Clamp field values (in-place) to the specified range.
dictionary dict
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Kinematic form of single-cell layer surface film model.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
virtual const volScalarField & alpha() const =0
Return the film coverage, 1 = covered, 0 = uncovered / [].
engineTime & runTime
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionSet dimless
Dimensionless.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
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.
virtual const volVectorField & Uw() const =0
Return the film wall velocity [m/s].
word timeName
Definition: getTimeIndex.H:3
Base class for surface film viscosity models.
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:41
Calculate the matrix for the first temporal derivative.
virtual const volScalarField & delta() const =0
Return the film thickness [m].
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:26
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
virtual void correct(const volScalarField &p, const volScalarField &T)
Correct.
const dimensionSet dimPressure
virtual const volVectorField & U() const =0
Return the film velocity [m/s].
dimensionedScalar muInf_
Limiting viscosity when lambda = 0.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Calculate the divergence of the given field.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
int debug
Static debugging option.
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1101
volScalarField & mu_
Reference to the viscosity field.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
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.
Represents 0/1 range or concept. Used for tagged dispatch or clamping.
Definition: pTraits.H:65
const dimensionedScalar mu
Atomic mass unit.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Calculate the matrix for the divergence of the given field and flux.
U
Definition: pEqn.H:72
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
void correctBoundaryConditions()
Correct boundary field.
dimensionedScalar c_
Model `c&#39; coefficient (read after d since dims depend on d value)
volScalarField & p
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
defineTypeNameAndDebug(kinematicSingleLayer, 0)
Calculate the finiteVolume matrix for implicit and explicit sources.
Namespace for OpenFOAM.
const dimensionSet & dimensions() const noexcept
Return dimensions.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127