mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
ENH: Added calcFrequency for steady cloud calcs
This commit is contained in:
@ -47,6 +47,7 @@ void Foam::KinematicCloud<ParcelType>::cloudSolution::read()
|
||||
|
||||
if (steadyState())
|
||||
{
|
||||
dict_.lookup("calcFrequency") >> calcFrequency_;
|
||||
dict_.lookup("maxTrackTime") >> maxTrackTime_;
|
||||
dict_.subDict("sourceTerms").lookup("resetOnStartup")
|
||||
>> resetSourcesOnStartup_;
|
||||
@ -65,6 +66,9 @@ Foam::KinematicCloud<ParcelType>::cloudSolution::cloudSolution
|
||||
dict_(dict),
|
||||
active_(dict.lookup("active")),
|
||||
transient_(false),
|
||||
calcFrequency_(1),
|
||||
iter_(1),
|
||||
deltaT_(0.0),
|
||||
coupled_(false),
|
||||
cellValueSourceCorrection_(false),
|
||||
maxTrackTime_(0.0),
|
||||
@ -87,6 +91,9 @@ Foam::KinematicCloud<ParcelType>::cloudSolution::cloudSolution
|
||||
dict_(cs.dict_),
|
||||
active_(cs.active_),
|
||||
transient_(cs.transient_),
|
||||
calcFrequency_(cs.calcFrequency_),
|
||||
iter_(cs.iter_),
|
||||
deltaT_(cs.deltaT_),
|
||||
coupled_(cs.coupled_),
|
||||
cellValueSourceCorrection_(cs.cellValueSourceCorrection_),
|
||||
maxTrackTime_(cs.maxTrackTime_),
|
||||
@ -104,6 +111,9 @@ Foam::KinematicCloud<ParcelType>::cloudSolution::cloudSolution
|
||||
dict_(dictionary::null),
|
||||
active_(false),
|
||||
transient_(false),
|
||||
calcFrequency_(0),
|
||||
iter_(0),
|
||||
deltaT_(0.0),
|
||||
coupled_(false),
|
||||
cellValueSourceCorrection_(false),
|
||||
maxTrackTime_(0.0),
|
||||
@ -137,14 +147,36 @@ bool Foam::KinematicCloud<ParcelType>::cloudSolution::semiImplicit
|
||||
|
||||
|
||||
template<class ParcelType>
|
||||
bool Foam::KinematicCloud<ParcelType>::cloudSolution::sourceActive() const
|
||||
bool Foam::KinematicCloud<ParcelType>::cloudSolution::solveThisStep() const
|
||||
{
|
||||
return coupled_ && (active_ || steadyState());
|
||||
return
|
||||
active_
|
||||
&& (
|
||||
mesh_.time().outputTime()
|
||||
|| (mesh_.time().timeIndex() % calcFrequency_ == 0)
|
||||
);
|
||||
}
|
||||
|
||||
|
||||
template<class ParcelType>
|
||||
bool Foam::KinematicCloud<ParcelType>::cloudSolution::writeThisStep() const
|
||||
bool Foam::KinematicCloud<ParcelType>::cloudSolution::canEvolve()
|
||||
{
|
||||
// Set the calculation time step
|
||||
if (transient_)
|
||||
{
|
||||
deltaT_ = mesh_.time().deltaTValue();
|
||||
}
|
||||
else
|
||||
{
|
||||
deltaT_ = maxTrackTime_;
|
||||
}
|
||||
|
||||
return solveThisStep();
|
||||
}
|
||||
|
||||
|
||||
template<class ParcelType>
|
||||
bool Foam::KinematicCloud<ParcelType>::cloudSolution::output() const
|
||||
{
|
||||
return active_ && mesh_.time().outputTime();
|
||||
}
|
||||
@ -278,10 +310,10 @@ void Foam::KinematicCloud<ParcelType>::evolveCloud
|
||||
{
|
||||
// this->surfaceFilm().injectStreadyState(td);
|
||||
|
||||
this->injection().injectSteadyState(td, solution_.maxTrackTime());
|
||||
this->injection().injectSteadyState(td, solution_.deltaT());
|
||||
|
||||
td.part() = ParcelType::trackData::tpLinearTrack;
|
||||
Cloud<ParcelType>::move(td, solution_.maxTrackTime());
|
||||
Cloud<ParcelType>::move(td, solution_.deltaT());
|
||||
}
|
||||
}
|
||||
|
||||
@ -331,10 +363,10 @@ void Foam::KinematicCloud<ParcelType>::moveCollide
|
||||
)
|
||||
{
|
||||
td.part() = ParcelType::trackData::tpVelocityHalfStep;
|
||||
Cloud<ParcelType>::move(td, this->db().time().deltaTValue());
|
||||
Cloud<ParcelType>::move(td, solution_.deltaT());
|
||||
|
||||
td.part() = ParcelType::trackData::tpLinearTrack;
|
||||
Cloud<ParcelType>::move(td, this->db().time().deltaTValue());
|
||||
Cloud<ParcelType>::move(td, solution_.deltaT());
|
||||
|
||||
// td.part() = ParcelType::trackData::tpRotationalTrack;
|
||||
// Cloud<ParcelType>::move(td);
|
||||
@ -344,7 +376,7 @@ void Foam::KinematicCloud<ParcelType>::moveCollide
|
||||
this->collision().collide();
|
||||
|
||||
td.part() = ParcelType::trackData::tpVelocityHalfStep;
|
||||
Cloud<ParcelType>::move(td, this->db().time().deltaTValue());
|
||||
Cloud<ParcelType>::move(td, solution_.deltaT());
|
||||
}
|
||||
|
||||
|
||||
@ -362,6 +394,8 @@ void Foam::KinematicCloud<ParcelType>::postEvolve()
|
||||
forces_.cacheFields(false, solution_.interpolationSchemes());
|
||||
|
||||
this->postProcessing().post();
|
||||
|
||||
solution_.nextIter();
|
||||
}
|
||||
|
||||
|
||||
@ -732,7 +766,7 @@ void Foam::KinematicCloud<ParcelType>::relaxSources
|
||||
template<class ParcelType>
|
||||
void Foam::KinematicCloud<ParcelType>::evolve()
|
||||
{
|
||||
if (solution_.active())
|
||||
if (solution_.canEvolve())
|
||||
{
|
||||
typename ParcelType::trackData td(*this);
|
||||
|
||||
|
||||
@ -132,6 +132,16 @@ public:
|
||||
//- Transient flag
|
||||
Switch transient_;
|
||||
|
||||
//- Calculation frequency - carrier steps per cloud step
|
||||
// NOTE: Steady operation only
|
||||
label calcFrequency_;
|
||||
|
||||
//- Current cloud iteration
|
||||
label iter_;
|
||||
|
||||
//- Cloud solution time step
|
||||
scalar deltaT_;
|
||||
|
||||
|
||||
// Run-time options
|
||||
|
||||
@ -205,6 +215,18 @@ public:
|
||||
//- Return const access to the steady flag
|
||||
inline const Switch steadyState() const;
|
||||
|
||||
//- Return const access to the calculation frequency
|
||||
inline label calcFrequency() const;
|
||||
|
||||
//- Return const access to the current cloud iteration
|
||||
inline label iter() const;
|
||||
|
||||
//- Increment and return iter counter
|
||||
inline label nextIter();
|
||||
|
||||
//- Return the time step
|
||||
inline scalar deltaT() const;
|
||||
|
||||
//- Return const access to the coupled flag
|
||||
inline const Switch coupled() const;
|
||||
|
||||
@ -229,11 +251,15 @@ public:
|
||||
|
||||
// Helper functions
|
||||
|
||||
//- Returns true if sources are active (at this time)
|
||||
bool sourceActive() const;
|
||||
//- Returns true if performing a cloud iteration this calc step
|
||||
bool solveThisStep() const;
|
||||
|
||||
//- Returns true if possible to evolve the cloud and sets timestep
|
||||
// parameters
|
||||
bool canEvolve();
|
||||
|
||||
//- Returns true if writing this step
|
||||
bool writeThisStep() const;
|
||||
bool output() const;
|
||||
};
|
||||
|
||||
|
||||
@ -439,6 +465,9 @@ public:
|
||||
//- Return const access to the solution properties
|
||||
inline const cloudSolution& solution() const;
|
||||
|
||||
//- Return access to the solution properties
|
||||
inline cloudSolution& solution();
|
||||
|
||||
//- Return the constant properties
|
||||
inline const typename ParcelType::constantProperties&
|
||||
constProps() const;
|
||||
|
||||
@ -91,6 +91,36 @@ Foam::KinematicCloud<ParcelType>::cloudSolution::steadyState() const
|
||||
}
|
||||
|
||||
|
||||
template<class ParcelType>
|
||||
inline Foam::label
|
||||
Foam::KinematicCloud<ParcelType>::cloudSolution::calcFrequency() const
|
||||
{
|
||||
return calcFrequency_;
|
||||
}
|
||||
|
||||
|
||||
template<class ParcelType>
|
||||
inline Foam::label Foam::KinematicCloud<ParcelType>::cloudSolution::iter() const
|
||||
{
|
||||
return iter_;
|
||||
}
|
||||
|
||||
|
||||
template<class ParcelType>
|
||||
inline Foam::label Foam::KinematicCloud<ParcelType>::cloudSolution::nextIter()
|
||||
{
|
||||
return ++iter_;
|
||||
}
|
||||
|
||||
|
||||
template<class ParcelType>
|
||||
inline Foam::scalar
|
||||
Foam::KinematicCloud<ParcelType>::cloudSolution::deltaT() const
|
||||
{
|
||||
return deltaT_;
|
||||
}
|
||||
|
||||
|
||||
template<class ParcelType>
|
||||
inline const Foam::Switch
|
||||
Foam::KinematicCloud<ParcelType>::cloudSolution::coupled() const
|
||||
@ -164,6 +194,14 @@ Foam::KinematicCloud<ParcelType>::solution() const
|
||||
}
|
||||
|
||||
|
||||
template<class ParcelType>
|
||||
inline typename Foam::KinematicCloud<ParcelType>::cloudSolution&
|
||||
Foam::KinematicCloud<ParcelType>::solution()
|
||||
{
|
||||
return solution_;
|
||||
}
|
||||
|
||||
|
||||
template<class ParcelType>
|
||||
inline const typename ParcelType::constantProperties&
|
||||
Foam::KinematicCloud<ParcelType>::constProps() const
|
||||
@ -445,7 +483,7 @@ Foam::KinematicCloud<ParcelType>::SU(volVectorField& U) const
|
||||
<< max(UCoeff()).value() << endl;
|
||||
}
|
||||
|
||||
if (solution_.sourceActive())
|
||||
if (solution_.coupled())
|
||||
{
|
||||
if (solution_.semiImplicit("U"))
|
||||
{
|
||||
|
||||
@ -297,7 +297,7 @@ void Foam::ReactingCloud<ParcelType>::relaxSources
|
||||
template<class ParcelType>
|
||||
void Foam::ReactingCloud<ParcelType>::evolve()
|
||||
{
|
||||
if (this->solution().active())
|
||||
if (this->solution().canEvolve())
|
||||
{
|
||||
typename ParcelType::trackData td(*this);
|
||||
|
||||
|
||||
@ -90,7 +90,7 @@ Foam::ReactingCloud<ParcelType>::SYi
|
||||
volScalarField& Yi
|
||||
) const
|
||||
{
|
||||
if (this->solution().sourceActive())
|
||||
if (this->solution().coupled())
|
||||
{
|
||||
if (this->solution().semiImplicit("Yi"))
|
||||
{
|
||||
@ -165,7 +165,7 @@ Foam::ReactingCloud<ParcelType>::Srho(const label i) const
|
||||
)
|
||||
);
|
||||
|
||||
if (this->solution().sourceActive())
|
||||
if (this->solution().coupled())
|
||||
{
|
||||
scalarField& rhoi = tRhoi();
|
||||
rhoi = rhoTrans_[i]/(this->db().time().deltaT()*this->mesh().V());
|
||||
@ -202,7 +202,7 @@ Foam::ReactingCloud<ParcelType>::Srho() const
|
||||
)
|
||||
);
|
||||
|
||||
if (this->solution().sourceActive())
|
||||
if (this->solution().coupled())
|
||||
{
|
||||
scalarField& sourceField = trhoTrans();
|
||||
forAll(rhoTrans_, i)
|
||||
@ -221,7 +221,7 @@ template<class ParcelType>
|
||||
inline Foam::tmp<Foam::fvScalarMatrix>
|
||||
Foam::ReactingCloud<ParcelType>::Srho(volScalarField& rho) const
|
||||
{
|
||||
if (this->solution().sourceActive())
|
||||
if (this->solution().coupled())
|
||||
{
|
||||
tmp<volScalarField> trhoTrans
|
||||
(
|
||||
|
||||
@ -220,7 +220,7 @@ void Foam::ReactingMultiphaseCloud<ParcelType>::resetSourceTerms()
|
||||
template<class ParcelType>
|
||||
void Foam::ReactingMultiphaseCloud<ParcelType>::evolve()
|
||||
{
|
||||
if (this->solution().active())
|
||||
if (this->solution().canEvolve())
|
||||
{
|
||||
typename ParcelType::trackData td(*this);
|
||||
|
||||
|
||||
@ -286,7 +286,7 @@ void Foam::ThermoCloud<ParcelType>::relaxSources
|
||||
template<class ParcelType>
|
||||
void Foam::ThermoCloud<ParcelType>::evolve()
|
||||
{
|
||||
if (this->solution().active())
|
||||
if (this->solution().canEvolve())
|
||||
{
|
||||
typename ParcelType::trackData td(*this);
|
||||
|
||||
|
||||
@ -133,7 +133,7 @@ Foam::ThermoCloud<ParcelType>::Sh(volScalarField& hs) const
|
||||
<< max(hsCoeff()).value() << endl;
|
||||
}
|
||||
|
||||
if (this->solution().sourceActive())
|
||||
if (this->solution().coupled())
|
||||
{
|
||||
if (this->solution().semiImplicit("hs"))
|
||||
{
|
||||
|
||||
@ -105,7 +105,7 @@ void Foam::ParticleTracks<CloudType>::postFace(const parcelType& p)
|
||||
{
|
||||
if
|
||||
(
|
||||
this->owner().solution().writeThisStep()
|
||||
this->owner().solution().output()
|
||||
|| this->owner().solution().transient()
|
||||
)
|
||||
{
|
||||
|
||||
Reference in New Issue
Block a user