CMS 3D CMS Logo

SiStripGainCosmicCalculator.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 // Package: SiStripChannelGain
3 // Class: SiStripGainCosmicCalculator
4 // Original Author: G. Bruno, D. Kcira
5 // Created: Mon May 20 10:04:31 CET 2007
18 #include "CLHEP/Random/RandFlat.h"
19 #include "CLHEP/Random/RandGauss.h"
34 //#include "DQM/SiStripCommon/interface/SiStripGenerateKey.h"
35 
36 //---------------------------------------------------------------------------------------------------------
38  : ConditionDBWriter<SiStripApvGain>(iConfig) {
39  edm::LogInfo("SiStripGainCosmicCalculator::SiStripGainCosmicCalculator");
41  edm::LogInfo("SiStripApvGainCalculator::SiStripApvGainCalculator")
42  << "ExpectedChargeDeposition=" << ExpectedChargeDeposition;
43 
44  TrackProducer = iConfig.getParameter<std::string>("TrackProducer");
45  TrackLabel = iConfig.getParameter<std::string>("TrackLabel");
46 
47  detModulesToBeExcluded.clear();
48  detModulesToBeExcluded = iConfig.getParameter<std::vector<unsigned> >("detModulesToBeExcluded");
49  MinNrEntries = iConfig.getUntrackedParameter<unsigned>("minNrEntries", 20);
50  MaxChi2OverNDF = iConfig.getUntrackedParameter<double>("maxChi2OverNDF", 5.);
51 
52  outputHistogramsInRootFile = iConfig.getParameter<bool>("OutputHistogramsInRootFile");
53  outputFileName = iConfig.getParameter<std::string>("OutputFileName");
54 
55  edm::LogInfo("SiStripApvGainCalculator")
56  << "Clusters from " << detModulesToBeExcluded.size() << " modules will be ignored in the calibration:";
57  edm::LogInfo("SiStripApvGainCalculator") << "The calibration for these DetIds will be set to a default value";
58  for (std::vector<uint32_t>::const_iterator imod = detModulesToBeExcluded.begin();
59  imod != detModulesToBeExcluded.end();
60  imod++) {
61  edm::LogInfo("SiStripApvGainCalculator") << "exclude detid = " << *imod;
62  }
63 
64  printdebug_ = iConfig.getUntrackedParameter<bool>("printDebug", false);
65  tTopo = nullptr;
66 }
67 
69  edm::LogInfo("SiStripGainCosmicCalculator::~SiStripGainCosmicCalculator");
70 }
71 
73 
75  //Retrieve tracker topology from geometry
77  iSetup.get<TrackerTopologyRcd>().get(tTopoHandle);
78  tTopo = tTopoHandle.product();
79 
80  edm::ESHandle<SiStripDetCabling> siStripDetCablingH;
81  iSetup.get<SiStripDetCablingRcd>().get(siStripDetCablingH);
82  siStripDetCabling = siStripDetCablingH.product();
83 
85  iSetup.get<TrackerDigiGeometryRecord>().get(tkGeomH);
86  tkGeom = tkGeomH.product();
87 
88  std::cout << "SiStripGainCosmicCalculator::algoBeginJob called" << std::endl;
90  HlistAPVPairs = new TObjArray();
91  HlistOtherHistos = new TObjArray();
92  //
93  HlistOtherHistos->Add(new TH1F(Form("APVPairCorrections"), Form("APVPairCorrections"), 50, -1., 4.));
94  HlistOtherHistos->Add(new TH1F(Form("APVPairCorrectionsTIB1mono"), Form("APVPairCorrectionsTIB1mono"), 50, -1., 4.));
95  HlistOtherHistos->Add(
96  new TH1F(Form("APVPairCorrectionsTIB1stereo"), Form("APVPairCorrectionsTIB1stereo"), 50, -1., 4.));
97  HlistOtherHistos->Add(new TH1F(Form("APVPairCorrectionsTIB2"), Form("APVPairCorrectionsTIB2"), 50, -1., 4.));
98  HlistOtherHistos->Add(new TH1F(Form("APVPairCorrectionsTOB1"), Form("APVPairCorrectionsTOB1"), 50, -1., 4.));
99  HlistOtherHistos->Add(new TH1F(Form("APVPairCorrectionsTOB2"), Form("APVPairCorrectionsTOB2"), 50, -1., 4.));
100  HlistOtherHistos->Add(new TH1F(Form("LocalAngle"), Form("LocalAngle"), 70, -0.1, 3.4));
101  HlistOtherHistos->Add(new TH1F(Form("LocalAngleAbsoluteCosine"), Form("LocalAngleAbsoluteCosine"), 48, -0.1, 1.1));
102  HlistOtherHistos->Add(new TH1F(Form("LocalPosition_cm"), Form("LocalPosition_cm"), 100, -5., 5.));
103  HlistOtherHistos->Add(new TH1F(Form("LocalPosition_normalized"), Form("LocalPosition_normalized"), 100, -1.1, 1.1));
104  TH1F* local_histo = new TH1F(Form("SiStripRecHitType"), Form("SiStripRecHitType"), 2, 0.5, 2.5);
105  HlistOtherHistos->Add(local_histo);
106  local_histo->GetXaxis()->SetBinLabel(1, "simple");
107  local_histo->GetXaxis()->SetBinLabel(2, "matched");
108 
109  // get cabling and find out list of active detectors
111  iSetup.get<SiStripDetCablingRcd>().get(siStripDetCabling);
112  std::vector<uint32_t> activeDets;
113  activeDets.clear();
114  SelectedDetIds.clear();
115  siStripDetCabling->addActiveDetectorsRawIds(activeDets);
116  // SelectedDetIds = activeDets; // all active detector modules
117  // use SiStripSubStructure for selecting certain regions
119  activeDets, SelectedDetIds, tTopo, 0, 0, 0, 0); // this adds rawDetIds to SelectedDetIds
121  activeDets, SelectedDetIds, tTopo, 0, 0, 0); // this adds rawDetIds to SelectedDetIds
122  // get tracker geometry and find nr. of apv pairs for each active detector
124  iSetup.get<TrackerDigiGeometryRecord>().get(tkGeom);
125  for (TrackerGeometry::DetContainer::const_iterator it = tkGeom->dets().begin(); it != tkGeom->dets().end();
126  it++) { // loop over detector modules
127  if (dynamic_cast<const StripGeomDetUnit*>((*it)) != nullptr) {
128  uint32_t detid = ((*it)->geographicalId()).rawId();
129  // get thickness for all detector modules, not just for active, this is strange
130  double module_thickness =
131  (*it)
132  ->surface()
133  .bounds()
134  .thickness(); // get thickness of detector from GeomDet (DetContainer == vector<GeomDet*>)
135  thickness_map.insert(std::make_pair(detid, module_thickness));
136  //
137  bool is_active_detector = false;
138  for (std::vector<uint32_t>::iterator iactive = SelectedDetIds.begin(); iactive != SelectedDetIds.end();
139  iactive++) {
140  if (*iactive == detid) {
141  is_active_detector = true;
142  break; // leave for loop if found matching detid
143  }
144  }
145  //
146  bool exclude_this_detid = false;
147  for (std::vector<uint32_t>::const_iterator imod = detModulesToBeExcluded.begin();
148  imod != detModulesToBeExcluded.end();
149  imod++) {
150  if (*imod == detid)
151  exclude_this_detid = true; // found in exclusion list
152  break;
153  }
154  //
155  if (is_active_detector &&
156  (!exclude_this_detid)) { // check whether is active detector and that should not be excluded
157  const StripTopology& p = dynamic_cast<const StripGeomDetUnit*>((*it))->specificTopology();
158  unsigned short NAPVPairs = p.nstrips() / 256;
159  if (NAPVPairs < 2 || NAPVPairs > 3) {
160  edm::LogError("SiStripGainCosmicCalculator")
161  << "Problem with Number of strips in detector: " << p.nstrips() << " Exiting program";
162  exit(1);
163  }
164  for (int iapp = 0; iapp < NAPVPairs; iapp++) {
165  TString hid = Form("ChargeAPVPair_%i_%i", detid, iapp);
166  HlistAPVPairs->Add(
167  new TH1F(hid, hid, 45, 0., 1350.)); // multiply by 3 to take into account division by width
168  }
169  }
170  }
171  }
172 }
173 
174 //---------------------------------------------------------------------------------------------------------
176  using namespace edm;
178 
179  //TO BE RESTORED
180  // anglefinder_->init(event,iSetup);
181 
182  // get seeds
183  // edm::Handle<TrajectorySeedCollection> seedcoll;
184  // event.getByType(seedcoll);
185  // get tracks
187  iEvent.getByLabel(TrackProducer, TrackLabel, trackCollection);
188  const reco::TrackCollection* tracks = trackCollection.product();
189 
190  // // get magnetic field
191  // edm::ESHandle<MagneticField> esmagfield;
192  // es.get<IdealMagneticFieldRecord>().get(esmagfield);
193  // magfield=&(*esmagfield);
194  // loop over tracks
195  for (reco::TrackCollection::const_iterator itr = tracks->begin(); itr != tracks->end();
196  itr++) { // looping over tracks
197 
198  //TO BE RESTORED
199  // std::vector<std::pair<const TrackingRecHit *,float> >hitangle =anglefinder_->findtrackangle((*(*seedcoll).begin()),*itr);
200  std::vector<std::pair<const TrackingRecHit*, float> >
201  hitangle; // =anglefinder_->findtrackangle((*(*seedcoll).begin()),*itr);
202 
203  for (std::vector<std::pair<const TrackingRecHit*, float> >::const_iterator hitangle_iter = hitangle.begin();
204  hitangle_iter != hitangle.end();
205  hitangle_iter++) {
206  const TrackingRecHit* trechit = hitangle_iter->first;
207  float local_angle = hitangle_iter->second;
208  LocalPoint local_position = trechit->localPosition();
209  const SiStripRecHit2D* sistripsimplehit = dynamic_cast<const SiStripRecHit2D*>(trechit);
210  const SiStripMatchedRecHit2D* sistripmatchedhit = dynamic_cast<const SiStripMatchedRecHit2D*>(trechit);
211  // std::cout<<" hit/matched "<<std::ios::hex<<sistripsimplehit<<" "<<sistripmatchedhit<<std::endl;
212  ((TH1F*)HlistOtherHistos->FindObject("LocalAngle"))->Fill(local_angle);
213  ((TH1F*)HlistOtherHistos->FindObject("LocalAngleAbsoluteCosine"))->Fill(fabs(cos(local_angle)));
214  if (sistripsimplehit) {
215  ((TH1F*)HlistOtherHistos->FindObject("SiStripRecHitType"))->Fill(1.);
216  const SiStripRecHit2D::ClusterRef& cluster = sistripsimplehit->cluster();
217  const auto& ampls = cluster->amplitudes();
218  uint32_t thedetid = 0; // is zero since long time cluster->geographicalId();
219  double module_width = moduleWidth(thedetid);
220  ((TH1F*)HlistOtherHistos->FindObject("LocalPosition_cm"))->Fill(local_position.x());
221  ((TH1F*)HlistOtherHistos->FindObject("LocalPosition_normalized"))->Fill(local_position.x() / module_width);
222  double module_thickness = moduleThickness(thedetid);
223  int ifirststrip = cluster->firstStrip();
224  int theapvpairid = int(float(ifirststrip) / 256.);
225  TH1F* histopointer = (TH1F*)HlistAPVPairs->FindObject(Form("ChargeAPVPair_%i_%i", thedetid, theapvpairid));
226  if (histopointer) {
227  short cCharge = 0;
228  for (unsigned int iampl = 0; iampl < ampls.size(); iampl++) {
229  cCharge += ampls[iampl];
230  }
231  double cluster_charge_over_path = ((double)cCharge) * fabs(cos(local_angle)) / (10. * module_thickness);
232  histopointer->Fill(cluster_charge_over_path);
233  }
234  } else {
235  if (sistripmatchedhit)
236  ((TH1F*)HlistOtherHistos->FindObject("SiStripRecHitType"))->Fill(2.);
237  }
238  }
239  }
240 }
241 
242 //---------------------------------------------------------------------------------------------------------
244  TH1F* inputHisto) { // automated fitting with finding of the appropriate nr. of ADCs
245  // set some default dummy value and return if no entries
246  double adcs = -0.5;
247  double error = 0.;
248  double nr_of_entries = inputHisto->GetEntries();
249  if (nr_of_entries < MinNrEntries) {
250  return std::make_pair(adcs, error);
251  }
252  //
253  // // fit with initial setting of parameter values
254  // double rms_of_histogram = inputHisto->GetRMS();
255  // TF1 *landaufit = new TF1("landaufit","landau",0.,450.);
256  // landaufit->SetParameters(nr_of_entries,mean_of_histogram,rms_of_histogram);
257  // inputHisto->Fit("landaufit","0Q+");
258  // delete landaufit;
259  //
260  // perform fit with standard landau
261  // make our own copy to avoid problems with threads
262  std::unique_ptr<TF1> fitfunction(new TF1("landaufit", "landau"));
263  inputHisto->Fit(fitfunction.get(), "0Q");
264  adcs = fitfunction->GetParameter("MPV");
265  error = fitfunction->GetParError(1); // MPV is parameter 1 (0=constant, 1=MPV, 2=Sigma)
266  double chi2 = fitfunction->GetChisquare();
267  double ndf = fitfunction->GetNDF();
268  double chi2overndf = chi2 / ndf;
269  // in case things went wrong, try to refit in smaller range
270  if (adcs < 2. || (error / adcs) > 1.8) {
271  inputHisto->Fit(fitfunction.get(), "0Q", nullptr, 0., 400.);
272  std::cout << "refitting landau for histogram " << inputHisto->GetTitle() << std::endl;
273  std::cout << "initial error/adcs =" << error << " / " << adcs << std::endl;
274  std::cout << "new error/adcs =" << fitfunction->GetParError(1) << " / " << fitfunction->GetParameter("MPV")
275  << std::endl;
276  adcs = fitfunction->GetParameter("MPV");
277  error = fitfunction->GetParError(1); // MPV is parameter 1 (0=constant, 1=MPV, 2=Sigma)
278  chi2 = fitfunction->GetChisquare();
279  ndf = fitfunction->GetNDF();
280  chi2overndf = chi2 / ndf;
281  }
282  // if still wrong, give up
283  if (adcs < 2. || chi2overndf > MaxChi2OverNDF) {
284  adcs = -0.5;
285  error = 0.;
286  }
287  return std::make_pair(adcs, error);
288 }
289 
290 //---------------------------------------------------------------------------------------------------------
291 double SiStripGainCosmicCalculator::moduleWidth(const uint32_t detid) // get width of the module detid
292 { //dk: copied from A. Giammanco and hacked, module_width values : 10.49 12.03 6.144 7.14 9.3696
293  double module_width = 0.;
294  const GeomDetUnit* it = tkGeom->idToDetUnit(DetId(detid));
295  if (dynamic_cast<const StripGeomDetUnit*>(it) == nullptr && dynamic_cast<const PixelGeomDetUnit*>(it) == nullptr) {
296  std::cout << "this detID doesn't seem to belong to the Tracker" << std::endl;
297  } else {
298  module_width = it->surface().bounds().width();
299  }
300  return module_width;
301 }
302 
303 //---------------------------------------------------------------------------------------------------------
304 double SiStripGainCosmicCalculator::moduleThickness(const uint32_t detid) // get thickness of the module detid
305 { //dk: copied from A. Giammanco and hacked
306  double module_thickness = 0.;
307  const GeomDetUnit* it = tkGeom->idToDetUnit(DetId(detid));
308  if (dynamic_cast<const StripGeomDetUnit*>(it) == nullptr && dynamic_cast<const PixelGeomDetUnit*>(it) == nullptr) {
309  std::cout << "this detID doesn't seem to belong to the Tracker" << std::endl;
310  } else {
311  module_thickness = it->surface().bounds().thickness();
312  }
313  return module_thickness;
314 }
315 
316 //---------------------------------------------------------------------------------------------------------
317 std::unique_ptr<SiStripApvGain> SiStripGainCosmicCalculator::getNewObject() {
318  std::cout << "SiStripGainCosmicCalculator::getNewObject called" << std::endl;
319 
320  std::cout << "total_nr_of_events=" << total_nr_of_events << std::endl;
321  // book some more histograms
322  TH1F* ChargeOfEachAPVPair = new TH1F("ChargeOfEachAPVPair", "ChargeOfEachAPVPair", 1, 0, 1);
323  ChargeOfEachAPVPair->SetCanExtend(TH1::kAllAxes);
324  TH1F* EntriesApvPairs = new TH1F("EntriesApvPairs", "EntriesApvPairs", 1, 0, 1);
325  EntriesApvPairs->SetCanExtend(TH1::kAllAxes);
326  TH1F* NrOfEntries =
327  new TH1F("NrOfEntries", "NrOfEntries", 351, -0.5, 350.5); // NrOfEntries->SetCanExtend(TH1::kAllAxes);
328  TH1F* ModuleThickness = new TH1F("ModuleThickness", "ModuleThickness", 2, 0.5, 2.5);
329  HlistOtherHistos->Add(ModuleThickness);
330  ModuleThickness->GetXaxis()->SetBinLabel(1, "320mu");
331  ModuleThickness->GetXaxis()->SetBinLabel(2, "500mu");
332  ModuleThickness->SetYTitle("Nr APVPairs");
333  TH1F* ModuleWidth = new TH1F("ModuleWidth", "ModuleWidth", 5, 0.5, 5.5);
334  HlistOtherHistos->Add(ModuleWidth);
335  ModuleWidth->GetXaxis()->SetBinLabel(1, "6.144cm");
336  ModuleWidth->GetXaxis()->SetBinLabel(2, "7.14cm");
337  ModuleWidth->GetXaxis()->SetBinLabel(3, "9.3696cm");
338  ModuleWidth->GetXaxis()->SetBinLabel(4, "10.49cm");
339  ModuleWidth->GetXaxis()->SetBinLabel(5, "12.03cm");
340  ModuleWidth->SetYTitle("Nr APVPairs");
341  // loop over single histograms and extract peak value of charge
342  HlistAPVPairs->Sort(); // sort alfabetically
343  TIter hiterator(HlistAPVPairs);
344  double MeanCharge = 0.;
345  double NrOfApvPairs = 0.;
346  TH1F* MyHisto = (TH1F*)hiterator();
347  while (MyHisto) {
348  TString histo_title = MyHisto->GetTitle();
349  if (histo_title.Contains("ChargeAPVPair_")) {
350  std::pair<double, double> two_values = getPeakOfLandau(MyHisto);
351  double local_nrofadcs = two_values.first;
352  double local_sigma = two_values.second;
353  ChargeOfEachAPVPair->Fill(histo_title, local_nrofadcs);
354  int ichbin = ChargeOfEachAPVPair->GetXaxis()->FindBin(histo_title.Data());
355  ChargeOfEachAPVPair->SetBinError(ichbin, local_sigma);
356  EntriesApvPairs->Fill(histo_title, MyHisto->GetEntries());
357  NrOfEntries->Fill(MyHisto->GetEntries());
358  if (local_nrofadcs > 0) { // if nr of adcs is negative, the fitting routine could not extract meaningfull numbers
359  MeanCharge += local_nrofadcs;
360  NrOfApvPairs += 1.; // count nr of apv pairs since do not know whether nr of bins of histogram is the same
361  }
362  }
363  MyHisto = (TH1F*)hiterator();
364  }
365  ChargeOfEachAPVPair->LabelsDeflate("X");
366  EntriesApvPairs->LabelsDeflate("X"); // trim nr. of bins to match active labels
367  HlistOtherHistos->Add(ChargeOfEachAPVPair);
368  HlistOtherHistos->Add(EntriesApvPairs);
369  HlistOtherHistos->Add(NrOfEntries);
370  MeanCharge = MeanCharge / NrOfApvPairs;
371  // calculate correction
372  TH1F* CorrectionOfEachAPVPair = (TH1F*)ChargeOfEachAPVPair->Clone("CorrectionOfEachAPVPair");
373  TH1F* ChargeOfEachAPVPairControlView =
374  new TH1F("ChargeOfEachAPVPairControlView", "ChargeOfEachAPVPairControlView", 1, 0, 1);
375  ChargeOfEachAPVPairControlView->SetCanExtend(TH1::kAllAxes);
376  TH1F* CorrectionOfEachAPVPairControlView =
377  new TH1F("CorrectionOfEachAPVPairControlView", "CorrectionOfEachAPVPairControlView", 1, 0, 1);
378  CorrectionOfEachAPVPairControlView->SetCanExtend(TH1::kAllAxes);
379  std::ofstream APVPairTextOutput("apvpair_corrections.txt");
380  APVPairTextOutput << "# MeanCharge = " << MeanCharge << std::endl;
381  APVPairTextOutput << "# Nr. of APVPairs = " << NrOfApvPairs << std::endl;
382  for (int ibin = 1; ibin <= ChargeOfEachAPVPair->GetNbinsX(); ibin++) {
383  TString local_bin_label = ChargeOfEachAPVPair->GetXaxis()->GetBinLabel(ibin);
384  double local_charge_over_path = ChargeOfEachAPVPair->GetBinContent(ibin);
385  if (local_bin_label.Contains("ChargeAPVPair_") &&
386  local_charge_over_path > 0.0000001) { // calculate correction only for meaningful numbers
387  uint32_t extracted_detid;
388  std::istringstream read_label((local_bin_label(14, 9)).Data());
389  read_label >> extracted_detid;
390  unsigned short extracted_apvpairid;
391  std::istringstream read_apvpair((local_bin_label(24, 1)).Data());
392  read_apvpair >> extracted_apvpairid;
393  double local_error_of_charge = ChargeOfEachAPVPair->GetBinError(ibin);
394  double local_correction = -0.5;
395  double local_error_correction = 0.;
396  local_correction =
397  MeanCharge / local_charge_over_path; // later use ExpectedChargeDeposition instead of MeanCharge
398  local_error_correction = local_correction * local_error_of_charge / local_charge_over_path;
399  if (local_error_correction > 1.8) { // understand why error too large sometimes
400  std::cout << "too large error " << local_error_correction << " for histogram " << local_bin_label << std::endl;
401  }
402  double nr_of_entries = EntriesApvPairs->GetBinContent(ibin);
403  APVPairTextOutput << local_bin_label << " " << local_correction << " " << local_charge_over_path << " "
404  << nr_of_entries << std::endl;
405  CorrectionOfEachAPVPair->SetBinContent(ibin, local_correction);
406  CorrectionOfEachAPVPair->SetBinError(ibin, local_error_correction);
407  ((TH1F*)HlistOtherHistos->FindObject("APVPairCorrections"))->Fill(local_correction);
408  DetId thedetId = DetId(extracted_detid);
409  unsigned int generalized_layer = 0;
410  // calculate generalized_layer: 31,32 = TIB1, 33 = TIB2, 33 = TIB3, 51 = TOB1, 52 = TOB2, 60 = TEC
411  if (thedetId.subdetId() == StripSubdetector::TIB) {
412  generalized_layer =
413  10 * thedetId.subdetId() + tTopo->tibLayer(thedetId.rawId()) + tTopo->tibStereo(thedetId.rawId());
414  if (tTopo->tibLayer(thedetId.rawId()) == 2) {
415  generalized_layer++;
416  if (tTopo->tibGlued(thedetId.rawId()))
417  edm::LogError("ClusterMTCCFilter") << "WRONGGGG" << std::endl;
418  }
419  } else {
420  generalized_layer = 10 * thedetId.subdetId();
421  if (thedetId.subdetId() == StripSubdetector::TOB) {
422  generalized_layer += tTopo->tobLayer(thedetId.rawId());
423  }
424  }
425  if (generalized_layer == 31) {
426  ((TH1F*)HlistOtherHistos->FindObject("APVPairCorrectionsTIB1mono"))->Fill(local_correction);
427  }
428  if (generalized_layer == 32) {
429  ((TH1F*)HlistOtherHistos->FindObject("APVPairCorrectionsTIB1stereo"))->Fill(local_correction);
430  }
431  if (generalized_layer == 33) {
432  ((TH1F*)HlistOtherHistos->FindObject("APVPairCorrectionsTIB2"))->Fill(local_correction);
433  }
434  if (generalized_layer == 51) {
435  ((TH1F*)HlistOtherHistos->FindObject("APVPairCorrectionsTOB1"))->Fill(local_correction);
436  }
437  if (generalized_layer == 52) {
438  ((TH1F*)HlistOtherHistos->FindObject("APVPairCorrectionsTOB2"))->Fill(local_correction);
439  }
440  // control view
441  const FedChannelConnection& fedchannelconnection =
442  siStripDetCabling->getConnection(extracted_detid, extracted_apvpairid);
443  std::ostringstream local_key;
444  // in S. Mersi's analysis the APVPair id seems to be used instead of the lldChannel, hence use the same here
445  local_key << "fecCrate" << fedchannelconnection.fecCrate() << "_fecSlot" << fedchannelconnection.fecSlot()
446  << "_fecRing" << fedchannelconnection.fecRing() << "_ccuAddr" << fedchannelconnection.ccuAddr()
447  << "_ccuChan" << fedchannelconnection.ccuChan() << "_apvPair" << extracted_apvpairid;
448  TString control_key = local_key.str();
449  ChargeOfEachAPVPairControlView->Fill(control_key, local_charge_over_path);
450  int ibin1 = ChargeOfEachAPVPairControlView->GetXaxis()->FindBin(control_key);
451  ChargeOfEachAPVPairControlView->SetBinError(ibin1, local_error_of_charge);
452  CorrectionOfEachAPVPairControlView->Fill(control_key, local_correction);
453  int ibin2 = CorrectionOfEachAPVPairControlView->GetXaxis()->FindBin(control_key);
454  CorrectionOfEachAPVPairControlView->SetBinError(ibin2, local_error_correction);
455  // thickness of each module
456  double module_thickness = moduleThickness(extracted_detid);
457  if (fabs(module_thickness - 0.032) < 0.001)
458  ModuleThickness->Fill(1);
459  if (fabs(module_thickness - 0.05) < 0.001)
460  ModuleThickness->Fill(2);
461  // width of each module
462  double module_width = moduleWidth(extracted_detid);
463  if (fabs(module_width - 6.144) < 0.01)
464  ModuleWidth->Fill(1);
465  if (fabs(module_width - 7.14) < 0.01)
466  ModuleWidth->Fill(2);
467  if (fabs(module_width - 9.3696) < 0.01)
468  ModuleWidth->Fill(3);
469  if (fabs(module_width - 10.49) < 0.01)
470  ModuleWidth->Fill(4);
471  if (fabs(module_width - 12.03) < 0.01)
472  ModuleWidth->Fill(5);
473  }
474  }
475  HlistOtherHistos->Add(CorrectionOfEachAPVPair);
476  ChargeOfEachAPVPairControlView->LabelsDeflate("X");
477  CorrectionOfEachAPVPairControlView->LabelsDeflate("X");
478  HlistOtherHistos->Add(ChargeOfEachAPVPairControlView);
479  HlistOtherHistos->Add(CorrectionOfEachAPVPairControlView);
480  // output histograms to file
481 
483  TFile* outputfile = new TFile(outputFileName, "RECREATE");
484  HlistAPVPairs->Write();
485  HlistOtherHistos->Write();
486  outputfile->Close();
487  }
488 
489  auto obj = std::make_unique<SiStripApvGain>();
490 
491  // for(std::map<uint32_t,OptoScanAnalysis*>::const_iterator it = analyses.begin(); it != analyses.end(); it++){
492  // //Generate Gain for det detid
493  // std::vector<float> theSiStripVector;
494  // for(unsigned short j=0; j<it->second; j++){
495  // float gain;
496 
497  // // if(sigmaGain_/meanGain_ < 0.00001) gain = meanGain_;
498  // // else{
499  // gain = CLHEP::RandGauss::shoot(meanGain_, sigmaGain_);
500  // if(gain<=minimumPosValue_) gain=minimumPosValue_;
501  // // }
502 
503  // if (printdebug_)
504  // edm::LogInfo("SiStripGainCalculator") << "detid " << it->first << " \t"
505  // << " apv " << j << " \t"
506  // << gain << " \t"
507  // << std::endl;
508  // theSiStripVector.push_back(gain);
509  // }
510  // SiStripApvGain::Range range(theSiStripVector.begin(),theSiStripVector.end());
511  // if ( ! obj->put(it->first,range) )
512  // edm::LogError("SiStripGainCalculator")<<"[SiStripGainCalculator::beginJob] detid already exists"<<std::endl;
513  // }
514 
515  return obj;
516 }
const uint16_t & fecSlot() const
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
std::vector< uint32_t > SelectedDetIds
const uint16_t & fecCrate() const
const FedChannelConnection & getConnection(uint32_t det_id, unsigned short apv_pair) const
void getTIBDetectors(const std::vector< uint32_t > &inputDetRawIds, std::vector< uint32_t > &tibDetRawIds, const TrackerTopology *trackerTopology, uint32_t layer=0, uint32_t bkw_frw=0, uint32_t int_ext=0, uint32_t string=0)
void addActiveDetectorsRawIds(std::vector< uint32_t > &) const
unsigned int tibLayer(const DetId &id) const
double moduleThickness(const uint32_t detid)
void algoBeginJob(const edm::EventSetup &) override
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
const Bounds & bounds() const
Definition: Surface.h:89
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
std::pair< double, double > getPeakOfLandau(TH1F *inputHisto)
virtual float width() const =0
const DetContainer & dets() const override
Returm a vector of all GeomDet (including all GeomDetUnits)
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
int iEvent
Definition: GenABIO.cc:224
const uint16_t & fecRing() const
Class containning control, module, detector and connection information, at the level of a FED channel...
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
const uint16_t & ccuChan() const
SiStripGainCosmicCalculator(const edm::ParameterSet &)
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
std::map< uint32_t, double > thickness_map
ClusterRef cluster() const
static constexpr auto TOB
SiStripDetCabling const * siStripDetCabling
virtual LocalPoint localPosition() const =0
const uint16_t & ccuAddr() const
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:488
std::unique_ptr< SiStripApvGain > getNewObject() override
Definition: DetId.h:17
static constexpr auto TIB
uint32_t tibGlued(const DetId &id) const
T const * product() const
Definition: Handle.h:69
virtual float thickness() const =0
virtual int nstrips() const =0
std::vector< uint32_t > detModulesToBeExcluded
HLT enums.
T get() const
Definition: EventSetup.h:73
void getTOBDetectors(const std::vector< uint32_t > &inputDetRawIds, std::vector< uint32_t > &tobDetRawIds, const TrackerTopology *trackerTopology, uint32_t layer=0, uint32_t bkw_frw=0, uint32_t rod=0)
uint32_t tibStereo(const DetId &id) const
double moduleWidth(const uint32_t detid)
T x() const
Definition: PV3DBase.h:59
T const * product() const
Definition: ESHandle.h:86
void algoAnalyze(const edm::Event &, const edm::EventSetup &) override
unsigned int tobLayer(const DetId &id) const
def exit(msg="")