initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / applications / solvers / multiphase / twoPhaseEulerFoam / interfacialModels / dragModels / Gibilaro / Gibilaro.C
blobfdb131324fd73d9221dbbd454c39477d08f74375
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
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.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "Gibilaro.H"
28 #include "addToRunTimeSelectionTable.H"
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 namespace Foam
34     defineTypeNameAndDebug(Gibilaro, 0);
36     addToRunTimeSelectionTable
37     (
38         dragModel,
39         Gibilaro,
40         dictionary
41     );
45 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
47 Foam::Gibilaro::Gibilaro
49     const dictionary& interfaceDict,
50     const volScalarField& alpha,
51     const phaseModel& phasea,
52     const phaseModel& phaseb
55     dragModel(interfaceDict, alpha, phasea, phaseb)
59 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
61 Foam::Gibilaro::~Gibilaro()
65 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
67 Foam::tmp<Foam::volScalarField> Foam::Gibilaro::K
69     const volScalarField& Ur
70 ) const
72     volScalarField beta = max(scalar(1) - alpha_, scalar(1.0e-6));
73     volScalarField bp = pow(beta, -2.8);
74     volScalarField Re = max(beta*Ur*phasea_.d()/phaseb_.nu(), scalar(1.0e-3));
76     return (17.3/Re + scalar(0.336))*phaseb_.rho()*Ur*bp/phasea_.d();
80 // ************************************************************************* //