Loading [MathJax]/extensions/tex2jax.js
Tomographerv4.1
Tomographer C++ Framework Documentation
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
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) 2016 ETH Zurich, Institute for Theoretical Physics, Philippe Faist
7  * Copyright (c) 2017 Caltech, Institute for Quantum Information and Matter, Philippe Faist
8  *
9  * Permission is hereby granted, free of charge, to any person obtaining a copy
10  * of this software and associated documentation files (the "Software"), to deal
11  * in the Software without restriction, including without limitation the rights
12  * to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
13  * copies of the Software, and to permit persons to whom the Software is
14  * furnished to do so, subject to the following conditions:
15  *
16  * The above copyright notice and this permission notice shall be included in
17  * all copies or substantial portions of the Software.
18  *
19  * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
20  * IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
21  * FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
22  * AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
23  * LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
24  * OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE
25  * SOFTWARE.
26  */
27 
28 #ifndef TOMOGRAPHER_MHRW_BIN_ERR_H
29 #define TOMOGRAPHER_MHRW_BIN_ERR_H
30 
31 
33 #include <tomographer/tools/eigenutil.h> // replicated(), powersOfTwo()
36 
37 #include <boost/math/constants/constants.hpp>
38 
39 
49 namespace Tomographer {
50 
51 namespace tomo_internal {
52 //
53 // internal helper to silence "left shift is negative" GCC warnings
54 //
55 template<int NumLevels, bool calculate>
56 struct helper_samples_size {
57  enum {
58  value = Eigen::Dynamic
59  };
60 };
61 template<int NumLevels>
62 struct helper_samples_size<NumLevels,true> {
63  enum {
64  value = (1 << NumLevels)
65  };
66 };
67 } // namespace tomo_internal
68 
69 
70 
108 template<typename ValueType_, int NumTrackValues_ = Eigen::Dynamic, int NumLevels_ = Eigen::Dynamic,
109  bool StoreBinSums_ = true, typename CountIntType_ = int>
110 TOMOGRAPHER_EXPORT struct BinningAnalysisParams
111 {
114  typedef ValueType_ ValueType;
117  typedef CountIntType_ CountIntType;
118 
122  static constexpr int NumTrackValuesCTime = NumTrackValues_;
126  static constexpr int NumLevelsCTime = NumLevels_;
130  static constexpr int NumLevelsPlusOneCTime = (NumLevelsCTime == Eigen::Dynamic
131  ? Eigen::Dynamic
132  : (NumLevelsCTime + 1));
155  static constexpr int SamplesSizeCTime =
156  tomo_internal::helper_samples_size<NumLevelsCTime, (NumLevelsCTime > 0 && NumLevelsCTime < 7)>::value;
157 
163  static constexpr bool StoreBinSums = StoreBinSums_;
164 
173 
189 
191  enum {
193  UNKNOWN_CONVERGENCE = 0,
197  NOT_CONVERGED
198  };
199 
200 };
201 
258 template<typename Params, typename LoggerType_>
259 TOMOGRAPHER_EXPORT class BinningAnalysis
260  // inheriting from this has some advantages over EIGEN_MAKE_ALIGNED_OPERATOR_NEW, such
261  // as not needing to explicitly declare the specialization
262  // NeedOwnOperatorNew<BinningAnalysis<...> >
264 {
265 public:
267  typedef typename Params::ValueType ValueType;
269  typedef typename Params::CountIntType CountIntType;
270 
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;
291 
295  typedef typename Params::SamplesArray SamplesArray;
299  typedef typename Params::BinSumArray BinSumArray;
303  typedef typename Params::BinSumSqArray BinSumSqArray;
304 
337 
339  enum {
341  UNKNOWN_CONVERGENCE = Params::UNKNOWN_CONVERGENCE,
343  CONVERGED = Params::CONVERGED,
345  NOT_CONVERGED = Params::NOT_CONVERGED
346  };
347 
349  typedef LoggerType_ LoggerType;
350 
351 private:
352 
361  SamplesArray samples;
362 
363  // where we store the flushed values
364 
369  CountIntType n_samples;
375  CountIntType n_flushes;
396  BinSumSqArray bin_sumsq;
397 
399  LoggerType & logger;
400 
401 
402 public:
403 
415  BinningAnalysis(int num_track_values_, int num_levels_, LoggerType & logger_)
416  : numTrackValues(num_track_values_),
417  numLevels(num_levels_),
418  samplesSize(1 << numLevels()),
419  samples(numTrackValues(), samplesSize()),
420  n_flushes(0),
421  bin_sum(BinSumArray::Zero(numTrackValues())),
422  bin_sumsq(BinSumSqArray::Zero(numTrackValues(), numLevels()+1)),
423  logger(logger_)
424  {
426  tomographer_assert(Tools::isPowerOfTwo(samplesSize()));
427  tomographer_assert( (1<<numLevels()) == samplesSize() );
428 
429  reset();
430  }
431 
437  inline void reset()
438  {
439  n_flushes = 0;
440  n_samples = 0;
441  helper_reset_bin_sum();
442  bin_sumsq = BinSumSqArray::Zero(numTrackValues(), numLevels()+1);
443  logger.longdebug("BinningAnalysis::reset()", "ready to go.");
444  }
445 
456  template<typename Derived>
457  inline void processNewValues(const Eigen::DenseBase<Derived> & vals)
458  {
459  const int ninbin = n_samples % samplesSize();
460 
461  ++n_samples;
462 
463  tomographer_assert(vals.rows() == numTrackValues());
464  tomographer_assert(vals.cols() == 1);
465 
466  // store the new values in the bins [also if ninbin == 0]
467  samples.col(ninbin) = vals;
468 
469  // add to our sum of values, if applicable.
470  helper_update_bin_sum(samples.col(ninbin));
471 
472  // see if we have to flush the bins (equivalent to `ninbin == samples_size()-1`)
473  if ( ninbin == samplesSize() - 1 ) {
474 
475  // we have filled all bins. Flush them. Re-use the beginning of the samples[] array
476  // to store the reduced bins while flushing them.
477  logger.longdebug("BinningAnalysis::processNewValues()", [&](std::ostream & str) {
478  str << "n_samples is now " << n_samples << "; flushing bins. samplesSize() = " << samplesSize();
479  });
480 
481  // the size of the samples at the current level of binning. Starts at samplesSize,
482  // and decreases by half at each higher level.
483 
484  for (int level = 0; level <= numLevels(); ++level) {
485 
486  const int binnedsize = 1 << (numLevels()-level);
487 
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);
492  });
493 
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));
499  }
500  }
501 
502  }
503 
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";
507  });
508 
509  ++n_flushes;
510  }
511 
512  }
513 
519  template<typename CalcValType, TOMOGRAPHER_ENABLED_IF_TMPL(NumTrackValuesCTime == 1)>
520  inline void processNewValue(const CalcValType val)
521  {
522  // for a single value
523  //processNewValues(Eigen::Map<const Eigen::Array<CalcValType,1,1> >(&val));
524  processNewValues(Eigen::Array<CalcValType,1,1>::Constant(val));
525  }
526 
527 
528 
529 private:
530 
531  TOMOGRAPHER_ENABLED_IF(StoreBinSums)
532  inline void helper_reset_bin_sum()
533  {
534  bin_sum.value = BinSumArray::Zero(numTrackValues());
535  }
536  TOMOGRAPHER_ENABLED_IF(!StoreBinSums)
537  inline void helper_reset_bin_sum() { }
538 
539  TOMOGRAPHER_ENABLED_IF(StoreBinSums)
540  inline void helper_update_bin_sum(const Eigen::Ref<const Eigen::Array<ValueType, NumTrackValuesCTime, 1> > &
541  new_samples)
542  {
543  bin_sum.value += new_samples;
544  }
545  template<typename Derived,
546  TOMOGRAPHER_ENABLED_IF_TMPL(!StoreBinSums)>
547  inline void helper_update_bin_sum(const Eigen::DenseBase<Derived> & ) { }
548 
549 
550 public:
551 
552  // retrieve results.
553 
559  inline CountIntType getNumFlushes() const { return n_flushes; }
560 
564  TOMOGRAPHER_ENABLED_IF(StoreBinSums)
565  inline auto getBinMeans() const
566 #ifndef TOMOGRAPHER_PARSED_BY_DOXYGEN
567  -> decltype(BinSumArray() / ValueType(n_samples))
568 #endif
569  {
570  return bin_sum.value / ValueType(n_samples);
571  }
572 
579  inline auto getBinSqmeans() const
580 #ifndef TOMOGRAPHER_PARSED_BY_DOXYGEN
581  -> decltype(
582  bin_sumsq.cwiseQuotient(n_flushes * Tools::replicated<NumTrackValuesCTime,1>(
584  .transpose().reverse(),
585  // replicated by:
586  numTrackValues(), 1
587  ))
588  )
589 #endif
590  {
591  return bin_sumsq.cwiseQuotient(n_flushes * Tools::replicated<NumTrackValuesCTime,1>(
593  .transpose().reverse(),
594  // replicated by:
595  numTrackValues(), 1
596  ));
597  }
598 
599 
604  TOMOGRAPHER_ENABLED_IF(StoreBinSums)
605  inline const BinSumArray & getBinSum() const { return bin_sum.value; }
606 
613  inline const BinSumSqArray & getBinSumsq() const {
614  return bin_sumsq;
615  }
616 
617 
628  template<typename Derived>
629  inline BinSumSqArray calcErrorLevels(const Eigen::ArrayBase<Derived> & means) const
630  {
631  tomographer_assert(means.rows() == numTrackValues());
632  tomographer_assert(means.cols() == 1);
633  const int n_levels_plus_one = numLevels()+1;
634  const int n_track_values = numTrackValues();
635 
640  return (
641  getBinSqmeans() - Tools::replicated<1,NumLevelsPlusOneCTime>(
642  means.cwiseProduct(means).template cast<ValueType>(),
643  // replicated by:
644  1, n_levels_plus_one
645  )
646  ).cwiseMax(0).cwiseQuotient(
647  // divide by the number of samples from which these bin-means were obtained, minus one.
648  Tools::replicated<NumTrackValuesCTime,1>(
650  .transpose().reverse(),
651  // replicated by:
652  n_track_values, 1
653  ) * n_flushes
655  numTrackValues(), n_levels_plus_one,
656  1 // the constant...
657  )
658  ).cwiseSqrt();
659  }
660 
661 
673  template<typename Derived>
674  inline BinSumArray calcErrorLastLevel(const Eigen::ArrayBase<Derived> & means) const {
675  tomographer_assert(means.rows() == numTrackValues());
676  tomographer_assert(means.cols() == 1);
677  return (
678  bin_sumsq.col(numLevels()) / ValueType(n_flushes) - means.cwiseProduct(means).template cast<ValueType>()
679  ).cwiseMax(0).cwiseSqrt() / std::sqrt(ValueType(n_flushes-1));
680  }
681 
693  TOMOGRAPHER_ENABLED_IF(StoreBinSums)
694  inline BinSumSqArray calcErrorLevels() const {
695  BinSumArray means = getBinMeans();
696  return calcErrorLevels(means);
697  }
698 
711  TOMOGRAPHER_ENABLED_IF(StoreBinSums)
712  inline BinSumArray calcErrorLastLevel() const {
713  BinSumArray means = getBinMeans();
714  return calcErrorLastLevel(means);
715  }
716 
729  inline Eigen::ArrayXi determineErrorConvergence(const Eigen::Ref<const BinSumSqArray> & error_levels) const
730  {
731  Eigen::ArrayXi converged_status(numTrackValues()); // RVO will help
732 
733  tomographer_assert(error_levels.rows() == numTrackValues());
734  tomographer_assert(error_levels.cols() == numLevels() + 1);
735 
736  logger.longdebug("BinningAnalysis::determineErrorConvergence", [&](std::ostream & str) {
737  str << "error_levels = \n" << error_levels << "\n";
738  });
739 
740  // verify that indeed the errors have converged. Inspired from ALPS code, see
741  // https://alps.comp-phys.org/svn/alps1/trunk/alps/src/alps/alea/simplebinning.h
742 
743  const int range = 4;
744  if (numLevels() < range-1) {
745 
746  converged_status = Eigen::ArrayXi::Constant(numTrackValues(), UNKNOWN_CONVERGENCE);
747 
748  } else {
749 
750  converged_status = Eigen::ArrayXi::Constant(numTrackValues(), CONVERGED);
751 
752  const auto & errors = error_levels.col(numLevels());
753 
754  for (int level = numLevels()+1 - range; level < numLevels(); ++level) {
755 
756  const auto & errors_thislevel = error_levels.col(level);
757 
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;
761  });
762 
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;
772  }
773  }
774 
775  }
776 
777  }
778 
779  logger.longdebug("BinningAnalysis::determineErrorConvergence", [&](std::ostream & str) {
780  str << "Done. converged_status [0=UNNOWN,1=CONVERGED,2=NOT CONVERGED] = \n"
781  << converged_status;
782  });
783 
784  return converged_status;
785  }
786 };
787 
788 
789 
790 } // namespace Tomographer
791 
792 
793 
794 
795 
796 
797 
798 
799 
800 
801 
802 
803 
804 
805 
806 #endif
void processNewValue(const CalcValType val)
Process a new value (if we&#39;re tracking a single function only)
Definition: mhrw_bin_err.h:520
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:303
auto getBinSqmeans() const
Get the raw average of the squared values observed, for each binning level.
Definition: mhrw_bin_err.h:579
Provides correct operator-new implementation for Eigen types via the NeedOwnOperatorNew mechanism...
Base namespace for the Tomographer project.
Definition: densellh.h:45
Group template parameters for BinningAnalysis.
Definition: mhrw_bin_err.h:110
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:729
Eigen::Array< ValueType, NumTrackValuesCTime, 1 > BinSumArray
Type used to store the sum of values.
Definition: mhrw_bin_err.h:180
Params::CountIntType CountIntType
Type used to count the number of samples. See BinningAnalysisParams::CountIntType.
Definition: mhrw_bin_err.h:269
void processNewValues(const Eigen::DenseBase< Derived > &vals)
Process new raw samples.
Definition: mhrw_bin_err.h:457
const Tools::StaticOrDynamic< int,(NumLevelsCTime==Eigen::Dynamic), NumLevelsCTime > numLevels
The number of levels in the binning analysis.
Definition: mhrw_bin_err.h:326
Params::ValueType ValueType
Type of the value(s) for which we are calculating error bars. See BinningAnalysisParams::ValueType.
Definition: mhrw_bin_err.h:267
constexpr bool isPowerOfTwo(IntType N)
Return true if the argument is a power of two, false otherwise.
Definition: cxxutil.h:282
ColXpr col(Index i)
CountIntType getNumFlushes() const
Return the number of times the collected samples were flushed.
Definition: mhrw_bin_err.h:559
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...
Definition: mhrw_bin_err.h:114
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:306
void reset()
Reset this object and start again as if freshly constructed.
Definition: mhrw_bin_err.h:437
BinSumSqArray calcErrorLevels(const Eigen::ArrayBase< Derived > &means) const
Calculate the error bars of samples at different binning levels.
Definition: mhrw_bin_err.h:629
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:321
The error bars appear to have converged.
Definition: mhrw_bin_err.h:195
Params::SamplesArray SamplesArray
Type of the internal buffer for the raw samples. See BinningAnalysisParams::SamplesArray.
Definition: mhrw_bin_err.h:295
Params::BinSumArray BinSumArray
Type used to store the sum of values. See BinningAnalysisParams::BinSumArray.
Definition: mhrw_bin_err.h:299
Some C++ utilities, with a tad of C++11 tricks.
Simple binning analysis for determining error bars.
Definition: mhrw_bin_err.h:259
LoggerType_ LoggerType
Type of the logger we will be logging debugging messages to.
Definition: mhrw_bin_err.h:349
BinningAnalysis(int num_track_values_, int num_levels_, LoggerType &logger_)
Constructor.
Definition: mhrw_bin_err.h:415
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.
Definition: mhrw_bin_err.h:188
VarValueDecoder< T >::RetType value(const Var &var)
Access the value of the given variable, as a C++ type.
Definition: ezmatio.h:878
const BinSumSqArray & getBinSumsq() const
Get the raw sums of the squared values observed, at each binning level.
Definition: mhrw_bin_err.h:613
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:336
T sqrt(T... args)
Eigen::Array< ValueType, NumTrackValuesCTime, SamplesSizeCTime > SamplesArray
Type of the internal buffer for the raw samples.
Definition: mhrw_bin_err.h:172
BinSumArray calcErrorLastLevel(const Eigen::ArrayBase< Derived > &means) const
Calculate the error bar of samples (from the last binning level).
Definition: mhrw_bin_err.h:674
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:117
#define tomographer_assert(...)
Assertion test macro.
Definition: cxxdefs.h:83
const Tools::StaticOrDynamic< int,(NumTrackValuesCTime==Eigen::Dynamic), NumTrackValuesCTime > numTrackValues
The number of functions being tracked/analyzed.
Definition: mhrw_bin_err.h:315
Utilities for logging messages.