JohnsonJacksonSchaefferFrictionalStress.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2016-2017 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
31 #include "unitConversion.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace kineticTheoryModels
38 {
39 namespace frictionalStressModels
40 {
41  defineTypeNameAndDebug(JohnsonJacksonSchaeffer, 0);
42 
44  (
45  frictionalStressModel,
46  JohnsonJacksonSchaeffer,
47  dictionary
48  );
49 }
50 }
51 }
52 
53 
54 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55 
58 (
59  const dictionary& dict
60 )
61 :
62  frictionalStressModel(dict),
63  coeffDict_(dict.optionalSubDict(typeName + "Coeffs")),
64  Fr_("Fr", dimensionSet(1, -1, -2, 0, 0), coeffDict_),
65  eta_("eta", dimless, coeffDict_),
66  p_("p", dimless, coeffDict_),
67  phi_("phi", dimless, coeffDict_),
68  alphaDeltaMin_("alphaDeltaMin", dimless, coeffDict_)
69 {
70  phi_ *= degToRad();
71 }
72 
73 
74 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
75 
78 {}
79 
80 
81 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
82 
86 (
87  const phaseModel& phase,
88  const dimensionedScalar& alphaMinFriction,
90 ) const
91 {
92  const volScalarField& alpha = phase;
93 
94  return
95  Fr_*pow(max(alpha - alphaMinFriction, scalar(0)), eta_)
96  /pow(max(alphaMax - alpha, alphaDeltaMin_), p_);
97 }
98 
99 
103 (
104  const phaseModel& phase,
105  const dimensionedScalar& alphaMinFriction,
107 ) const
108 {
109  const volScalarField& alpha = phase;
110 
111  return Fr_*
112  (
113  eta_*pow(max(alpha - alphaMinFriction, scalar(0)), eta_ - 1.0)
114  *(alphaMax-alpha)
115  + p_*pow(max(alpha - alphaMinFriction, scalar(0)), eta_)
116  )/pow(max(alphaMax - alpha, alphaDeltaMin_), p_ + 1.0);
117 }
118 
119 
123 (
124  const phaseModel& phase,
125  const dimensionedScalar& alphaMinFriction,
127  const volScalarField& pf,
128  const volSymmTensorField& D
129 ) const
130 {
131  const volScalarField& alpha = phase;
132 
133  tmp<volScalarField> tnu
134  (
135  new volScalarField
136  (
137  IOobject
138  (
139  "JohnsonJacksonSchaeffer:nu",
140  phase.mesh().time().timeName(),
141  phase.mesh(),
144  false
145  ),
146  phase.mesh(),
147  dimensionedScalar(dimensionSet(0, 2, -1, 0, 0), Zero)
148  )
149  );
150 
151  volScalarField& nuf = tnu.ref();
152 
153  forAll(D, celli)
154  {
155  if (alpha[celli] > alphaMinFriction.value())
156  {
157  nuf[celli] =
158  0.5*pf[celli]*sin(phi_.value())
159  /(
160  sqrt((1.0/3.0)*sqr(tr(D[celli])) - invariantII(D[celli]))
161  + SMALL
162  );
163  }
164  }
165 
166  const fvPatchList& patches = phase.mesh().boundary();
167  const volVectorField& U = phase.U();
168 
169  volScalarField::Boundary& nufBf = nuf.boundaryFieldRef();
170 
171  forAll(patches, patchi)
172  {
173  if (!patches[patchi].coupled())
174  {
175  nufBf[patchi] =
176  (
177  pf.boundaryField()[patchi]*sin(phi_.value())
178  /(
179  mag(U.boundaryField()[patchi].snGrad())
180  + SMALL
181  )
182  );
183  }
184  }
185 
186  // Correct coupled BCs
187  nuf.correctBoundaryConditions();
188 
189  return tnu;
190 }
191 
192 
195 {
196  coeffDict_ <<= dict_.optionalSubDict(typeName + "Coeffs");
197 
198  Fr_.read(coeffDict_);
199  eta_.read(coeffDict_);
200  p_.read(coeffDict_);
201 
202  phi_.read(coeffDict_);
203  phi_ *= degToRad();
204 
205  alphaDeltaMin_.read(coeffDict_);
206 
207  return true;
208 }
209 
210 
211 // ************************************************************************* //
Cmpt invariantII(const SymmTensor< Cmpt > &st)
Return the 2nd invariant of a SymmTensor.
Definition: SymmTensorI.H:375
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:88
const Type & value() const noexcept
Return const reference to value.
dictionary dict
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Unit conversion functions.
bool coupled(solutionDict.getOrDefault("coupledEnergyField", false))
dimensionedScalar sqrt(const dimensionedScalar &ds)
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
Type of boundary fields.
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:85
addToRunTimeSelectionTable(frictionalStressModel, JohnsonJackson, dictionary)
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:84
dimensionedScalar alphaMax("alphaMax", dimless/dimTime, laminarTransport)
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)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
U
Definition: pEqn.H:72
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const polyBoundaryMesh & patches
Nothing to be read.
virtual tmp< volScalarField > frictionalPressurePrime(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
virtual tmp< volScalarField > frictionalPressure(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const
const dimensionedScalar & D
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
Namespace for OpenFOAM.
PtrList< fvPatch > fvPatchList
Store lists of fvPatch as a PtrList.
Definition: fvPatch.H:57
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157