Loading...
Searching...
No Matches
Beetstra.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-------------------------------------------------------------------------------
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
28#include "Beetstra.H"
29#include "phasePair.H"
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:
54 residualRe_("residualRe", dimless, dict)
55{}
56
57
58// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
61{}
62
63
64// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
65
67{
69 (
70 max(pair_.dispersed(), pair_.continuous().residualAlpha())
71 );
72
74 (
75 max(scalar(1) - pair_.dispersed(), pair_.continuous().residualAlpha())
76 );
77
78 volScalarField Res(alpha2*pair_.Re());
79
80 volScalarField ResLim
81 (
82 "ReLim",
83 max(Res, residualRe_)
84 );
85
87 (
88 "F0",
89 10*alpha1/sqr(alpha2) + sqr(alpha2)*(1 + 1.5*sqrt(alpha1))
90 );
91
93 (
94 "F1",
95 0.413*Res/(24*sqr(alpha2))*(1.0/alpha2
96 + 3*alpha1*alpha2 + 8.4*pow(ResLim, -0.343))
97 /(1 + pow(10.0, 3*alpha1)*pow(ResLim, -(1 + 4*alpha1)/2.0))
98 );
99
100 return 24*alpha2*(F0 + F1);
101}
102
103
104// ************************************************************************* //
#define F1(B, C, D)
Definition SHA1.C:149
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 & alpha1
const volScalarField & alpha2
bool registerObject() const noexcept
Should objects created with this IOobject be registered?
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
dragModel(const phasePair &pair, const bool registerObject)
Definition dragModel.C:43
const phasePair & pair_
Phase pair.
Definition dragModel.H:63
Drag model of Beetstra et al. for monodisperse gas-particle flows obtained with direct numerical simu...
Definition Beetstra.H:64
virtual ~Beetstra()
Destructor.
Definition Beetstra.C:53
Beetstra(const dictionary &dict, const phasePair &pair, const bool registerObject)
Construct from a dictionary and a phase pair.
Definition Beetstra.C:40
virtual tmp< volScalarField > CdRe() const
Drag coefficient.
Definition Beetstra.C:59
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition phasePair.H:52
A class for managing temporary objects.
Definition tmp.H:75
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
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
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dictionary dict