JohnsonJacksonSchaefferFrictionalStress.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) 2016-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 
52 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
53 
54 bool Foam::kineticTheoryModels::frictionalStressModels::
55 JohnsonJacksonSchaeffer::readCoeffs
56 (
57  const dictionary& coeffDict
58 )
59 {
60  Fr_.read(coeffDict);
61  eta_.read(coeffDict);
62  p_.read(coeffDict);
63 
64  phi_.read(coeffDict);
65  phi_ *= constant::mathematical::pi/180.0;
66 
67  alphaDeltaMin_.read(coeffDict);
68 
69  return true;
70 }
71 
72 
73 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
74 
77 (
78  const dictionary& coeffDict
79 )
80 :
81  frictionalStressModel(coeffDict),
82  Fr_("Fr", dimensionSet(1, -1, -2, 0, 0), coeffDict),
83  eta_("eta", dimless, coeffDict),
84  p_("p", dimless, coeffDict),
85  phi_("phi", dimless, coeffDict),
86  alphaDeltaMin_("alphaDeltaMin", dimless, coeffDict)
87 {
88  phi_ *= constant::mathematical::pi/180.0;
89 }
90 
91 
92 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
93 
96 {}
97 
98 
99 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
100 
104 (
105  const phaseModel& phase,
106  const dimensionedScalar& alphaMinFriction,
108 ) const
109 {
110  const volScalarField& alpha = phase;
111 
112  return
113  Fr_*pow(max(alpha - alphaMinFriction, scalar(0)), eta_)
114  /pow(max(alphaMax - alpha, alphaDeltaMin_), p_);
115 }
116 
117 
121 (
122  const phaseModel& phase,
123  const dimensionedScalar& alphaMinFriction,
125 ) const
126 {
127  const volScalarField& alpha = phase;
128 
129  return Fr_*
130  (
131  eta_*pow(max(alpha - alphaMinFriction, scalar(0)), eta_ - 1)
132  *(alphaMax - alpha)
133  + p_*pow(max(alpha - alphaMinFriction, scalar(0)), eta_)
134  )/pow(max(alphaMax - alpha, alphaDeltaMin_), p_ + 1);
135 }
136 
137 
141 (
142  const phaseModel& phase,
143  const dimensionedScalar& alphaMinFriction,
145  const volScalarField& pf,
146  const volSymmTensorField& D
147 ) const
148 {
149  const volScalarField& alpha = phase;
150 
152  (
154  (
156  (
157  Foam::typedName<frictionalStressModel>("nu"),
158  phase.group()
159  ),
160  phase.mesh(),
161  dimensionedScalar(dimensionSet(0, 2, -1, 0, 0), 0)
162  )
163  );
164 
165  volScalarField& nuf = tnu.ref();
166 
167  forAll(D, celli)
168  {
169  if (alpha[celli] > alphaMinFriction.value())
170  {
171  nuf[celli] =
172  0.5*pf[celli]*sin(phi_.value())
173  /(
174  sqrt((1.0/3.0)*sqr(tr(D[celli])) - invariantII(D[celli]))
175  + small
176  );
177  }
178  }
179 
180  const fvPatchList& patches = phase.mesh().boundary();
182 
184  {
185  if (!patches[patchi].coupled())
186  {
187  nufBf[patchi] =
188  (
189  pf.boundaryField()[patchi]*sin(phi_.value())
190  /(
191  mag(phase.U()().boundaryField()[patchi].snGrad())
192  + small
193  )
194  );
195  }
196  }
197 
198  // Correct coupled BCs
200 
201  return tnu;
202 }
203 
204 
205 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Macros for easy insertion into run-time selection tables.
const Mesh & mesh() const
Return mesh.
Generic GeometricBoundaryField class.
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
void correctBoundaryConditions()
Correct boundary field.
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 templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
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>
const Type & value() const
Return const reference to value.
JohnsonJacksonSchaeffer(const dictionary &coeffDict)
Construct from the coefficients dictionary.
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
virtual tmp< volVectorField > U() const =0
Return the velocity.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
label patchi
const fvPatchList & patches
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)
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
Cmpt invariantII(const SymmTensor< Cmpt > &st)
Return the 2nd invariant of a symmetric tensor.
Definition: SymmTensorI.H:421
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
void tr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< tensor > &f1)
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)