lambdaThixotropic.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2020-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "lambdaThixotropic.H"
27 #include "fvmDiv.H"
28 #include "fvmSup.H"
29 #include "fvModels.H"
30 #include "fvConstraints.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace laminarModels
38 {
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
42 template<class BasicMomentumTransportModel>
44 (
45  const alphaField& alpha,
46  const rhoField& rho,
47  const volVectorField& U,
48  const surfaceScalarField& alphaRhoPhi,
49  const surfaceScalarField& phi,
50  const viscosity& viscosity
51 )
52 :
53  linearViscousStress<laminarModel<BasicMomentumTransportModel>>
54  (
55  typeName,
56  alpha,
57  rho,
58  U,
59  alphaRhoPhi,
60  phi,
61  viscosity
62  ),
63 
64  a_("a", dimless/dimTime, this->coeffDict_),
65  b_("b", dimless, this->coeffDict_),
66  d_("d", dimless, this->coeffDict_),
67  c_("c", pow(dimTime, d_.value() - scalar(1)), this->coeffDict_),
68  nu0_("nu0", dimKinematicViscosity, this->coeffDict_),
69  nuInf_("nuInf", dimKinematicViscosity, this->coeffDict_),
70  K_(1 - sqrt(nuInf_/nu0_)),
71  BinghamPlastic_(this->coeffDict_.found("sigmay")),
72  sigmay_
73  (
74  BinghamPlastic_
75  ? dimensionedScalar("sigmay", dimPressure/dimDensity, this->coeffDict_)
77  ),
78  residualAlpha_
79  (
80  dimensioned<scalar>::lookupOrAddToDict
81  (
82  "residualAlpha",
83  this->coeffDict_,
84  dimless,
85  1e-6
86  )
87  ),
88  lambda_
89  (
90  IOobject
91  (
92  this->groupName(typedName("lambda")),
93  this->runTime_.name(),
94  this->mesh_,
95  IOobject::MUST_READ,
96  IOobject::AUTO_WRITE
97  ),
98  this->mesh_
99  ),
100 
101  nu_
102  (
103  IOobject
104  (
105  this->groupName(typedName("nu")),
106  this->runTime_.name(),
107  this->mesh_,
108  IOobject::NO_READ,
109  IOobject::AUTO_WRITE
110  ),
111  calcNu(this->strainRate())
112  )
113 {}
114 
115 
116 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
117 
118 template<class BasicMomentumTransportModel>
121 (
122  const volScalarField& strainRate
123 ) const
124 {
126  (
127  nuInf_/(sqr(1 - K_*lambda_) + rootVSmall)
128  );
129 
130  // Add optional yield stress contribution to the viscosity
131  if (BinghamPlastic_)
132  {
133  dimensionedScalar sigmaySmall
134  (
135  "sigmaySmall",
136  sigmay_.dimensions(),
137  small
138  );
139 
140  // Limit the Bingham viscosity to 100x the thixotropic viscosity
141  // for numerical stability
142  dimensionedScalar nuMax_("nuMax", 100*nu0_);
143 
144  nu.ref() = min
145  (
146  sigmay_/(strainRate + 1e-4*(sigmay_ + sigmaySmall)/nu0_) + nu(),
147  nuMax_
148  );
149  }
150 
151  return nu;
152 }
153 
154 
155 template<class BasicMomentumTransportModel>
156 tmp<volScalarField>
157 lambdaThixotropic<BasicMomentumTransportModel>::strainRate() const
158 {
159  return sqrt(2.0)*mag(symm(fvc::grad(this->U())));
160 }
161 
162 
163 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
164 
165 template<class BasicMomentumTransportModel>
167 {
169  {
170  a_.read(this->coeffDict());
171  b_.read(this->coeffDict());
172  d_.read(this->coeffDict());
173 
174  c_ = dimensionedScalar
175  (
176  "c",
177  pow(dimTime, d_.value() - scalar(1)),
178  this->coeffDict_
179  );
180 
181  nu0_.read(this->coeffDict());
182  nuInf_.read(this->coeffDict());
183 
184  K_ = (1 - sqrt(nuInf_/nu0_));
185 
186  return true;
187  }
188  else
189  {
190  return false;
191  }
192 }
193 
194 
195 template<class BasicMomentumTransportModel>
198 {
199  return volScalarField::New
200  (
201  this->groupName("nuEff"),
202  nu_
203  );
204 }
205 
206 
207 template<class BasicMomentumTransportModel>
210 (
211  const label patchi
212 ) const
213 {
214  return nu_.boundaryField()[patchi];
215 }
216 
217 
218 template<class BasicMomentumTransportModel>
220 {
221  // Local references
222  const alphaField& alpha = this->alpha_;
223  const rhoField& rho = this->rho_;
224  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
225  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
227  (
228  Foam::fvConstraints::New(this->mesh_)
229  );
230 
231  const volScalarField strainRate(this->strainRate());
232 
233  tmp<fvScalarMatrix> lambdaEqn
234  (
235  fvm::ddt(alpha, rho, lambda_)
236  + fvm::div(alphaRhoPhi, lambda_)
237  ==
238  alpha()*rho()*a_*pow(1 - lambda_(), b_)
239  - fvm::Sp
240  (
241  rho()
242  *(
243  alpha()*c_*pow(strainRate(), d_)
244  + max(residualAlpha_ - alpha(), dimensionedScalar(dimless, 0))*a_
245  ),
246  lambda_
247  )
248  + fvModels.source(alpha, rho, lambda_)
249  );
250 
251  lambdaEqn.ref().relax();
252  fvConstraints.constrain(lambdaEqn.ref());
253  solve(lambdaEqn);
254  fvConstraints.constrain(lambda_);
255 
256  lambda_.maxMin(scalar(0), scalar(1));
257 
258  nu_ = calcNu(strainRate);
259 
261 }
262 
263 
264 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
265 
266 } // End namespace laminarModels
267 } // End namespace Foam
268 
269 // ************************************************************************* //
bool found
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:100
Finite volume constraints.
Definition: fvConstraints.H:67
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
Finite volume models.
Definition: fvModels.H:65
tmp< fvMatrix< Type > > source(const VolField< Type > &field) const
Return source for an equation.
Templated abstract base class for laminar transport models.
Definition: laminarModel.H:52
virtual void correct()
Predict the laminar viscosity.
Definition: laminarModel.C:267
Thixotropic viscosity momentum transport model based on the evolution of the structural parameter :
BasicMomentumTransportModel::alphaField alphaField
virtual void correct()
Correct the lambdaThixotropic viscosity.
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
lambdaThixotropic(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
Construct from components.
virtual bool read()
Read momentumTransport dictionary.
BasicMomentumTransportModel::rhoField rhoField
Linear viscous stress turbulence model base class.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
Foam::fvConstraints & fvConstraints(Foam::fvConstraints::New(mesh))
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
Calculate the matrix for the divergence of the given field and flux.
Calculate the matrix for implicit and explicit sources.
label patchi
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
Namespace for OpenFOAM.
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
const doubleScalar e
Definition: doubleScalar.H:106
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const dimensionSet dimPressure
dimensionedSymmTensor sqr(const dimensionedVector &dv)
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimKinematicViscosity
const dimensionSet dimless
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const dimensionSet dimTime
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionSet dimDensity
dimensioned< scalar > mag(const dimensioned< Type > &)
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:176
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.