All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
ana::EnsembleSpectrum Class Reference

#include <EnsembleSpectrum.h>

Public Member Functions

 EnsembleSpectrum (SpectrumLoaderBase &loader, const HistAxis &axis, const SpillCut &spillcut, const Cut &cut, const std::vector< SystShifts > &univ_shifts, const Var &cv_wei=kUnweighted)
 
 EnsembleSpectrum (SpectrumLoaderBase &loader, const HistAxis &axis, const SpillCut &spillcut, const Cut &cut, const std::vector< Var > &univ_weis, const Var &cv_wei=kUnweighted)
 
Spectrum Nominal () const
 
unsigned int NUniverses () const
 
Spectrum Universe (unsigned int i) const
 
double POT () const
 
double Livetime () const
 
TGraphAsymmErrors * ErrorBand (double exposure, EExposureType expotype=kPOT, EBinType bintype=kBinContent) const
 Result can be painted prettily with DrawErrorBand. More...
 
void Scale (double c)
 
EnsembleSpectrumoperator+= (const EnsembleSpectrum &rhs)
 
EnsembleSpectrum operator+ (const EnsembleSpectrum &rhs) const
 
EnsembleSpectrumoperator-= (const EnsembleSpectrum &rhs)
 
EnsembleSpectrum operator- (const EnsembleSpectrum &rhs) const
 
EnsembleSpectrumoperator*= (const EnsembleRatio &rhs)
 
EnsembleSpectrum operator* (const EnsembleRatio &rhs) const
 
EnsembleSpectrumoperator/= (const EnsembleRatio &rhs)
 
EnsembleSpectrum operator/ (const EnsembleRatio &rhs) const
 
void SaveTo (TDirectory *dir) const
 
unsigned int NDimensions () const
 
std::vector< std::string > GetLabels () const
 
std::vector< BinningGetBinnings () const
 

Static Public Member Functions

static std::unique_ptr
< EnsembleSpectrum
LoadFrom (TDirectory *dir)
 

Protected Member Functions

 EnsembleSpectrum (const Spectrum &nom)
 

Protected Attributes

Spectrum fNom
 
std::vector< SpectrumfUnivs
 

Detailed Description

Definition at line 16 of file EnsembleSpectrum.h.

Constructor & Destructor Documentation

ana::EnsembleSpectrum::EnsembleSpectrum ( SpectrumLoaderBase loader,
const HistAxis axis,
const SpillCut spillcut,
const Cut cut,
const std::vector< SystShifts > &  univ_shifts,
const Var cv_wei = kUnweighted 
)

Definition at line 14 of file EnsembleSpectrum.cxx.

20  : fNom(loader, axis, spillcut, cut, kNoShift, cv_wei)
21  {
22  fUnivs.reserve(univ_shifts.size());
23  for(const SystShifts& ss: univ_shifts){
24  fUnivs.emplace_back(loader, axis, spillcut, cut, ss, cv_wei);
25  }
26  }
tuple loader
Definition: demo.py:7
const SystShifts kNoShift
Definition: SystShifts.h:61
std::vector< Spectrum > fUnivs
ana::EnsembleSpectrum::EnsembleSpectrum ( SpectrumLoaderBase loader,
const HistAxis axis,
const SpillCut spillcut,
const Cut cut,
const std::vector< Var > &  univ_weis,
const Var cv_wei = kUnweighted 
)

Definition at line 29 of file EnsembleSpectrum.cxx.

35  : fNom(loader, axis, spillcut, cut, kNoShift, cv_wei)
36  {
37  fUnivs.reserve(univ_weis.size());
38  for(const Var& w: univ_weis){
39  fUnivs.emplace_back(loader, axis, spillcut, cut, kNoShift, cv_wei * w);
40  }
41  }
tuple loader
Definition: demo.py:7
_Var< caf::SRSliceProxy > Var
Representation of a variable to be retrieved from a caf::StandardRecord object.
Definition: Var.h:73
const SystShifts kNoShift
Definition: SystShifts.h:61
std::vector< Spectrum > fUnivs
ana::EnsembleSpectrum::EnsembleSpectrum ( const Spectrum nom)
inlineprotected

Definition at line 72 of file EnsembleSpectrum.h.

72 : fNom(nom) {}

Member Function Documentation

TGraphAsymmErrors * ana::EnsembleSpectrum::ErrorBand ( double  exposure,
EExposureType  expotype = kPOT,
EBinType  bintype = kBinContent 
) const

Result can be painted prettily with DrawErrorBand.

Definition at line 44 of file EnsembleSpectrum.cxx.

47  {
48  std::unique_ptr<TH1D> hnom(fNom.ToTH1(exposure, expotype, bintype));
49 
50  std::vector<std::unique_ptr<TH1D>> hunivs;
51  hunivs.reserve(fUnivs.size());
52  for(const Spectrum& u: fUnivs) hunivs.emplace_back(u.ToTH1(exposure, expotype, bintype));
53 
54  TGraphAsymmErrors* g = new TGraphAsymmErrors;
55 
56  for(int binIdx = 0; binIdx < hnom->GetNbinsX()+2; ++binIdx){
57  const double xnom = hnom->GetXaxis()->GetBinCenter(binIdx);
58  const double ynom = hnom->GetBinContent(binIdx);
59  g->SetPoint(binIdx, xnom, ynom);
60 
61  const double dx = hnom->GetXaxis()->GetBinWidth(binIdx);
62 
63  std::vector<double> ys;
64  ys.reserve(hunivs.size());
65  for(const std::unique_ptr<TH1D>& hu: hunivs){
66  ys.push_back(hu->GetBinContent(binIdx));
67  }
68 
69  // 1 sigma
70  const double y0 = FindQuantile(.5-0.6827/2, ys);
71  const double y1 = FindQuantile(.5+0.6827/2, ys);
72 
73  // It's theoretically possible for the central value to be outside the
74  // error bands - clamp to zero in that case
75  g->SetPointError(binIdx, dx/2, dx/2,
76  std::max(ynom-y0, 0.),
77  std::max(y1-ynom, 0.));
78  } // end for binIdx
79 
80  return g;
81  }
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.
double FindQuantile(double frac, std::vector< double > &xs)
BEGIN_PROLOG g
std::vector< Spectrum > fUnivs
std::vector<Binning> ana::EnsembleSpectrum::GetBinnings ( ) const
inline

Definition at line 69 of file EnsembleSpectrum.h.

69 {return fNom.GetBinnings();}
std::vector< Binning > GetBinnings() const
Definition: Spectrum.h:325
std::vector<std::string> ana::EnsembleSpectrum::GetLabels ( ) const
inline

Definition at line 68 of file EnsembleSpectrum.h.

68 {return fNom.GetLabels();}
std::vector< std::string > GetLabels() const
Definition: Spectrum.h:324
double ana::EnsembleSpectrum::Livetime ( ) const
inline

Definition at line 43 of file EnsembleSpectrum.h.

43 {return fNom.Livetime();}
double Livetime() const
Seconds. For informational purposes only. No calculations use this.
Definition: Spectrum.h:292
std::unique_ptr< EnsembleSpectrum > ana::EnsembleSpectrum::LoadFrom ( TDirectory *  dir)
static

Definition at line 175 of file EnsembleSpectrum.cxx.

176  {
177  DontAddDirectory guard;
178 
179  TObjString* tag = (TObjString*)dir->Get("type");
180  assert(tag);
181  assert(tag->GetString() == "EnsembleSpectrum");
182  delete tag;
183 
184  std::unique_ptr<EnsembleSpectrum> ret(new EnsembleSpectrum(*Spectrum::LoadFrom(dir->GetDirectory("nom"))));
185 
186  for(unsigned int i = 0; ; ++i){
187  TDirectory* d = dir->GetDirectory(("univ"+std::to_string(i)).c_str());
188  if(!d) break;
189  ret->fUnivs.push_back(*Spectrum::LoadFrom(d));
190  }
191 
192  return ret;
193  }
static std::unique_ptr< Spectrum > LoadFrom(TDirectory *dir)
EnsembleSpectrum(SpectrumLoaderBase &loader, const HistAxis &axis, const SpillCut &spillcut, const Cut &cut, const std::vector< SystShifts > &univ_shifts, const Var &cv_wei=kUnweighted)
tuple dir
Definition: dropbox.py:28
std::string to_string(WindowPattern const &pattern)
unsigned int ana::EnsembleSpectrum::NDimensions ( ) const
inline

Definition at line 67 of file EnsembleSpectrum.h.

67 {return fNom.NDimensions();}
unsigned int NDimensions() const
Definition: Spectrum.h:323
Spectrum ana::EnsembleSpectrum::Nominal ( ) const
inline

Definition at line 33 of file EnsembleSpectrum.h.

33 {return fNom;}
unsigned int ana::EnsembleSpectrum::NUniverses ( ) const
inline

Definition at line 34 of file EnsembleSpectrum.h.

34 {return fUnivs.size();}
std::vector< Spectrum > fUnivs
EnsembleSpectrum ana::EnsembleSpectrum::operator* ( const EnsembleRatio rhs) const

Definition at line 134 of file EnsembleSpectrum.cxx.

135  {
136  EnsembleSpectrum ret = *this;
137  ret *= rhs;
138  return ret;
139  }
EnsembleSpectrum(SpectrumLoaderBase &loader, const HistAxis &axis, const SpillCut &spillcut, const Cut &cut, const std::vector< SystShifts > &univ_shifts, const Var &cv_wei=kUnweighted)
EnsembleSpectrum & ana::EnsembleSpectrum::operator*= ( const EnsembleRatio rhs)

Definition at line 125 of file EnsembleSpectrum.cxx.

126  {
127  fNom *= rhs.Nominal();
128  assert(rhs.NUniverses() == fUnivs.size());
129  for(unsigned int i = 0; i < fUnivs.size(); ++i) fUnivs[i] *= rhs.Universe(i);
130  return *this;
131  }
std::vector< Spectrum > fUnivs
EnsembleSpectrum ana::EnsembleSpectrum::operator+ ( const EnsembleSpectrum rhs) const

Definition at line 100 of file EnsembleSpectrum.cxx.

101  {
102  EnsembleSpectrum ret = *this;
103  ret += rhs;
104  return ret;
105  }
EnsembleSpectrum(SpectrumLoaderBase &loader, const HistAxis &axis, const SpillCut &spillcut, const Cut &cut, const std::vector< SystShifts > &univ_shifts, const Var &cv_wei=kUnweighted)
EnsembleSpectrum & ana::EnsembleSpectrum::operator+= ( const EnsembleSpectrum rhs)

Definition at line 91 of file EnsembleSpectrum.cxx.

92  {
93  fNom += rhs.fNom;
94  assert(fUnivs.size() == rhs.fUnivs.size());
95  for(unsigned int i = 0; i < fUnivs.size(); ++i) fUnivs[i] += rhs.fUnivs[i];
96  return *this;
97  }
std::vector< Spectrum > fUnivs
EnsembleSpectrum ana::EnsembleSpectrum::operator- ( const EnsembleSpectrum rhs) const

Definition at line 117 of file EnsembleSpectrum.cxx.

118  {
119  EnsembleSpectrum ret = *this;
120  ret -= rhs;
121  return ret;
122  }
EnsembleSpectrum(SpectrumLoaderBase &loader, const HistAxis &axis, const SpillCut &spillcut, const Cut &cut, const std::vector< SystShifts > &univ_shifts, const Var &cv_wei=kUnweighted)
EnsembleSpectrum & ana::EnsembleSpectrum::operator-= ( const EnsembleSpectrum rhs)

Definition at line 108 of file EnsembleSpectrum.cxx.

109  {
110  fNom -= rhs.fNom;
111  assert(fUnivs.size() == rhs.fUnivs.size());
112  for(unsigned int i = 0; i < fUnivs.size(); ++i) fUnivs[i] -= rhs.fUnivs[i];
113  return *this;
114  }
std::vector< Spectrum > fUnivs
EnsembleSpectrum ana::EnsembleSpectrum::operator/ ( const EnsembleRatio rhs) const

Definition at line 151 of file EnsembleSpectrum.cxx.

152  {
153  EnsembleSpectrum ret = *this;
154  ret /= rhs;
155  return ret;
156  }
EnsembleSpectrum(SpectrumLoaderBase &loader, const HistAxis &axis, const SpillCut &spillcut, const Cut &cut, const std::vector< SystShifts > &univ_shifts, const Var &cv_wei=kUnweighted)
EnsembleSpectrum & ana::EnsembleSpectrum::operator/= ( const EnsembleRatio rhs)

Definition at line 142 of file EnsembleSpectrum.cxx.

143  {
144  fNom /= rhs.Nominal();
145  assert(rhs.NUniverses() == fUnivs.size());
146  for(unsigned int i = 0; i < fUnivs.size(); ++i) fUnivs[i] /= rhs.Universe(i);
147  return *this;
148  }
std::vector< Spectrum > fUnivs
double ana::EnsembleSpectrum::POT ( ) const
inline

Definition at line 41 of file EnsembleSpectrum.h.

41 {return fNom.POT();}
double POT() const
Definition: Spectrum.h:289
void ana::EnsembleSpectrum::SaveTo ( TDirectory *  dir) const

Definition at line 159 of file EnsembleSpectrum.cxx.

160  {
161  TDirectory* tmp = gDirectory;
162  dir->cd();
163 
164  TObjString("EnsembleSpectrum").Write("type");
165 
166  fNom.SaveTo(dir->mkdir("nom"));
167  for(unsigned int i = 0; i < fUnivs.size(); ++i){
168  fUnivs[i].SaveTo(dir->mkdir(("univ"+std::to_string(i)).c_str()));
169  }
170 
171  tmp->cd();
172  }
tuple dir
Definition: dropbox.py:28
std::vector< Spectrum > fUnivs
std::string to_string(WindowPattern const &pattern)
void SaveTo(TDirectory *dir) const
void ana::EnsembleSpectrum::Scale ( double  c)

Definition at line 84 of file EnsembleSpectrum.cxx.

85  {
86  fNom.Scale(c);
87  for(Spectrum& s: fUnivs) s.Scale(c);
88  }
void Scale(double c)
Multiply this spectrum by a constant c.
std::vector< Spectrum > fUnivs
then echo File list $list not found else cat $list while read file do echo $file sed s
Definition: file_to_url.sh:60
Spectrum ana::EnsembleSpectrum::Universe ( unsigned int  i) const
inline

Definition at line 35 of file EnsembleSpectrum.h.

36  {
37  assert(i < fUnivs.size());
38  return fUnivs[i];
39  }
std::vector< Spectrum > fUnivs

Member Data Documentation

Spectrum ana::EnsembleSpectrum::fNom
protected

Definition at line 74 of file EnsembleSpectrum.h.

std::vector<Spectrum> ana::EnsembleSpectrum::fUnivs
protected

Definition at line 75 of file EnsembleSpectrum.h.


The documentation for this class was generated from the following files: