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-2020 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 {
38  defineTypeNameAndDebug(JohnsonJackson, 0);
39 
41  (
42  frictionalStressModel,
43  JohnsonJackson,
44  dictionary
45  );
46 }
47 }
48 }
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
55 (
56  const dictionary& dict
57 )
58 :
59  frictionalStressModel(dict),
60  coeffDict_(dict.optionalSubDict(typeName + "Coeffs")),
61  Fr_("Fr", dimensionSet(1, -1, -2, 0, 0), coeffDict_),
62  eta_("eta", dimless, coeffDict_),
63  p_("p", dimless, coeffDict_),
64  phi_("phi", dimless, coeffDict_),
65  alphaDeltaMin_("alphaDeltaMin", dimless, coeffDict_)
66 {
67  phi_ *= constant::mathematical::pi/180.0;
68 }
69 
70 
71 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
72 
75 {}
76 
77 
78 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
79 
83 (
84  const phaseModel& phase,
85  const dimensionedScalar& alphaMinFriction,
86  const dimensionedScalar& alphaMax
87 ) const
88 {
89  const volScalarField& alpha = phase;
90 
91  return
92  Fr_*pow(max(alpha - alphaMinFriction, scalar(0)), eta_)
93  /pow(max(alphaMax - alpha, alphaDeltaMin_), p_);
94 }
95 
96 
100 (
101  const phaseModel& phase,
102  const dimensionedScalar& alphaMinFriction,
103  const dimensionedScalar& alphaMax
104 ) const
105 {
106  const volScalarField& alpha = phase;
107 
108  return Fr_*
109  (
110  eta_*pow(max(alpha - alphaMinFriction, scalar(0)), eta_ - 1)
111  *(alphaMax - alpha)
112  + p_*pow(max(alpha - alphaMinFriction, scalar(0)), eta_)
113  )/pow(max(alphaMax - alpha, alphaDeltaMin_), p_ + 1);
114 }
115 
116 
119 (
120  const phaseModel& phase,
121  const dimensionedScalar& alphaMinFriction,
122  const dimensionedScalar& alphaMax,
123  const volScalarField& pf,
124  const volSymmTensorField& D
125 ) const
126 {
127  return volScalarField::New
128  (
130  (
131  IOobject::modelName("nu", frictionalStressModel::typeName),
132  phase.group()
133  ),
134  dimensionedScalar(dimTime, 0.5)*pf*sin(phi_)
135  );
136 }
137 
138 
140 {
141  coeffDict_ <<= dict_.optionalSubDict(typeName + "Coeffs");
142 
143  Fr_.read(coeffDict_);
144  eta_.read(coeffDict_);
145  p_.read(coeffDict_);
146 
147  phi_.read(coeffDict_);
148  phi_ *= constant::mathematical::pi/180.0;
149 
150  alphaDeltaMin_.read(coeffDict_);
151 
152  return true;
153 }
154 
155 
156 // ************************************************************************* //
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:62
dictionary dict
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
virtual tmp< volScalarField > nu(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax, const volScalarField &pf, const volSymmTensorField &D) const
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
const dimensionSet dimless
bool read(Istream &, const bool keepHeader=false)
Read dictionary from Istream, optionally keeping the header.
Definition: dictionaryIO.C:104
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
const dimensionSet dimTime
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:1060
static word groupName(Name name, const word &group)
dimensionedScalar sin(const dimensionedScalar &ds)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
const dictionary & dict_
Reference to higher-level dictionary for re-read.
defineTypeNameAndDebug(combustionModel, 0)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
static word modelName(Name name, const word &model)
Return the name of the object within the given model.
A class for managing temporary objects.
Definition: PtrList.H:53
virtual tmp< volScalarField > frictionalPressurePrime(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const
JohnsonJackson(const dictionary &dict)
Construct from components.
Namespace for OpenFOAM.
virtual tmp< volScalarField > frictionalPressure(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const