Updating InteractingKinematicParcel with changes in KinematicParcel.

This commit is contained in:
graham
2009-09-03 12:25:03 +01:00
parent 0e49835fd7
commit ad9c4d5b5c
3 changed files with 42 additions and 6 deletions

View File

@ -102,6 +102,12 @@ void Foam::InteractingKinematicParcel<ParcelType>::calc
const scalar rho0 = rho_;
const scalar mass0 = mass();
// Reynolds number
const scalar Re = this->Re(U0, d0, rhoc_, muc_);
// Sources
//~~~~~~~~
// Explicit momentum source for particle
vector Su = vector::zero;
@ -115,7 +121,8 @@ void Foam::InteractingKinematicParcel<ParcelType>::calc
// ~~~~~~
// Calculate new particle velocity
vector U1 = calcVelocity(td, dt, cellI, d0, U0, rho0, mass0, Su, dUTrans);
vector U1 =
calcVelocity(td, dt, cellI, Re, muc_, d0, U0, rho0, mass0, Su, dUTrans);
// Accumulate carrier phase source terms
@ -140,6 +147,8 @@ const Foam::vector Foam::InteractingKinematicParcel<ParcelType>::calcVelocity
TrackData& td,
const scalar dt,
const label cellI,
const scalar Re,
const scalar mu,
const scalar d,
const vector& U,
const scalar rho,
@ -151,8 +160,7 @@ const Foam::vector Foam::InteractingKinematicParcel<ParcelType>::calcVelocity
const polyMesh& mesh = this->cloud().pMesh();
// Momentum transfer coefficient
const scalar utc =
td.cloud().drag().utc(U - Uc_, d, rhoc_, muc_) + ROOTVSMALL;
const scalar utc = td.cloud().drag().utc(Re, d, mu) + ROOTVSMALL;
// Momentum source due to particle forces
const vector FCoupled =

View File

@ -49,7 +49,6 @@ SourceFiles
#include "autoPtr.H"
#include "interpolationCellPoint.H"
#include "contiguous.H"
#include "mathConstants.H"
#include "InteractingKinematicCloud.H"
@ -243,6 +242,8 @@ protected:
TrackData& td,
const scalar dt, // timestep
const label cellI, // owner cell
const scalar Re, // Reynolds number
const scalar mu, // local carrier viscosity
const scalar d, // diameter
const vector& U, // velocity
const scalar rho, // density
@ -401,6 +402,15 @@ public:
//- Surface area for given diameter
inline scalar areaS(const scalar d) const;
//- Reynolds number
inline scalar Re
(
const vector& U, // particle velocity
const scalar d, // particle diameter
const scalar rhoc, // carrier density
const scalar muc // carrier dynamic viscosity
) const;
// Main calculation loop

View File

@ -24,6 +24,10 @@ License
\*---------------------------------------------------------------------------*/
#include "mathConstants.H"
using namespace Foam::constant;
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
template <class ParcelType>
@ -374,7 +378,7 @@ template <class ParcelType>
inline Foam::scalar
Foam::InteractingKinematicParcel<ParcelType>::volume(const scalar d) const
{
return constant::math::pi/6.0*pow3(d);
return math::pi/6.0*pow3(d);
}
@ -404,7 +408,21 @@ template <class ParcelType>
inline Foam::scalar
Foam::InteractingKinematicParcel<ParcelType>::areaS(const scalar d) const
{
return constant::math::pi*d*d;
return math::pi*d*d;
}
template <class ParcelType>
inline Foam::scalar
Foam::InteractingKinematicParcel<ParcelType>::Re
(
const vector& U,
const scalar d,
const scalar rhoc,
const scalar muc
) const
{
return rhoc*mag(U - Uc_)*d/muc;
}