SchaefferFrictionalStress.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-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 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace kineticTheoryModels
34 {
35 namespace frictionalStressModels
36 {
38 
40  (
42  Schaeffer,
44  );
45 }
46 }
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
53 (
54  const dictionary& dict
55 )
56 :
58  coeffDict_(dict.optionalSubDict(typeName + "Coeffs")),
59  phi_("phi", dimless, coeffDict_)
60 {
61  phi_ *= constant::mathematical::pi/180.0;
62 }
63 
64 
65 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
66 
68 {}
69 
70 
71 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
72 
76 (
77  const phaseModel& phase,
78  const dimensionedScalar& alphaMinFriction,
80 ) const
81 {
82  const volScalarField& alpha = phase;
83 
84  return
85  dimensionedScalar(dimensionSet(1, -1, -2, 0, 0), 1e24)
86  *pow(Foam::max(alpha - alphaMinFriction, scalar(0)), 10.0);
87 }
88 
89 
93 (
94  const phaseModel& phase,
95  const dimensionedScalar& alphaMinFriction,
97 ) const
98 {
99  const volScalarField& alpha = phase;
100 
101  return
102  dimensionedScalar(dimensionSet(1, -1, -2, 0, 0), 1e25)
103  *pow(Foam::max(alpha - alphaMinFriction, scalar(0)), 9.0);
104 }
105 
106 
109 (
110  const phaseModel& phase,
111  const dimensionedScalar& alphaMinFriction,
113  const volScalarField& pf,
114  const volSymmTensorField& D
115 ) const
116 {
117  const volScalarField& alpha = phase;
118 
120  (
122  (
124  (
125  Foam::typedName<frictionalStressModel>("nu"),
126  phase.group()
127  ),
128  phase.mesh(),
129  dimensionedScalar(dimensionSet(0, 2, -1, 0, 0), 0)
130  )
131  );
132 
133  volScalarField& nuf = tnu.ref();
134 
135  forAll(D, celli)
136  {
137  if (alpha[celli] > alphaMinFriction.value())
138  {
139  nuf[celli] =
140  0.5*pf[celli]*sin(phi_.value())
141  /(
142  sqrt((1.0/3.0)*sqr(tr(D[celli])) - invariantII(D[celli]))
143  + small
144  );
145  }
146  }
147 
148  const fvPatchList& patches = phase.mesh().boundary();
150 
152  {
153  if (!patches[patchi].coupled())
154  {
155  nufBf[patchi] =
156  (
157  pf.boundaryField()[patchi]*sin(phi_.value())
158  /(
159  mag(phase.U()().boundaryField()[patchi].snGrad())
160  + small
161  )
162  );
163  }
164  }
165 
166  // Correct coupled BCs
168 
169  return tnu;
170 }
171 
172 
174 {
175  coeffDict_ <<= dict_.optionalSubDict(typeName + "Coeffs");
176 
177  phi_.read(coeffDict_);
178  phi_ *= constant::mathematical::pi/180.0;
179 
180  return true;
181 }
182 
183 
184 // ************************************************************************* //
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.
Schaeffer(const dictionary &dict)
Construct from components.
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)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict