SyamlalOBrien.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 -------------------------------------------------------------------------------
10 License
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 "SyamlalOBrien.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace dragModels
36 {
37  defineTypeNameAndDebug(SyamlalOBrien, 0);
38 
40  (
41  dragModel,
42  SyamlalOBrien,
43  dictionary
44  );
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
52 (
53  const dictionary& interfaceDict,
54  const phaseModel& phase1,
55  const phaseModel& phase2
56 )
57 :
58  dragModel(interfaceDict, phase1, phase2)
59 {}
60 
61 
62 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
63 
65 {}
66 
67 
68 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
69 
71 (
72  const volScalarField& Ur
73 ) const
74 {
75  volScalarField alpha2(max(phase2_, scalar(1.0e-6)));
76  volScalarField A(pow(alpha2, 4.14));
78  (
79  neg(alpha2 - 0.85)*(0.8*pow(alpha2, 1.28))
80  + pos0(alpha2 - 0.85)*(pow(alpha2, 2.65))
81  );
82 
83  volScalarField Re(max(Ur*phase1_.d()/phase2_.nu(), scalar(1.0e-3)));
84 
86  (
87  0.5*
88  (
89  A - 0.06*Re + sqrt(sqr(0.06*Re) + 0.12*Re*(2.0*B - A) + sqr(A))
90  )
91  );
92 
93  volScalarField Cds(sqr(0.63 + 4.8*sqrt(Vr/Re)));
94 
95  return 0.75*Cds*phase2_.rho()*Ur/(phase1_.d()*sqr(Vr));
96 }
97 
98 
99 // ************************************************************************* //
Foam::phaseModel
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:54
Foam::dragModels::defineTypeNameAndDebug
defineTypeNameAndDebug(blended, 0)
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
alpha2
const volScalarField & alpha2
Definition: setRegionFluidFields.H:9
Foam::dragModel::K
virtual tmp< volScalarField > K() const
Return the drag coefficient K.
Definition: dragModel.C:153
B
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
A
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
Foam::pos0
dimensionedScalar pos0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:188
Foam::dragModels::SyamlalOBrien::~SyamlalOBrien
virtual ~SyamlalOBrien()
Destructor.
Definition: SyamlalOBrien.C:64
Foam::dragModels::SyamlalOBrien::SyamlalOBrien
SyamlalOBrien(const dictionary &interfaceDict, const phaseModel &phase1, const phaseModel &phase2)
Construct from components.
Definition: SyamlalOBrien.C:52
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::dragModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(dragModel, blended, dictionary)
phase2
phaseModel & phase2
Definition: setRegionFluidFields.H:6
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::Re
scalarField Re(const UList< complex > &cf)
Extract real component.
Definition: complexField.C:159
SyamlalOBrien.H
Foam::dragModel
Definition: dragModel.H:51
phase1
phaseModel & phase1
Definition: setRegionFluidFields.H:5
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::neg
dimensionedScalar neg(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:199