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-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 
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 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
52 bool Foam::kineticTheoryModels::frictionalStressModels::Schaeffer::readCoeffs
53 (
54  const dictionary& coeffDict
55 )
56 {
57  phi_.read(coeffDict);
58  phi_ *= constant::mathematical::pi/180.0;
59 
60  return true;
61 }
62 
63 
64 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
65 
67 (
68  const dictionary& coeffDict
69 )
70 :
71  frictionalStressModel(coeffDict),
72  phi_("phi", dimless, coeffDict)
73 {
74  phi_ *= constant::mathematical::pi/180.0;
75 }
76 
77 
78 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
79 
81 {}
82 
83 
84 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
85 
89 (
90  const phaseModel& phase,
91  const dimensionedScalar& alphaMinFriction,
93 ) const
94 {
95  const volScalarField& alpha = phase;
96 
97  return
98  dimensionedScalar(dimensionSet(1, -1, -2, 0, 0), 1e24)
99  *pow(Foam::max(alpha - alphaMinFriction, scalar(0)), 10.0);
100 }
101 
102 
106 (
107  const phaseModel& phase,
108  const dimensionedScalar& alphaMinFriction,
110 ) const
111 {
112  const volScalarField& alpha = phase;
113 
114  return
115  dimensionedScalar(dimensionSet(1, -1, -2, 0, 0), 1e25)
116  *pow(Foam::max(alpha - alphaMinFriction, scalar(0)), 9.0);
117 }
118 
119 
122 (
123  const phaseModel& phase,
124  const dimensionedScalar& alphaMinFriction,
126  const volScalarField& pf,
127  const volSymmTensorField& D
128 ) const
129 {
130  const volScalarField& alpha = phase;
131 
133  (
135  (
137  (
138  Foam::typedName<frictionalStressModel>("nu"),
139  phase.group()
140  ),
141  phase.mesh(),
142  dimensionedScalar(dimensionSet(0, 2, -1, 0, 0), 0)
143  )
144  );
145 
146  volScalarField& nuf = tnu.ref();
147 
148  forAll(D, celli)
149  {
150  if (alpha[celli] > alphaMinFriction.value())
151  {
152  nuf[celli] =
153  0.5*pf[celli]*sin(phi_.value())
154  /(
155  sqrt((1.0/3.0)*sqr(tr(D[celli])) - invariantII(D[celli]))
156  + small
157  );
158  }
159  }
160 
161  const fvPatchList& patches = phase.mesh().boundary();
163 
165  {
166  if (!patches[patchi].coupled())
167  {
168  nufBf[patchi] =
169  (
170  pf.boundaryField()[patchi]*sin(phi_.value())
171  /(
172  mag(phase.U()().boundaryField()[patchi].snGrad())
173  + small
174  )
175  );
176  }
177  }
178 
179  // Correct coupled BCs
181 
182  return tnu;
183 }
184 
185 
186 // ************************************************************************* //
#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.
Schaeffer(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)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.