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

#include <MedianSurface.h>

Inheritance diagram for ana::MedianSurface:
ana::Surface

Public Member Functions

 MedianSurface (const std::vector< Surface > &throws)
 
void DrawEnsemble (TH2 *fc, Color_t color=kGray)
 
void DrawBand (TH2 *fc)
 
void SaveTo (TDirectory *dir) const
 
- Public Member Functions inherited from ana::Surface
 Surface (const IExperiment *expt, osc::IOscCalcAdjustable *calc, const IFitVar *xvar, int nbinsx, double xmin, double xmax, const IFitVar *yvar, int nbinsy, double ymin, double ymax, const std::vector< const IFitVar * > &profVars={}, const std::vector< const ISyst * > &profSysts={}, const std::map< const IFitVar *, std::vector< double >> &seedPts={}, const std::vector< SystShifts > &systSeedPts={}, bool parallel=false, Fitter::Precision prec=Fitter::kNormal)
 
 Surface (const IExperiment *expt, osc::IOscCalcAdjustable *calc, const FitAxis &xax, const FitAxis &yax, const std::vector< const IFitVar * > &profVars={}, const std::vector< const ISyst * > &profSysts={}, const std::map< const IFitVar *, std::vector< double >> &seedPts={}, const std::vector< SystShifts > &systSeedPts={}, bool parallel=false, Fitter::Precision prec=Fitter::kNormal)
 
void Draw () const
 Draw the surface itself. More...
 
void DrawBestFit (Color_t color, Int_t marker=kFullCircle) const
 Draw the best fit point. More...
 
double MinChi () const
 
double GetMinX () const
 
double GetMinY () const
 
void DrawContour (TH2 *fc, Style_t style, Color_t color, double minchi=-1)
 
TH2F * ToTH2 (double minchi=-1) const
 
void SetTitle (const char *str)
 
std::vector< TH2 * > GetProfiledHists ()
 Maps of the values taken on by the profiled parameters. More...
 
std::vector< TH2 * > GetMarginalizedHists ()
 Deprecated. Retained for backwards compatibility. More...
 
std::vector< TGraph * > GetGraphs (TH2 *fc, double minchi=-1)
 For expert use, custom painting of contours. More...
 
void SaveTo (TDirectory *dir) const
 

Static Public Member Functions

static std::unique_ptr
< MedianSurface
LoadFrom (TDirectory *dir)
 
- Static Public Member Functions inherited from ana::Surface
static std::unique_ptr< SurfaceLoadFrom (TDirectory *dir)
 
static std::unique_ptr< SurfaceLoadFromMulti (const std::vector< TFile * > &files, const std::string &label)
 
static std::unique_ptr< SurfaceLoadFromMulti (const std::string &wildcard, const std::string &label)
 

Protected Attributes

std::vector< SurfacefThrows
 
TH2F * fHistUp1
 
TH2F * fHistDn1
 
TH2F * fHistUp2
 
TH2F * fHistDn2
 
- Protected Attributes inherited from ana::Surface
bool fParallel
 
Fitter::Precision fPrec
 
double fMinChi
 
double fMinX
 
double fMinY
 
TH2F * fHist
 
bool fLogX
 
bool fLogY
 
std::vector< TH2 * > fProfHists
 
std::vector< double > fSeedValues
 
std::vector< int > fBinMask
 

Additional Inherited Members

- Protected Member Functions inherited from ana::Surface
 Surface ()
 
void EnsureAxes () const
 
void CheckMask (const std::string &func) const
 
void FillSurface (const std::string &progTitle, const IExperiment *expt, osc::IOscCalcAdjustable *calc, const FitAxis &xax, const FitAxis &yax, const std::vector< const IFitVar * > &profVars, const std::vector< const ISyst * > &profSysts, const std::map< const IFitVar *, std::vector< double >> &seedPts, const std::vector< SystShifts > &systSeedPts)
 
void FillSurfacePoint (const IExperiment *expt, osc::IOscCalcAdjustable *calc, const FitAxis &xax, double x, const FitAxis &yax, double y, const std::vector< const IFitVar * > &profVars, const std::vector< const ISyst * > &profSysts, const std::map< const IFitVar *, std::vector< double >> &seedPts, const std::vector< SystShifts > &systSeedPts)
 

Detailed Description

Definition at line 9 of file MedianSurface.h.

Constructor & Destructor Documentation

ana::MedianSurface::MedianSurface ( const std::vector< Surface > &  throws)

Definition at line 15 of file MedianSurface.cxx.

16  : fThrows(throws)
17  {
18  assert(!throws.empty());
19 
20  for(const Surface& s: throws){
21  assert(s.fHist->GetNbinsX() == throws[0].fHist->GetNbinsX());
22  assert(s.fHist->GetNbinsY() == throws[0].fHist->GetNbinsY());
23  // TODO check min and max match at least
24  }
25 
26  fLogX = throws[0].fLogX;
27  fLogY = throws[0].fLogY;
28 
29  // Think about what to do with these
30  fMinX = -1;
31  fMinY = -1;
32  // Though this is right
33  fMinChi = 0;
34 
35  fHist = new TH2F(*throws[0].fHist);
36  fHist->Reset();
37  fHistUp1 = new TH2F(*fHist);
38  fHistDn1 = new TH2F(*fHist);
39  fHistUp2 = new TH2F(*fHist);
40  fHistDn2 = new TH2F(*fHist);
41 
42  for(int ix = 0; ix < fHist->GetNbinsX()+2; ++ix){
43  for(int iy = 0; iy < fHist->GetNbinsY()+2; ++iy){
44  std::vector<float> chis;
45  chis.reserve(throws.size());
46  for(const Surface& s: throws){
47  chis.push_back(s.fHist->GetBinContent(ix, iy));
48  }
49  std::sort(chis.begin(), chis.end());
50 
51  // TODO think about off-by-one errors and interpolation here
52  // Median
53  fHist->SetBinContent(ix, iy, chis[throws.size()/2]);
54  // One sigma bounds
55  const double tail1 = (1-0.6827)/2;
56  fHistUp1->SetBinContent(ix, iy, chis[throws.size()*tail1]);
57  fHistDn1->SetBinContent(ix, iy, chis[throws.size()*(1-tail1)]);
58 
59  // Two sigma bounds
60  const double tail2 = (1-0.9545)/2;
61  fHistUp2->SetBinContent(ix, iy, chis[throws.size()*tail2]);
62  fHistDn2->SetBinContent(ix, iy, chis[(throws.size()-1)*(1-tail2)]);
63  } // end for iy
64  } // end for ix
65  }
std::vector< Surface > fThrows
Definition: MedianSurface.h:20
std::size_t size(FixedBins< T, C > const &) noexcept
Definition: FixedBins.h:561
then echo File list $list not found else cat $list while read file do echo $file sed s
Definition: file_to_url.sh:60
bool empty(FixedBins< T, C > const &) noexcept
Definition: FixedBins.h:555

Member Function Documentation

void ana::MedianSurface::DrawBand ( TH2 *  fc)

Definition at line 86 of file MedianSurface.cxx.

87  {
88  EnsureAxes();
89 
90  TH2F surf1(*fHist);
91  TH2F surf2(*fHist);
92 
93  for(int ix = 0; ix < fHist->GetNbinsX()+2; ++ix){
94  for(int iy = 0; iy < fHist->GetNbinsY()+2; ++iy){
95  const double c = fc->GetBinContent(ix, iy);
96 
97  const double u1 = fHistUp1->GetBinContent(ix, iy);
98  const double d1 = fHistDn1->GetBinContent(ix, iy);
99  surf1.SetBinContent(ix, iy, std::max(u1-c, c-d1));
100 
101  const double u2 = fHistUp2->GetBinContent(ix, iy);
102  const double d2 = fHistDn2->GetBinContent(ix, iy);
103  surf2.SetBinContent(ix, iy, std::max(u2-c, c-d2));
104  }
105  }
106 
107  const double level = 0;
108 
109  // I haven't been able to figure out how to draw these filled properly. cont0 does it, but can't handle countours that reach the edge of the space
110  surf2.SetContour(1, &level);
111  surf2.SetLineColor(kYellow);
112  surf2.DrawCopy("cont3 same");
113 
114  surf1.SetContour(1, &level);
115  surf1.SetLineColor(kGreen);
116  surf1.SetFillColor(kGreen+2);
117  surf1.DrawCopy("cont3 same");
118 
119  gPad->Update();
120  }
void ana::MedianSurface::DrawEnsemble ( TH2 *  fc,
Color_t  color = kGray 
)

Definition at line 68 of file MedianSurface.cxx.

69  {
70  EnsureAxes();
71 
72  for(Surface& s: fThrows){
73  std::vector<TGraph*> gs = s.GetGraphs(fc, -1);
74 
75  for(TGraph* g: gs){
76  g->SetLineWidth(1);
77  g->SetLineColor(color);
78  g->Draw("l");
79  }
80  }
81 
82  gPad->Update();
83  }
std::vector< Surface > fThrows
Definition: MedianSurface.h:20
BEGIN_PROLOG g
then echo File list $list not found else cat $list while read file do echo $file sed s
Definition: file_to_url.sh:60
std::unique_ptr< MedianSurface > ana::MedianSurface::LoadFrom ( TDirectory *  dir)
static

Definition at line 137 of file MedianSurface.cxx.

138  {
139  DontAddDirectory guard;
140 
141  TObjString* tag = (TObjString*)dir->Get("type");
142  assert(tag);
143  assert(tag->GetString() == "MedianSurface");
144 
145  std::vector<Surface> surfs;
146  for(unsigned int i = 0; ; ++i){
147  TDirectory* surfdir = dir->GetDirectory(TString::Format("surf%d", i));
148  if(!surfdir) break; // we got all of them
149  surfs.push_back(*ana::LoadFrom<Surface>(surfdir));
150  }
151 
152  return std::make_unique<MedianSurface>(surfs);
153  }
tuple dir
Definition: dropbox.py:28
void ana::MedianSurface::SaveTo ( TDirectory *  dir) const

Definition at line 123 of file MedianSurface.cxx.

124  {
125  TDirectory* tmp = gDirectory;
126  dir->cd();
127  TObjString("MedianSurface").Write("type");
128 
129  for(unsigned int i = 0; i < fThrows.size(); ++i){
130  fThrows[i].SaveTo(dir->mkdir(TString::Format("surf%d", i)));
131  }
132 
133  tmp->cd();
134  }
std::vector< Surface > fThrows
Definition: MedianSurface.h:20
tuple dir
Definition: dropbox.py:28

Member Data Documentation

TH2F * ana::MedianSurface::fHistDn1
protected

Definition at line 22 of file MedianSurface.h.

TH2F * ana::MedianSurface::fHistDn2
protected

Definition at line 22 of file MedianSurface.h.

TH2F* ana::MedianSurface::fHistUp1
protected

Definition at line 22 of file MedianSurface.h.

TH2F * ana::MedianSurface::fHistUp2
protected

Definition at line 22 of file MedianSurface.h.

std::vector<Surface> ana::MedianSurface::fThrows
protected

Definition at line 20 of file MedianSurface.h.


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