Loading...
Searching...
No Matches
SchaefferFrictionalStress.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) 2011-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 * * * * * * * * * * * * * //
35namespace Foam
37namespace kineticTheoryModels
38{
40{
42
44 (
48 );
50}
51}
52
53
54// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55
57(
58 const dictionary& dict
59)
60:
61 frictionalStressModel(dict),
62 coeffDict_(dict.optionalSubDict(typeName + "Coeffs")),
63 phi_("phi", dimless, coeffDict_)
65 phi_ *= degToRad();
66}
67
68
69// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
70
72{}
73
74
75// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
76
80(
81 const phaseModel& phase,
82 const dimensionedScalar& alphaMinFriction,
84) const
85{
86 const volScalarField& alpha = phase;
87
88 return
89 dimensionedScalar("", dimensionSet(1, -1, -2, 0, 0), 1e24)
90 *pow(Foam::max(alpha - alphaMinFriction, scalar(0)), 10.0);
91}
92
93
97(
98 const phaseModel& phase,
99 const dimensionedScalar& alphaMinFriction,
101) const
102{
103 const volScalarField& alpha = phase;
104
105 return
106 dimensionedScalar("", dimensionSet(1, -1, -2, 0, 0), 1e25)
107 *pow(Foam::max(alpha - alphaMinFriction, scalar(0)), 9.0);
108}
109
110
113(
114 const phaseModel& phase,
115 const dimensionedScalar& alphaMinFriction,
117 const volScalarField& pf,
118 const volSymmTensorField& D
119) const
120{
121 const volScalarField& alpha = phase;
122
123 auto tnu = volScalarField::New
124 (
125 IOobject::scopedName("Schaeffer", "nu"),
127 phase.mesh(),
128 dimensionedScalar(dimensionSet(0, 2, -1, 0, 0), Zero)
129 );
130 auto& nuf = tnu.ref();
131
132 forAll(D, celli)
133 {
134 if (alpha[celli] > alphaMinFriction.value())
135 {
136 nuf[celli] =
137 0.5*pf[celli]*sin(phi_.value())
138 /(
139 sqrt((1.0/3.0)*sqr(tr(D[celli])) - invariantII(D[celli]))
140 + SMALL
141 );
142 }
143 }
144
145 const fvPatchList& patches = phase.mesh().boundary();
146 const tmp<volVectorField> tU(phase.U());
147 const volVectorField& U = tU();
148
149 volScalarField::Boundary& nufBf = nuf.boundaryFieldRef();
150
151 forAll(patches, patchi)
152 {
153 if (!patches[patchi].coupled())
154 {
155 nufBf[patchi] =
156 (
157 pf.boundaryField()[patchi]*sin(phi_.value())
158 /(
159 mag(U.boundaryField()[patchi].snGrad())
160 + SMALL
161 )
162 );
163 }
164 }
165
166 // Correct coupled BCs
167 nuf.correctBoundaryConditions();
168
169 return tnu;
170}
171
172
174{
175 coeffDict_ <<= dict_.optionalSubDict(typeName + "Coeffs");
176
177 phi_.read(coeffDict_);
178 phi_ *= degToRad();
179
180 return true;
181}
182
183
184// ************************************************************************* //
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 noexcept
Return const reference to mesh.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ NO_REGISTER
Do not request registration (bool: false).
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
const Type & value() const noexcept
Return const reference to value.
frictionalStressModel(const dictionary &dict)
Construct from components.
const dictionary & dict_
Reference to higher-level dictionary for re-read.
Schaeffer(const dictionary &dict)
Construct from components.
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
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phaseModel.H:56
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phase.H:53
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.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
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)
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
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.