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-2018 OpenFOAM Foundation
9 Copyright (C) 2019 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
35namespace Foam
36{
37namespace kineticTheoryModels
38{
39namespace frictionalStressModels
40{
42
44 (
48 );
49}
50}
51}
52
53
54// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55
56Foam::kineticTheoryModels::frictionalStressModels::
57JohnsonJacksonSchaeffer::JohnsonJacksonSchaeffer
58(
59 const dictionary& dict
60)
61:
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)
114 *(alphaMax - alpha)
115 + p_*pow(max(alpha - alphaMinFriction, scalar(0)), eta_)
116 )/pow(max(alphaMax - alpha, alphaDeltaMin_), p_ + 1);
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
134 (
136 (
137 "JohnsonJacksonSchaeffer:nu",
138 phase.mesh(),
139 dimensionedScalar(dimensionSet(0, 2, -1, 0, 0))
140 )
141 );
142
143 volScalarField& nuf = tnu.ref();
144
145 forAll(D, celli)
146 {
147 if (alpha[celli] > alphaMinFriction.value())
148 {
149 nuf[celli] =
150 0.5*pf[celli]*sin(phi_.value())
151 /(
152 sqrt((1.0/3.0)*sqr(tr(D[celli])) - invariantII(D[celli]))
153 + SMALL
154 );
155 }
156 }
157
158 const fvPatchList& patches = phase.mesh().boundary();
159 const volVectorField& U = phase.U();
160
162
163 forAll(patches, patchi)
164 {
165 if (!patches[patchi].coupled())
166 {
167 nufBf[patchi] =
168 (
169 pf.boundaryField()[patchi]*sin(phi_.value())
170 /(
171 mag(U.boundaryField()[patchi].snGrad())
172 + SMALL
173 )
174 );
175 }
176 }
177
178 // Correct coupled BCs
180
181 return tnu;
182}
183
184
187{
188 coeffDict_ <<= dict_.optionalSubDict(typeName + "Coeffs");
189
190 Fr_.read(coeffDict_);
191 eta_.read(coeffDict_);
192 p_.read(coeffDict_);
193
194 phi_.read(coeffDict_);
195 phi_ *= degToRad();
196
197 alphaDeltaMin_.read(coeffDict_);
198
199 return true;
200}
201
202
203// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const Mesh & mesh() const
Return mesh.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
void correctBoundaryConditions()
Correct boundary field.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:109
const Type & value() const
Return const reference to value.
virtual tmp< volScalarField > frictionalPressure(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const
virtual tmp< volScalarField > frictionalPressurePrime(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:61
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:57
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
U
Definition: pEqn.H:72
const polyBoundaryMesh & patches
bool coupled(solutionDict.getOrDefault("coupledEnergyField", false))
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
const dimensionSet dimless
Dimensionless.
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Cmpt invariantII(const SymmTensor< Cmpt > &st)
Return the 2nd invariant of a SymmTensor.
Definition: SymmTensorI.H:468
dimensionedScalar sin(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
constexpr scalar degToRad() noexcept
Multiplication factor for degrees to radians conversion.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
volScalarField & nu
volScalarField & alpha
dictionary dict
const dimensionedScalar & D
dimensionedScalar alphaMax("alphaMax", dimless/dimTime, laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Unit conversion functions.