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