GidaspowSchillerNaumann.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-2017 OpenFOAM Foundation
9 Copyright (C) 2022 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
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
36namespace multiphaseEuler
37{
38namespace dragModels
39{
41
43 (
47 );
48}
49}
50}
51
52
53// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54
55Foam::multiphaseEuler::dragModels::GidaspowSchillerNaumann::
56GidaspowSchillerNaumann
57(
58 const dictionary& interfaceDict,
59 const phaseModel& phase1,
60 const phaseModel& phase2
61)
62:
63 dragModel(interfaceDict, phase1, phase2)
64{}
65
66
67// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
68
71{}
72
73
74// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
75
78(
79 const volScalarField& Ur
80) const
81{
82 volScalarField alpha2(max(phase2_, scalar(1e-6)));
83 volScalarField bp(pow(alpha2, -2.65));
84
85 volScalarField Re(max(alpha2*Ur*phase1_.d()/phase2_.nu(), scalar(1.0e-3)));
87 (
88 neg(Re - 1000)*(24.0*(1.0 + 0.15*pow(Re, 0.687))/Re)
89 + pos0(Re - 1000)*0.44
90 );
91
92 return 0.75*Cds*phase2_.rho()*Ur*bp/phase1_.d();
93}
94
95
96// ************************************************************************* //
CGAL::Exact_predicates_exact_constructions_kernel K
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
phaseModel & phase1
const volScalarField & alpha2
phaseModel & phase2
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
H, Enwald, E. Peirano, A-E Almstedt 'Eulerian Two-Phase Flow Theory Applied to Fluidization' Int....
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:61
A class for managing temporary objects.
Definition: tmp.H:65
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
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
scalarField Re(const UList< complex > &cf)
Extract real component.
Definition: complexField.C:159
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar neg(const dimensionedScalar &ds)
volScalarField & e
Definition: createFields.H:11