multiphaseEulerFoam/.../aerosolDrag: Improvements

Expanded the documentation and updated the mean free path calculation

Patch contributed by Institute of Fluid Dynamics,
Helmholtz-Zentrum Dresden - Rossendorf (HZDR)
This commit is contained in:
Will Bainbridge
2020-08-26 14:19:51 +01:00
parent 1844051d91
commit 36ce8b31ae
4 changed files with 59 additions and 31 deletions

View File

@ -29,8 +29,6 @@ License
#include "addToRunTimeSelectionTable.H"
#include "constants.H"
using namespace Foam::constant;
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
@ -42,6 +40,9 @@ namespace dragModels
}
}
using Foam::constant::physicoChemical::k;
using Foam::constant::mathematical::pi;
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
@ -56,7 +57,7 @@ Foam::dragModels::aerosolDrag::aerosolDrag
A1_(dict.lookupOrDefault<scalar>("A1", 2.514)),
A2_(dict.lookupOrDefault<scalar>("A2", 0.8)),
A3_(dict.lookupOrDefault<scalar>("A3", 0.55)),
dm_(dimensionedScalar::lookupOrDefault("dm", dict, dimLength, 364e-12))
sigma_("sigma", dimLength, dict)
{}
@ -70,24 +71,14 @@ Foam::dragModels::aerosolDrag::~aerosolDrag()
Foam::tmp<Foam::volScalarField> Foam::dragModels::aerosolDrag::CdRe() const
{
const volScalarField lambda
(
physicoChemical::k
*pair_.continuous().thermo().T()
/sqrt(2.0)
/mathematical::pi
/pair_.continuous().thermo().p()
/sqr(dm_)
);
const volScalarField& T = pair_.continuous().thermo().T();
const volScalarField& p = pair_.continuous().thermo().p();
tmp<volScalarField> td(pair_.dispersed().d());
const volScalarField& d = td();
const volScalarField dp(pair_.dispersed().d());
const volScalarField lambda(k*T/(sqrt(2.0)*pi*p*sqr(sigma_)));
const volScalarField Cc
(
1 + lambda/dp*(A1_ + A2_*exp(-A3_*dp/lambda))
);
return 24/Cc;
return 24/(1 + lambda/d*(A1_ + A2_*exp(-A3_*d/lambda)));
}

View File

@ -25,9 +25,42 @@ Class
Foam::dragModels::aerosolDrag
Description
This drag model is intended for simulations involving dispersed phase
with a very small diameter. It combines the Stokes drag with a
Cunningham correction factor.
Stokes drag with Cunningham slip correction. The slip correction coefficient
is implemented in the following form:
\f[
C_c = 1 + \lambda [A_1 + A_2 \exp(-A_3 d/\lambda)]/d\,.
\f]
The coefficients default to the values proposed by Davis (1945). The mean
free path is computed by
\f[
\lambda = \frac{kT}{\sqrt{2} \pi p \sigma^{2}}\,.
\f]
\vartable
A_1 | Coefficient [-]
A_2 | Coefficient [-]
A_3 | Coefficient [-]
\sigma | Lennard-Jones parameter [m]
\endvartable
Reference:
\verbatim
Davies, C. N. (1945).
Definitive equations for the fluid resistance of spheres.
Proceedings of the Physical Society, 57(4), 259.
\endverbatim
Usage
\table
Property | Description | Required | Default value
A1 | Coefficient A1 | no | 2.514
A2 | Coefficient A2 | no | 0.8
A3 | Coefficient A2 | no | 0.55
sigma | Lennard-Jones parameter | yes | none
\endtable
SourceFiles
aerosolDrag.C
@ -47,7 +80,7 @@ namespace dragModels
{
/*---------------------------------------------------------------------------*\
Class aerosolDrag Declaration
Class aerosolDrag Declaration
\*---------------------------------------------------------------------------*/
class aerosolDrag
@ -56,17 +89,17 @@ class aerosolDrag
{
// Private Data
//- Model coefficient A1
const scalar A1_;
//- Cunningham slip correction coefficient A1
const dimensionedScalar A1_;
//- Model coefficient A2
const scalar A2_;
//- Cunningham slip correction coefficient A2
const dimensionedScalar A2_;
//- Model coefficient A3
const scalar A3_;
//- Cunningham slip correction coefficient A3
const dimensionedScalar A3_;
//- Diameter of the continuous phase molecules, defaults to value for N2
const dimensionedScalar dm_;
//- Lennard-Jones sigma parameter
const dimensionedScalar sigma_;
public:

View File

@ -159,6 +159,8 @@ drag
{
type aerosolDrag;
sigma 340e-12;
swarmCorrection
{
type none;

View File

@ -173,6 +173,8 @@ drag
{
type aerosolDrag;
sigma 340e-12;
swarmCorrection
{
type none;