CMS 3D CMS Logo

MeasureLA.cc
Go to the documentation of this file.
5 
6 #include <boost/lexical_cast.hpp>
7 #include <TChain.h>
8 #include <TFile.h>
9 
10 
11 namespace sistrip {
12 
13 void MeasureLA::
15  for(auto const& p : vpset) {
16  methods |= p.getParameter<int32_t>("Method");
17  byModule = byModule || p.getParameter<int32_t>("Granularity");
18  byLayer = byLayer || !p.getParameter<int32_t>("Granularity");
19  }
20 }
21 
23  inputFiles( conf.getParameter<std::vector<std::string> >("InputFiles") ),
24  inFileLocation( conf.getParameter<std::string>("InFileLocation")),
25  fp_(conf.getParameter<edm::FileInPath>("SiStripDetInfo") ),
26  reports( conf.getParameter<edm::VParameterSet>("Reports")),
27  measurementPreferences( conf.getParameter<edm::VParameterSet>("MeasurementPreferences")),
28  calibrations(conf.getParameter<edm::VParameterSet>("Calibrations")),
30  localybin(conf.getUntrackedParameter<double>("LocalYBin",0.0)),
31  stripsperbin(conf.getUntrackedParameter<unsigned>("StripsPerBin",0)),
32  maxEvents( conf.getUntrackedParameter<unsigned>("MaxEvents",0)),
33  tTopo_(StandaloneTrackerTopology::fromTrackerParametersXMLFile(conf.getParameter<edm::FileInPath>("TrackerParameters").fullPath()))
34 {
38 
39  TChain*const chain = new TChain("la_data");
40  for(auto const& file : inputFiles) chain->Add((file+inFileLocation).c_str());
41 
43  laff.fill(chain, book);
44  laff.fit(book);
47 
49 }
50 
51 
52 std::unique_ptr<SiStripLorentzAngle> MeasureLA::
54  auto lorentzAngle = std::make_unique<SiStripLorentzAngle>();
55  /*
56  std::map<uint32_t,LA_Filler_Fitter::Result>
57  module_results = LA_Filler_Fitter::module_results(book, LA_Filler_Fitter::SQRTVAR);
58 
59  BOOST_FOREACH(const uint32_t& detid, SiStripDetInfoFileReader(fp_.fullPath()).getAllDetIds()) {
60  float la = module_results[detid].measure / module_results[detid].field ;
61  lorentzAngle->putLorentzAngle( detid, la );
62  }
63  */
64  return lorentzAngle;
65 }
66 
67 
68 
69 
70 void MeasureLA::
74  for(auto& result : LA_Filler_Fitter::module_results(book, method)) {
75 
76  calibrate( calibration_key(result.first,method), result.second);
78  label = boost::regex_replace(label,boost::regex("layer"),"");
79 
80  const double mu_H = -result.second.calMeasured.first / result.second.field;
81  const double sigma_mu_H = result.second.calMeasured.second / result.second.field;
82  const double weight = pow(1./sigma_mu_H, 2);
83 
84  book.fill(mu_H, label, 150,-0.05,0.1, weight);
85  }
86  for(Book::iterator it = book.begin(".*"+granularity(MODULESUMMARY)+".*"); it!=book.end(); ++it) {
87  if(it->second->GetEntries()) it->second->Fit("gaus","LLQ");
88  }
89  }
90 }
91 
92 void MeasureLA::
93 process_reports() const {
94  for(auto const& p : reports) {
95  const GRANULARITY gran = (GRANULARITY) p.getParameter<int32_t>("Granularity");
96  const std::string name = p.getParameter<std::string>("ReportName");
97  const LA_Filler_Fitter::Method method = (LA_Filler_Fitter::Method) p.getParameter<int32_t>("Method");
98 
99  write_report_plots( name, method, gran);
100  switch(gran) {
101  case LAYER: write_report_text( name, method, LA_Filler_Fitter::layer_results(book, method) ); break;
102  case MODULE: write_report_text( name, method, LA_Filler_Fitter::module_results(book, method) ); break;
103  case MODULESUMMARY: write_report_text_ms( name, method); break;
104  }
105  }
106 
107  { TFile widthsFile("widths.root","RECREATE");
108  for(Book::const_iterator it = book.begin(".*_width"); it!=book.end(); it++) if(it->second) it->second->Write();
109  widthsFile.Close();}
110 }
111 
112 void MeasureLA::
114  TFile file((name+".root").c_str(),"RECREATE");
115  const std::string key = ".*" + granularity(gran) + ".*("+LA_Filler_Fitter::method(method)+"|"+LA_Filler_Fitter::method(method,false)+".*)";
116  for(Book::const_iterator hist = book.begin(key); hist!=book.end(); ++hist)
117  if(hist->second) hist->second->Write();
118  file.Close();
119 }
120 
121 template <class T>
122 void MeasureLA::
123 write_report_text(std::string name, const LA_Filler_Fitter::Method& _method, const std::map<T,LA_Filler_Fitter::Result>& _results) const {
125  std::map<T,LA_Filler_Fitter::Result>results = _results;
126  std::fstream file((name+".dat").c_str(),std::ios::out);
127  for(auto& result : results) {
128  calibrate( calibration_key(result.first,method), result.second);
129  file << result.first << "\t" << result.second << std::endl;
130  }
131  file.close();
132 }
133 
134 void MeasureLA::
136  std::fstream file((name+".dat").c_str(),std::ios::out);
138  for(Book::const_iterator it = book.begin(key); it!=book.end(); ++it) {
139  const TF1*const f = it->second->GetFunction("gaus");
140  if(f) {
141  file << it->first << "\t"
142  << f->GetParameter(1) << "\t"
143  << f->GetParError(1) << "\t"
144  << f->GetParameter(2) << "\t"
145  << f->GetParError(2) << std::endl;
146  }
147  }
148  file.close();
149 }
150 
151 void MeasureLA::
153  for(auto const& p : calibrations) {
154  LA_Filler_Fitter::Method method = (LA_Filler_Fitter::Method) p.getParameter<int32_t>("Method");
155  std::vector<double> slopes(p.getParameter<std::vector<double> >("Slopes")); assert(slopes.size()==14);
156  std::vector<double> offsets(p.getParameter<std::vector<double> >("Offsets")); assert(offsets.size()==14);
157  std::vector<double> pulls(p.getParameter<std::vector<double> >("Pulls")); assert(pulls.size()==14);
158 
159  for(unsigned i=0; i<14; i++) {
160  const std::pair<unsigned,LA_Filler_Fitter::Method> key( i, method);
161  offset[key] = offsets[i];
162  slope[key] = slopes[i];
163  error_scaling[key] = pulls[i];
164  }
165  }
166 }
167 
168 inline
169 void MeasureLA::
170 calibrate(const std::pair<unsigned,LA_Filler_Fitter::Method> key, LA_Filler_Fitter::Result& result) const {
171  result.calMeasured = std::make_pair<float,float>( ( result.measured.first - offset.find(key)->second ) / slope.find(key)->second ,
172  result.measured.second * error_scaling.find(key)->second / slope.find(key)->second );
173 }
174 
175 std::pair<uint32_t,LA_Filler_Fitter::Method> MeasureLA::
177  boost::regex format(".*(T[IO]B)_layer(\\d)([as]).*");
178  const bool isTIB = "TIB" == boost::regex_replace(layer, format, "\\1");
179  const bool stereo = "s" == boost::regex_replace(layer, format, "\\3");
180  const unsigned layerNum = boost::lexical_cast<unsigned>(boost::regex_replace(layer, format, "\\2"));
181  return std::make_pair(LA_Filler_Fitter::layer_index(isTIB,stereo,layerNum),method);
182 }
183 
184 std::pair<uint32_t,LA_Filler_Fitter::Method> MeasureLA::
185 calibration_key(const uint32_t detid, const LA_Filler_Fitter::Method method) const {
186  const bool isTIB = SiStripDetId(detid).subDetector() == SiStripDetId::TIB;
187  const bool stereo = isTIB ? tTopo_.tibStereo(detid) : tTopo_.tobStereo(detid);
188  const unsigned layer = isTIB ? tTopo_.tibLayer(detid) : tTopo_.tobStereo(detid);
189 
190  return std::make_pair(LA_Filler_Fitter::layer_index(isTIB,stereo,layer),method);
191 }
192 
193 }
auto setWhatProduced(T *iThis, const es::Label &iLabel={})
Definition: ESProducer.h:124
const edm::FileInPath fp_
Definition: MeasureLA.h:49
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:62
void calibrate(const std::pair< unsigned, LA_Filler_Fitter::Method >, LA_Filler_Fitter::Result &) const
Definition: MeasureLA.cc:170
void summarize_module_muH_byLayer(const LA_Filler_Fitter &)
Definition: MeasureLA.cc:71
static unsigned layer_index(bool TIB, bool stereo, unsigned layer)
static void fit(Book &book)
unsigned int tibLayer(const DetId &id) const
Definition: chain.py:1
const unsigned stripsperbin
Definition: MeasureLA.h:55
const edm::VParameterSet measurementPreferences
Definition: MeasureLA.h:50
std::vector< ParameterSet > VParameterSet
Definition: ParameterSet.h:33
const std::vector< std::string > inputFiles
Definition: MeasureLA.h:47
const unsigned maxEvents
Definition: MeasureLA.h:55
uint32_t tobStereo(const DetId &id) const
const float localybin
Definition: MeasureLA.h:54
void process_reports() const
Definition: MeasureLA.cc:93
Definition: weight.py:1
std::map< std::pair< uint32_t, LA_Filler_Fitter::Method >, float > slope
Definition: MeasureLA.h:51
static std::string granularity(int32_t g)
Definition: MeasureLA.h:27
void fill(TTree *, Book &) const
Definition: LA_Filler.cc:9
std::pair< float, float > measured
sistrip classes
void store_methods_and_granularity(const edm::VParameterSet &)
Definition: MeasureLA.cc:14
char const * label
const edm::VParameterSet calibrations
Definition: MeasureLA.h:50
MeasureLA(const edm::ParameterSet &)
Definition: MeasureLA.cc:22
const edm::VParameterSet reports
Definition: MeasureLA.h:50
double f[11][100]
iterator begin(string_t re=".*")
Definition: Book.h:48
format
Some error handling for the usage.
const std::string inFileLocation
Definition: MeasureLA.h:48
static std::string method(Method m, bool fit=true)
std::unique_ptr< SiStripLorentzAngle > produce(const SiStripLorentzAngleRcd &)
Definition: MeasureLA.cc:53
Detector identifier class for the strip tracker.
Definition: SiStripDetId.h:17
SubDetector subDetector() const
Definition: SiStripDetId.h:102
void write_report_text_ms(const std::string, const LA_Filler_Fitter::Method) const
Definition: MeasureLA.cc:135
void store_calibrations()
Definition: MeasureLA.cc:152
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:50
std::map< std::pair< uint32_t, LA_Filler_Fitter::Method >, float > offset
Definition: MeasureLA.h:51
void write_report_plots(const std::string, const LA_Filler_Fitter::Method, const GRANULARITY) const
Definition: MeasureLA.cc:113
boost::filter_iterator< match_name, book_t::const_iterator > const_iterator
Definition: Book.h:47
HLT enums.
void write_report_text(const std::string, const LA_Filler_Fitter::Method &, const std::map< T, LA_Filler_Fitter::Result > &) const
Definition: MeasureLA.cc:123
std::map< std::pair< uint32_t, LA_Filler_Fitter::Method >, float > error_scaling
Definition: MeasureLA.h:51
TrackerTopology fromTrackerParametersXMLFile(const std::string &xmlFileName)
std::pair< float, float > calMeasured
uint32_t tibStereo(const DetId &id) const
boost::filter_iterator< match_name, book_t::iterator > iterator
Definition: Book.h:46
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
std::pair< unsigned, LA_Filler_Fitter::Method > calibration_key(const std::string layer, const LA_Filler_Fitter::Method) const
Definition: MeasureLA.cc:176
TrackerTopology tTopo_
Definition: MeasureLA.h:58
iterator end(string_t re=".*")
Definition: Book.h:50