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-------------------------------------------------------------------------------
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
31
32#include "fvmDdt.H"
33#include "fvmDiv.H"
34#include "fvcDiv.H"
35#include "fvmSup.H"
36
37// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38
39namespace Foam
40{
41namespace regionModels
42{
43namespace surfaceFilmModels
44{
45
46// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
47
49
51(
55);
56
57
58// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
59
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 (
78 (
79 typeName + ":lambda",
80 film.regionMesh().time().timeName(),
81 film.regionMesh(),
82 IOobject::MUST_READ,
83 IOobject::AUTO_WRITE
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_;
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 (
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);
174}
175
176
177// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
178
179} // End namespace surfaceFilmModels
180} // End namespace regionModels
181} // End namespace Foam
182
183// ************************************************************************* //
scalar delta
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
surfaceScalarField & phi
const dimensionSet & dimensions() const
Return dimensions.
void max(const dimensioned< Type > &dt)
Use the maximum of the field and specified value.
void min(const dimensioned< Type > &dt)
Use the minimum of the field and specified value.
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.
bool writeTime() const noexcept
True if this is a write time.
Definition: TimeStateI.H:67
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:1092
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
virtual bool write(const bool valid=true) const
Write using setting from DB.
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.
virtual const volVectorField & U() const =0
Return the film velocity [m/s].
virtual const volScalarField & alpha() const =0
Return the film coverage, 1 = covered, 0 = uncovered / [].
virtual const volScalarField & delta() const =0
Return the film thickness [m].
virtual const volVectorField & Uw() const =0
Return the film wall velocity [m/s].
Thixotropic viscosity model based on the evolution of the structural parameter :
dimensionedScalar c_
Model `c' coefficient (read after d since dims depend on d value)
dimensionedScalar muInf_
Limiting viscosity when lambda = 0.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
U
Definition: pEqn.H:72
volScalarField & p
const volScalarField & T
const volScalarField & mu
engineTime & runTime
Calculate the divergence of the given field.
Calculate the matrix for the first temporal derivative.
Calculate the matrix for the divergence of the given field and flux.
Calculate the finiteVolume matrix for implicit and explicit sources.
word timeName
Definition: getTimeIndex.H:3
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.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Namespace for OpenFOAM.
const dimensionSet dimPressure
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
volScalarField & alpha
dictionary dict