BinghamPlastic.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) 2014-2023 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 "BinghamPlastic.H"
28 #include "fvcGrad.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace mixtureViscosityModels
36 {
38 
40  (
44  );
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
52 (
54 )
55 :
57  yieldStressCoeff_
58  (
59  "BinghamCoeff",
60  dimensionSet(1, -1, -2, 0, 0),
61  plasticCoeffs_
62  ),
63  yieldStressExponent_
64  (
65  "BinghamExponent",
66  dimless,
67  plasticCoeffs_
68  ),
69  yieldStressOffset_
70  (
71  "BinghamOffset",
72  dimless,
73  plasticCoeffs_
74  )
75 {}
76 
77 
78 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
79 
82 (
83  const volScalarField& muc,
84  const volVectorField& U
85 ) const
86 {
87  volScalarField tauy
88  (
89  yieldStressCoeff_
90  *(
91  pow
92  (
93  scalar(10),
94  min
95  (
96  log10(vGreat),
97  yieldStressExponent_
98  *(max(mixture_.alphad(), scalar(0)) + yieldStressOffset_)
99  )
100  )
101  - pow
102  (
103  scalar(10),
104  yieldStressExponent_*yieldStressOffset_
105  )
106  )
107  );
108 
109  volScalarField mup(plastic::mu(muc, U));
110 
111  dimensionedScalar tauySmall("tauySmall", tauy.dimensions(), small);
112 
113  return min
114  (
115  tauy
116  /(
117  sqrt(2.0)*mag(symm(fvc::grad(U)))
118  + 1.0e-4*(tauy + tauySmall)/mup
119  )
120  + mup,
121  muMax_
122  );
123 }
124 
125 
127 {
128  if (plastic::read())
129  {
130  plasticCoeffs_.lookup("yieldStressCoeff") >> yieldStressCoeff_;
131  plasticCoeffs_.lookup("yieldStressExponent") >> yieldStressExponent_;
132  plasticCoeffs_.lookup("yieldStressOffset") >> yieldStressOffset_;
133 
134  return true;
135  }
136  else
137  {
138  return false;
139  }
140 }
141 
142 
143 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
const dimensionSet & dimensions() const
Return dimensions.
Generic GeometricField class.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Dimension set for the base types.
Definition: dimensionSet.H:122
Class to represent a mixture of two constant density phases.
An abstract base class for incompressible mixtureViscosityModels.
Viscosity correction model for Bingham plastics.
virtual tmp< volScalarField > mu(const volScalarField &muc, const volVectorField &U) const
Return the mixture viscosity.
virtual bool read()
Read phaseProperties dictionary.
BinghamPlastic(const incompressibleDriftFluxMixture &mixture)
Construct from mixture.
Viscosity correction model for a generic power-law plastic.
Definition: plastic.H:58
virtual tmp< volScalarField > mu(const volScalarField &muc, const volVectorField &U) const
Return the mixture viscosity.
Definition: plastic.C:77
virtual bool read()
Read phaseProperties dictionary.
Definition: plastic.C:98
A class for managing temporary objects.
Definition: tmp.H:55
Calculate the gradient of the given field.
U
Definition: pEqn.H:72
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
defineTypeNameAndDebug(BinghamPlastic, 0)
addToRunTimeSelectionTable(mixtureViscosityModel, BinghamPlastic, dictionary)
Namespace for OpenFOAM.
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
const dimensionSet dimless
dimensionedScalar log10(const dimensionedScalar &ds)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)