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-2018 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.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,
79  const dimensionedScalar& alphaMax
80 ) const
81 {
82  const volScalarField& alpha = phase;
83 
84  return
85  dimensionedScalar("1e24", 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,
96  const dimensionedScalar& alphaMax
97 ) const
98 {
99  const volScalarField& alpha = phase;
100 
101  return
102  dimensionedScalar("1e25", 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,
112  const dimensionedScalar& alphaMax,
113  const volScalarField& pf,
114  const volSymmTensorField& D
115 ) const
116 {
117  const volScalarField& alpha = phase;
118 
119  tmp<volScalarField> tnu
120  (
121  new volScalarField
122  (
123  IOobject
124  (
125  "Schaeffer:nu",
126  phase.mesh().time().timeName(),
127  phase.mesh(),
130  false
131  ),
132  phase.mesh(),
133  dimensionedScalar("nu", dimensionSet(0, 2, -1, 0, 0), 0.0)
134  )
135  );
136 
137  volScalarField& nuf = tnu.ref();
138 
139  forAll(D, celli)
140  {
141  if (alpha[celli] > alphaMinFriction.value())
142  {
143  nuf[celli] =
144  0.5*pf[celli]*sin(phi_.value())
145  /(
146  sqrt((1.0/3.0)*sqr(tr(D[celli])) - invariantII(D[celli]))
147  + small
148  );
149  }
150  }
151 
152  const fvPatchList& patches = phase.mesh().boundary();
153  const volVectorField& U = phase.U();
154 
155  volScalarField::Boundary& nufBf = nuf.boundaryFieldRef();
156 
157  forAll(patches, patchi)
158  {
159  if (!patches[patchi].coupled())
160  {
161  nufBf[patchi] =
162  (
163  pf.boundaryField()[patchi]*sin(phi_.value())
164  /(
165  mag(U.boundaryField()[patchi].snGrad())
166  + small
167  )
168  );
169  }
170  }
171 
172  // Correct coupled BCs
173  nuf.correctBoundaryConditions();
174 
175  return tnu;
176 }
177 
178 
180 {
181  coeffDict_ <<= dict_.optionalSubDict(typeName + "Coeffs");
182 
183  phi_.read(coeffDict_);
184  phi_ *= constant::mathematical::pi/180.0;
185 
186  return true;
187 }
188 
189 
190 // ************************************************************************* //
Cmpt invariantII(const SymmTensor< Cmpt > &st)
Return the 2nd invariant of a symmetric tensor.
Definition: SymmTensorI.H:403
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:58
virtual tmp< volScalarField > frictionalPressurePrime(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const
dictionary dict
Schaeffer(const dictionary &dict)
Construct from components.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
void read(const dictionary &)
Update the value of dimensioned<Type>
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:759
const Type & value() const
Return const reference to value.
virtual tmp< volScalarField > nu(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax, const volScalarField &pf, const volSymmTensorField &D) const
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
dimensionedScalar sin(const dimensionedScalar &ds)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
virtual tmp< volScalarField > frictionalPressure(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const
defineTypeNameAndDebug(combustionModel, 0)
Internal & ref()
Return a reference to the dimensioned internal field.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
PtrList< fvPatch > fvPatchList
container classes for fvPatch
Definition: fvPatchList.H:45
label patchi
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
const dictionary & dict_
Reference to higher-level dictionary for re-read.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.