Tomographer  v2.0
Tomographer C++ Framework Documentation
mhrw_bin_err.h
Go to the documentation of this file.
1 /* This file is part of the Tomographer project, which is distributed under the
2  * terms of the MIT license.
3  *
4  * The MIT License (MIT)
5  *
6  * Copyright (c) 2015 ETH Zurich, Institute for Theoretical Physics, Philippe Faist
7  *
8  * Permission is hereby granted, free of charge, to any person obtaining a copy
9  * of this software and associated documentation files (the "Software"), to deal
10  * in the Software without restriction, including without limitation the rights
11  * to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
12  * copies of the Software, and to permit persons to whom the Software is
13  * furnished to do so, subject to the following conditions:
14  *
15  * The above copyright notice and this permission notice shall be included in
16  * all copies or substantial portions of the Software.
17  *
18  * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
19  * IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
20  * FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
21  * AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
22  * LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
23  * OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE
24  * SOFTWARE.
25  */
26 
27 #ifndef TOMOGRAPHER_MHRW_BIN_ERR_H
28 #define TOMOGRAPHER_MHRW_BIN_ERR_H
29 
30 
32 #include <tomographer2/tools/eigenutil.h> // replicated(), powersOfTwo()
35 
36 #include <boost/math/constants/constants.hpp>
37 
38 
48 namespace Tomographer {
49 
50 namespace tomo_internal {
51 //
52 // internal helper to silence "left shift is negative" GCC warnings
53 //
54 template<int NumLevels, bool calculate>
55 struct helper_samples_size {
56  enum {
57  value = Eigen::Dynamic
58  };
59 };
60 template<int NumLevels>
61 struct helper_samples_size<NumLevels,true> {
62  enum {
63  value = (1 << NumLevels)
64  };
65 };
66 } // namespace tomo_internal
67 
68 
69 
107 template<typename ValueType_, int NumTrackValues_ = Eigen::Dynamic, int NumLevels_ = Eigen::Dynamic,
108  bool StoreBinSums_ = true, typename CountIntType_ = int>
110 {
113  typedef ValueType_ ValueType;
116  typedef CountIntType_ CountIntType;
117 
121  static constexpr int NumTrackValuesCTime = NumTrackValues_;
125  static constexpr int NumLevelsCTime = NumLevels_;
129  static constexpr int NumLevelsPlusOneCTime = (NumLevelsCTime == Eigen::Dynamic
130  ? Eigen::Dynamic
131  : (NumLevelsCTime + 1));
154  static constexpr int SamplesSizeCTime =
155  tomo_internal::helper_samples_size<NumLevelsCTime, (NumLevelsCTime > 0 && NumLevelsCTime < 7)>::value;
156 
162  static constexpr bool StoreBinSums = StoreBinSums_;
163 
172 
188 
190  enum {
192  UNKNOWN_CONVERGENCE = 0,
196  NOT_CONVERGED
197  };
198 
199 };
200 
257 template<typename Params, typename LoggerType_>
259  // inheriting from this has some advantages over EIGEN_MAKE_ALIGNED_OPERATOR_NEW, such
260  // as not needing to explicitly declare the specialization
261  // NeedOwnOperatorNew<BinningAnalysis<...> >
263 {
264 public:
266  typedef typename Params::ValueType ValueType;
268  typedef typename Params::CountIntType CountIntType;
269 
273  static constexpr int NumTrackValuesCTime = Params::NumTrackValuesCTime;
277  static constexpr int NumLevelsCTime = Params::NumLevelsCTime;
281  static constexpr int NumLevelsPlusOneCTime = Params::NumLevelsPlusOneCTime;
285  static constexpr int SamplesSizeCTime = Params::SamplesSizeCTime;
289  static constexpr bool StoreBinSums = Params::StoreBinSums;
290 
294  typedef typename Params::SamplesArray SamplesArray;
298  typedef typename Params::BinSumArray BinSumArray;
302  typedef typename Params::BinSumSqArray BinSumSqArray;
303 
336 
338  enum {
340  UNKNOWN_CONVERGENCE = Params::UNKNOWN_CONVERGENCE,
342  CONVERGED = Params::CONVERGED,
344  NOT_CONVERGED = Params::NOT_CONVERGED
345  };
346 
348  typedef LoggerType_ LoggerType;
349 
350 private:
351 
360  SamplesArray samples;
361 
362  // where we store the flushed values
363 
368  CountIntType n_samples;
374  CountIntType n_flushes;
395  BinSumSqArray bin_sumsq;
396 
398  LoggerType & logger;
399 
400 
401 public:
402 
414  BinningAnalysis(int num_track_values_, int num_levels_, LoggerType & logger_)
415  : numTrackValues(num_track_values_),
416  numLevels(num_levels_),
417  samplesSize(1 << numLevels()),
418  samples(numTrackValues(), samplesSize()),
419  n_flushes(0),
420  bin_sum(BinSumArray::Zero(numTrackValues())),
421  bin_sumsq(BinSumSqArray::Zero(numTrackValues(), numLevels()+1)),
422  logger(logger_)
423  {
424  tomographer_assert(Tools::isPositive(numLevels()));
425  tomographer_assert(Tools::isPowerOfTwo(samplesSize()));
426  tomographer_assert( (1<<numLevels()) == samplesSize() );
427 
428  reset();
429  }
430 
436  inline void reset()
437  {
438  n_flushes = 0;
439  n_samples = 0;
440  helper_reset_bin_sum();
441  bin_sumsq = BinSumSqArray::Zero(numTrackValues(), numLevels()+1);
442  logger.longdebug("BinningAnalysis::reset()", "ready to go.");
443  }
444 
455  template<typename Derived>
456  inline void processNewValues(const Eigen::DenseBase<Derived> & vals)
457  {
458  const int ninbin = n_samples % samplesSize();
459 
460  ++n_samples;
461 
462  tomographer_assert(vals.rows() == numTrackValues());
463  tomographer_assert(vals.cols() == 1);
464 
465  // store the new values in the bins [also if ninbin == 0]
466  samples.col(ninbin) = vals;
467 
468  // add to our sum of values, if applicable.
469  helper_update_bin_sum(samples.col(ninbin));
470 
471  // see if we have to flush the bins (equivalent to `ninbin == samples_size()-1`)
472  if ( ninbin == samplesSize() - 1 ) {
473 
474  // we have filled all bins. Flush them. Re-use the beginning of the samples[] array
475  // to store the reduced bins while flushing them.
476  logger.longdebug("BinningAnalysis::processNewValues()", [&](std::ostream & str) {
477  str << "n_samples is now " << n_samples << "; flushing bins. samplesSize() = " << samplesSize();
478  });
479 
480  // the size of the samples at the current level of binning. Starts at samplesSize,
481  // and decreases by half at each higher level.
482 
483  for (int level = 0; level <= numLevels(); ++level) {
484 
485  const int binnedsize = 1 << (numLevels()-level);
486 
487  logger.longdebug("BinningAnalysis::processNewValues()", [&](std::ostream & str) {
488  str << "Processing binning level = " << level << ": binnedsize="<<binnedsize
489  << "; n_flushes=" << n_flushes << "\n";
490  str << "\tbinned samples = \n" << samples.block(0,0,numTrackValues(),binnedsize);
491  });
492 
493  for (int ksample = 0; ksample < binnedsize; ++ksample) {
494  bin_sumsq.col(level) += samples.col(ksample).cwiseProduct(samples.col(ksample));
495  if (ksample % 2 == 0 && binnedsize > 1) {
496  samples.col(ksample/2) = boost::math::constants::half<ValueType>() *
497  (samples.col(ksample) + samples.col(ksample+1));
498  }
499  }
500 
501  }
502 
503  logger.longdebug("BinningAnalysis::processNewValues()", [&](std::ostream & str) {
504  str << "Flushing #" << n_flushes << " done. bin_sum is = \n" << bin_sum << "\n"
505  << "\tbin_sumsq is = \n" << bin_sumsq << "\n";
506  });
507 
508  ++n_flushes;
509  }
510 
511  }
512 
518  template<typename CalcValType, TOMOGRAPHER_ENABLED_IF_TMPL(NumTrackValuesCTime == 1)>
519  inline void processNewValue(const CalcValType val)
520  {
521  // for a single value
522  //processNewValues(Eigen::Map<const Eigen::Array<CalcValType,1,1> >(&val));
523  processNewValues(Eigen::Array<CalcValType,1,1>::Constant(val));
524  }
525 
526 
527 
528 private:
529 
530  TOMOGRAPHER_ENABLED_IF(StoreBinSums)
531  inline void helper_reset_bin_sum()
532  {
533  bin_sum.value = BinSumArray::Zero(numTrackValues());
534  }
535  TOMOGRAPHER_ENABLED_IF(!StoreBinSums)
536  inline void helper_reset_bin_sum() { }
537 
538  TOMOGRAPHER_ENABLED_IF(StoreBinSums)
539  inline void helper_update_bin_sum(const Eigen::Ref<const Eigen::Array<ValueType, NumTrackValuesCTime, 1> > &
540  new_samples)
541  {
542  bin_sum.value += new_samples;
543  }
544  template<typename Derived,
545  TOMOGRAPHER_ENABLED_IF_TMPL(!StoreBinSums)>
546  inline void helper_update_bin_sum(const Eigen::DenseBase<Derived> & ) { }
547 
548 
549 public:
550 
551  // retrieve results.
552 
558  inline CountIntType getNumFlushes() const { return n_flushes; }
559 
563  TOMOGRAPHER_ENABLED_IF(StoreBinSums)
564  inline auto getBinMeans() const
565 #ifndef TOMOGRAPHER_PARSED_BY_DOXYGEN
566  -> decltype(BinSumArray() / ValueType(n_samples))
567 #endif
568  {
569  return bin_sum.value / ValueType(n_samples);
570  }
571 
578  inline auto getBinSqmeans() const
579 #ifndef TOMOGRAPHER_PARSED_BY_DOXYGEN
580  -> decltype(
581  bin_sumsq.cwiseQuotient(n_flushes * Tools::replicated<NumTrackValuesCTime,1>(
583  .transpose().reverse(),
584  // replicated by:
585  numTrackValues(), 1
586  ))
587  )
588 #endif
589  {
590  return bin_sumsq.cwiseQuotient(n_flushes * Tools::replicated<NumTrackValuesCTime,1>(
592  .transpose().reverse(),
593  // replicated by:
594  numTrackValues(), 1
595  ));
596  }
597 
598 
603  TOMOGRAPHER_ENABLED_IF(StoreBinSums)
604  inline const BinSumArray & getBinSum() const { return bin_sum.value; }
605 
612  inline const BinSumSqArray & getBinSumsq() const {
613  return bin_sumsq;
614  }
615 
616 
627  template<typename Derived>
628  inline BinSumSqArray calcErrorLevels(const Eigen::ArrayBase<Derived> & means) const
629  {
630  tomographer_assert(means.rows() == numTrackValues());
631  tomographer_assert(means.cols() == 1);
632  const int n_levels_plus_one = numLevels()+1;
633  const int n_track_values = numTrackValues();
634 
639  return (
640  getBinSqmeans() - Tools::replicated<1,NumLevelsPlusOneCTime>(
641  means.cwiseProduct(means).template cast<ValueType>(),
642  // replicated by:
643  1, n_levels_plus_one
644  )
645  ).cwiseMax(0).cwiseQuotient(
646  // divide by the number of samples from which these bin-means were obtained, minus one.
647  Tools::replicated<NumTrackValuesCTime,1>(
649  .transpose().reverse(),
650  // replicated by:
651  n_track_values, 1
652  ) * n_flushes
654  numTrackValues(), n_levels_plus_one,
655  1 // the constant...
656  )
657  ).cwiseSqrt();
658  }
659 
660 
672  template<typename Derived>
673  inline BinSumArray calcErrorLastLevel(const Eigen::ArrayBase<Derived> & means) const {
674  tomographer_assert(means.rows() == numTrackValues());
675  tomographer_assert(means.cols() == 1);
676  return (
677  bin_sumsq.col(numLevels()) / ValueType(n_flushes) - means.cwiseProduct(means).template cast<ValueType>()
678  ).cwiseMax(0).cwiseSqrt() / std::sqrt(ValueType(n_flushes-1));
679  }
680 
692  TOMOGRAPHER_ENABLED_IF(StoreBinSums)
693  inline BinSumSqArray calcErrorLevels() const {
694  BinSumArray means = getBinMeans();
695  return calcErrorLevels(means);
696  }
697 
710  TOMOGRAPHER_ENABLED_IF(StoreBinSums)
711  inline BinSumArray calcErrorLastLevel() const {
712  BinSumArray means = getBinMeans();
713  return calcErrorLastLevel(means);
714  }
715 
728  inline Eigen::ArrayXi determineErrorConvergence(const Eigen::Ref<const BinSumSqArray> & error_levels) const
729  {
730  Eigen::ArrayXi converged_status(numTrackValues()); // RVO will help
731 
732  tomographer_assert(error_levels.rows() == numTrackValues());
733  tomographer_assert(error_levels.cols() == numLevels() + 1);
734 
735  logger.longdebug("BinningAnalysis::determineErrorConvergence", [&](std::ostream & str) {
736  str << "error_levels = \n" << error_levels << "\n";
737  });
738 
739  // verify that indeed the errors have converged. Inspired from ALPS code, see
740  // https://alps.comp-phys.org/svn/alps1/trunk/alps/src/alps/alea/simplebinning.h
741 
742  const int range = 4;
743  if (numLevels() < range-1) {
744 
745  converged_status = Eigen::ArrayXi::Constant(numTrackValues(), UNKNOWN_CONVERGENCE);
746 
747  } else {
748 
749  converged_status = Eigen::ArrayXi::Constant(numTrackValues(), CONVERGED);
750 
751  const auto & errors = error_levels.col(numLevels());
752 
753  for (int level = numLevels()+1 - range; level < numLevels(); ++level) {
754 
755  const auto & errors_thislevel = error_levels.col(level);
756 
757  logger.longdebug("BinningAnalysis::determineErrorConvergence", [&](std::ostream & str) {
758  str << "About to study level " << level << ": at this point, converged_status = \n"
759  << converged_status << "\nand errors_thislevel = \n" << errors_thislevel;
760  });
761 
762  for (int val_it = 0; val_it < numTrackValues(); ++val_it) {
763  if (errors_thislevel(val_it) >= errors(val_it) &&
764  converged_status(val_it) != NOT_CONVERGED) {
765  converged_status(val_it) = CONVERGED;
766  } else if (errors_thislevel(val_it) < 0.824 * errors(val_it)) {
767  converged_status(val_it) = NOT_CONVERGED;
768  } else if ((errors_thislevel(val_it) < 0.9 * errors(val_it)) &&
769  converged_status(val_it) != NOT_CONVERGED) {
770  converged_status(val_it) = UNKNOWN_CONVERGENCE;
771  }
772  }
773 
774  }
775 
776  }
777 
778  logger.longdebug("BinningAnalysis::determineErrorConvergence", [&](std::ostream & str) {
779  str << "Done. converged_status [0=UNNOWN,1=CONVERGED,2=NOT CONVERGED] = \n"
780  << converged_status;
781  });
782 
783  return converged_status;
784  }
785 };
786 
787 
788 
789 } // namespace Tomographer
790 
791 
792 
793 
794 
795 
796 
797 
798 
799 
800 
801 
802 
803 
804 
805 #endif
void processNewValue(const CalcValType val)
Process a new value (if we&#39;re tracking a single function only)
Definition: mhrw_bin_err.h:519
Params::BinSumSqArray BinSumSqArray
Type used to store the sum of squares of values at each binning level. See BinningAnalysisParams::Bin...
Definition: mhrw_bin_err.h:302
Provides correct operator-new implementation for Eigen types via the NeedOwnOperatorNew mechanism...
Base namespace for the Tomographer project.
Definition: densellh.h:44
Group template parameters for BinningAnalysis.
Definition: mhrw_bin_err.h:109
BinSumSqArray calcErrorLevels(const Eigen::ArrayBase< Derived > &means) const
Calculate the error bars of samples at different binning levels.
Definition: mhrw_bin_err.h:628
Eigen::Array< ValueType, NumTrackValuesCTime, 1 > BinSumArray
Type used to store the sum of values.
Definition: mhrw_bin_err.h:179
Params::CountIntType CountIntType
Type used to count the number of samples. See BinningAnalysisParams::CountIntType.
Definition: mhrw_bin_err.h:268
Eigen::ArrayXi determineErrorConvergence(const Eigen::Ref< const BinSumSqArray > &error_levels) const
Attempt to determine if the error bars have converged.
Definition: mhrw_bin_err.h:728
void processNewValues(const Eigen::DenseBase< Derived > &vals)
Process new raw samples.
Definition: mhrw_bin_err.h:456
const Tools::StaticOrDynamic< int,(NumLevelsCTime==Eigen::Dynamic), NumLevelsCTime > numLevels
The number of levels in the binning analysis.
Definition: mhrw_bin_err.h:325
auto powersOfTwo(IndexTypes...sizes) -> const Eigen::CwiseNullaryOp< tomo_internal::powers_of_two_generator< typename Eigen::internal::traits< Der >::Scalar >, Der >
Expression for a 1-D expression of powers of two.
Definition: eigenutil.h:305
Params::ValueType ValueType
Type of the value(s) for which we are calculating error bars. See BinningAnalysisParams::ValueType.
Definition: mhrw_bin_err.h:266
constexpr bool isPowerOfTwo(IntType N)
Return true if the argument is a power of two, false otherwise.
Definition: cxxutil.h:325
ColXpr col(Index i)
const CwiseBinaryOp< internal::scalar_product_op< typename Derived::Scalar, typename OtherDerived::Scalar >, const Derived, const OtherDerived > cwiseProduct(const Eigen::ArrayBase< OtherDerived > &other) const
BinSumArray calcErrorLastLevel(const Eigen::ArrayBase< Derived > &means) const
Calculate the error bar of samples (from the last binning level).
Definition: mhrw_bin_err.h:673
ValueType_ ValueType
Type of the value which we are calculating error bars for. Also the type of the error bars themselves...
Definition: mhrw_bin_err.h:113
void reset()
Reset this object and start again as if freshly constructed.
Definition: mhrw_bin_err.h:436
Basic utilities for dealing with Eigen matrices and other types.
std::enable_if< std::is_unsigned< X >::value, bool >::type isPositive(const X)
Test whether the given value is positive or zero.
Definition: cxxutil.h:376
The error bars appear to have converged.
Definition: mhrw_bin_err.h:194
Params::SamplesArray SamplesArray
Type of the internal buffer for the raw samples. See BinningAnalysisParams::SamplesArray.
Definition: mhrw_bin_err.h:294
Params::BinSumArray BinSumArray
Type used to store the sum of values. See BinningAnalysisParams::BinSumArray.
Definition: mhrw_bin_err.h:298
Some C++ utilities, with a tad of C++11 tricks.
Simple binning analysis for determining error bars.
Definition: mhrw_bin_err.h:258
LoggerType_ LoggerType
Type of the logger we will be logging debugging messages to.
Definition: mhrw_bin_err.h:348
BinningAnalysis(int num_track_values_, int num_levels_, LoggerType &logger_)
Constructor.
Definition: mhrw_bin_err.h:414
Managing the need for specific overrides to operator new() for some types (especially Eigen types) ...
const BinSumSqArray & getBinSumsq() const
Get the raw sums of the squared values observed, at each binning level.
Definition: mhrw_bin_err.h:612
Eigen::Array< ValueType, NumTrackValuesCTime, NumLevelsPlusOneCTime > BinSumSqArray
Type used to store the sum of squares of values at each binning level.
Definition: mhrw_bin_err.h:187
VarValueDecoder< T >::RetType value(const Var &var)
Access the value of the given variable, as a C++ type.
Definition: ezmatio.h:854
const Tools::StaticOrDynamic< CountIntType,(SamplesSizeCTime==Eigen::Dynamic), SamplesSizeCTime > samplesSize
The size of our samples buffer. (avoid using, might change in the future)
Definition: mhrw_bin_err.h:335
T sqrt(T...args)
Eigen::Array< ValueType, NumTrackValuesCTime, SamplesSizeCTime > SamplesArray
Type of the internal buffer for the raw samples.
Definition: mhrw_bin_err.h:171
STL class.
CountIntType_ CountIntType
Type used to count the number of samples. Usually int is fine, except if you are really taking a LOT ...
Definition: mhrw_bin_err.h:116
CountIntType getNumFlushes() const
Return the number of times the collected samples were flushed.
Definition: mhrw_bin_err.h:558
const Tools::StaticOrDynamic< int,(NumTrackValuesCTime==Eigen::Dynamic), NumTrackValuesCTime > numTrackValues
The number of functions being tracked/analyzed.
Definition: mhrw_bin_err.h:314
Utilities for logging messages.
auto getBinSqmeans() const
Get the raw average of the squared values observed, for each binning level.
Definition: mhrw_bin_err.h:578