mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
ENH: Improvements to particle integration methods
This commit is contained in:
@ -61,7 +61,7 @@ Foam::Analytical<Type>::integrate
|
|||||||
(
|
(
|
||||||
const Type& phi,
|
const Type& phi,
|
||||||
const scalar dt,
|
const scalar dt,
|
||||||
const Type& alpha,
|
const Type& alphaBeta,
|
||||||
const scalar beta
|
const scalar beta
|
||||||
) const
|
) const
|
||||||
{
|
{
|
||||||
@ -69,9 +69,18 @@ Foam::Analytical<Type>::integrate
|
|||||||
|
|
||||||
const scalar expTerm = exp(min(50, -beta*dt));
|
const scalar expTerm = exp(min(50, -beta*dt));
|
||||||
|
|
||||||
retValue.average() =
|
if (beta > ROOTVSMALL)
|
||||||
alpha + (phi - alpha)*(1 - expTerm)/(beta*dt + ROOTVSMALL);
|
{
|
||||||
retValue.value() = alpha + (phi - alpha)*expTerm;
|
const Type alpha = alphaBeta/beta;
|
||||||
|
retValue.average() = alpha + (phi - alpha)*(1 - expTerm)/(beta*dt);
|
||||||
|
retValue.value() = alpha + (phi - alpha)*expTerm;
|
||||||
|
}
|
||||||
|
else
|
||||||
|
{
|
||||||
|
retValue.average() = phi;
|
||||||
|
retValue.value() = phi;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
return retValue;
|
return retValue;
|
||||||
}
|
}
|
||||||
|
|||||||
@ -83,7 +83,7 @@ public:
|
|||||||
(
|
(
|
||||||
const Type& phi,
|
const Type& phi,
|
||||||
const scalar dt,
|
const scalar dt,
|
||||||
const Type& alpha,
|
const Type& alphaBeta,
|
||||||
const scalar beta
|
const scalar beta
|
||||||
) const;
|
) const;
|
||||||
};
|
};
|
||||||
|
|||||||
@ -61,12 +61,12 @@ Foam::Euler<Type>::integrate
|
|||||||
(
|
(
|
||||||
const Type& phi,
|
const Type& phi,
|
||||||
const scalar dt,
|
const scalar dt,
|
||||||
const Type& alpha,
|
const Type& alphaBeta,
|
||||||
const scalar beta
|
const scalar beta
|
||||||
) const
|
) const
|
||||||
{
|
{
|
||||||
typename IntegrationScheme<Type>::integrationResult retValue;
|
typename IntegrationScheme<Type>::integrationResult retValue;
|
||||||
retValue.value() = (phi + beta*dt*alpha)/(1.0 + beta*dt);
|
retValue.value() = (phi + alphaBeta*dt)/(1.0 + beta*dt);
|
||||||
retValue.average() = 0.5*(phi + retValue.value());
|
retValue.average() = 0.5*(phi + retValue.value());
|
||||||
|
|
||||||
return retValue;
|
return retValue;
|
||||||
|
|||||||
@ -80,7 +80,7 @@ public:
|
|||||||
(
|
(
|
||||||
const Type& phi,
|
const Type& phi,
|
||||||
const scalar dt,
|
const scalar dt,
|
||||||
const Type& alpha,
|
const Type& alphaBeta,
|
||||||
const scalar beta
|
const scalar beta
|
||||||
) const;
|
) const;
|
||||||
};
|
};
|
||||||
|
|||||||
@ -62,7 +62,7 @@ Foam::IntegrationScheme<Type>::integrate
|
|||||||
(
|
(
|
||||||
const Type& phi,
|
const Type& phi,
|
||||||
const scalar dt,
|
const scalar dt,
|
||||||
const Type& alpha,
|
const Type& alphaBeta,
|
||||||
const scalar beta
|
const scalar beta
|
||||||
) const
|
) const
|
||||||
{
|
{
|
||||||
|
|||||||
@ -185,7 +185,7 @@ public:
|
|||||||
(
|
(
|
||||||
const Type& phi,
|
const Type& phi,
|
||||||
const scalar dt,
|
const scalar dt,
|
||||||
const Type& alpha,
|
const Type& alphaBeta,
|
||||||
const scalar beta
|
const scalar beta
|
||||||
) const;
|
) const;
|
||||||
};
|
};
|
||||||
|
|||||||
@ -192,21 +192,20 @@ const Foam::vector Foam::KinematicParcel<ParcelType>::calcVelocity
|
|||||||
const parcelType& p = static_cast<const parcelType&>(*this);
|
const parcelType& p = static_cast<const parcelType&>(*this);
|
||||||
const forceSuSp Fcp = forces.calcCoupled(p, dt, mass, Re, mu);
|
const forceSuSp Fcp = forces.calcCoupled(p, dt, mass, Re, mu);
|
||||||
const forceSuSp Fncp = forces.calcNonCoupled(p, dt, mass, Re, mu);
|
const forceSuSp Fncp = forces.calcNonCoupled(p, dt, mass, Re, mu);
|
||||||
forceSuSp Feff = Fcp + Fncp;
|
const forceSuSp Feff = Fcp + Fncp;
|
||||||
Feff.Sp() += ROOTVSMALL;
|
|
||||||
|
|
||||||
|
|
||||||
// New particle velocity
|
// New particle velocity
|
||||||
//~~~~~~~~~~~~~~~~~~~~~~
|
//~~~~~~~~~~~~~~~~~~~~~~
|
||||||
|
|
||||||
// Update velocity - treat as 3-D
|
// Update velocity - treat as 3-D
|
||||||
const vector ap = Uc_ + (Feff.Su() + Su)/Feff.Sp();
|
const vector abp = (Feff.Sp()*Uc_ + (Feff.Su() + Su))/mass;
|
||||||
const scalar bp = Feff.Sp()/mass;
|
const scalar bp = Feff.Sp()/mass;
|
||||||
|
|
||||||
Spu = Feff.Sp();
|
Spu = Feff.Sp();
|
||||||
|
|
||||||
IntegrationScheme<vector>::integrationResult Ures =
|
IntegrationScheme<vector>::integrationResult Ures =
|
||||||
td.cloud().UIntegrator().integrate(U, dt, ap, bp);
|
td.cloud().UIntegrator().integrate(U, dt, abp, bp);
|
||||||
|
|
||||||
vector Unew = Ures.value();
|
vector Unew = Ures.value();
|
||||||
|
|
||||||
|
|||||||
@ -333,7 +333,7 @@ Foam::scalar Foam::ThermoParcel<ParcelType>::calcHeatTransfer
|
|||||||
|
|
||||||
// Integrate to find the new parcel temperature
|
// Integrate to find the new parcel temperature
|
||||||
IntegrationScheme<scalar>::integrationResult Tres =
|
IntegrationScheme<scalar>::integrationResult Tres =
|
||||||
td.cloud().TIntegrator().integrate(T, dt, ap, bp);
|
td.cloud().TIntegrator().integrate(T, dt, ap*bp, bp);
|
||||||
|
|
||||||
scalar Tnew = max(Tres.value(), td.cloud().constProps().TMin());
|
scalar Tnew = max(Tres.value(), td.cloud().constProps().TMin());
|
||||||
|
|
||||||
|
|||||||
Reference in New Issue
Block a user