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-2022 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"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace mixtureViscosityModels
35 {
36  defineTypeNameAndDebug(BinghamPlastic, 0);
37 
39  (
40  mixtureViscosityModel,
41  BinghamPlastic,
42  dictionary
43  );
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
51 (
52  const incompressibleTwoPhaseInteractingMixture& mixture
53 )
54 :
55  plastic(mixture),
56  yieldStressCoeff_
57  (
58  "BinghamCoeff",
59  dimensionSet(1, -1, -2, 0, 0),
60  plasticCoeffs_
61  ),
62  yieldStressExponent_
63  (
64  "BinghamExponent",
65  dimless,
66  plasticCoeffs_
67  ),
68  yieldStressOffset_
69  (
70  "BinghamOffset",
71  dimless,
72  plasticCoeffs_
73  )
74 {}
75 
76 
77 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
78 
81 (
82  const volScalarField& muc,
83  const volVectorField& U
84 ) const
85 {
86  volScalarField tauy
87  (
88  yieldStressCoeff_
89  *(
90  pow
91  (
92  scalar(10),
93  min
94  (
95  log10(vGreat),
96  yieldStressExponent_
97  *(max(mixture_.alphad(), scalar(0)) + yieldStressOffset_)
98  )
99  )
100  - pow
101  (
102  scalar(10),
103  yieldStressExponent_*yieldStressOffset_
104  )
105  )
106  );
107 
108  volScalarField mup(plastic::mu(muc, U));
109 
110  dimensionedScalar tauySmall("tauySmall", tauy.dimensions(), small);
111 
112  return min
113  (
114  tauy
115  /(
116  sqrt(2.0)*mag(symm(fvc::grad(U)))
117  + 1.0e-4*(tauy + tauySmall)/mup
118  )
119  + mup,
120  muMax_
121  );
122 }
123 
124 
126 {
127  if (plastic::read())
128  {
129  plasticCoeffs_.lookup("yieldStressCoeff") >> yieldStressCoeff_;
130  plasticCoeffs_.lookup("yieldStressExponent") >> yieldStressExponent_;
131  plasticCoeffs_.lookup("yieldStressOffset") >> yieldStressOffset_;
132 
133  return true;
134  }
135  else
136  {
137  return false;
138  }
139 }
140 
141 
142 // ************************************************************************* //
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
virtual tmp< volScalarField > mu(const volScalarField &muc, const volVectorField &U) const
Return the mixture viscosity.
virtual bool read()
Read phaseProperties dictionary.
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionSet dimless
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:59
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
Calculate the gradient of the given field.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
const dimensionedScalar mu
Atomic mass unit.
defineTypeNameAndDebug(combustionModel, 0)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
BinghamPlastic(const incompressibleTwoPhaseInteractingMixture &mixture)
Construct from mixture.
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensioned< scalar > mag(const dimensioned< Type > &)
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
A class for managing temporary objects.
Definition: PtrList.H:53
dimensionedScalar log10(const dimensionedScalar &ds)
Namespace for OpenFOAM.