Tomographer  v5.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) 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 
464  SamplesArray samples;
465 
466  // where we store the flushed values
467 
472  CountIntType n_samples;
478  CountIntType n_flushes;
499  BinSumSqArray bin_sumsq;
500 
502  LoggerType & logger;
503 
504 
505 public:
506 
518  BinningAnalysis(int num_track_values_, int num_levels_, LoggerType & logger_)
519  : numTrackValues(num_track_values_),
520  numLevels(num_levels_),
521  samplesSize(1 << numLevels()),
522  samples((Eigen::Index)numTrackValues(), (Eigen::Index)samplesSize()),
523  n_flushes(0),
524  bin_sum(BinSumArray::Zero((Eigen::Index)numTrackValues())),
525  bin_sumsq(BinSumSqArray::Zero((Eigen::Index)numTrackValues(), (Eigen::Index)numLevels()+1)),
526  logger(logger_)
527  {
529  tomographer_assert(Tools::isPowerOfTwo(samplesSize()));
530  tomographer_assert( (1<<numLevels()) == samplesSize() );
531 
532  reset();
533  }
534 
540  inline void reset()
541  {
542  n_flushes = 0;
543  n_samples = 0;
544  helper_reset_bin_sum();
545  bin_sumsq = BinSumSqArray::Zero((Eigen::Index)numTrackValues(), (Eigen::Index)numLevels()+1);
546  logger.longdebug("BinningAnalysis::reset()", "ready to go.");
547  }
548 
559  template<typename Derived>
560  inline void processNewValues(const Eigen::DenseBase<Derived> & vals)
561  {
562  const CountIntType ninbin = n_samples % samplesSize();
563 
564  ++n_samples;
565 
566  tomographer_assert(vals.rows() == numTrackValues());
567  tomographer_assert(vals.cols() == 1);
568 
569  // store the new values in the bins [also if ninbin == 0]
570  samples.col(ninbin) = vals;
571 
572  // add to our sum of values, if applicable.
573  helper_update_bin_sum(samples.col(ninbin));
574 
575  // see if we have to flush the bins (equivalent to `ninbin == samples_size()-1`)
576  if ( ninbin == samplesSize() - 1 ) {
577 
578  // we have filled all bins. Flush them. Re-use the beginning of the samples[] array
579  // to store the reduced bins while flushing them.
580  logger.longdebug("BinningAnalysis::processNewValues()", [&](std::ostream & str) {
581  str << "n_samples is now " << n_samples << "; flushing bins. samplesSize() = " << samplesSize();
582  });
583 
584  // the size of the samples at the current level of binning. Starts at samplesSize,
585  // and decreases by half at each higher level.
586 
587  for (int level = 0; level <= numLevels(); ++level) {
588 
589  const CountIntType binnedsize = CountIntType(1) << (numLevels()-level);
590 
591  logger.longdebug("BinningAnalysis::processNewValues()", [&](std::ostream & str) {
592  str << "Processing binning level = " << level << ": binnedsize="<<binnedsize
593  << "; n_flushes=" << n_flushes << "\n";
594  str << "\tbinned samples = \n"
595  << samples.block(0,0,(Eigen::Index)numTrackValues(),(Eigen::Index)binnedsize);
596  });
597 
598  for (CountIntType ksample = 0; ksample < binnedsize; ++ksample) {
599  bin_sumsq.col(level) += samples.col(ksample).cwiseProduct(samples.col(ksample));
600  if (ksample % 2 == 0 && binnedsize > 1) {
601  samples.col(ksample/2) = boost::math::constants::half<ValueType>() *
602  (samples.col(ksample) + samples.col(ksample+1));
603  }
604  }
605 
606  }
607 
608  logger.longdebug("BinningAnalysis::processNewValues()", [&](std::ostream & str) {
609  str << "Flushing #" << n_flushes << " done. bin_sum is = \n" << bin_sum << "\n"
610  << "\tbin_sumsq is = \n" << bin_sumsq << "\n";
611  });
612 
613  ++n_flushes;
614  }
615 
616  }
617 
623  template<typename CalcValType, TOMOGRAPHER_ENABLED_IF_TMPL(NumTrackValuesCTime == 1)>
624  inline void processNewValue(const CalcValType val)
625  {
626  // for a single value
627  //processNewValues(Eigen::Map<const Eigen::Array<CalcValType,1,1> >(&val));
628  processNewValues(Eigen::Array<CalcValType,1,1>::Constant(val));
629  }
630 
631 
632 
633 private:
634 
635  TOMOGRAPHER_ENABLED_IF(StoreBinSums)
636  inline void helper_reset_bin_sum()
637  {
638  bin_sum.value = BinSumArray::Zero((Eigen::Index)numTrackValues());
639  }
640  TOMOGRAPHER_ENABLED_IF(!StoreBinSums)
641  inline void helper_reset_bin_sum() { }
642 
643  TOMOGRAPHER_ENABLED_IF(StoreBinSums)
644  inline void helper_update_bin_sum(const Eigen::Ref<const Eigen::Array<ValueType, NumTrackValuesCTime, 1> > &
645  new_samples)
646  {
647  bin_sum.value += new_samples;
648  }
649  template<typename Derived,
650  TOMOGRAPHER_ENABLED_IF_TMPL(!StoreBinSums)>
651  inline void helper_update_bin_sum(const Eigen::DenseBase<Derived> & ) { }
652 
653 
654 public:
655 
656  // retrieve results.
657 
663  inline CountIntType getNumFlushes() const { return n_flushes; }
664 
668  TOMOGRAPHER_ENABLED_IF(StoreBinSums)
669  inline auto getBinMeans() const
670 #ifndef TOMOGRAPHER_PARSED_BY_DOXYGEN
671  -> decltype(BinSumArray() / ValueType(n_samples))
672 #endif
673  {
674  return bin_sum.value / ValueType(n_samples);
675  }
676 
683  inline auto getBinSqmeans() const
684 #ifndef TOMOGRAPHER_PARSED_BY_DOXYGEN
685  -> decltype(
686  bin_sumsq.cwiseQuotient(n_flushes * Tools::replicated<NumTrackValuesCTime,1>(
688  .transpose().reverse(),
689  // replicated by:
690  numTrackValues(), 1
691  ))
692  )
693 #endif
694  {
695  return bin_sumsq.cwiseQuotient(n_flushes * Tools::replicated<NumTrackValuesCTime,1>(
697  .transpose().reverse(),
698  // replicated by:
699  numTrackValues(), 1
700  ));
701  }
702 
703 
708  TOMOGRAPHER_ENABLED_IF(StoreBinSums)
709  inline const BinSumArray & getBinSum() const { return bin_sum.value; }
710 
717  inline const BinSumSqArray & getBinSumsq() const {
718  return bin_sumsq;
719  }
720 
721 
732  template<typename Derived>
733  inline BinSumSqArray calcErrorLevels(const Eigen::ArrayBase<Derived> & means) const
734  {
735  tomographer_assert(means.rows() == (Eigen::Index)numTrackValues());
736  tomographer_assert(means.cols() == 1);
737  const int n_levels_plus_one = numLevels()+1;
738  const int n_track_values = numTrackValues();
739 
744  return (
745  getBinSqmeans() - Tools::replicated<1,NumLevelsPlusOneCTime>(
746  means.cwiseProduct(means).template cast<ValueType>(),
747  // replicated by:
748  1, n_levels_plus_one
749  )
750  ).cwiseMax(0).cwiseQuotient(
751  // divide by the number of samples from which these bin-means were obtained, minus one.
752  Tools::replicated<NumTrackValuesCTime,1>(
754  .transpose().reverse(),
755  // replicated by:
756  n_track_values, 1
757  ) * n_flushes
759  (Eigen::Index)numTrackValues(), (Eigen::Index)n_levels_plus_one,
760  1 // the constant...
761  )
762  ).cwiseSqrt();
763  }
764 
765 
777  template<typename Derived>
778  inline BinSumArray calcErrorLastLevel(const Eigen::ArrayBase<Derived> & means) const {
779  tomographer_assert(means.rows() == (Eigen::Index)numTrackValues());
780  tomographer_assert(means.cols() == 1);
781  return (
782  bin_sumsq.col((Eigen::Index)numLevels()) / ValueType(n_flushes)
783  - means.cwiseProduct(means).template cast<ValueType>()
784  ).cwiseMax(0).cwiseSqrt() / std::sqrt(ValueType(n_flushes-1));
785  }
786 
798  TOMOGRAPHER_ENABLED_IF(StoreBinSums)
799  inline BinSumSqArray calcErrorLevels() const {
800  BinSumArray means = getBinMeans();
801  return calcErrorLevels(std::move(means));
802  }
803 
816  TOMOGRAPHER_ENABLED_IF(StoreBinSums)
817  inline BinSumArray calcErrorLastLevel() const {
818  BinSumArray means = getBinMeans();
819  return calcErrorLastLevel(std::move(means));
820  }
821 
834  inline Eigen::ArrayXi determineErrorConvergence(const Eigen::Ref<const BinSumSqArray> & error_levels) const
835  {
836  Eigen::ArrayXi converged_status(numTrackValues()); // RVO will help
837 
838  tomographer_assert(error_levels.rows() == (Eigen::Index)numTrackValues());
839  tomographer_assert(error_levels.cols() == (Eigen::Index)numLevels() + 1);
840 
841  logger.longdebug("BinningAnalysis::determineErrorConvergence", [&](std::ostream & str) {
842  str << "error_levels = \n" << error_levels << "\n";
843  });
844 
845  // verify that indeed the errors have converged. Inspired from ALPS code, see
846  // https://alps.comp-phys.org/svn/alps1/trunk/alps/src/alps/alea/simplebinning.h
847 
848  const int range = 4;
849  if (numLevels() < range-1) {
850 
851  converged_status = Eigen::ArrayXi::Constant((Eigen::Index)numTrackValues(), UNKNOWN_CONVERGENCE);
852 
853  } else {
854 
855  converged_status = Eigen::ArrayXi::Constant((Eigen::Index)numTrackValues(), CONVERGED);
856 
857  const auto & errors = error_levels.col((Eigen::Index)numLevels());
858 
859  for (int level = numLevels()+1 - range; level < numLevels(); ++level) {
860 
861  const auto & errors_thislevel = error_levels.col((Eigen::Index)level);
862 
863  logger.longdebug("BinningAnalysis::determineErrorConvergence", [&](std::ostream & str) {
864  str << "About to study level " << level << ": at this point, converged_status = \n"
865  << converged_status << "\nand errors_thislevel = \n" << errors_thislevel;
866  });
867 
868  for (Eigen::Index val_it = 0; val_it < (Eigen::Index)numTrackValues(); ++val_it) {
869  if (errors_thislevel(val_it) >= errors(val_it) &&
870  converged_status(val_it) != NOT_CONVERGED) {
871  converged_status(val_it) = CONVERGED;
872  } else if (errors_thislevel(val_it) < 0.824 * errors(val_it)) {
873  converged_status(val_it) = NOT_CONVERGED;
874  } else if ((errors_thislevel(val_it) < 0.9 * errors(val_it)) &&
875  converged_status(val_it) != NOT_CONVERGED) {
876  converged_status(val_it) = UNKNOWN_CONVERGENCE;
877  }
878  }
879 
880  }
881 
882  }
883 
884  logger.longdebug("BinningAnalysis::determineErrorConvergence", [&](std::ostream & str) {
885  str << "Done. converged_status [0=UNNOWN,1=CONVERGED,2=NOT CONVERGED] = \n"
886  << converged_status;
887  });
888 
889  return converged_status;
890  }
891 
892 
893 };
894 
895 
896 
897 } // namespace Tomographer
898 
899 
900 
901 
902 
903 
904 
905 
906 
907 
908 
909 
910 
911 
912 
913 #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:624
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:683
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:834
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:560
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:311
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:663
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
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:540
BinSumSqArray calcErrorLevels(const Eigen::ArrayBase< Derived > &means) const
Calculate the error bars of samples at different binning levels.
Definition: mhrw_bin_err.h:733
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:350
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:518
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:717
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:778
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.