plastic.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-2017 OpenFOAM Foundation
9  Copyright (C) 2020 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 "plastic.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace mixtureViscosityModels
37 {
38  defineTypeNameAndDebug(plastic, 0);
39 
41  (
42  mixtureViscosityModel,
43  plastic,
44  dictionary
45  );
46 }
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
53 (
54  const word& name,
55  const dictionary& viscosityProperties,
56  const volVectorField& U,
57  const surfaceScalarField& phi,
58  const word modelName
59 )
60 :
61  mixtureViscosityModel(name, viscosityProperties, U, phi),
62  plasticCoeffs_(viscosityProperties.optionalSubDict(modelName + "Coeffs")),
63  plasticViscosityCoeff_("coeff", dimDynamicViscosity, plasticCoeffs_),
64  plasticViscosityExponent_("exponent", dimless, plasticCoeffs_),
65  muMax_("muMax", dimDynamicViscosity, plasticCoeffs_),
66  alpha_
67  (
68  U.mesh().lookupObject<volScalarField>
69  (
70  IOobject::groupName
71  (
72  viscosityProperties.getOrDefault<word>("alpha", "alpha"),
73  viscosityProperties.dictName()
74  )
75  )
76  )
77 {}
78 
79 
80 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
81 
84 {
85  return min
86  (
87  muc
88  + plasticViscosityCoeff_
89  *(
90  pow
91  (
92  scalar(10),
93  plasticViscosityExponent_*alpha_
94  ) - scalar(1)
95  ),
96  muMax_
97  );
98 }
99 
100 
102 (
103  const dictionary& viscosityProperties
104 )
105 {
106  mixtureViscosityModel::read(viscosityProperties);
107 
108  plasticCoeffs_ = viscosityProperties.optionalSubDict(typeName + "Coeffs");
109 
110  plasticCoeffs_.readEntry("k", plasticViscosityCoeff_);
111  plasticCoeffs_.readEntry("n", plasticViscosityExponent_);
112  plasticCoeffs_.readEntry("muMax", muMax_);
113 
114  return true;
115 }
116 
117 
118 // ************************************************************************* //
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.
const word dictName("faMeshDefinition")
const dimensionSet dimless
Dimensionless.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
const dimensionSet dimDynamicViscosity
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
virtual bool read(const dictionary &viscosityProperties)=0
Read transportProperties dictionary.
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
defineTypeNameAndDebug(combustionModel, 0)
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
U
Definition: pEqn.H:72
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.