All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
UnfoldSVD.cxx
Go to the documentation of this file.
2 
4 
5 #include "TSVDUnfold.h"
6 
7 #include "TH2.h"
8 
9 namespace ana
10 {
11  //----------------------------------------------------------------------
13  const ReweightableSpectrum& recoVsTrue,
14  unsigned int reg)
15  {
16  const double pot = reco.POT();
17 
18  DontAddDirectory guard;
19 
20  std::unique_ptr<TH1D> hReco(reco.ToTH1(pot));
21  std::unique_ptr<TH1D> hMCReco(recoVsTrue.UnWeighted().ToTH1(pot));
22  std::unique_ptr<TH1D> hMCTrue(recoVsTrue.WeightingVariable().ToTH1(pot));
23  std::unique_ptr<TH2D> hRT(recoVsTrue.ToTH2(pot));
24 
25  TSVDUnfold uf(hReco.get(), hMCReco.get(), hMCTrue.get(), hRT.get());
26 
27  std::unique_ptr<TH1D> h_unf(uf.Unfold(reg));
28 
29  // Enforce matching normalization
30  h_unf->Scale(reco.Integral(pot)/h_unf->Integral(0, -1));
31 
32  // TODO in principle these should be the true labels and bins. Will be
33  // easier with cafanacore
34  return Spectrum(std::move(h_unf), reco.GetLabels(), reco.GetBinnings(),
35  pot, reco.Livetime());
36  }
37 }
TH1D * ToTH1(double exposure, Color_t col=kBlack, Style_t style=kSolid, EExposureType expotype=kPOT, EBinType bintype=kBinContent) const
Histogram made from this Spectrum, scaled to some exposure.
Spectrum with the value of a second variable, allowing for reweighting
process_name opflashCryoW ana
double Integral(double exposure, double *err=0, EExposureType expotype=kPOT) const
Return total number of events scaled to pot.
Representation of a spectrum in any variable, with associated POT.
Definition: Spectrum.h:30
Spectrum UnfoldSVD(const Spectrum &reco, const ReweightableSpectrum &recoVsTrue, unsigned int reg)
Singular value decomposition unfolding using ROOT&#39;s TSVDDecomp.
Definition: UnfoldSVD.cxx:12
process_name standard_reco_uboone reco
double POT() const
Definition: Spectrum.h:289
std::vector< std::string > GetLabels() const
Definition: Spectrum.h:324
Prevent histograms being added to the current directory.
double Livetime() const
Seconds. For informational purposes only. No calculations use this.
Definition: Spectrum.h:292
std::vector< Binning > GetBinnings() const
Definition: Spectrum.h:325
TH2D * ToTH2(double pot) const