CMS 3D CMS Logo

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