All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2018 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"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace mixtureViscosityModels
34 {
35  defineTypeNameAndDebug(plastic, 0);
36 
38  (
39  mixtureViscosityModel,
40  plastic,
41  dictionary
42  );
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
50 (
51  const word& name,
52  const dictionary& viscosityProperties,
53  const volVectorField& U,
54  const surfaceScalarField& phi,
55  const word modelName
56 )
57 :
58  mixtureViscosityModel(name, viscosityProperties, U, phi),
59  plasticCoeffs_(viscosityProperties.optionalSubDict(modelName + "Coeffs")),
60  plasticViscosityCoeff_
61  (
62  "coeff",
63  dimensionSet(1, -1, -1, 0, 0),
64  plasticCoeffs_.lookup("coeff")
65  ),
66  plasticViscosityExponent_
67  (
68  "exponent",
69  dimless,
70  plasticCoeffs_.lookup("exponent")
71  ),
72  muMax_
73  (
74  "muMax",
75  dimensionSet(1, -1, -1, 0, 0),
76  plasticCoeffs_.lookup("muMax")
77  ),
78  alpha_
79  (
80  U.mesh().lookupObject<volScalarField>
81  (
82  IOobject::groupName
83  (
84  viscosityProperties.lookupOrDefault<word>("alpha", "alpha"),
85  viscosityProperties.dictName()
86  )
87  )
88  )
89 {}
90 
91 
92 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
93 
96 {
97  return min
98  (
99  muc
100  + plasticViscosityCoeff_
101  *(
102  pow
103  (
104  scalar(10),
105  plasticViscosityExponent_*alpha_
106  ) - scalar(1)
107  ),
108  muMax_
109  );
110 }
111 
112 
114 (
115  const dictionary& viscosityProperties
116 )
117 {
118  mixtureViscosityModel::read(viscosityProperties);
119 
120  plasticCoeffs_ = viscosityProperties.optionalSubDict(typeName + "Coeffs");
121 
122  plasticCoeffs_.lookup("k") >> plasticViscosityCoeff_;
123  plasticCoeffs_.lookup("n") >> plasticViscosityExponent_;
124  plasticCoeffs_.lookup("muMax") >> muMax_;
125 
126  return true;
127 }
128 
129 
130 // ************************************************************************* //
plastic(const word &name, const dictionary &viscosityProperties, const volVectorField &U, const surfaceScalarField &phi, const word modelName=typeName)
Construct from components.
tmp< volScalarField > mu(const volScalarField &muc) const
Return the mixture viscosity.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
Macros for easy insertion into run-time selection tables.
phi
Definition: pEqn.H:104
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
virtual bool read(const dictionary &viscosityProperties)=0
Read transportProperties dictionary.
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
const word dictName("particleTrackDict")
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
U
Definition: pEqn.H:72
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
bool read(const dictionary &viscosityProperties)
Read transportProperties dictionary.
Namespace for OpenFOAM.