BinghamPlastic.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) 2014-2015 OpenFOAM Foundation
9  Copyright (C) 2021 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 "BinghamPlastic.H"
31 #include "fvcGrad.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace mixtureViscosityModels
38 {
39  defineTypeNameAndDebug(BinghamPlastic, 0);
40 
42  (
43  mixtureViscosityModel,
44  BinghamPlastic,
45  dictionary
46  );
47 }
48 }
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
54 (
55  const word& name,
56  const dictionary& viscosityProperties,
57  const volVectorField& U,
58  const surfaceScalarField& phi
59 )
60 :
61  plastic(name, viscosityProperties, U, phi, typeName),
62  yieldStressCoeff_
63  (
64  "BinghamCoeff",
65  dimensionSet(1, -1, -2, 0, 0),
66  plasticCoeffs_
67  ),
68  yieldStressExponent_
69  (
70  "BinghamExponent",
71  dimless,
72  plasticCoeffs_
73  ),
74  yieldStressOffset_
75  (
76  "BinghamOffset",
77  dimless,
78  plasticCoeffs_
79  ),
80  U_(U)
81 {}
82 
83 
84 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
85 
88 (
89  const volScalarField& muc
90 ) const
91 {
92  volScalarField tauy
93  (
94  yieldStressCoeff_
95  *(
96  pow
97  (
98  scalar(10),
99  yieldStressExponent_
100  *(max(alpha_, scalar(0)) + yieldStressOffset_)
101  )
102  - pow
103  (
104  scalar(10),
105  yieldStressExponent_*yieldStressOffset_
106  )
107  )
108  );
109 
110  volScalarField mup(plastic::mu(muc));
111 
112  dimensionedScalar tauySmall("tauySmall", tauy.dimensions(), SMALL);
113 
114  return min
115  (
116  tauy
117  /(
118  sqrt(2.0)*mag(symm(fvc::grad(U_)))
119  + 1.0e-4*(tauy + tauySmall)/mup
120  )
121  + mup,
122  muMax_
123  );
124 }
125 
126 
128 (
129  const dictionary& viscosityProperties
130 )
131 {
132  plastic::read(viscosityProperties);
133 
134  plasticCoeffs_.readEntry("BinghamCoeff", yieldStressCoeff_);
135  plasticCoeffs_.readEntry("BinghamExponent", yieldStressExponent_);
136  plasticCoeffs_.readEntry("BinghamOffset", yieldStressOffset_);
137 
138  return true;
139 }
140 
141 
142 // ************************************************************************* //
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionSet dimless
Dimensionless.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:76
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
Calculate the gradient of the given field.
tmp< volScalarField > mu(const volScalarField &muc) const
Return the mixture viscosity.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
BinghamPlastic(const word &name, const dictionary &viscosityProperties, const volVectorField &U, const surfaceScalarField &phi)
Construct from components.
void read(Istream &, label &val, const dictionary &)
In-place read with dictionary lookup.
defineTypeNameAndDebug(combustionModel, 0)
const dimensionedScalar mu
Atomic mass unit.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
U
Definition: pEqn.H:72
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
bool read(const dictionary &viscosityProperties)
Read transportProperties dictionary.
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)