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  Copyright (C) 2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 "thixotropicViscosity.H"
30 #include "kinematicSingleLayer.H"
32 
33 #include "fvmDdt.H"
34 #include "fvmDiv.H"
35 #include "fvcDiv.H"
36 #include "fvmSup.H"
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42 namespace regionModels
43 {
44 namespace surfaceFilmModels
45 {
46 
47 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
48 
49 defineTypeNameAndDebug(thixotropicViscosity, 0);
50 
52 (
53  filmViscosityModel,
54  thixotropicViscosity,
56 );
57 
58 
59 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
60 
61 thixotropicViscosity::thixotropicViscosity
62 (
64  const dictionary& dict,
66 )
67 :
68  filmViscosityModel(typeName, film, dict, mu),
69  a_("a", dimless/dimTime, coeffDict_),
70  b_("b", dimless, coeffDict_),
71  d_("d", dimless, coeffDict_),
72  c_("c", pow(dimTime, d_.value() - scalar(1)), coeffDict_),
73  mu0_("mu0", dimPressure*dimTime, coeffDict_),
74  muInf_("muInf", mu0_.dimensions(), coeffDict_),
75  K_(1 - sqrt(muInf_/mu0_)),
76  lambda_
77  (
78  IOobject
79  (
80  IOobject::scopedName(typeName, "lambda"),
81  film.regionMesh().time().timeName(),
82  film.regionMesh().thisDb(),
83  IOobject::MUST_READ,
84  IOobject::AUTO_WRITE,
85  IOobject::REGISTER
86  ),
87  film.regionMesh()
88  )
89 {
91 
92  // Initialise viscosity to inf value because it cannot be evaluated yet
95 }
96 
97 
98 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
99 
101 {}
102 
103 
104 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
105 
107 (
108  const volScalarField& p,
109  const volScalarField& T
110 )
111 {
112  const kinematicSingleLayer& film = filmType<kinematicSingleLayer>();
113 
114  const volVectorField& U = film.U();
115  const volVectorField& Uw = film.Uw();
116  const volScalarField& delta = film.delta();
117  const volScalarField& deltaRho = film.deltaRho();
118  const surfaceScalarField& phi = film.phi();
119  const volScalarField& alpha = film.alpha();
120  const Time& runTime = this->film().regionMesh().time();
121 
122  // Shear rate
123  const volScalarField gDot
124  (
125  "gDot",
126  alpha*mag(U - Uw)/(delta + film.deltaSmall())
127  );
128 
129  if (debug && runTime.writeTime())
130  {
131  gDot.write();
132  }
133 
134  const dimensionedScalar deltaRho0
135  (
136  "deltaRho0",
137  deltaRho.dimensions(),
138  ROOTVSMALL
139  );
140 
141  const surfaceScalarField phiU(phi/fvc::interpolate(deltaRho + deltaRho0));
142 
143  const dimensionedScalar c0("c0", dimless/dimTime, ROOTVSMALL);
144  const volScalarField coeff("coeff", -c_*pow(gDot, d_) + c0);
145 
146  fvScalarMatrix lambdaEqn
147  (
149  + fvm::div(phiU, lambda_)
150  - fvm::Sp(fvc::div(phiU), lambda_)
151  ==
152  a_*pow((1 - lambda_), b_)
153  + fvm::SuSp(coeff, lambda_)
154 
155  // Include the effect of the impinging droplets added with lambda = 0
156  - fvm::Sp
157  (
158  max
159  (
160  -film.rhoSp(),
161  dimensionedScalar(film.rhoSp().dimensions(), Zero)
162  )/(deltaRho + deltaRho0),
163  lambda_
164  )
165  );
166 
167  lambdaEqn.relax();
168  lambdaEqn.solve();
169 
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
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:1095
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:180
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