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-2018 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
29#include "phasePair.H"
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
36namespace dragModels
37{
40}
41}
42
43
44// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45
47(
48 const dictionary& dict,
49 const phasePair& pair,
50 const bool registerObject
51)
52:
53 dragModel(dict, pair, registerObject),
54 residualRe_("residualRe", dimless, dict)
55{}
56
57
58// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
59
61{}
62
63
64// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
65
68{
70 (
71 max(scalar(1) - pair_.dispersed(), pair_.continuous().residualAlpha())
72 );
73
74 volScalarField Re(alpha2*pair_.Re());
75
76 volScalarField CdsRe
77 (
78 neg(Re - 1000)*24*(1.0 + 0.15*pow(Re, 0.687))/alpha2
79 + pos0(Re - 1000)*0.44*max(Re, residualRe_)
80 );
81
82 return
83 CdsRe
84 *pow(alpha2, -2.65)
85 *max(pair_.continuous(), pair_.continuous().residualAlpha());
86}
87
88
89// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const volScalarField & alpha2
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Gidaspow, Schiller and Naumann drag model.
virtual tmp< volScalarField > CdRe() const
Drag coefficient.
H, Enwald, E. Peirano, A-E Almstedt 'Eulerian Two-Phase Flow Theory Applied to Fluidization' Int....
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:56
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
const dimensionSet dimless
Dimensionless.
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)
dictionary dict