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