CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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) :
20  ensembleSize_(M),
23  localYbin_(0),
24  stripsPerBin_(0),
25  maxEvents_(max),
26  methods_(methods) {};
27 
28  //measurement constructor
29  LA_Filler_Fitter(int methods, bool layer, bool module, float localybin, unsigned stripbin, unsigned max) :
30  ensembleSize_(0),
32  byLayer_(layer),byModule_(module),
33  localYbin_(localybin),
34  stripsPerBin_(stripbin),
35  maxEvents_(max),
36  methods_(methods) {};
37 
38  enum Method { WIDTH =1<<0, FIRST_METHOD=1<<0,
39  PROB1 =1<<1,
40  AVGV2 =1<<2,
41  AVGV3 =1<<3,
42  RMSV2 =1<<4,
43  RMSV3 =1<<5, LAST_METHOD=1<<5};
44 
45  static std::string method(Method m,bool fit=true) {
46  switch(m) {
47  case WIDTH: return "_width_prof";
48  case PROB1: return fit? SymmetryFit::name(method(m,0)): "_prob_w1" ;
49  case AVGV2: return fit? SymmetryFit::name(method(m,0)): "_avg_var_w2";
50  case AVGV3: return fit? SymmetryFit::name(method(m,0)): "_avg_var_w3";
51  case RMSV2: return fit? SymmetryFit::name(method(m,0)): "_rms_var_w2";
52  case RMSV3: return fit? SymmetryFit::name(method(m,0)): "_rms_var_w3";
53  default: return "_UNKNOWN";
54  }
55  }
56 
57  struct Result {
58  std::pair<float,float> reco,measured,calMeasured;
59  float field,chi2;
60  unsigned ndof,entries;
61  Result() : reco(1000,0), measured(1000,0), calMeasured(1000,0),
62  field(0), chi2(0), ndof(0), entries(0) {}
63  };
64 
65  struct EnsembleSummary {
66  std::pair<float,float> meanMeasured, sigmaMeasured, meanUncertainty, pull;
67  float truth;
68  unsigned samples;
70  truth(0), samples(0) {}
71  };
72 
73  //Located in src/LA_Filler.cc
74  void fill(TTree*, Book&) const;
75  void fill_one_cluster(Book&,
76  const poly<std::string>&,
77  const unsigned, const float, const float, const float, const float) const;
78  poly<std::string> allAndOne(const unsigned width) const;
79  poly<std::string> varWidth(const unsigned width) const;
80  poly<std::string> granularity(const SiStripDetId, const float, const Long64_t, const float, const unsigned) const;
81  static std::string subdetLabel(const SiStripDetId);
82  static std::string moduleLabel(const SiStripDetId);
83  static std::string layerLabel(const SiStripDetId);
84  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 );}
85 
86  //Located in src/LA_Fitter.cc
87  static void fit(Book& book) { make_and_fit_symmchi2(book); fit_width_profile(book); }
88  static void make_and_fit_symmchi2(Book&);
89  static void fit_width_profile(Book&);
90  static TH1* rms_profile(const std::string, const TProfile* const);
91  static TH1* subset_probability(const std::string name, const TH1* const , const TH1* const );
92  static unsigned find_rebin(const TH1* const);
93 
94  //Located in src/LA_Results.cc
95  void summarize_ensembles(Book&) const;
96  static Result result(Method, const std::string name, const Book&);
97  static std::map< std::string, Result > layer_results(const Book&, const Method);
98  static std::map< uint32_t, Result > module_results(const Book&, const Method);
99  static std::map< std::string, std::vector<Result> > ensemble_results(const Book&, const Method );
100  static std::map< std::string, std::vector<EnsembleSummary> > ensemble_summary(const Book& );
101  static std::pair<std::pair<float,float>, std::pair<float,float> > offset_slope(const std::vector<EnsembleSummary>&);
102  static float pull(const std::vector<EnsembleSummary>&);
103 
104  private:
105 
107  const double ensembleLow_, ensembleUp_;
108  const bool byLayer_, byModule_;
109  const float localYbin_;
110  const unsigned stripsPerBin_;
111  const Long64_t maxEvents_;
112  const int methods_;
113 };
114 
115 std::ostream& operator<<(std::ostream&, const LA_Filler_Fitter::Result&);
116 std::ostream& operator<<(std::ostream&, const LA_Filler_Fitter::EnsembleSummary&);
117 
118 #endif
LA_Filler_Fitter(int methods, int M, int N, double low, double up, unsigned max)
std::pair< float, float > pull
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
std::pair< float, float > meanMeasured
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:72
static std::pair< std::pair< float, float >, std::pair< float, float > > offset_slope(const std::vector< EnsembleSummary > &)
Definition: LA_Results.cc:130
static std::map< std::string, std::vector< EnsembleSummary > > ensemble_summary(const Book &)
Definition: LA_Results.cc:105
std::ostream & operator<<(std::ostream &out, const ALILine &li)
Definition: ALILine.cc:187
void fill(TTree *, Book &) const
Definition: LA_Filler.cc:8
static std::string moduleLabel(const SiStripDetId)
Definition: LA_Filler.cc:90
std::pair< float, float > measured
LA_Filler_Fitter(int methods, bool layer, bool module, float localybin, unsigned stripbin, unsigned max)
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
static std::string layerLabel(const SiStripDetId)
Definition: LA_Filler.cc:92
const float localYbin_
poly< std::string > varWidth(const unsigned width) const
Definition: LA_Filler.cc:64
Definition: poly.h:11
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
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
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
std::pair< float, float > meanUncertainty
static std::string subdetLabel(const SiStripDetId)
Definition: LA_Filler.cc:88
std::pair< float, float > reco
static float pull(const std::vector< EnsembleSummary > &)
Definition: LA_Results.cc:152
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:45
std::pair< float, float > calMeasured
volatile std::atomic< bool > shutdown_flag false
const double ensembleLow_
Definition: vlib.h:208
poly< std::string > allAndOne(const unsigned width) const
Definition: LA_Filler.cc:60