CMS 3D CMS Logo

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