CMS 3D CMS Logo

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