GidaspowErgunWenYu.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 | Copyright (C) 2011-2012 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "GidaspowErgunWenYu.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace dragModels
34 {
35  defineTypeNameAndDebug(GidaspowErgunWenYu, 0);
36 
38  (
39  dragModel,
40  GidaspowErgunWenYu,
41  dictionary
42  );
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
50 (
51  const dictionary& interfaceDict,
52  const phaseModel& phase1,
53  const phaseModel& phase2
54 )
55 :
56  dragModel(interfaceDict, phase1, phase2)
57 {}
58 
59 
60 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
61 
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
67 
69 (
70  const volScalarField& Ur
71 ) const
72 {
73  volScalarField alpha2(max(phase2_, scalar(1.0e-6)));
75  volScalarField bp(pow(alpha2, -2.65));
76  volScalarField Re(max(Ur*d/phase2_.nu(), scalar(1.0e-3)));
77 
78  volScalarField Cds
79  (
80  neg(Re - 1000)*(24.0*(1.0 + 0.15*pow(Re, 0.687))/Re)
81  + pos(Re - 1000)*0.44
82  );
83 
84  // Wen and Yu (1966)
85  return
86  (
87  pos(alpha2 - 0.8)
88  *(0.75*Cds*phase2_.rho()*Ur*bp/d)
89  + neg(alpha2 - 0.8)
90  *(
91  150.0*phase1_*phase2_.nu()*phase2_.rho()/(sqr(alpha2*d))
92  + 1.75*phase2_.rho()*Ur/(alpha2*d)
93  )
94  );
95 }
96 
97 
98 // ************************************************************************* //
const double e
Elementary charge.
Definition: doubleFloat.H:78
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const phaseModel & phase1_
Definition: dragModel.H:57
virtual ~GidaspowErgunWenYu()
Destructor.
dimensionedScalar neg(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dimensionedScalar pos(const dimensionedScalar &ds)
const phaseModel & phase2_
Definition: dragModel.H:58
virtual tmp< volScalarField > K() const
Return the drag coefficient K.
const dimensionedScalar & rho() const
Definition: phaseModel.H:165
const dimensionedScalar & nu() const
Return the laminar viscosity.
Definition: phaseModel.H:150
alpha2
Definition: alphaEqn.H:112
phaseModel & phase1
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
tmp< volScalarField > d() const
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
GidaspowErgunWenYu(const dictionary &interfaceDict, const phaseModel &phase1, const phaseModel &phase2)
Construct from components.
phaseModel & phase2
A class for managing temporary objects.
Definition: PtrList.H:54
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
Namespace for OpenFOAM.