Tomographer  v5.2
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) 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 
84 };
85 
86 
87 
88 
89 
90 
107 struct TOMOGRAPHER_EXPORT BinningErrorBarConvergenceSummary {
108  BinningErrorBarConvergenceSummary(Eigen::Index n_bins_ = 0,
109  Eigen::Index n_converged_ = 0,
110  Eigen::Index n_unknown_ = 0,
111  Eigen::Index n_unknown_isolated_ = 0,
112  Eigen::Index n_not_converged_ = 0)
113  : n_bins(n_bins_),
114  n_converged(n_converged_),
115  n_unknown(n_unknown_),
116  n_unknown_isolated(n_unknown_isolated_),
117  n_not_converged(n_not_converged_)
118  {
119  }
120 
121  Eigen::Index n_bins;
122  Eigen::Index n_converged;
123  Eigen::Index n_unknown;
124  Eigen::Index n_unknown_isolated;
125  Eigen::Index n_not_converged;
126 
130  static BinningErrorBarConvergenceSummary
132  {
133  const Eigen::Index n_bins = converged_status.size();
134 
135  const Eigen::ArrayXi unkn_arr =
136  converged_status.cwiseEqual(BINNING_UNKNOWN_CONVERGENCE).cast<int>();
137  const Eigen::ArrayXi conv_arr =
138  converged_status.cwiseEqual(BINNING_CONVERGED).cast<int>();
139 
140  const Eigen::Index n_unknown = unkn_arr.cast<Eigen::Index>().count();
141  const Eigen::Index n_converged = conv_arr.cast<Eigen::Index>().count();
142  const Eigen::Index n_not_converged =
143  converged_status.cwiseEqual(BINNING_NOT_CONVERGED).cast<Eigen::Index>().count();
144 
145  // Little heuristic to see whether the "unknown" converged error bars are isolated or
146  // not. Use conv_arr shifted by one in each direction as a mask -- an unconverged
147  // error bar is isolated if it is surrounded by converged error bars
148  const Eigen::Index n_unknown_isol = (
149  ( unkn_arr.segment(1,n_bins-2)
150  .cwiseProduct(conv_arr.segment(0,n_bins-2))
151  .cwiseProduct(conv_arr.segment(2,n_bins-2))
152  ).count()
153  // edges
154  + unkn_arr(0)*conv_arr(1) + unkn_arr(n_bins-1)*conv_arr(n_bins-2)
155  );
156 
157  return BinningErrorBarConvergenceSummary(n_bins, n_converged, n_unknown, n_unknown_isol, n_not_converged);
158  }
159 };
160 
161 
162 inline std::ostream & operator<<(std::ostream & stream, const BinningErrorBarConvergenceSummary & s)
163 {
164  return stream << s.n_converged << " converged / "
165  << s.n_unknown << " maybe (" << s.n_unknown_isolated << " isolated) / "
166  << s.n_not_converged << " not converged";
167 }
168 
169 
170 
171 
172 
210 template<typename ValueType_, int NumTrackValues_ = Eigen::Dynamic, int NumLevels_ = Eigen::Dynamic,
211  bool StoreBinSums_ = true, typename CountIntType_ = int>
212 struct TOMOGRAPHER_EXPORT BinningAnalysisParams
213 {
216  typedef ValueType_ ValueType;
219  typedef CountIntType_ CountIntType;
220 
224  static constexpr int NumTrackValuesCTime = NumTrackValues_;
228  static constexpr int NumLevelsCTime = NumLevels_;
232  static constexpr int NumLevelsPlusOneCTime = (NumLevelsCTime == Eigen::Dynamic
233  ? Eigen::Dynamic
234  : (NumLevelsCTime + 1));
257  static constexpr int SamplesSizeCTime =
258  tomo_internal::helper_samples_size<NumLevelsCTime, (NumLevelsCTime > 0 && NumLevelsCTime < 7)>::value;
259 
265  static constexpr bool StoreBinSums = StoreBinSums_;
266 
275 
291 
292  // //! Constants for error bar convergence analysis.
293  // enum {
294  // //! Unable to determine whether the error bars have converged.
295  // UNKNOWN_CONVERGENCE = BINNING_UNKNOWN_CONVERGENCE,
296  // //! The error bars appear to have converged.
297  // CONVERGED = BINNING_CONVERGED,
298  // //! The error bars don't seem to have converged.
299  // NOT_CONVERGED = BINNING_NOT_CONVERGED
300  // };
301 
302 };
303 
360 template<typename Params, typename LoggerType_>
361 class TOMOGRAPHER_EXPORT BinningAnalysis
362  // inheriting from this has some advantages over EIGEN_MAKE_ALIGNED_OPERATOR_NEW, such
363  // as not needing to explicitly declare the specialization
364  // NeedOwnOperatorNew<BinningAnalysis<...> >
366 {
367 public:
369  typedef typename Params::ValueType ValueType;
371  typedef typename Params::CountIntType CountIntType;
372 
376  static constexpr int NumTrackValuesCTime = Params::NumTrackValuesCTime;
380  static constexpr int NumLevelsCTime = Params::NumLevelsCTime;
384  static constexpr int NumLevelsPlusOneCTime = Params::NumLevelsPlusOneCTime;
388  static constexpr int SamplesSizeCTime = Params::SamplesSizeCTime;
392  static constexpr bool StoreBinSums = Params::StoreBinSums;
393 
397  typedef typename Params::SamplesArray SamplesArray;
401  typedef typename Params::BinSumArray BinSumArray;
405  typedef typename Params::BinSumSqArray BinSumSqArray;
406 
438  const Tools::StaticOrDynamic<CountIntType, (SamplesSizeCTime==Eigen::Dynamic),
439  (SamplesSizeCTime < 0 ? 0 : SamplesSizeCTime)> samplesSize;
440 
442  enum {
444  UNKNOWN_CONVERGENCE = BINNING_UNKNOWN_CONVERGENCE,
446  CONVERGED = BINNING_CONVERGED,
448  NOT_CONVERGED = BINNING_NOT_CONVERGED
449  };
450 
452  typedef LoggerType_ LoggerType;
453 
454 private:
455 
456  // private not picked up by doxygen anyway
457 
458  /* * \brief The array in which we store samples that arrive from the simulation.
459  *
460  * This array has size \a samplesSize() (for each tracking value). Once this array is
461  * filled, it is <em>flushed</em>, i.e. the values are processed and stored as
462  * appropriate in \ref bin_sum and \ref bin_sumsq.
463  *
464  * This array has \a numTrackingValues() rows and \a samplesSize() columns.
465  */
466  SamplesArray samples;
467 
468  // where we store the flushed values
469 
470  /* * \brief Number of samples seen.
471  *
472  * This is equal to the number of times \ref processNewValues() was called.
473  */
474  CountIntType n_samples;
475  /* * \brief Number of flushes.
476  *
477  * A flush corresponds to having filled all the samples in the sample vector (of size \a
478  * samplesSize()), and pushing new values into \ref bin_sum and \ref bin_sumsq.
479  */
480  CountIntType n_flushes;
481  /* * \brief Sum of all values seen.
482  *
483  * This is a column vector of \a numTrackingValues() entries.
484  *
485  * \note This member is only available if the template parameter \a StoreBinSums is set
486  * to \c true.
487  *
488  * \note values are added to this array as soon as they are seen, not when the samples
489  * array is flushed. In particular, if the total number of values is not a multiple of
490  * \a samplesSize(), then there will be samples counted into \a bin_sum but not into \a
491  * bin_sumsq.
492  */
494  /* * \brief Sum of the squares of all flushed & processed values, at different binning
495  * levels.
496  *
497  * This is a matrix of \a numTrackingValues() rows and <em>numLevels()+1</em>
498  * columns.
499  *
500  */
501  BinSumSqArray bin_sumsq;
502 
503  // ! Just a boring logger...
504  LoggerType & logger;
505 
506 
507 public:
508 
520  BinningAnalysis(int num_track_values_, int num_levels_, LoggerType & logger_)
521  : numTrackValues(num_track_values_),
522  numLevels(num_levels_),
523  samplesSize(1 << numLevels()),
524  samples((Eigen::Index)numTrackValues(), (Eigen::Index)samplesSize()),
525  n_flushes(0),
526  bin_sum(BinSumArray::Zero((Eigen::Index)numTrackValues())),
527  bin_sumsq(BinSumSqArray::Zero((Eigen::Index)numTrackValues(), (Eigen::Index)numLevels()+1)),
528  logger(logger_)
529  {
531  tomographer_assert(Tools::isPowerOfTwo(samplesSize()));
532  tomographer_assert( (1<<numLevels()) == samplesSize() );
533 
534  reset();
535  }
536 
542  inline void reset()
543  {
544  n_flushes = 0;
545  n_samples = 0;
546  helper_reset_bin_sum();
547  bin_sumsq = BinSumSqArray::Zero((Eigen::Index)numTrackValues(), (Eigen::Index)numLevels()+1);
548  logger.longdebug("BinningAnalysis::reset()", "ready to go.");
549  }
550 
561  template<typename Derived>
562  inline void processNewValues(const Eigen::DenseBase<Derived> & vals)
563  {
564  const CountIntType ninbin = n_samples % samplesSize();
565 
566  ++n_samples;
567 
568  tomographer_assert(vals.rows() == numTrackValues());
569  tomographer_assert(vals.cols() == 1);
570 
571  // store the new values in the bins [also if ninbin == 0]
572  samples.col(ninbin) = vals;
573 
574  // add to our sum of values, if applicable.
575  helper_update_bin_sum(samples.col(ninbin));
576 
577  // see if we have to flush the bins (equivalent to `ninbin == samples_size()-1`)
578  if ( ninbin == samplesSize() - 1 ) {
579 
580  // we have filled all bins. Flush them. Re-use the beginning of the samples[] array
581  // to store the reduced bins while flushing them.
582  logger.longdebug("BinningAnalysis::processNewValues()", [&](std::ostream & str) {
583  str << "n_samples is now " << n_samples << "; flushing bins. samplesSize() = " << samplesSize();
584  });
585 
586  // the size of the samples at the current level of binning. Starts at samplesSize,
587  // and decreases by half at each higher level.
588 
589  for (int level = 0; level <= numLevels(); ++level) {
590 
591  const CountIntType binnedsize = CountIntType(1) << (numLevels()-level);
592 
593  logger.longdebug("BinningAnalysis::processNewValues()", [&](std::ostream & str) {
594  str << "Processing binning level = " << level << ": binnedsize="<<binnedsize
595  << "; n_flushes=" << n_flushes << "\n";
596  str << "\tbinned samples = \n"
597  << samples.block(0,0,(Eigen::Index)numTrackValues(),(Eigen::Index)binnedsize);
598  });
599 
600  for (CountIntType ksample = 0; ksample < binnedsize; ++ksample) {
601  bin_sumsq.col(level) += samples.col(ksample).cwiseProduct(samples.col(ksample));
602  if (ksample % 2 == 0 && binnedsize > 1) {
603  samples.col(ksample/2) = boost::math::constants::half<ValueType>() *
604  (samples.col(ksample) + samples.col(ksample+1));
605  }
606  }
607 
608  }
609 
610  logger.longdebug("BinningAnalysis::processNewValues()", [&](std::ostream & str) {
611  str << "Flushing #" << n_flushes << " done. bin_sum is = \n" << bin_sum << "\n"
612  << "\tbin_sumsq is = \n" << bin_sumsq << "\n";
613  });
614 
615  ++n_flushes;
616  }
617 
618  }
619 
625  template<typename CalcValType, TOMOGRAPHER_ENABLED_IF_TMPL(NumTrackValuesCTime == 1)>
626  inline void processNewValue(const CalcValType val)
627  {
628  // for a single value
629  //processNewValues(Eigen::Map<const Eigen::Array<CalcValType,1,1> >(&val));
630  processNewValues(Eigen::Array<CalcValType,1,1>::Constant(val));
631  }
632 
633 
634 
635 private:
636 
637  TOMOGRAPHER_ENABLED_IF(StoreBinSums)
638  inline void helper_reset_bin_sum()
639  {
640  bin_sum.value = BinSumArray::Zero((Eigen::Index)numTrackValues());
641  }
642  TOMOGRAPHER_ENABLED_IF(!StoreBinSums)
643  inline void helper_reset_bin_sum() { }
644 
645  TOMOGRAPHER_ENABLED_IF(StoreBinSums)
646  inline void helper_update_bin_sum(const Eigen::Ref<const Eigen::Array<ValueType, NumTrackValuesCTime, 1> > &
647  new_samples)
648  {
649  bin_sum.value += new_samples;
650  }
651  template<typename Derived,
652  TOMOGRAPHER_ENABLED_IF_TMPL(!StoreBinSums)>
653  inline void helper_update_bin_sum(const Eigen::DenseBase<Derived> & ) { }
654 
655 
656 public:
657 
658  // retrieve results.
659 
665  inline CountIntType getNumFlushes() const { return n_flushes; }
666 
670  TOMOGRAPHER_ENABLED_IF(StoreBinSums)
671  inline auto getBinMeans() const
672 #ifndef TOMOGRAPHER_PARSED_BY_DOXYGEN
673  -> decltype(BinSumArray() / ValueType(n_samples))
674 #endif
675  {
676  return bin_sum.value / ValueType(n_samples);
677  }
678 
685  inline auto getBinSqmeans() const
686 #ifndef TOMOGRAPHER_PARSED_BY_DOXYGEN
687  -> decltype(
688  bin_sumsq.cwiseQuotient(n_flushes * Tools::replicated<NumTrackValuesCTime,1>(
690  .transpose().reverse(),
691  // replicated by:
692  numTrackValues(), 1
693  ))
694  )
695 #endif
696  {
697  return bin_sumsq.cwiseQuotient(n_flushes * Tools::replicated<NumTrackValuesCTime,1>(
699  .transpose().reverse(),
700  // replicated by:
701  numTrackValues(), 1
702  ));
703  }
704 
705 
710  TOMOGRAPHER_ENABLED_IF(StoreBinSums)
711  inline const BinSumArray & getBinSum() const { return bin_sum.value; }
712 
719  inline const BinSumSqArray & getBinSumsq() const {
720  return bin_sumsq;
721  }
722 
723 
734  template<typename Derived>
735  inline BinSumSqArray calcErrorLevels(const Eigen::ArrayBase<Derived> & means) const
736  {
737  tomographer_assert(means.rows() == (Eigen::Index)numTrackValues());
738  tomographer_assert(means.cols() == 1);
739  const int n_levels_plus_one = numLevels()+1;
740  const int n_track_values = numTrackValues();
741 
746  return (
747  getBinSqmeans() - Tools::replicated<1,NumLevelsPlusOneCTime>(
748  means.cwiseProduct(means).template cast<ValueType>(),
749  // replicated by:
750  1, n_levels_plus_one
751  )
752  ).cwiseMax(0).cwiseQuotient(
753  // divide by the number of samples from which these bin-means were obtained, minus one.
754  Tools::replicated<NumTrackValuesCTime,1>(
756  .transpose().reverse(),
757  // replicated by:
758  n_track_values, 1
759  ) * n_flushes
761  (Eigen::Index)numTrackValues(), (Eigen::Index)n_levels_plus_one,
762  1 // the constant...
763  )
764  ).cwiseSqrt();
765  }
766 
767 
779  template<typename Derived>
780  inline BinSumArray calcErrorLastLevel(const Eigen::ArrayBase<Derived> & means) const {
781  tomographer_assert(means.rows() == (Eigen::Index)numTrackValues());
782  tomographer_assert(means.cols() == 1);
783  return (
784  bin_sumsq.col((Eigen::Index)numLevels()) / ValueType(n_flushes)
785  - means.cwiseProduct(means).template cast<ValueType>()
786  ).cwiseMax(0).cwiseSqrt() / std::sqrt(ValueType(n_flushes-1));
787  }
788 
800  TOMOGRAPHER_ENABLED_IF(StoreBinSums)
801  inline BinSumSqArray calcErrorLevels() const {
802  BinSumArray means = getBinMeans();
803  return calcErrorLevels(std::move(means));
804  }
805 
818  TOMOGRAPHER_ENABLED_IF(StoreBinSums)
819  inline BinSumArray calcErrorLastLevel() const {
820  BinSumArray means = getBinMeans();
821  return calcErrorLastLevel(std::move(means));
822  }
823 
836  inline Eigen::ArrayXi determineErrorConvergence(const Eigen::Ref<const BinSumSqArray> & error_levels) const
837  {
838  Eigen::ArrayXi converged_status(numTrackValues()); // RVO will help
839 
840  tomographer_assert(error_levels.rows() == (Eigen::Index)numTrackValues());
841  tomographer_assert(error_levels.cols() == (Eigen::Index)numLevels() + 1);
842 
843  logger.longdebug("BinningAnalysis::determineErrorConvergence", [&](std::ostream & str) {
844  str << "error_levels = \n" << error_levels << "\n";
845  });
846 
847  // verify that indeed the errors have converged. Inspired from ALPS code, see
848  // https://alps.comp-phys.org/svn/alps1/trunk/alps/src/alps/alea/simplebinning.h
849 
850  const int range = 4;
851  if (numLevels() < range-1) {
852 
853  converged_status = Eigen::ArrayXi::Constant((Eigen::Index)numTrackValues(), UNKNOWN_CONVERGENCE);
854 
855  } else {
856 
857  converged_status = Eigen::ArrayXi::Constant((Eigen::Index)numTrackValues(), CONVERGED);
858 
859  const auto & errors = error_levels.col((Eigen::Index)numLevels());
860 
861  for (int level = numLevels()+1 - range; level < numLevels(); ++level) {
862 
863  const auto & errors_thislevel = error_levels.col((Eigen::Index)level);
864 
865  logger.longdebug("BinningAnalysis::determineErrorConvergence", [&](std::ostream & str) {
866  str << "About to study level " << level << ": at this point, converged_status = \n"
867  << converged_status << "\nand errors_thislevel = \n" << errors_thislevel;
868  });
869 
870  for (Eigen::Index val_it = 0; val_it < (Eigen::Index)numTrackValues(); ++val_it) {
871  if (errors_thislevel(val_it) >= errors(val_it) &&
872  converged_status(val_it) != NOT_CONVERGED) {
873  converged_status(val_it) = CONVERGED;
874  } else if (errors_thislevel(val_it) < 0.824 * errors(val_it)) {
875  converged_status(val_it) = NOT_CONVERGED;
876  } else if ((errors_thislevel(val_it) < 0.9 * errors(val_it)) &&
877  converged_status(val_it) != NOT_CONVERGED) {
878  converged_status(val_it) = UNKNOWN_CONVERGENCE;
879  }
880  }
881 
882  }
883 
884  }
885 
886  logger.longdebug("BinningAnalysis::determineErrorConvergence", [&](std::ostream & str) {
887  str << "Done. converged_status [0=UNNOWN,1=CONVERGED,2=NOT CONVERGED] = \n"
888  << converged_status;
889  });
890 
891  return converged_status;
892  }
893 
894 
895 };
896 
897 
898 
915 template<typename IterCountIntType, typename LocalLoggerType>
916 inline int sanitizeBinningLevels(int binning_num_levels,
917  IterCountIntType n_run,
918  IterCountIntType samples_last_level,
919  LocalLoggerType & logger)
920 {
921  if (binning_num_levels > 0) {
922  // provided manual value
923  if (binning_num_levels < 4) {
924  logger.warning([&](std::ostream & stream) {
925  stream << "You are using few binning levels = " << binning_num_levels
926  << " and the resulting error bars might not be reliable." ;
927  });
928  }
929  } else {
930  // choose automatically. Make sure that the last level has
931  // ~samples_last_level samples to calculate std deviation.
932  binning_num_levels = (int)(
933  std::floor(std::log(n_run/samples_last_level)
934  / std::log(2)) + 1e-3
935  ) ;
936  if (binning_num_levels < 1) {
937  binning_num_levels = 1;
938  }
939  if (binning_num_levels < 4) {
940  logger.warning([&](std::ostream & stream) {
941  stream << "Because n_run is small, you will be using few binning levels = "
942  << binning_num_levels << ", and the resulting error bars "
943  << "might not be reliable." ;
944  });
945  }
946  }
947  const IterCountIntType binning_last_level_num_samples =
948  (IterCountIntType) std::ldexp((double)n_run, - binning_num_levels);
949 
950  logger.debug([&](std::ostream & stream) {
951  stream << "Binning analysis: " << binning_num_levels << " levels, with "
952  << binning_last_level_num_samples << " samples at last level";
953  });
954 
955  // warn if number of samples @ last level is below recommended value
956  if (binning_last_level_num_samples < samples_last_level) {
957  logger.warning([&](std::ostream & stream) {
958  stream << "The number of samples (" << binning_last_level_num_samples
959  << ") at the last binning level is below the recommended value ("
960  << samples_last_level << "). Consider increasing n_run "
961  << "or decreasing binning_num_levels.";
962  });
963  }
964  return binning_num_levels;
965 }
966 
967 
968 
969 
970 
971 
972 } // namespace Tomographer
973 
974 
975 
976 
977 
978 
979 
980 
981 
982 
983 
984 
985 
986 
987 
988 #endif
It is unknown whether the error bar has converged or not.
Definition: mhrw_bin_err.h:79
void processNewValue(const CalcValType val)
Process a new value (if we&#39;re tracking a single function only)
Definition: mhrw_bin_err.h:626
int sanitizeBinningLevels(int binning_num_levels, IterCountIntType n_run, IterCountIntType samples_last_level, LocalLoggerType &logger)
Utility for sanitizing/determining number of binning levels.
Definition: mhrw_bin_err.h:916
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:405
static BinningErrorBarConvergenceSummary fromConvergedStatus(const Eigen::Ref< const Eigen::ArrayXi > &converged_status)
Construct a summary object from a list of converged status obtained from BinningAnalysis::determineEr...
Definition: mhrw_bin_err.h:131
auto getBinSqmeans() const
Get the raw average of the squared values observed, for each binning level.
Definition: mhrw_bin_err.h:685
Provides correct operator-new implementation for Eigen types via the NeedOwnOperatorNew mechanism...
Base namespace for the Tomographer project.
Definition: densellh.h:45
The error bar appears to have converged.
Definition: mhrw_bin_err.h:81
Group template parameters for BinningAnalysis.
Definition: mhrw_bin_err.h:212
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:836
T log(T... args)
Eigen::Array< ValueType, NumTrackValuesCTime, 1 > BinSumArray
Type used to store the sum of values.
Definition: mhrw_bin_err.h:282
Params::CountIntType CountIntType
Type used to count the number of samples. See BinningAnalysisParams::CountIntType.
Definition: mhrw_bin_err.h:371
void processNewValues(const Eigen::DenseBase< Derived > &vals)
Process new raw samples.
Definition: mhrw_bin_err.h:562
const Tools::StaticOrDynamic< int,(NumLevelsCTime==Eigen::Dynamic), NumLevelsCTime > numLevels
The number of levels in the binning analysis.
Definition: mhrw_bin_err.h:428
Params::ValueType ValueType
Type of the value(s) for which we are calculating error bars. See BinningAnalysisParams::ValueType.
Definition: mhrw_bin_err.h:369
constexpr bool isPowerOfTwo(IntType N)
Return true if the argument is a power of two, false otherwise.
Definition: cxxutil.h:321
T floor(T... args)
const Tools::StaticOrDynamic< CountIntType,(SamplesSizeCTime==Eigen::Dynamic),(SamplesSizeCTime< 0 ? 0 :SamplesSizeCTime)> samplesSize
The size of our samples buffer. (avoid using, might change in the future)
Definition: mhrw_bin_err.h:439
CountIntType getNumFlushes() const
Return the number of times the collected samples were flushed.
Definition: mhrw_bin_err.h:665
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:216
T ldexp(T... args)
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:542
BinSumSqArray calcErrorLevels(const Eigen::ArrayBase< Derived > &means) const
Calculate the error bars of samples at different binning levels.
Definition: mhrw_bin_err.h:735
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:360
Params::SamplesArray SamplesArray
Type of the internal buffer for the raw samples. See BinningAnalysisParams::SamplesArray.
Definition: mhrw_bin_err.h:397
Params::BinSumArray BinSumArray
Type used to store the sum of values. See BinningAnalysisParams::BinSumArray.
Definition: mhrw_bin_err.h:401
Some C++ utilities, with a tad of C++11 tricks.
T move(T... args)
Simple binning analysis for determining error bars.
Definition: mhrw_bin_err.h:361
LoggerType_ LoggerType
Type of the logger we will be logging debugging messages to.
Definition: mhrw_bin_err.h:452
BinningAnalysis(int num_track_values_, int num_levels_, LoggerType &logger_)
Constructor.
Definition: mhrw_bin_err.h:520
Managing the need for specific overrides to operator new() for some types (especially Eigen types) ...
The error bar appears NOT to have converged.
Definition: mhrw_bin_err.h:83
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:290
VarValueDecoder< T >::RetType value(const Var &var)
Access the value of the given variable, as a C++ type.
Definition: ezmatio.h:878
Brief summary of convergence of error bars from a binning analysis.
Definition: mhrw_bin_err.h:107
const BinSumSqArray & getBinSumsq() const
Get the raw sums of the squared values observed, at each binning level.
Definition: mhrw_bin_err.h:719
T sqrt(T... args)
Eigen::Array< ValueType, NumTrackValuesCTime, SamplesSizeCTime > SamplesArray
Type of the internal buffer for the raw samples.
Definition: mhrw_bin_err.h:274
BinSumArray calcErrorLastLevel(const Eigen::ArrayBase< Derived > &means) const
Calculate the error bar of samples (from the last binning level).
Definition: mhrw_bin_err.h:780
BinningConvergence
Convergence status of binning error bar.
Definition: mhrw_bin_err.h:77
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:219
#define tomographer_assert(...)
Assertion test macro.
Definition: cxxdefs.h:84
const Tools::StaticOrDynamic< int,(NumTrackValuesCTime==Eigen::Dynamic), NumTrackValuesCTime > numTrackValues
The number of functions being tracked/analyzed.
Definition: mhrw_bin_err.h:417
Utilities for logging messages.