All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
ReweightableSpectrum.h
Go to the documentation of this file.
1 #pragma once
2 
4 
5 #include <string>
6 
7 class TDirectory;
8 class TH2;
9 class TH2D;
10 
11 namespace ana
12 {
13  /// %Spectrum with the value of a second variable, allowing for reweighting
15  {
16  public:
17  friend class SpectrumLoaderBase;
18  friend class SpectrumLoader;
19  friend class NullLoader;
20  friend class MRCCLoader;
21 
23  const HistAxis& recoAxis,
24  const HistAxis& trueAxis,
25  const Cut& cut,
26  const SystShifts& shift = kNoShift,
27  const Var& wei = kUnweighted);
28 
30  const HistAxis& recoAxis,
31  const HistAxis& trueAxis,
32  const SpillCut& spillcut,
33  const SliceCut& slicecut,
34  const SystShifts& shift = kNoShift,
35  const Var& wei = kUnweighted);
36 
37  ReweightableSpectrum(const Var& rwVar,
38  const std::string& xlabel, const std::string& ylabel,
39  double pot,
40  int nbinsx, double xmin, double xmax,
41  int nbinsy, double ymin, double ymax);
42 
43  ReweightableSpectrum(const Var& rwVar,
44  TH2* h,
45  const std::vector<std::string>& labels,
46  const std::vector<Binning>& bins,
47  double pot, double livetime);
48 
49  ReweightableSpectrum(const Var& rwVar,
50  std::unique_ptr<TH2D> h,
51  const std::vector<std::string>& labels,
52  const std::vector<Binning>& bins,
53  double pot, double livetime);
54 
55  virtual ~ReweightableSpectrum();
56 
59 
60  /// \brief The variable that will be used to fill the y-axis
61  ///
62  /// By convention, return zero if the information can't be obtained, and
63  /// this event will be skipped.
64  const Var& ReweightVar() const {return fRWVar;}
65 
66  void Fill(double x, double y, double w = 1);
67 
68  TH2D* ToTH2(double pot) const;
69 
70  Spectrum UnWeighted() const;
71 
73 
74  Spectrum WeightedBy(const TH1* weights) const;
75 
76  /// Rescale bins so that \ref WeightingVariable will return \a target
77  void ReweightToTrueSpectrum(const Spectrum& target);
78  /// Recale bins so that \ref Unweighted will return \a target
79  void ReweightToRecoSpectrum(const Spectrum& target);
80 
81  void Clear();
82 
83  /// Function to save a ReweightableSpectrum to file
84  /// the fRWVar member is not written to file, so when
85  /// the spectrum is loaded back from file, ReweightVar
86  /// should not be accessed, but reweighting still works
87  void SaveTo(TDirectory* dir) const;
88 
89  static std::unique_ptr<ReweightableSpectrum> LoadFrom(TDirectory* dir);
90 
91  unsigned int NDimensions() const{return fLabels.size();}
92  std::vector<std::string> GetLabels() const {return fLabels;}
93  std::vector<Binning> GetBinnings() const {return fBins;}
94 
95  protected:
96  // Derived classes can be trusted take care of their own construction
97  ReweightableSpectrum(const std::vector<std::string>& labels,
98  const std::vector<Binning>& bins,
99  const Var& rwVar)
100  : fRWVar(rwVar),
101  fHist(0), fPOT(0), fLivetime(0),
102  fLabels(labels), fBins(bins)
103  {
104  }
105 
106  ReweightableSpectrum(const std::string& label,
107  const Binning& bins,
108  const Var& rwVar)
109  : fRWVar(rwVar),
110  fHist(0), fPOT(0), fLivetime(0),
111  fLabels(1, label), fBins(1, bins)
112  {
113  }
114 
115  /// Constructor needed by LoadFrom. Since there's no good
116  /// way to store a Var, ReweightVar will return nonsense
117  /// for ReweightableSpectrum that are loaded from a file
119  const std::vector<std::string>& labels,
120  const std::vector<Binning>& bins,
121  double pot, double livetime)
122  : ReweightableSpectrum(kUnweighted, h, labels, bins, pot, livetime)
123  {
124  }
125 
128 
129  Var fRWVar; ///< What goes on the y axis?
130 
131  TH2D* fHist;
132  double fPOT;
133  double fLivetime;
134 
135  std::vector<std::string> fLabels;
136  std::vector<Binning> fBins;
137 
138  std::string fTrueLabel;
139 
140  /// This count is maintained by SpectrumLoader, as a sanity check
141  std::set<SpectrumLoaderBase*> fLoaderCount;
142  };
143 }
const Var & ReweightVar() const
The variable that will be used to fill the y-axis.
* labels
Represent the binning of a Spectrum&#39;s x-axis.
Definition: Binning.h:18
process_name opflash particleana ie x
void Fill(double x, double y, double w=1)
tuple loader
Definition: demo.py:7
Simple record of shifts applied to systematic parameters.
Definition: SystShifts.h:16
Var fRWVar
What goes on the y axis?
Spectrum with the value of a second variable, allowing for reweighting
void ReweightToRecoSpectrum(const Spectrum &target)
Recale bins so that Unweighted will return target.
process_name opflashCryoW ana
ReweightableSpectrum(const std::vector< std::string > &labels, const std::vector< Binning > &bins, const Var &rwVar)
shift
Definition: fcl_checks.sh:26
Representation of a spectrum in any variable, with associated POT.
Definition: Spectrum.h:30
void RemoveLoader(SpectrumLoaderBase *)
process_name pandoraGausCryo1 vertexChargeCryo1 vertexStubCryo1 xmin
while getopts h
process_name opflash particleana ie ie y
std::vector< Binning > GetBinnings() const
std::vector< std::string > GetLabels() const
void SaveTo(TDirectory *dir) const
std::set< SpectrumLoaderBase * > fLoaderCount
This count is maintained by SpectrumLoader, as a sanity check.
std::vector< std::string > fLabels
const SystShifts kNoShift
Definition: SystShifts.h:61
tuple dir
Definition: dropbox.py:28
Base class for the various types of spectrum loader.
ReweightableSpectrum & operator=(const ReweightableSpectrum &rhs)
Collaborates with Spectrum and OscillatableSpectrum to fill spectra from CAF files.
const Var kUnweighted([](const caf::SRSliceProxy *){return 1;})
The simplest possible Var, always 1. Used as a default weight.
Spectrum WeightedBy(const TH1 *weights) const
static std::unique_ptr< ReweightableSpectrum > LoadFrom(TDirectory *dir)
std::vector< Binning > fBins
Dummy loader that doesn&#39;t load any files.
void ReweightToTrueSpectrum(const Spectrum &target)
Rescale bins so that WeightingVariable will return target.
ReweightableSpectrum(SpectrumLoaderBase &loader, const HistAxis &recoAxis, const HistAxis &trueAxis, const Cut &cut, const SystShifts &shift=kNoShift, const Var &wei=kUnweighted)
void AddLoader(SpectrumLoaderBase *)
ReweightableSpectrum(TH2 *h, const std::vector< std::string > &labels, const std::vector< Binning > &bins, double pot, double livetime)
ReweightableSpectrum(const std::string &label, const Binning &bins, const Var &rwVar)
unsigned int NDimensions() const
TH2D * ToTH2(double pot) const