JohnsonJacksonFrictionalStress.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) 2011-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 
28 #include "mathematicalConstants.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace kineticTheoryModels
35 {
36 namespace frictionalStressModels
37 {
39 
41  (
45  );
46 }
47 }
48 }
49 
50 
51 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52 
53 bool Foam::kineticTheoryModels::frictionalStressModels::JohnsonJackson::
54 readCoeffs
55 (
56  const dictionary& coeffDict
57 )
58 {
59  Fr_.read(coeffDict);
60  eta_.read(coeffDict);
61  p_.read(coeffDict);
62 
63  phi_.read(coeffDict);
64  phi_ *= constant::mathematical::pi/180.0;
65 
66  alphaDeltaMin_.read(coeffDict);
67 
68  return true;
69 }
70 
71 
72 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
73 
76 (
77  const dictionary& coeffDict
78 )
79 :
80  frictionalStressModel(coeffDict),
81  Fr_("Fr", dimensionSet(1, -1, -2, 0, 0), coeffDict),
82  eta_("eta", dimless, coeffDict),
83  p_("p", dimless, coeffDict),
84  phi_("phi", dimless, coeffDict),
85  alphaDeltaMin_("alphaDeltaMin", dimless, coeffDict)
86 {
87  phi_ *= constant::mathematical::pi/180.0;
88 }
89 
90 
91 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
92 
95 {}
96 
97 
98 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
99 
103 (
104  const phaseModel& phase,
105  const dimensionedScalar& alphaMinFriction,
107 ) const
108 {
109  const volScalarField& alpha = phase;
110 
111  return
112  Fr_*pow(max(alpha - alphaMinFriction, scalar(0)), eta_)
113  /pow(max(alphaMax - alpha, alphaDeltaMin_), p_);
114 }
115 
116 
120 (
121  const phaseModel& phase,
122  const dimensionedScalar& alphaMinFriction,
124 ) const
125 {
126  const volScalarField& alpha = phase;
127 
128  return Fr_*
129  (
130  eta_*pow(max(alpha - alphaMinFriction, scalar(0)), eta_ - 1)
131  *(alphaMax - alpha)
132  + p_*pow(max(alpha - alphaMinFriction, scalar(0)), eta_)
133  )/pow(max(alphaMax - alpha, alphaDeltaMin_), p_ + 1);
134 }
135 
136 
139 (
140  const phaseModel& phase,
141  const dimensionedScalar& alphaMinFriction,
143  const volScalarField& pf,
144  const volSymmTensorField& D
145 ) const
146 {
147  return volScalarField::New
148  (
150  (
151  Foam::typedName<frictionalStressModel>("nu"),
152  phase.group()
153  ),
154  dimensionedScalar(dimTime, 0.5)*pf*sin(phi_)
155  );
156 }
157 
158 
159 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
static word group(const word &name)
Return group (extension part of name)
Definition: IOobject.C:131
static word groupName(Name name, const word &group)
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Dimension set for the base types.
Definition: dimensionSet.H:125
void read(const dictionary &, const unitConversion &defaultUnits=NullObjectRef< unitConversion >())
Update the value of dimensioned<Type>
virtual tmp< volScalarField > frictionalPressure(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const
virtual tmp< volScalarField > nu(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax, const volScalarField &pf, const volSymmTensorField &D) const
virtual tmp< volScalarField > frictionalPressurePrime(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const
JohnsonJackson(const dictionary &coeffDict)
Construct from the coefficients dictionary.
A class for managing temporary objects.
Definition: tmp.H:55
dimensionedScalar alphaMax(viscosity->lookup("alphaMax"))
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
addToRunTimeSelectionTable(frictionalStressModel, JohnsonJackson, dictionary)
static const coefficient D("D", dimTemperature, 257.14)
Namespace for OpenFOAM.
const dimensionSet dimless
dimensionedScalar sin(const dimensionedScalar &ds)
const dimensionSet dimTime
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)