28 #ifndef TOMOGRAPHER_MHRW_BIN_ERR_H 29 #define TOMOGRAPHER_MHRW_BIN_ERR_H 37 #include <boost/math/constants/constants.hpp> 51 namespace tomo_internal {
55 template<
int NumLevels,
bool calculate>
56 struct helper_samples_size {
58 value = Eigen::Dynamic
61 template<
int NumLevels>
62 struct helper_samples_size<NumLevels,true> {
64 value = (1 << NumLevels)
108 template<
typename ValueType_,
int NumTrackValues_ = Eigen::Dynamic,
int NumLevels_ = Eigen::Dynamic,
109 bool StoreBinSums_ =
true,
typename CountIntType_ =
int>
122 static constexpr
int NumTrackValuesCTime = NumTrackValues_;
126 static constexpr
int NumLevelsCTime = NumLevels_;
130 static constexpr
int NumLevelsPlusOneCTime = (NumLevelsCTime == Eigen::Dynamic
132 : (NumLevelsCTime + 1));
155 static constexpr
int SamplesSizeCTime =
163 static constexpr
bool StoreBinSums = StoreBinSums_;
193 UNKNOWN_CONVERGENCE = 0,
258 template<
typename Params,
typename LoggerType_>
274 static constexpr
int NumTrackValuesCTime = Params::NumTrackValuesCTime;
278 static constexpr
int NumLevelsCTime = Params::NumLevelsCTime;
282 static constexpr
int NumLevelsPlusOneCTime = Params::NumLevelsPlusOneCTime;
286 static constexpr
int SamplesSizeCTime = Params::SamplesSizeCTime;
290 static constexpr
bool StoreBinSums = Params::StoreBinSums;
341 UNKNOWN_CONVERGENCE = Params::UNKNOWN_CONVERGENCE,
343 CONVERGED = Params::CONVERGED,
345 NOT_CONVERGED = Params::NOT_CONVERGED
361 SamplesArray samples;
369 CountIntType n_samples;
375 CountIntType n_flushes;
396 BinSumSqArray bin_sumsq;
416 : numTrackValues(num_track_values_),
417 numLevels(num_levels_),
418 samplesSize(1 << numLevels()),
419 samples(numTrackValues(), samplesSize()),
421 bin_sum(BinSumArray::Zero(numTrackValues())),
422 bin_sumsq(BinSumSqArray::Zero(numTrackValues(), numLevels()+1)),
441 helper_reset_bin_sum();
442 bin_sumsq = BinSumSqArray::Zero(numTrackValues(), numLevels()+1);
443 logger.longdebug(
"BinningAnalysis::reset()",
"ready to go.");
456 template<
typename Derived>
459 const int ninbin = n_samples % samplesSize();
467 samples.col(ninbin) = vals;
470 helper_update_bin_sum(samples.col(ninbin));
473 if ( ninbin == samplesSize() - 1 ) {
477 logger.longdebug(
"BinningAnalysis::processNewValues()", [&](
std::ostream & str) {
478 str <<
"n_samples is now " << n_samples <<
"; flushing bins. samplesSize() = " << samplesSize();
484 for (
int level = 0; level <= numLevels(); ++level) {
486 const int binnedsize = 1 << (numLevels()-level);
488 logger.longdebug(
"BinningAnalysis::processNewValues()", [&](
std::ostream & str) {
489 str <<
"Processing binning level = " << level <<
": binnedsize="<<binnedsize
490 <<
"; n_flushes=" << n_flushes <<
"\n";
491 str <<
"\tbinned samples = \n" << samples.block(0,0,numTrackValues(),binnedsize);
494 for (
int ksample = 0; ksample < binnedsize; ++ksample) {
495 bin_sumsq.
col(level) += samples.col(ksample).cwiseProduct(samples.col(ksample));
496 if (ksample % 2 == 0 && binnedsize > 1) {
497 samples.col(ksample/2) = boost::math::constants::half<ValueType>() *
498 (samples.col(ksample) + samples.col(ksample+1));
504 logger.longdebug(
"BinningAnalysis::processNewValues()", [&](
std::ostream & str) {
505 str <<
"Flushing #" << n_flushes <<
" done. bin_sum is = \n" << bin_sum <<
"\n" 506 <<
"\tbin_sumsq is = \n" << bin_sumsq <<
"\n";
519 template<
typename CalcValType, TOMOGRAPHER_ENABLED_IF_TMPL(NumTrackValuesCTime == 1)>
531 TOMOGRAPHER_ENABLED_IF(StoreBinSums)
532 inline void helper_reset_bin_sum()
534 bin_sum.value = BinSumArray::Zero(numTrackValues());
536 TOMOGRAPHER_ENABLED_IF(!StoreBinSums)
537 inline void helper_reset_bin_sum() { }
539 TOMOGRAPHER_ENABLED_IF(StoreBinSums)
543 bin_sum.value += new_samples;
545 template<
typename Derived,
546 TOMOGRAPHER_ENABLED_IF_TMPL(!StoreBinSums)>
564 TOMOGRAPHER_ENABLED_IF(StoreBinSums)
565 inline auto getBinMeans() const
566 #ifndef TOMOGRAPHER_PARSED_BY_DOXYGEN 567 -> decltype(BinSumArray() / ValueType(n_samples))
570 return bin_sum.value / ValueType(n_samples);
580 #ifndef TOMOGRAPHER_PARSED_BY_DOXYGEN 582 bin_sumsq.cwiseQuotient(n_flushes * Tools::replicated<NumTrackValuesCTime,1>(
584 .transpose().reverse(),
591 return bin_sumsq.cwiseQuotient(n_flushes * Tools::replicated<NumTrackValuesCTime,1>(
593 .transpose().reverse(),
604 TOMOGRAPHER_ENABLED_IF(StoreBinSums)
605 inline const BinSumArray & getBinSum()
const {
return bin_sum.value; }
628 template<
typename Derived>
633 const int n_levels_plus_one = numLevels()+1;
634 const int n_track_values = numTrackValues();
641 getBinSqmeans() - Tools::replicated<1,NumLevelsPlusOneCTime>(
646 ).cwiseMax(0).cwiseQuotient(
648 Tools::replicated<NumTrackValuesCTime,1>(
650 .transpose().reverse(),
655 numTrackValues(), n_levels_plus_one,
673 template<
typename Derived>
678 bin_sumsq.col(numLevels()) / ValueType(n_flushes) - means.
cwiseProduct(means).template cast<ValueType>()
679 ).cwiseMax(0).cwiseSqrt() /
std::sqrt(ValueType(n_flushes-1));
693 TOMOGRAPHER_ENABLED_IF(StoreBinSums)
694 inline BinSumSqArray calcErrorLevels()
const {
695 BinSumArray means = getBinMeans();
696 return calcErrorLevels(means);
711 TOMOGRAPHER_ENABLED_IF(StoreBinSums)
712 inline BinSumArray calcErrorLastLevel()
const {
713 BinSumArray means = getBinMeans();
714 return calcErrorLastLevel(means);
731 Eigen::ArrayXi converged_status(numTrackValues());
736 logger.longdebug(
"BinningAnalysis::determineErrorConvergence", [&](
std::ostream & str) {
737 str <<
"error_levels = \n" << error_levels <<
"\n";
744 if (numLevels() < range-1) {
746 converged_status = Eigen::ArrayXi::Constant(numTrackValues(), UNKNOWN_CONVERGENCE);
750 converged_status = Eigen::ArrayXi::Constant(numTrackValues(), CONVERGED);
752 const auto & errors = error_levels.col(numLevels());
754 for (
int level = numLevels()+1 - range; level < numLevels(); ++level) {
756 const auto & errors_thislevel = error_levels.col(level);
758 logger.longdebug(
"BinningAnalysis::determineErrorConvergence", [&](
std::ostream & str) {
759 str <<
"About to study level " << level <<
": at this point, converged_status = \n" 760 << converged_status <<
"\nand errors_thislevel = \n" << errors_thislevel;
763 for (
int val_it = 0; val_it < numTrackValues(); ++val_it) {
764 if (errors_thislevel(val_it) >= errors(val_it) &&
765 converged_status(val_it) != NOT_CONVERGED) {
766 converged_status(val_it) = CONVERGED;
767 }
else if (errors_thislevel(val_it) < 0.824 * errors(val_it)) {
768 converged_status(val_it) = NOT_CONVERGED;
769 }
else if ((errors_thislevel(val_it) < 0.9 * errors(val_it)) &&
770 converged_status(val_it) != NOT_CONVERGED) {
771 converged_status(val_it) = UNKNOWN_CONVERGENCE;
779 logger.longdebug(
"BinningAnalysis::determineErrorConvergence", [&](
std::ostream & str) {
780 str <<
"Done. converged_status [0=UNNOWN,1=CONVERGED,2=NOT CONVERGED] = \n" 784 return converged_status;
void processNewValue(const CalcValType val)
Process a new value (if we're tracking a single function only)
Params::BinSumSqArray BinSumSqArray
Type used to store the sum of squares of values at each binning level. See BinningAnalysisParams::Bin...
auto getBinSqmeans() const
Get the raw average of the squared values observed, for each binning level.
Base namespace for the Tomographer project.
Group template parameters for BinningAnalysis.
Eigen::ArrayXi determineErrorConvergence(const Eigen::Ref< const BinSumSqArray > &error_levels) const
Attempt to determine if the error bars have converged.
Eigen::Array< ValueType, NumTrackValuesCTime, 1 > BinSumArray
Type used to store the sum of values.
Params::CountIntType CountIntType
Type used to count the number of samples. See BinningAnalysisParams::CountIntType.
void processNewValues(const Eigen::DenseBase< Derived > &vals)
Process new raw samples.
const Tools::StaticOrDynamic< int,(NumLevelsCTime==Eigen::Dynamic), NumLevelsCTime > numLevels
The number of levels in the binning analysis.
Params::ValueType ValueType
Type of the value(s) for which we are calculating error bars. See BinningAnalysisParams::ValueType.
CountIntType getNumFlushes() const
Return the number of times the collected samples were flushed.
const CwiseBinaryOp< internal::scalar_product_op< typename Derived::Scalar, typename OtherDerived::Scalar >, const Derived, const OtherDerived > cwiseProduct(const Eigen::ArrayBase< OtherDerived > &other) const
ValueType_ ValueType
Type of the value which we are calculating error bars for. Also the type of the error bars themselves...
void reset()
Reset this object and start again as if freshly constructed.
BinSumSqArray calcErrorLevels(const Eigen::ArrayBase< Derived > &means) const
Calculate the error bars of samples at different binning levels.
Basic utilities for dealing with Eigen matrices and other types.
The error bars appear to have converged.
Params::SamplesArray SamplesArray
Type of the internal buffer for the raw samples. See BinningAnalysisParams::SamplesArray.
Params::BinSumArray BinSumArray
Type used to store the sum of values. See BinningAnalysisParams::BinSumArray.
Some C++ utilities, with a tad of C++11 tricks.
Simple binning analysis for determining error bars.
LoggerType_ LoggerType
Type of the logger we will be logging debugging messages to.
BinningAnalysis(int num_track_values_, int num_levels_, LoggerType &logger_)
Constructor.
Managing the need for specific overrides to operator new() for some types (especially Eigen types) ...
Eigen::Array< ValueType, NumTrackValuesCTime, NumLevelsPlusOneCTime > BinSumSqArray
Type used to store the sum of squares of values at each binning level.
VarValueDecoder< T >::RetType value(const Var &var)
Access the value of the given variable, as a C++ type.
const BinSumSqArray & getBinSumsq() const
Get the raw sums of the squared values observed, at each binning level.
const Tools::StaticOrDynamic< CountIntType,(SamplesSizeCTime==Eigen::Dynamic), SamplesSizeCTime > samplesSize
The size of our samples buffer. (avoid using, might change in the future)
Eigen::Array< ValueType, NumTrackValuesCTime, SamplesSizeCTime > SamplesArray
Type of the internal buffer for the raw samples.
BinSumArray calcErrorLastLevel(const Eigen::ArrayBase< Derived > &means) const
Calculate the error bar of samples (from the last binning level).
CountIntType_ CountIntType
Type used to count the number of samples. Usually int is fine, except if you are really taking a LOT ...
#define tomographer_assert(...)
Assertion test macro.
const Tools::StaticOrDynamic< int,(NumTrackValuesCTime==Eigen::Dynamic), NumTrackValuesCTime > numTrackValues
The number of functions being tracked/analyzed.
Utilities for logging messages.