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 | Copyright (C) 2011-2015 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 "GidaspowSchillerNaumann.H"
27 #include "phasePair.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace dragModels
35 {
36  defineTypeNameAndDebug(GidaspowSchillerNaumann, 0);
37  addToRunTimeSelectionTable(dragModel, GidaspowSchillerNaumann, dictionary);
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
45 (
46  const dictionary& dict,
47  const phasePair& pair,
48  const bool registerObject
49 )
50 :
51  dragModel(dict, pair, registerObject),
52  residualRe_("residualRe", dimless, dict)
53 {}
54 
55 
56 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
57 
59 {}
60 
61 
62 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
63 
66 {
68  (
69  max(scalar(1) - pair_.dispersed(), pair_.continuous().residualAlpha())
70  );
71 
73 
74  volScalarField CdsRe
75  (
76  neg(Re - 1000)*24.0*(1.0 + 0.15*pow(Re, 0.687))/alpha2
77  + pos(Re - 1000)*0.44*max(Re, residualRe_)
78  );
79 
80  return
81  CdsRe
82  *pow(alpha2, -2.65)
84 }
85 
86 
87 // ************************************************************************* //
dictionary dict
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
virtual const phaseModel & continuous() const
Continuous phase.
const dimensionedScalar & residualAlpha() const
Return the residual phase-fraction for given phase.
dimensionedScalar neg(const dimensionedScalar &ds)
const phasePair & pair_
Phase pair.
Definition: dragModel.H:62
Macros for easy insertion into run-time selection tables.
virtual ~GidaspowSchillerNaumann()
Destructor.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dimensionedScalar pos(const dimensionedScalar &ds)
virtual tmp< volScalarField > CdRe() const
Drag coefficient.
GidaspowSchillerNaumann(const dictionary &interfaceDict, const phaseModel &phase1, const phaseModel &phase2)
Construct from components.
alpha2
Definition: alphaEqn.H:112
tmp< volScalarField > Re() const
Reynolds number.
virtual const phaseModel & dispersed() const
Dispersed phase.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
A class for managing temporary objects.
Definition: PtrList.H:54
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
Namespace for OpenFOAM.