plastic.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-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 "plastic.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace mixtureViscosityModels
35 {
37 
39  (
41  plastic,
43  );
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
51 (
53 )
54 :
56  plasticViscosityCoeff_
57  (
58  "coeff",
60  coeffDict()
61  ),
62  plasticViscosityExponent_
63  (
64  "exponent",
65  dimless,
66  coeffDict()
67  ),
68  muMax_
69  (
70  "muMax",
72  coeffDict()
73  )
74 {}
75 
76 
77 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
78 
81 (
82  const volScalarField& muc,
83  const volVectorField& U
84 ) const
85 {
86  return min
87  (
88  muc
89  + plasticViscosityCoeff_
90  *(
91  pow
92  (
93  scalar(10),
94  plasticViscosityExponent_*mixture_.alphad()
95  ) - scalar(1)
96  ),
97  muMax_
98  );
99 }
100 
101 
103 {
105  {
106  const dictionary& plasticCoeffs = coeffDict();
107 
108  plasticViscosityCoeff_.read(plasticCoeffs);
109  plasticViscosityExponent_.read(plasticCoeffs);
110  muMax_.read(plasticCoeffs);
111 
112  return true;
113  }
114  else
115  {
116  return false;
117  }
118 }
119 
120 
121 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Class to represent a mixture of two constant density phases.
An abstract base class for incompressible mixtureViscosityModels.
Viscosity correction model for a generic power-law plastic.
Definition: plastic.H:57
virtual tmp< volScalarField > mu(const volScalarField &muc, const volVectorField &U) const
Return the mixture viscosity.
Definition: plastic.C:81
plastic(const incompressibleDriftFluxMixture &mixture)
Construct from mixture.
Definition: plastic.C:51
virtual bool read()
Read phaseProperties dictionary.
Definition: plastic.C:102
A class for managing temporary objects.
Definition: tmp.H:55
virtual bool read()=0
Read physicalProperties dictionary.
U
Definition: pEqn.H:72
defineTypeNameAndDebug(BinghamPlastic, 0)
addToRunTimeSelectionTable(mixtureViscosityModel, BinghamPlastic, dictionary)
Namespace for OpenFOAM.
const dimensionSet dimDynamicViscosity
const dimensionSet dimless
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)