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-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 "JohnsonJacksonSchaefferFrictionalStress.H"
28 #include "mathematicalConstants.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace kineticTheoryModels
35 {
36 namespace frictionalStressModels
37 {
38  defineTypeNameAndDebug(JohnsonJacksonSchaeffer, 0);
39 
41  (
42  frictionalStressModel,
43  JohnsonJacksonSchaeffer,
44  dictionary
45  );
46 }
47 }
48 }
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
55 (
56  const dictionary& dict
57 )
58 :
59  frictionalStressModel(dict),
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,
86  const dimensionedScalar& alphaMax
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,
103  const dimensionedScalar& alphaMax
104 ) const
105 {
106  const volScalarField& alpha = phase;
107 
108  return Fr_*
109  (
110  eta_*pow(max(alpha - alphaMinFriction, scalar(0)), eta_ - 1.0)
111  *(alphaMax-alpha)
112  + p_*pow(max(alpha - alphaMinFriction, scalar(0)), eta_)
113  )/pow(max(alphaMax - alpha, alphaDeltaMin_), p_ + 1.0);
114 }
115 
116 
120 (
121  const phaseModel& phase,
122  const dimensionedScalar& alphaMinFriction,
123  const dimensionedScalar& alphaMax,
124  const volScalarField& pf,
125  const volSymmTensorField& D
126 ) const
127 {
128  const volScalarField& alpha = phase;
129 
130  tmp<volScalarField> tnu
131  (
133  (
134  "JohnsonJacksonSchaeffer:nu",
135  phase.mesh(),
136  dimensionedScalar(dimensionSet(0, 2, -1, 0, 0), 0)
137  )
138  );
139 
140  volScalarField& nuf = tnu.ref();
141 
142  forAll(D, celli)
143  {
144  if (alpha[celli] > alphaMinFriction.value())
145  {
146  nuf[celli] =
147  0.5*pf[celli]*sin(phi_.value())
148  /(
149  sqrt((1.0/3.0)*sqr(tr(D[celli])) - invariantII(D[celli]))
150  + small
151  );
152  }
153  }
154 
155  const fvPatchList& patches = phase.mesh().boundary();
156  const volVectorField& U = phase.U();
157 
158  volScalarField::Boundary& nufBf = nuf.boundaryFieldRef();
159 
160  forAll(patches, patchi)
161  {
162  if (!patches[patchi].coupled())
163  {
164  nufBf[patchi] =
165  (
166  pf.boundaryField()[patchi]*sin(phi_.value())
167  /(
168  mag(U.boundaryField()[patchi].snGrad())
169  + small
170  )
171  );
172  }
173  }
174 
175  // Correct coupled BCs
176  nuf.correctBoundaryConditions();
177 
178  return tnu;
179 }
180 
181 
184 {
185  coeffDict_ <<= dict_.optionalSubDict(typeName + "Coeffs");
186 
187  Fr_.read(coeffDict_);
188  eta_.read(coeffDict_);
189  p_.read(coeffDict_);
190 
191  phi_.read(coeffDict_);
192  phi_ *= constant::mathematical::pi/180.0;
193 
194  alphaDeltaMin_.read(coeffDict_);
195 
196  return true;
197 }
198 
199 
200 // ************************************************************************* //
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
dictionary dict
virtual tmp< volScalarField > nu(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax, const volScalarField &pf, const volSymmTensorField &D) const
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Mesh &, const dimensionSet &, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return a temporary field constructed from name, mesh, dimensionSet.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
virtual tmp< volScalarField > frictionalPressurePrime(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const
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:766
const Type & value() const
Return const reference to value.
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
dimensionedScalar sin(const dimensionedScalar &ds)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
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 > &)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
A class for managing temporary objects.
Definition: PtrList.H:53
JohnsonJacksonSchaeffer(const dictionary &dict)
Construct from components.
Namespace for OpenFOAM.
virtual tmp< volScalarField > frictionalPressure(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const