SchaefferFrictionalStress.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 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 
26 #include "SchaefferFrictionalStress.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace kineticTheoryModels
34 {
35 namespace frictionalStressModels
36 {
37  defineTypeNameAndDebug(Schaeffer, 0);
38 
40  (
41  frictionalStressModel,
42  Schaeffer,
43  dictionary
44  );
45 }
46 }
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
53 (
54  const dictionary& dict
55 )
56 :
57  frictionalStressModel(dict),
58  coeffDict_(dict.subDict(typeName + "Coeffs")),
59  phi_("phi", dimless, coeffDict_)
60 {
62 }
63 
64 
65 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
66 
68 {}
69 
70 
71 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
72 
76 (
77  const volScalarField& alpha1,
78  const dimensionedScalar& alphaMinFriction,
80 ) const
81 {
82  return
83  dimensionedScalar("1e24", dimensionSet(1, -1, -2, 0, 0), 1e24)
84  *pow(Foam::max(alpha1 - alphaMinFriction, scalar(0)), 10.0);
85 }
86 
87 
91 (
92  const volScalarField& alpha1,
93  const dimensionedScalar& alphaMinFriction,
95 ) const
96 {
97  return
98  dimensionedScalar("1e25", dimensionSet(1, -1, -2, 0, 0), 1e25)
99  *pow(Foam::max(alpha1 - alphaMinFriction, scalar(0)), 9.0);
100 }
101 
102 
105 (
106  const volScalarField& alpha1,
107  const dimensionedScalar& alphaMinFriction,
109  const volScalarField& pf,
110  const volSymmTensorField& D
111 ) const
112 {
113  const scalar I2Dsmall = 1.0e-15;
114 
115  // Creating nu assuming it should be 0 on the boundary which may not be
116  // true
117  tmp<volScalarField> tnu
118  (
119  new volScalarField
120  (
121  IOobject
122  (
123  "Schaeffer:nu",
124  alpha1.mesh().time().timeName(),
125  alpha1.mesh(),
128  false
129  ),
130  alpha1.mesh(),
131  dimensionedScalar("nu", dimensionSet(0, 2, -1, 0, 0), 0.0)
132  )
133  );
134 
135  volScalarField& nuf = tnu();
136 
137  forAll(D, celli)
138  {
139  if (alpha1[celli] > alphaMinFriction.value())
140  {
141  nuf[celli] =
142  0.5*pf[celli]*sin(phi_.value())
143  /(
144  sqrt(1.0/6.0*(sqr(D[celli].xx() - D[celli].yy())
145  + sqr(D[celli].yy() - D[celli].zz())
146  + sqr(D[celli].zz() - D[celli].xx()))
147  + sqr(D[celli].xy()) + sqr(D[celli].xz())
148  + sqr(D[celli].yz())) + I2Dsmall
149  );
150  }
151  }
152 
153  // Correct coupled BCs
154  nuf.correctBoundaryConditions();
155 
156  return tnu;
157 }
158 
159 
161 {
162  coeffDict_ <<= dict_.subDict(typeName + "Coeffs");
163 
164  phi_.read(coeffDict_);
165  phi_ *= constant::mathematical::pi/180.0;
166 
167  return true;
168 }
169 
170 
171 // ************************************************************************* //
Foam::kineticTheoryModels::frictionalStressModels::JohnsonJackson::phi_
dimensionedScalar phi_
Angle of internal friction.
Definition: JohnsonJacksonFrictionalStress.H:69
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::kineticTheoryModels::frictionalStressModels::Schaeffer::nu
virtual tmp< volScalarField > nu(const volScalarField &alpha1, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax, const volScalarField &pf, const volSymmTensorField &D) const
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:255
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::volSymmTensorField
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:58
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::kineticTheoryModels::frictionalStressModels::Schaeffer::read
virtual bool read()
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
dict
dictionary dict
Definition: searchingEngine.H:14
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
alphaMax
dimensionedScalar alphaMax(laminarTransport.lookup("alphaMax"))
alpha1
volScalarField & alpha1
Definition: createFields.H:15
Foam::kineticTheoryModels::frictionalStressModels::Schaeffer::frictionalPressurePrime
virtual tmp< volScalarField > frictionalPressurePrime(const volScalarField &alpha1, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const
Foam::kineticTheoryModels::frictionalStressModels::Schaeffer::Schaeffer
Schaeffer(const dictionary &dict)
Construct from components.
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::constant::mathematical::pi
const scalar pi(M_PI)
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::kineticTheoryModels::frictionalStressModels::Schaeffer::~Schaeffer
virtual ~Schaeffer()
Destructor.
Foam::kineticTheoryModels::frictionalStressModels::Schaeffer::frictionalPressure
virtual tmp< volScalarField > frictionalPressure(const volScalarField &alpha1, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const