Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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{
40{
42
44 (
45 frictionalStressModel,
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,
89 const dimensionedScalar& alphaMax
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,
106 const dimensionedScalar& alphaMax
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,
126 const dimensionedScalar& alphaMax,
127 const volScalarField& pf,
128 const volSymmTensorField& D
129) const
130{
131 const volScalarField& alpha = phase;
132
133 auto tnu = volScalarField::New
134 (
135 IOobject::scopedName("JohnsonJacksonSchaeffer", "nu"),
136 IOobject::NO_REGISTER,
137 phase.mesh(),
138 dimensionedScalar(dimensionSet(0, 2, -1, 0, 0), Zero)
139 );
140 auto& 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_ *= degToRad();
193
194 alphaDeltaMin_.read(coeffDict_);
195
196 return true;
197}
198
199
200// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
virtual tmp< volScalarField > frictionalPressure(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const
virtual tmp< volScalarField > nu(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax, const volScalarField &pf, const volSymmTensorField &D) const
virtual tmp< volScalarField > frictionalPressurePrime(const phaseModel &phase, const dimensionedScalar &alphaMinFriction, const dimensionedScalar &alphaMax) const
A class for managing temporary objects.
Definition tmp.H:75
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
const polyBoundaryMesh & patches
bool coupled
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:40
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
PtrList< fvPatch > fvPatchList
Store lists of fvPatch as a PtrList.
Definition fvPatch.H:64
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar sin(const dimensionedScalar &ds)
Cmpt invariantII(const SymmTensor2D< Cmpt > &st)
Return the 2nd invariant of a SymmTensor2D.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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:299
Unit conversion functions.