mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
Removing the cached deltaT from the cloud, using the new deltaTValue
access memeber instead - no longer expensive to obtain.
This commit is contained in:
@ -297,7 +297,7 @@ void Foam::DsmcCloud<ParcelType>::collisions()
|
|||||||
// Temporary storage for subCells
|
// Temporary storage for subCells
|
||||||
List<DynamicList<label> > subCells(8);
|
List<DynamicList<label> > subCells(8);
|
||||||
|
|
||||||
scalar deltaT = cachedDeltaT();
|
scalar deltaT = mesh().time().deltaTValue();
|
||||||
|
|
||||||
label collisionCandidates = 0;
|
label collisionCandidates = 0;
|
||||||
|
|
||||||
@ -1042,9 +1042,6 @@ Foam::DsmcCloud<ParcelType>::~DsmcCloud()
|
|||||||
template<class ParcelType>
|
template<class ParcelType>
|
||||||
void Foam::DsmcCloud<ParcelType>::evolve()
|
void Foam::DsmcCloud<ParcelType>::evolve()
|
||||||
{
|
{
|
||||||
// cache the value of deltaT for this timestep
|
|
||||||
storeDeltaT();
|
|
||||||
|
|
||||||
typename ParcelType::trackData td(*this);
|
typename ParcelType::trackData td(*this);
|
||||||
|
|
||||||
// Reset the data collection fields
|
// Reset the data collection fields
|
||||||
|
|||||||
@ -137,10 +137,6 @@ class DsmcCloud
|
|||||||
//- Random number generator
|
//- Random number generator
|
||||||
Random rndGen_;
|
Random rndGen_;
|
||||||
|
|
||||||
//- In-cloud cache of deltaT, lookup in submodels and parcel is
|
|
||||||
// expensive
|
|
||||||
scalar cachedDeltaT_;
|
|
||||||
|
|
||||||
|
|
||||||
// boundary value fields
|
// boundary value fields
|
||||||
|
|
||||||
@ -267,12 +263,6 @@ public:
|
|||||||
//- Return refernce to the random object
|
//- Return refernce to the random object
|
||||||
inline Random& rndGen();
|
inline Random& rndGen();
|
||||||
|
|
||||||
//- Store (cache) the current value of deltaT
|
|
||||||
inline void storeDeltaT();
|
|
||||||
|
|
||||||
//- Return the cached value of deltaT
|
|
||||||
inline scalar cachedDeltaT() const;
|
|
||||||
|
|
||||||
|
|
||||||
// References to the boundary fields for surface data collection
|
// References to the boundary fields for surface data collection
|
||||||
|
|
||||||
|
|||||||
@ -125,20 +125,6 @@ inline Foam::Random& Foam::DsmcCloud<ParcelType>::rndGen()
|
|||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
template<class ParcelType>
|
|
||||||
inline void Foam::DsmcCloud<ParcelType>::storeDeltaT()
|
|
||||||
{
|
|
||||||
cachedDeltaT_ = mesh().time().deltaT().value();
|
|
||||||
}
|
|
||||||
|
|
||||||
|
|
||||||
template<class ParcelType>
|
|
||||||
inline Foam::scalar Foam::DsmcCloud<ParcelType>::cachedDeltaT() const
|
|
||||||
{
|
|
||||||
return cachedDeltaT_;
|
|
||||||
}
|
|
||||||
|
|
||||||
|
|
||||||
template<class ParcelType>
|
template<class ParcelType>
|
||||||
inline Foam::volScalarField::GeometricBoundaryField&
|
inline Foam::volScalarField::GeometricBoundaryField&
|
||||||
Foam::DsmcCloud<ParcelType>::qBF()
|
Foam::DsmcCloud<ParcelType>::qBF()
|
||||||
|
|||||||
@ -44,7 +44,7 @@ bool Foam::DsmcParcel<ParcelType>::move
|
|||||||
const polyMesh& mesh = td.cloud().pMesh();
|
const polyMesh& mesh = td.cloud().pMesh();
|
||||||
const polyBoundaryMesh& pbMesh = mesh.boundaryMesh();
|
const polyBoundaryMesh& pbMesh = mesh.boundaryMesh();
|
||||||
|
|
||||||
const scalar deltaT = td.cloud().cachedDeltaT();
|
const scalar deltaT = mesh.time().deltaTValue();
|
||||||
scalar tEnd = (1.0 - p.stepFraction())*deltaT;
|
scalar tEnd = (1.0 - p.stepFraction())*deltaT;
|
||||||
const scalar dtMax = tEnd;
|
const scalar dtMax = tEnd;
|
||||||
|
|
||||||
@ -135,7 +135,7 @@ void Foam::DsmcParcel<ParcelType>::hitWallPatch
|
|||||||
|
|
||||||
const scalar fA = mag(wpp.faceAreas()[wppLocalFace]);
|
const scalar fA = mag(wpp.faceAreas()[wppLocalFace]);
|
||||||
|
|
||||||
const scalar deltaT = td.cloud().cachedDeltaT();
|
const scalar deltaT = td.cloud().pMesh().time().deltaTValue();
|
||||||
|
|
||||||
const constantProperties& constProps(td.cloud().constProps(typeId_));
|
const constantProperties& constProps(td.cloud().constProps(typeId_));
|
||||||
|
|
||||||
|
|||||||
@ -129,7 +129,7 @@ void Foam::FreeStream<CloudType>::inflow()
|
|||||||
|
|
||||||
const polyMesh& mesh(cloud.mesh());
|
const polyMesh& mesh(cloud.mesh());
|
||||||
|
|
||||||
const scalar deltaT = cloud.cachedDeltaT();
|
const scalar deltaT = mesh.time().deltaTValue();
|
||||||
|
|
||||||
Random& rndGen(cloud.rndGen());
|
Random& rndGen(cloud.rndGen());
|
||||||
|
|
||||||
|
|||||||
Reference in New Issue
Block a user