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