CMS 3D CMS Logo

LA_Filler_Fitter.h
Go to the documentation of this file.
1 #ifndef LA_FILLER_FITTER_H
2 #define LA_FILLER_FITTER_H
3 
4 #include <string>
5 #include <vector>
6 #include <map>
10 #include <TTree.h>
11 #include "SymmetryFit.h"
12 class TProfile;
13 
15 
16  public:
17 
18  //ensemble constructor
19  LA_Filler_Fitter(int methods, int M, int N, double low, double up, unsigned max, const TrackerTopology *tTopo) :
20  ensembleSize_(M),
23  localYbin_(0),
24  stripsPerBin_(0),
25  maxEvents_(max),
26  methods_(methods),
27  tTopo_(tTopo) {}
28 
29  //measurement constructor
30  LA_Filler_Fitter(int methods, bool layer, bool module, float localybin, unsigned stripbin, unsigned max, const TrackerTopology *tTopo) :
31  ensembleSize_(0),
33  byLayer_(layer),byModule_(module),
34  localYbin_(localybin),
35  stripsPerBin_(stripbin),
36  maxEvents_(max),
37  methods_(methods),
38  tTopo_(tTopo) {}
39 
40  enum Method { WIDTH =1<<0, FIRST_METHOD=1<<0,
41  PROB1 =1<<1,
42  AVGV2 =1<<2,
43  AVGV3 =1<<3,
44  RMSV2 =1<<4,
45  RMSV3 =1<<5, LAST_METHOD=1<<5};
46 
47  static std::string method(Method m,bool fit=true) {
48  switch(m) {
49  case WIDTH: return "_width_prof";
50  case PROB1: return fit? SymmetryFit::name(method(m,false)): "_prob_w1" ;
51  case AVGV2: return fit? SymmetryFit::name(method(m,false)): "_avg_var_w2";
52  case AVGV3: return fit? SymmetryFit::name(method(m,false)): "_avg_var_w3";
53  case RMSV2: return fit? SymmetryFit::name(method(m,false)): "_rms_var_w2";
54  case RMSV3: return fit? SymmetryFit::name(method(m,false)): "_rms_var_w3";
55  default: return "_UNKNOWN";
56  }
57  }
58 
59  struct Result {
60  std::pair<float,float> reco,measured,calMeasured;
61  float field,chi2;
62  unsigned ndof,entries;
63  Result() : reco(1000,0), measured(1000,0), calMeasured(1000,0),
64  field(0), chi2(0), ndof(0), entries(0) {}
65  };
66 
67  struct EnsembleSummary {
68  std::pair<float,float> meanMeasured, sigmaMeasured, meanUncertainty, pull;
69  float truth;
70  unsigned samples;
71  EnsembleSummary() : meanMeasured(0,0), sigmaMeasured(0,0), meanUncertainty(0,0), pull(0,0),
72  truth(0), samples(0) {}
73  };
74 
75  //Located in src/LA_Filler.cc
76  void fill(TTree*, Book&) const;
77  void fill_one_cluster(Book&,
78  const poly<std::string>&,
79  const unsigned, const float, const float, const float, const float) const;
80  poly<std::string> allAndOne(const unsigned width) const;
81  poly<std::string> varWidth(const unsigned width) const;
82  poly<std::string> granularity(const SiStripDetId, const float, const Long64_t, const float, const unsigned) const;
83  static std::string subdetLabel(const SiStripDetId);
84  static std::string moduleLabel(const SiStripDetId);
85  std::string layerLabel(const SiStripDetId) const;
86  static unsigned layer_index(bool TIB, bool stereo, unsigned layer) { return layer + (TIB?0:6) +(stereo?1:0) + ( (layer>2)?1:(layer==1)?-1:0 );}
87 
88  //Located in src/LA_Fitter.cc
89  static void fit(Book& book) { make_and_fit_symmchi2(book); fit_width_profile(book); }
90  static void make_and_fit_symmchi2(Book&);
91  static void fit_width_profile(Book&);
92  static TH1* rms_profile(const std::string, const TProfile* const);
93  static TH1* subset_probability(const std::string name, const TH1* const , const TH1* const );
94  static unsigned find_rebin(const TH1* const);
95 
96  //Located in src/LA_Results.cc
97  void summarize_ensembles(Book&) const;
98  static Result result(Method, const std::string name, const Book&);
99  static std::map< std::string, Result > layer_results(const Book&, const Method);
100  static std::map< uint32_t, Result > module_results(const Book&, const Method);
101  static std::map< std::string, std::vector<Result> > ensemble_results(const Book&, const Method );
102  static std::map< std::string, std::vector<EnsembleSummary> > ensemble_summary(const Book& );
103  static std::pair<std::pair<float,float>, std::pair<float,float> > offset_slope(const std::vector<EnsembleSummary>&);
104  static float pull(const std::vector<EnsembleSummary>&);
105 
106  private:
107 
109  const double ensembleLow_, ensembleUp_;
110  const bool byLayer_, byModule_;
111  const float localYbin_;
112  const unsigned stripsPerBin_;
113  const Long64_t maxEvents_;
114  const int methods_;
116 };
117 
118 std::ostream& operator<<(std::ostream&, const LA_Filler_Fitter::Result&);
119 std::ostream& operator<<(std::ostream&, const LA_Filler_Fitter::EnsembleSummary&);
120 
121 #endif
Definition: BitonicSort.h:8
std::string layerLabel(const SiStripDetId) const
Definition: LA_Filler.cc:93
static std::map< std::string, Result > layer_results(const Book &, const Method)
Definition: LA_Results.cc:63
static unsigned layer_index(bool TIB, bool stereo, unsigned layer)
static void fit(Book &book)
static std::string name(std::string base)
Definition: SymmetryFit.h:14
static TH1 * subset_probability(const std::string name, const TH1 *const , const TH1 *const )
Definition: LA_Fitter.cc:108
std::pair< float, float > sigmaMeasured
static void fit_width_profile(Book &)
Definition: LA_Fitter.cc:10
const Long64_t maxEvents_
poly< std::string > granularity(const SiStripDetId, const float, const Long64_t, const float, const unsigned) const
Definition: LA_Filler.cc:73
static std::pair< std::pair< float, float >, std::pair< float, float > > offset_slope(const std::vector< EnsembleSummary > &)
Definition: LA_Results.cc:132
static std::map< std::string, std::vector< EnsembleSummary > > ensemble_summary(const Book &)
Definition: LA_Results.cc:107
void fill(TTree *, Book &) const
Definition: LA_Filler.cc:9
static std::string moduleLabel(const SiStripDetId)
Definition: LA_Filler.cc:91
std::pair< float, float > measured
void summarize_ensembles(Book &) const
Definition: LA_Results.cc:83
static Result result(Method, const std::string name, const Book &)
Definition: LA_Results.cc:14
const float localYbin_
poly< std::string > varWidth(const unsigned width) const
Definition: LA_Filler.cc:65
Definition: poly.h:11
std::ostream & operator<<(std::ostream &, const LA_Filler_Fitter::Result &)
Definition: LA_Results.cc:167
static std::string method(Method m, bool fit=true)
static unsigned find_rebin(const TH1 *const)
Definition: LA_Fitter.cc:81
Detector identifier class for the strip tracker.
Definition: SiStripDetId.h:17
LA_Filler_Fitter(int methods, int M, int N, double low, double up, unsigned max, const TrackerTopology *tTopo)
static std::map< std::string, std::vector< Result > > ensemble_results(const Book &, const Method)
Definition: LA_Results.cc:73
#define N
Definition: blowfish.cc:9
static TH1 * rms_profile(const std::string, const TProfile *const)
Definition: LA_Fitter.cc:95
const TrackerTopology * tTopo_
static std::map< uint32_t, Result > module_results(const Book &, const Method)
Definition: LA_Results.cc:51
const double ensembleUp_
static void make_and_fit_symmchi2(Book &)
Definition: LA_Fitter.cc:34
static std::string subdetLabel(const SiStripDetId)
Definition: LA_Filler.cc:89
std::pair< float, float > reco
static float pull(const std::vector< EnsembleSummary > &)
Definition: LA_Results.cc:156
const unsigned stripsPerBin_
Definition: Book.h:16
void fill_one_cluster(Book &, const poly< std::string > &, const unsigned, const float, const float, const float, const float) const
Definition: LA_Filler.cc:46
std::pair< float, float > calMeasured
const double ensembleLow_
LA_Filler_Fitter(int methods, bool layer, bool module, float localybin, unsigned stripbin, unsigned max, const TrackerTopology *tTopo)
Definition: vlib.h:208
poly< std::string > allAndOne(const unsigned width) const
Definition: LA_Filler.cc:61