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-2023 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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
55 (
56  const dictionary& dict
57 )
58 :
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,
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,
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 
120 (
121  const phaseModel& phase,
122  const dimensionedScalar& alphaMinFriction,
124  const volScalarField& pf,
125  const volSymmTensorField& D
126 ) const
127 {
128  const volScalarField& alpha = phase;
129 
131  (
133  (
135  (
136  Foam::typedName<frictionalStressModel>("nu"),
137  phase.group()
138  ),
139  phase.mesh(),
140  dimensionedScalar(dimensionSet(0, 2, -1, 0, 0), 0)
141  )
142  );
143 
144  volScalarField& nuf = tnu.ref();
145 
146  forAll(D, celli)
147  {
148  if (alpha[celli] > alphaMinFriction.value())
149  {
150  nuf[celli] =
151  0.5*pf[celli]*sin(phi_.value())
152  /(
153  sqrt((1.0/3.0)*sqr(tr(D[celli])) - invariantII(D[celli]))
154  + small
155  );
156  }
157  }
158 
159  const fvPatchList& patches = phase.mesh().boundary();
161 
163  {
164  if (!patches[patchi].coupled())
165  {
166  nufBf[patchi] =
167  (
168  pf.boundaryField()[patchi]*sin(phi_.value())
169  /(
170  mag(phase.U()().boundaryField()[patchi].snGrad())
171  + small
172  )
173  );
174  }
175  }
176 
177  // Correct coupled BCs
179 
180  return tnu;
181 }
182 
183 
186 {
187  coeffDict_ <<= dict_.optionalSubDict(typeName + "Coeffs");
188 
189  Fr_.read(coeffDict_);
190  eta_.read(coeffDict_);
191  p_.read(coeffDict_);
192 
193  phi_.read(coeffDict_);
194  phi_ *= constant::mathematical::pi/180.0;
195 
196  alphaDeltaMin_.read(coeffDict_);
197 
198  return true;
199 }
200 
201 
202 // ************************************************************************* //
static const Foam::dimensionedScalar D("D", Foam::dimTemperature, 257.14)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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.
void correctBoundaryConditions()
Correct boundary field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, 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:134
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 keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Dimension set for the base types.
Definition: dimensionSet.H:125
const Type & value() const
Return const reference to value.
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:181
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)
Namespace for OpenFOAM.
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimless
dimensionedScalar sin(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
Cmpt invariantII(const SymmTensor< Cmpt > &st)
Return the 2nd invariant of a symmetric tensor.
Definition: SymmTensorI.H:403
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dictionary dict