CMS 3D CMS Logo

EnsembleCalibrationLA.cc
Go to the documentation of this file.
6 #include <TChain.h>
7 #include <TFile.h>
8 #include <boost/lexical_cast.hpp>
9 #include <fstream>
10 
11 namespace sistrip {
12 
14  inputFiles( conf.getParameter<std::vector<std::string> >("InputFiles") ),
15  inFileLocation( conf.getParameter<std::string>("InFileLocation")),
16  Prefix( conf.getUntrackedParameter<std::string>("Prefix","")),
17  maxEvents( conf.getUntrackedParameter<unsigned>("MaxEvents",0)),
18  samples( conf.getParameter<unsigned>("Samples")),
19  nbins( conf.getParameter<unsigned>("NBins")),
20  lowBin( conf.getParameter<double>("LowBin")),
21  highBin( conf.getParameter<double>("HighBin")),
22  vMethods( conf.getParameter<std::vector<int> >("Methods"))
23 {}
24 
26 {
27  Book book("la_ensemble");
28  TChain*const chain = new TChain("la_ensemble");
29  for(auto const& file : inputFiles) chain->Add((file+inFileLocation).c_str());
30 
31  int methods = 0;
32  for(unsigned int method : vMethods) methods|=method;
33 
35  laff.fill(chain,book);
36  laff.fit(book);
37  laff.summarize_ensembles(book);
38 
41  write_samples_plots(book);
43 }
44 
46 {
48  eSetup.get<TrackerTopologyRcd>().get(tTopoHandle);
49  tTopo_ = tTopoHandle.product();
50 }
51 
53 {
54  for(auto const& ensemble : LA_Filler_Fitter::ensemble_summary(book)) {
55  std::fstream file((Prefix+ensemble.first+".dat").c_str(),std::ios::out);
56  for(auto const& summary : ensemble.second)
57  file << summary << std::endl;
58 
59  const std::pair<std::pair<float,float>,std::pair<float,float> > line = LA_Filler_Fitter::offset_slope(ensemble.second);
60  const float pull = LA_Filler_Fitter::pull(ensemble.second);
61 
62  unsigned index = 15;
64  {
65  std::cout << ensemble.first << std::endl;
66  boost::regex format(".*(T[IO]B)_layer(\\d)([as])_(.*)");
67  if(boost::regex_match(ensemble.first,format)) {
68  const bool TIB = "TIB" == boost::regex_replace(ensemble.first, format, "\\1");
69  const bool stereo = "s" == boost::regex_replace(ensemble.first, format, "\\3");
70  const unsigned layer = boost::lexical_cast<unsigned>(boost::regex_replace(ensemble.first, format, "\\2"));
71  label = boost::regex_replace(ensemble.first, format, "\\4");
72  index = LA_Filler_Fitter::layer_index(TIB,stereo,layer);
73 
74  calibrations[label].slopes[index]=line.second.first;
75  calibrations[label].offsets[index]=line.first.first;
76  calibrations[label].pulls[index]=pull;
77  }
78  }
79 
80  file << std::endl << std::endl
81  << "# Best Fit Line: "
82  << line.first.first <<"("<< line.first.second<<") + x* "
83  << line.second.first<<"("<< line.second.second<<")"
84  << std::endl
85  << "# Pull (average sigma of (x_measure-x_truth)/e_measure): " << pull
86  << std::endl
87  << "LA_Calibration( METHOD_XXXXX , xxx, " << line.second.first << ", " << line.first.first << ", " << pull << ")," << std::endl;
88  file.close();
89  }
90 }
91 
93 write_ensembles_plots(const Book& book) const {
94  TFile file((Prefix+"sampleFits.root").c_str(),"RECREATE");
95  for(Book::const_iterator hist = book.begin(".*(profile|ratio|reconstruction|symm|symmchi2|_w\\d)"); hist!=book.end(); ++hist)
96  hist->second->Write();
97  file.Close();
98 }
99 
101 write_samples_plots(const Book& book) const {
102  TFile file((Prefix+"ensembleFits.root").c_str(),"RECREATE");
103  for(Book::const_iterator hist = book.begin(".*(measure|merr|ensembleReco|pull)"); hist!=book.end(); ++hist)
104  hist->second->Write();
105  file.Close();
106 }
107 
110  std::fstream file((Prefix+"calibrations.dat").c_str(),std::ios::out);
111  for(auto const& cal : calibrations) {
112  file << cal.first << std::endl
113  << "\t slopes("; for(float i : cal.second.slopes) file << i<< ","; file << ")" << std::endl
114  << "\t offsets("; for(float i : cal.second.offsets) file << i<< ","; file << ")" << std::endl
115  << "\t pulls("; for(float i : cal.second.pulls) file << i<< ","; file << ")" << std::endl;
116  }
117  file.close();
118 }
119 
120 }
121 
std::map< std::string, MethodCalibrations > calibrations
static unsigned layer_index(bool TIB, bool stereo, unsigned layer)
static void fit(Book &book)
Definition: chain.py:1
void write_ensembles_plots(const Book &) const
static std::pair< std::pair< float, float >, std::pair< float, float > > offset_slope(const std::vector< EnsembleSummary > &)
Definition: LA_Results.cc:131
static std::map< std::string, std::vector< EnsembleSummary > > ensemble_summary(const Book &)
Definition: LA_Results.cc:106
void fill(TTree *, Book &) const
Definition: LA_Filler.cc:9
sistrip classes
void summarize_ensembles(Book &) const
Definition: LA_Results.cc:82
char const * label
const std::vector< std::string > inputFiles
iterator begin(string_t re=".*")
Definition: Book.h:48
format
Some error handling for the usage.
void write_samples_plots(const Book &) const
EnsembleCalibrationLA(const edm::ParameterSet &)
inputFiles
Definition: merge.py:6
boost::filter_iterator< match_name, book_t::const_iterator > const_iterator
Definition: Book.h:47
static float pull(const std::vector< EnsembleSummary > &)
Definition: LA_Results.cc:155
T get() const
Definition: EventSetup.h:71
Definition: Book.h:16
T const * product() const
Definition: ESHandle.h:86
void endRun(const edm::Run &, const edm::EventSetup &) override
Definition: Run.h:45
iterator end(string_t re=".*")
Definition: Book.h:50