Skip to content
Snippets Groups Projects
Commit c282a02a authored by Henry's avatar Henry
Browse files

Added LopezDeBertodano turbulentDispersionModel

parent 8fa9f894
Branches
Tags
No related merge requests found
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2014 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "LopezDeBertodano.H"
#include "phasePair.H"
#include "fvc.H"
#include "PhaseCompressibleTurbulenceModel.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
namespace turbulentDispersionModels
{
defineTypeNameAndDebug(LopezDeBertodano, 0);
addToRunTimeSelectionTable
(
turbulentDispersionModel,
LopezDeBertodano,
dictionary
);
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::turbulentDispersionModels::LopezDeBertodano::LopezDeBertodano
(
const dictionary& dict,
const phasePair& pair
)
:
turbulentDispersionModel(dict, pair),
Ctd_("Ctd", dimless, dict.lookup("Ctd"))
{}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
Foam::turbulentDispersionModels::LopezDeBertodano::~LopezDeBertodano()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
Foam::tmp<Foam::volVectorField>
Foam::turbulentDispersionModels::LopezDeBertodano::F() const
{
return
Ctd_
*pair_.continuous().rho()
*pair_.continuous().turbulence().k()
*fvc::grad(pair_.dispersed());
}
// ************************************************************************* //
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2014 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Class
Foam::turbulentDispersionModels::LopezDeBertodano
Description
Lopez de Bertodano (1992) turbulent dispersion model.
\verbatim
"Turbulent bubbly two-phase flow in a triangular
duct"
Lopez de Bertodano, M.
Ph.D. Thesis, Rensselaer Polytechnic Institution, New York, USA, 1992.
\endverbatim
\verbatim
"The Favre averaged drag model for turbulent dispersion in Eulerian
multi-phase flows"
Burns, A.D., Frank, T., Hamill, I., Shi, J.M.,
5th international conference on multiphase flow
Volume 4, Paper 392, May 2004
\endverbatim
SourceFiles
LopezDeBertodano.C
\*---------------------------------------------------------------------------*/
#ifndef LopezDeBertodano_H
#define LopezDeBertodano_H
#include "turbulentDispersionModel.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
class phasePair;
namespace turbulentDispersionModels
{
/*---------------------------------------------------------------------------*\
Class LopezDeBertodano Declaration
\*---------------------------------------------------------------------------*/
class LopezDeBertodano
:
public turbulentDispersionModel
{
// Private data
//- Constant turbulent dispersion coefficient
const dimensionedScalar Ctd_;
public:
//- Runtime type information
TypeName("LopezDeBertodano");
// Constructors
//- Construct from a dictionary and a phase pair
LopezDeBertodano
(
const dictionary& dict,
const phasePair& pair
);
//- Destructor
virtual ~LopezDeBertodano();
// Member Functions
//- Turbulent dispersion force
virtual tmp<volVectorField> F() const;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace turbulentDispersionModels
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //
0% or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment