CMS 3D CMS Logo

LA_Results.cc
Go to the documentation of this file.
3 
4 #include <cmath>
5 #include <boost/foreach.hpp>
6 #include <boost/regex.hpp>
7 #include <boost/lexical_cast.hpp>
8 #include <boost/algorithm/string/erase.hpp>
9 #include <TF1.h>
10 #include <TGraphErrors.h>
11 #include <TProfile.h>
12 
14 result(Method m, const std::string name, const Book& book) {
15  Result p;
16  const std::string base = boost::erase_all_copy(name,method(m));
17 
18  const TH1* const h = book[name];
19  const TH1* const reco = book[base+"_reconstruction"];
20  const TH1* const field = book[base+"_field"];
21 
22  if(reco) p.reco = std::make_pair<float,float>( reco->GetMean(), reco->GetMeanError());
23  if(field) p.field = field->GetMean();
24  if(h) { switch(m) {
25  case WIDTH: {
26  const TF1*const f = h->GetFunction("LA_profile_fit"); if(!f) break;
27  p.measured = std::make_pair<float,float>(f->GetParameter(0),f->GetParError(0));
28  p.chi2 = f->GetChisquare();
29  p.ndof = f->GetNDF();
30  p.entries = (unsigned)(h->GetEntries());
31  break;
32  }
33  case PROB1: case AVGV2: case AVGV3: case RMSV2: case RMSV3: /*case MULTI:*/ {
34  const TF1*const f = h->GetFunction("SymmetryFit"); if(!f) break;
35  p.measured = std::make_pair<float,float>( p.reco.first + f->GetParameter(0), f->GetParameter(1) );
36  p.chi2 = f->GetParameter(2);
37  p.ndof = (unsigned) (f->GetParameter(3));
38  p.entries =
39  (m&PROB1) ? (unsigned) book[base+"_w1"]->GetEntries() :
40  (m&(AVGV2|RMSV2)) ? (unsigned) book[base+method(AVGV2,false)]->GetEntries() :
41  (m&(AVGV3|RMSV3)) ? (unsigned) book[base+method(AVGV3,false)]->GetEntries() : 0 ;
42  break;
43  }
44  default: break;
45  }
46  }
47  return p;
48 }
49 
50 std::map<uint32_t,LA_Filler_Fitter::Result> LA_Filler_Fitter::
51 module_results( const Book& book, const Method m) {
52  std::map<uint32_t,Result> results;
53  for(Book::const_iterator it = book.begin(".*_module\\d*"+method(m)); it!=book.end(); ++it ) {
54  const uint32_t detid = boost::lexical_cast<uint32_t>( boost::regex_replace( it->first,
55  boost::regex(".*_module(\\d*)_.*"),
56  std::string("\\1")));
57  results[detid] = result(m,it->first,book);
58  }
59  return results;
60 }
61 
62 std::map<std::string,LA_Filler_Fitter::Result> LA_Filler_Fitter::
63 layer_results( const Book& book, const Method m) {
64  std::map<std::string,Result> results;
65  for(Book::const_iterator it = book.begin(".*layer\\d.*"+method(m)); it!=book.end(); ++it ) {
66  const std::string name = boost::erase_all_copy(it->first,method(m));
67  results[name] = result(m,it->first,book);
68  }
69  return results;
70 }
71 
72 std::map<std::string, std::vector<LA_Filler_Fitter::Result> > LA_Filler_Fitter::
73 ensemble_results( const Book& book, const Method m) {
74  std::map<std::string, std::vector<Result> > results;
75  for(Book::const_iterator it = book.begin(".*_sample.*"+method(m)); it!=book.end(); ++it ) {
76  const std::string name = boost::regex_replace(it->first,boost::regex("sample\\d*_"),"");
77  results[name].push_back(result(m,it->first,book));
78  }
79  return results;
80 }
81 
83 summarize_ensembles(Book& book) const {
84  typedef std::map<std::string, std::vector<Result> > results_t;
85  results_t results;
86  for(int m = FIRST_METHOD; m <= LAST_METHOD; m<<=1)
87  if(methods_ & m) { results_t g = ensemble_results(book,(Method)(m)); results.insert(g.begin(),g.end());}
88 
89  BOOST_FOREACH(const results_t::value_type group, results) {
90  const std::string name = group.first;
91  BOOST_FOREACH(const Result p, group.second) {
92  const float pad = (ensembleUp_-ensembleLow_)/10;
93  book.fill( p.reco.first, name+"_ensembleReco", 12*ensembleBins_, ensembleLow_-pad, ensembleUp_+pad );
94  book.fill( p.measured.first, name+"_measure", 12*ensembleBins_, ensembleLow_-pad, ensembleUp_+pad );
95  book.fill( p.measured.second, name+"_merr", 500, 0, 0.01);
96  book.fill( (p.measured.first-p.reco.first)/p.measured.second, name+"_pull", 500, -10,10);
97  }
98  //Need our own copy for thread safety
99  TF1 gaus("mygaus","gaus");
100  book[name+"_measure"]->Fit(&gaus,"LLQ");
101  book[name+"_merr"]->Fit(&gaus,"LLQ");
102  book[name+"_pull"]->Fit(&gaus,"LLQ");
103  }
104 }
105 
106 std::map<std::string, std::vector<LA_Filler_Fitter::EnsembleSummary> > LA_Filler_Fitter::
107 ensemble_summary(const Book& book) {
108  std::map<std::string, std::vector<EnsembleSummary> > summary;
109  for(Book::const_iterator it = book.begin(".*_ensembleReco"); it!=book.end(); ++it) {
110  const std::string base = boost::erase_all_copy(it->first,"_ensembleReco");
111 
112  const TH1*const reco = it->second;
113  const TH1*const measure = book[base+"_measure"];
114  const TH1*const merr = book[base+"_merr"];
115  const TH1*const pull = book[base+"_pull"];
116 
118  s.samples = reco->GetEntries();
119  s.truth = reco->GetMean();
120  s.meanMeasured = std::make_pair<float,float>( measure->GetFunction("gaus")->GetParameter(1), measure->GetFunction("gaus")->GetParError(1) );
121  s.sigmaMeasured = std::make_pair<float,float>( measure->GetFunction("gaus")->GetParameter(2), measure->GetFunction("gaus")->GetParError(2) );
122  s.meanUncertainty = std::make_pair<float,float>( merr->GetFunction("gaus")->GetParameter(1), merr->GetFunction("gaus")->GetParError(1) );
123  s.pull = std::make_pair<float,float>( pull->GetFunction("gaus")->GetParameter(2), pull->GetFunction("gaus")->GetParError(2) );
124 
125  const std::string name = boost::regex_replace(base,boost::regex("ensembleBin\\d*_"),"");
126  summary[name].push_back(s);
127  }
128  return summary;
129 }
130 
131 std::pair<std::pair<float,float>, std::pair<float,float> > LA_Filler_Fitter::
132 offset_slope(const std::vector<LA_Filler_Fitter::EnsembleSummary>& ensembles) {
133  try {
134  std::vector<float> x,y,xerr,yerr;
135  BOOST_FOREACH(EnsembleSummary ensemble, ensembles) {
136  x.push_back(ensemble.truth);
137  xerr.push_back(0);
138  y.push_back(ensemble.meanMeasured.first);
139  yerr.push_back(ensemble.meanMeasured.second);
140  }
141  TGraphErrors graph(x.size(),&(x[0]),&(y[0]),&(xerr[0]),&(yerr[0]));
142  //Need our own copy for thread safety
143  TF1 pol1("mypol1","pol1");
144  graph.Fit(&pol1);
145  const TF1*const fit = graph.GetFunction("pol1");
146 
147  return std::make_pair( std::make_pair(fit->GetParameter(0), fit->GetParError(0)),
148  std::make_pair(fit->GetParameter(1), fit->GetParError(1)) );
149  } catch(edm::Exception e) {
150  std::cerr << "Fitting Line Failed " << std::endl << e << std::endl;
151  return std::make_pair( std::make_pair(0,0), std::make_pair(0,0));
152  }
153 }
154 
155 float LA_Filler_Fitter::
156 pull(const std::vector<LA_Filler_Fitter::EnsembleSummary>& ensembles) {
157  float p(0),w(0);
158  BOOST_FOREACH(EnsembleSummary ensemble, ensembles) {
159  const float unc2 = pow(ensemble.pull.second,2);
160  p+= ensemble.pull.first / unc2;
161  w+= 1/unc2;
162  }
163  return p/w;
164 }
165 
166 
167 std::ostream& operator<<(std::ostream& strm, const LA_Filler_Fitter::Result& r) {
168  return strm << r.reco.first <<"\t"<< r.reco.second <<"\t"
169  << r.measured.first <<"\t"<< r.measured.second <<"\t"
170  << r.calMeasured.first <<"\t"<< r.calMeasured.second <<"\t"
171  << r.field <<"\t"
172  << r.chi2 <<"\t"
173  << r.ndof <<"\t"
174  << r.entries;
175 }
176 
177 std::ostream& operator<<(std::ostream& strm, const LA_Filler_Fitter::EnsembleSummary& e) {
178  return strm << e.truth <<"\t"
179  << e.meanMeasured.first <<"\t"<< e.meanMeasured.second <<"\t"
180  << e.sigmaMeasured.first <<"\t"<< e.sigmaMeasured.second <<"\t"
181  << e.meanUncertainty.first <<"\t"<< e.meanUncertainty.second << "\t"
182  << e.pull.first <<"\t"<< e.pull.second << "\t"
183  << e.samples;
184 }
185 
186 
187 
std::pair< float, float > pull
static std::map< std::string, Result > layer_results(const Book &, const Method)
Definition: LA_Results.cc:63
const double w
Definition: UKUtility.cc:23
std::pair< float, float > meanMeasured
std::pair< float, float > sigmaMeasured
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
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
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
double f[11][100]
iterator begin(string_t re=".*")
Definition: Book.h:48
base
Make Sure CMSSW is Setup ##.
static std::string method(Method m, bool fit=true)
static std::map< std::string, std::vector< Result > > ensemble_results(const Book &, const Method)
Definition: LA_Results.cc:73
void fill(double_t X, const char *name, uint_t NbinsX, double_t Xlow, double_t Xup, double_t W=1)
Definition: Book.h:60
static std::map< uint32_t, Result > module_results(const Book &, const Method)
Definition: LA_Results.cc:51
const double ensembleUp_
std::pair< float, float > meanUncertainty
boost::filter_iterator< match_name, book_t::const_iterator > const_iterator
Definition: Book.h:47
fixed size matrix
std::pair< float, float > reco
static float pull(const std::vector< EnsembleSummary > &)
Definition: LA_Results.cc:156
Definition: Book.h:16
std::pair< float, float > calMeasured
std::ostream & operator<<(std::ostream &strm, const LA_Filler_Fitter::Result &r)
Definition: LA_Results.cc:167
const double ensembleLow_
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
iterator end(string_t re=".*")
Definition: Book.h:50