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_("residualAlpha", dimless, this->coeffDict(), 1e-6),
79  lambda_
80  (
81  IOobject
82  (
83  this->groupName(typedName("lambda")),
84  this->runTime_.name(),
85  this->mesh_,
86  IOobject::MUST_READ,
87  IOobject::AUTO_WRITE
88  ),
89  this->mesh_
90  ),
91 
92  nu_
93  (
94  IOobject
95  (
96  this->groupName(typedName("nu")),
97  this->runTime_.name(),
98  this->mesh_,
99  IOobject::NO_READ,
100  IOobject::AUTO_WRITE
101  ),
102  calcNu(this->strainRate())
103  )
104 {}
105 
106 
107 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
108 
109 template<class BasicMomentumTransportModel>
112 (
113  const volScalarField& strainRate
114 ) const
115 {
117  (
118  nuInf_/(sqr(1 - K_*lambda_) + rootVSmall)
119  );
120 
121  // Add optional yield stress contribution to the viscosity
122  if (BinghamPlastic_)
123  {
124  dimensionedScalar sigmaySmall
125  (
126  "sigmaySmall",
127  sigmay_.dimensions(),
128  small
129  );
130 
131  // Limit the Bingham viscosity to 100x the thixotropic viscosity
132  // for numerical stability
133  dimensionedScalar nuMax_("nuMax", 100*nu0_);
134 
135  nu.ref() = min
136  (
137  sigmay_/(strainRate + 1e-4*(sigmay_ + sigmaySmall)/nu0_) + nu(),
138  nuMax_
139  );
140  }
141 
142  return nu;
143 }
144 
145 
146 template<class BasicMomentumTransportModel>
147 tmp<volScalarField>
148 lambdaThixotropic<BasicMomentumTransportModel>::strainRate() const
149 {
150  return sqrt(2.0)*mag(symm(fvc::grad(this->U())));
151 }
152 
153 
154 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
155 
156 template<class BasicMomentumTransportModel>
158 {
160  {
161  a_.read(this->coeffDict());
162  b_.read(this->coeffDict());
163  d_.read(this->coeffDict());
164 
165  c_ = dimensionedScalar
166  (
167  "c",
168  pow(dimTime, d_.value() - scalar(1)),
169  this->coeffDict()
170  );
171 
172  nu0_.read(this->coeffDict());
173  nuInf_.read(this->coeffDict());
174 
175  K_ = (1 - sqrt(nuInf_/nu0_));
176 
177  return true;
178  }
179  else
180  {
181  return false;
182  }
183 }
184 
185 
186 template<class BasicMomentumTransportModel>
189 {
190  return volScalarField::New
191  (
192  this->groupName("nuEff"),
193  nu_
194  );
195 }
196 
197 
198 template<class BasicMomentumTransportModel>
201 (
202  const label patchi
203 ) const
204 {
205  return nu_.boundaryField()[patchi];
206 }
207 
208 
209 template<class BasicMomentumTransportModel>
211 {
212  // Local references
213  const alphaField& alpha = this->alpha_;
214  const rhoField& rho = this->rho_;
215  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
216  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
218  (
219  Foam::fvConstraints::New(this->mesh_)
220  );
221 
222  const volScalarField strainRate(this->strainRate());
223 
224  tmp<fvScalarMatrix> lambdaEqn
225  (
226  fvm::ddt(alpha, rho, lambda_)
227  + fvm::div(alphaRhoPhi, lambda_)
228  ==
229  alpha()*rho()*a_*pow(1 - lambda_(), b_)
230  - fvm::Sp
231  (
232  rho()
233  *(
234  alpha()*c_*pow(strainRate(), d_)
235  + max(residualAlpha_ - alpha(), dimensionedScalar(dimless, 0))*a_
236  ),
237  lambda_
238  )
239  + fvModels.source(alpha, rho, lambda_)
240  );
241 
242  lambdaEqn.ref().relax();
243  fvConstraints.constrain(lambdaEqn.ref());
244  solve(lambdaEqn);
245  fvConstraints.constrain(lambda_);
246 
247  lambda_.maxMin(scalar(0), scalar(1));
248 
249  nu_ = calcNu(strainRate);
250 
252 }
253 
254 
255 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
256 
257 } // End namespace laminarModels
258 } // End namespace Foam
259 
260 // ************************************************************************* //
bool found
Generic GeometricField class.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, 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:103
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:274
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:197
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.
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
const dimensionSet dimKinematicViscosity
const dimensionSet dimless
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
void symm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimDensity
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:181
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
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.