CMS 3D CMS Logo

ValidIsoTrkCalib.cc
Go to the documentation of this file.
1 //
2 // Package: ValidIsoTrkCalib
3 // Class: ValidIsoTrkCalib
4 //
5 
6 /*
7  Description: Validation for Isolated tracks Calibration
8 
9  Implementation:
10 See the twiki page for details:
11 https://twiki.cern.ch/twiki/bin/view/CMS/ValidIsoTrkCalib
12 
13 */
14 
15 //
16 // Original Author: Andrey Pozdnyakov
17 // Created: Tue Nov 4 01:16:05 CET 2008
18 //
19 
20 // system include files
21 
22 #include <memory>
23 
24 // user include files
28 
29 //#include "FWCore/Framework/interface/Event.h"
30 //#include "FWCore/Framework/interface/MakerMacros.h"
31 //#include "FWCore/ParameterSet/interface/ParameterSet.h"
32 
41 
51 
53 
54 #include "TROOT.h"
55 #include "TFile.h"
56 #include "TTree.h"
57 #include "TH1F.h"
58 
59 
60 #include <fstream>
61 #include <map>
62 
63 using namespace edm;
64 using namespace std;
65 using namespace reco;
66 
68 public:
69  explicit ValidIsoTrkCalib(const edm::ParameterSet&);
70  ~ValidIsoTrkCalib() override;
71 
72  // double getDistInPlaneSimple(const GlobalPoint caloPoint, const GlobalPoint rechitPoint);
73 
74 private:
75 
76  void beginJob() override ;
77  void analyze(const edm::Event&, const edm::EventSetup&) override;
78  void endJob() override ;
79 
80 
81 
82  // ----------member data ---------------------------
83 
84  //Variables from HcalIsotrackAnalyzer
85 
88  double taECALCone_;
89  double taHCALCone_;
90 
91  const CaloGeometry* geo;
92  // nothing is done with these tags, so I leave it - cowden
95  std::vector<edm::InputTag> genecalLabel_;
96 
102 
103  //std::string m_inputTrackLabel;
104  //std::string m_hcalLabel;
105 
107  string AxB_;
109 
112  //string calibFactorsFileName_;
113  // string corrfile;
114 
116  //bool takeAllRecHits_, takeGenTracks_;
117  int gen, iso, pix;
118  float genPt[500], genPhi[500], genEta[500];
119  float isoPt[500], isoPhi[500], isoEta[500];
120  float pixPt[500], pixPhi[500], pixEta[500];
121  //int hbheiEta[5000],hbheiPhi[5000],hbheDepth[5000];
122  //float hbheEnergy[5000];
123 
124 
125  int NisoTrk;
126  float trackPt, trackE, trackEta, trackPhi;
127  float ptNear;
128  float ptrack, rvert;
129  //float eecal, ehcal;
130 
131  int MinNTrackHitsBarrel, MinNTECHitsEndcap;
132  double energyECALmip, maxPNear;
133  double energyMinIso, energyMaxIso;
134 
135 
136  Float_t emEnergy;
137  Float_t targetE;
138 
139  //TFile* rootFile;
140  // TTree* tree;
141  TTree *tTree, *fTree;
142 
143  Float_t xTrkEcal;
144  Float_t yTrkEcal;
145  Float_t zTrkEcal;
146 
147  Float_t xTrkHcal;
148  Float_t yTrkHcal;
149  Float_t zTrkHcal;
150 
151  int Nhits;
152  float eClustBefore; //Calo energy before calibration
153  float eClustAfter; //After calibration
154  float eTrack; //Track energy
155  float etaTrack;
156  float phiTrack;
157  float eECAL; // Energy deposited in ECAL
158  int numHits; //number of rechits
159  //float e3x3;
160 
168  int iEta;
169  int iPhi;
170  int iEtaTr;
171  int iPhiTr;
172  float iDr, delR;
173  int dietatr;
174  int diphitr;
175 
176  float iTime;
177  float HTime[100];
178  float e3x3Before;
179  float e3x3After;
180  float e5x5Before;
181  float e5x5After;
184  float PtNearBy;
185  float numVH, numVS, numValidTrkHits, numValidTrkStrips;
186 
188  // map<UInt_t, Float_t> CalibFactors;
189  // Bool_t ReadCalibFactors(string);
190 
191  TH1F *nTracks;
192 
194  // int Lumi_n;
195 
196 
197 };
198 
199 
201 
202 {
203 
204  //takeAllRecHits_=iConfig.getUntrackedParameter<bool>("takeAllRecHits");
205  takeGenTracks_=iConfig.getUntrackedParameter<bool>("takeGenTracks");
206 
207  tok_genTrack_ = consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("genTracksLabel"));
208  genhbheLabel_= iConfig.getParameter<edm::InputTag>("genHBHE");
209  //genhoLabel_=iConfig.getParameter<edm::InputTag>("genHO");
210  //genecalLabel_=iConfig.getParameter<std::vector<edm::InputTag> >("genECAL");
211 
212  // m_hcalLabel = iConfig.getUntrackedParameter<std::string> ("hcalRecHitsLabel","hbhereco");
213 
214  tok_hbhe_ = consumes<HBHERecHitCollection>(iConfig.getParameter<edm::InputTag>("hbheInput"));
215  tok_ho_ = consumes<HORecHitCollection>(iConfig.getParameter<edm::InputTag>("hoInput"));
216  //eLabel_=iConfig.getParameter<edm::InputTag>("eInput");
217  tok_track_ = consumes<reco::IsolatedPixelTrackCandidateCollection>(iConfig.getParameter<edm::InputTag>("HcalIsolTrackInput"));
218  tok_track1_ = consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("trackInput"));
219 
220  associationConeSize_=iConfig.getParameter<double>("associationConeSize");
221  allowMissingInputs_=iConfig.getUntrackedParameter<bool>("allowMissingInputs", true);
222 // outputFileName_=iConfig.getParameter<std::string>("outputFileName");
223  // calibFactorsFileName_=iConfig.getParameter<std::string>("calibFactorsFileName");
224 
225 
226  AxB_=iConfig.getParameter<std::string>("AxB");
227  calibrationConeSize_=iConfig.getParameter<double>("calibrationConeSize");
228 
229  MinNTrackHitsBarrel = iConfig.getParameter<int>("MinNTrackHitsBarrel");
230  MinNTECHitsEndcap = iConfig.getParameter<int>("MinNTECHitsEndcap");
231  energyECALmip = iConfig.getParameter<double>("energyECALmip");
232  energyMinIso = iConfig.getParameter<double>("energyMinIso");
233  energyMaxIso = iConfig.getParameter<double>("energyMaxIso");
234  maxPNear = iConfig.getParameter<double>("maxPNear");
235 
236  edm::ParameterSet parameters = iConfig.getParameter<edm::ParameterSet>("TrackAssociatorParameters");
237  edm::ConsumesCollector iC = consumesCollector();
238  parameters_.loadParameters( parameters, iC );
239  trackAssociator_.useDefaultPropagator();
240 
241  // taECALCone_=iConfig.getUntrackedParameter<double>("TrackAssociatorECALCone",0.5);
242  //taHCALCone_=iConfig.getUntrackedParameter<double>("TrackAssociatorHCALCone",0.6);
243 
244 }
245 
246 
248 {
249 
250  // do anything here that needs to be done at desctruction time
251  // (e.g. close files, deallocate resources etc.)
252 
253 }
254 
255 
256 // ------------ method called to for each event ------------
257 void
259 {
260  using namespace edm;
261 
262  try{
263  edm::ESHandle <HcalRespCorrs> recalibCorrs;
264  iSetup.get<HcalRespCorrsRcd>().get("recalibrate",recalibCorrs);
265  respRecalib = recalibCorrs.product();
266 
267  LogInfo("CalibConstants")<<" Loaded: OK ";
268 
269  }catch(const cms::Exception & e) {
270  LogWarning("CalibConstants")<<" Not Found!! ";
271  }
272 
274  iEvent.getByToken(tok_genTrack_, generalTracks);
275 
277  iEvent.getByToken(tok_track1_,isoProdTracks);
278 
279 
281  //edm::Handle<reco::TrackCollection> isoPixelTracks;
282  iEvent.getByToken(tok_track_,isoPixelTracks);
283 
284  /*
285  edm::Handle<EcalRecHitCollection> ecal;
286  iEvent.getByLabel(eLabel_,ecal);
287  const EcalRecHitCollection Hitecal = *(ecal.product());
288  */
289 
291  iEvent.getByToken(tok_hbhe_,hbhe);
292  const HBHERecHitCollection Hithbhe = *(hbhe.product());
293 
295  iSetup.get<CaloGeometryRecord>().get(pG);
296  geo = pG.product();
297 
298  const HcalGeometry* gHcal = static_cast<const HcalGeometry*>(geo->getSubdetectorGeometry(DetId::Hcal,HcalBarrel));
299  //Note: even though it says HcalBarrel, we actually get the whole Hcal detector geometry!
300 
301  // Lumi_n=iEvent.luminosityBlock();
302  parameters_.useEcal = true;
303  parameters_.useHcal = true;
304  parameters_.useCalo = false;
305  parameters_.useMuon = false;
306  //parameters_.dREcal = taECALCone_;
307  //parameters_.dRHcal = taHCALCone_;
308 
309  //cout<<"Hello World. TrackCollectionSize: "<< isoPixelTracks->size()<<endl;
310 
311  if (isoPixelTracks->empty()) return;
312 
313 
314 for (reco::TrackCollection::const_iterator trit=isoProdTracks->begin(); trit!=isoProdTracks->end(); trit++)
315  {
316 
317  reco::IsolatedPixelTrackCandidateCollection::const_iterator isoMatched=isoPixelTracks->begin();
318  //reco::TrackCollection::const_iterator isoMatched=isoPixelTracks->begin();
319  bool matched=false;
320 
321  //for (reco::IsolatedPixelTrackCandidateCollection::const_iterator trit = isoPixelTracks->begin(); trit!=isoPixelTracks->end(); trit++)
322  for (reco::IsolatedPixelTrackCandidateCollection::const_iterator it = isoPixelTracks->begin(); it!=isoPixelTracks->end(); it++)
323  //for (reco::TrackCollection::const_iterator it = isoPixelTracks->begin(); it!=isoPixelTracks->end(); it++)
324  {
325 
326  if (abs((trit->pt() - it->pt())/it->pt()) < 0.005 && abs(trit->eta() - it->eta()) < 0.01)
327  {
328  isoMatched=it;
329  matched=true;
330  break;
331  }
332  }
333  // CUT
334 
335  if (!matched) continue;
336  if(isoMatched->maxPtPxl() > maxPNear) continue;
337 
338  ptNear = isoMatched->maxPtPxl();
339  //cout<<"Point 0.1 isoMatch. ptnear: "<<ptNear<<endl;
340 
341 
342  // CUT
343  if (trit->hitPattern().numberOfValidHits()<MinNTrackHitsBarrel) continue;
344  if (fabs(trit->eta())>1.47&&trit->hitPattern().numberOfValidStripTECHits()<MinNTECHitsEndcap) continue;
345 
346  //cout<<"Point 0.2.1 after numofvalidhits HB: "<<trit->hitPattern().numberOfValidHits()<<endl;
347  //cout<<"Point 0.2.2 after numofvalidstrips HE: "<<trit->hitPattern().numberOfValidStripTECHits()<<endl;
348 
349  numVH = trit->hitPattern().numberOfValidHits();
350  numVS = trit->hitPattern().numberOfValidStripTECHits();
351 
352 
353 
354  trackE = sqrt(trit->px()*trit->px()+trit->py()*trit->py()+trit->pz()*trit->pz()+0.14*0.14);
355  trackPt = trit->pt();
356  trackEta = trit->eta();
357  trackPhi = trit->phi();
358 
359  emEnergy = isoMatched->energyIn();
360 
361  //cout<<"Point 0.3. Matched :: pt: "<<trit->pt()<<" wholeEnergy: "<<trackE<<" emEnergy: "<<emEnergy<<" eta: "<<etahcal<<" phi: "<<phihcal<<endl;
362  //cout<<"Point 0.4. EM energy in cone: "<<emEnergy<<" EtaHcal: "<<etahcal<<" PhiHcal: "<<phihcal<<endl;
363 
364  TrackDetMatchInfo info = trackAssociator_.associate(iEvent, iSetup,trackAssociator_.getFreeTrajectoryState(iSetup, *trit), parameters_);
365 
366  //float etaecal=info.trkGlobPosAtEcal.eta();
367  //float phiecal=info.trkGlobPosAtEcal.phi();
368  // float etahcal=info.trkGlobPosAtHcal.eta();
369  // float phihcal=info.trkGlobPosAtHcal.phi();
370 
371 
372  xTrkEcal=info.trkGlobPosAtEcal.x();
373  yTrkEcal=info.trkGlobPosAtEcal.y();
374  zTrkEcal=info.trkGlobPosAtEcal.z();
375 
376  xTrkHcal=info.trkGlobPosAtHcal.x();
377  yTrkHcal=info.trkGlobPosAtHcal.y();
378  zTrkHcal=info.trkGlobPosAtHcal.z();
379 
380  if (xTrkEcal==0 && yTrkEcal==0&& zTrkEcal==0) {cout<<"zero point at Ecal"<<endl; continue;}
381  if (xTrkHcal==0 && yTrkHcal==0&& zTrkHcal==0) {cout<<"zero point at Hcal"<<endl; continue;}
382 
383  /*GlobalVector trackMomAtEcal = info.trkMomAtEcal;
384  GlobalVector trackMomAtHcal = info.trkMomAtHcal;
385 
386  PxTrkHcal = trackMomAtHcal.x();
387  PyTrkHcal = trackMomAtHcal.y();
388  PzTrkHcal = trackMomAtHcal.z();
389  */
390 
391  GlobalPoint gPointEcal(xTrkEcal,yTrkEcal,zTrkEcal);
392  GlobalPoint gPointHcal(xTrkHcal,yTrkHcal,zTrkHcal);
393 
394  int iphitrue = -10;
395  int ietatrue = 100;
396  const HcalDetId tempId = gHcal->getClosestCell(gPointHcal);
397  ietatrue = tempId.ieta();
398  iphitrue = tempId.iphi();
399 
400 
401  MaxHit_struct MaxHit;
402 
403  MaxHit.hitenergy=-100.;
404 
405  //container for used recHits
406  std::vector<DetId> usedHits;
407  //
408  usedHits.clear();
409  //cout <<"Point 1. Entrance to HBHECollection"<<endl;
410 
411  //float dddeta = 1000.;
412  //float dddphi = 1000.;
413  //int iphitrue = 1234;
414  //int ietatrue = 1234;
415 
416  GlobalPoint gPhot;
417 
418 
419  for (HBHERecHitCollection::const_iterator hhit=Hithbhe.begin(); hhit!=Hithbhe.end(); hhit++)
420  {
421 
422  //check that this hit was not considered before and push it into usedHits
423  bool hitIsUsed=false;
424  for (uint32_t i=0; i<usedHits.size(); i++)
425  {
426  if (usedHits[i]==hhit->id()) hitIsUsed=true;
427  }
428  if (hitIsUsed) continue;
429  usedHits.push_back(hhit->id());
430  //
431 
432  // rof 16.05.2008 start: include the possibility for recalibration
433  float recal = 1;
434  // rof end
435 
436  GlobalPoint pos = geo->getPosition(hhit->detid());
437  //float phihit = pos.phi();
438  //float etahit = pos.eta();
439 
440  int iphihitm = (hhit->id()).iphi();
441  int ietahitm = (hhit->id()).ieta();
442  int depthhit = (hhit->id()).depth();
443  float enehit = hhit->energy() * recal;
444 
445  if (depthhit!=1) continue;
446 
447  /*
448  float dphi = fabs(phihcal - phihit);
449  if(dphi > 4.*atan(1.)) dphi = 8.*atan(1.) - dphi;
450  float deta = fabs(etahcal - etahit);
451  float dr = sqrt(dphi*dphi + deta*deta);
452  */
453 
454 
455  //double distAtHcal = getDistInPlaneTrackDir(gPointHcal, trackMomAtHcal, pos);
456  double distAtHcal = getDistInPlaneSimple(gPointHcal, pos);
457 
458  if(distAtHcal < associationConeSize_)
459  {
460 
461  for (HBHERecHitCollection::const_iterator hhit2=Hithbhe.begin(); hhit2!=Hithbhe.end(); hhit2++)
462  {
463  int iphihitm2 = (hhit2->id()).iphi();
464  int ietahitm2 = (hhit2->id()).ieta();
465  int depthhit2 = (hhit2->id()).depth();
466  float enehit2 = hhit2->energy() * recal;
467 
468  if (iphihitm==iphihitm2 && ietahitm==ietahitm2 && depthhit!=depthhit2) enehit = enehit+enehit2;
469 
470  }
471 
472 
473  //cout<<"IN CONE ieta: "<<ietahitm<<" iphi: "<<iphihitm<<" depthhit: "<<depthhit<<" dr: "<<dr<<" energy: "<<enehit<<endl;
474 
475  //Find a Hit with Maximum Energy
476 
477  if(enehit > MaxHit.hitenergy)
478  {
479  MaxHit.hitenergy = enehit;
480  MaxHit.ietahitm = (hhit->id()).ieta();
481  MaxHit.iphihitm = (hhit->id()).iphi();
482  MaxHit.dr = distAtHcal;
483  //MaxHit.depthhit = (hhit->id()).depth();
484  MaxHit.depthhit = 1;
485 
486  //gPhot = geo->getPosition(hhit->detid());
487  }
488 
489  }
490  } //end of all HBHE hits cycle
491 
492  usedHits.clear();
493 
494  //cout<<"Hottest ieta: "<<MaxHit.ietahitm<<" iphi: "<<MaxHit.iphihitm<<" dr: "<<MaxHit.dr<<endl;
495  //cout<<"Track ieta: "<<ietatrue<<" iphi: "<<iphitrue<<endl;
496 
497  //cout<<"Point 3. MaxHit :::En "<<MaxHit.hitenergy<<" ieta: "<<MaxHit.ietahitm<<" iphi: "<<MaxHit.iphihitm<<endl;
498 
499  Bool_t passCuts = kFALSE;
500  if(trackE > energyMinIso && trackE < energyMaxIso && emEnergy < energyECALmip && MaxHit.hitenergy > 0. && abs(MaxHit.ietahitm)<29) passCuts = kTRUE;
501 
502  //cout<<"Pont 0.1.1. trackE:"<<trackE<<" emEn: "<<emEnergy<<endl;
503 
504 
505  numHits=0;
506 
507  eClustBefore = 0.0;
508  eClustAfter = 0.0;
509  eBeforeDepth1 = 0.0;
510  eAfterDepth1 = 0.0;
511  eBeforeDepth2 = 0.0;
512  eAfterDepth2 = 0.0;
513  CentHitFactor = 0.0;
514  e3x3After = 0.0;
515  e3x3Before = 0.0;
516  e5x5After = 0.0;
517  e5x5Before = 0.0;
518 
519 
520  for (HBHERecHitCollection::const_iterator hhit=Hithbhe.begin(); hhit!=Hithbhe.end(); hhit++)
521  {
522 
523  //check that this hit was not considered before and push it into usedHits
524  bool hitIsUsed=false;
525  for (uint32_t i=0; i<usedHits.size(); i++)
526  {
527  if (usedHits[i]==hhit->id()) hitIsUsed=true;
528  }
529  if (hitIsUsed) continue;
530  usedHits.push_back(hhit->id());
531 
532  int DIETA = 100;
533  if(MaxHit.ietahitm*(hhit->id()).ieta()>0)
534  {
535  DIETA = MaxHit.ietahitm - (hhit->id()).ieta();
536  }
537  if(MaxHit.ietahitm*(hhit->id()).ieta()<0)
538  {
539  DIETA = MaxHit.ietahitm - (hhit->id()).ieta();
540  DIETA = DIETA>0 ? DIETA-1 : DIETA+1;
541  }
542 
543  int DIPHI = abs(MaxHit.iphihitm - (hhit->id()).iphi());
544  DIPHI = DIPHI>36 ? 72-DIPHI : DIPHI;
545 
546 
547 
548  int numbercell=100; //always collect Wide clastor!
549 
550  //if(AxB_=="5x5" || AxB_=="3x3" || AxB_=="7x7"|| AxB_=="Cone")
551 
552  //if(AxB_=="3x3") numbercell = 1;
553  //if(AxB_=="5x5") numbercell = 2;
554  //if(AxB_=="Cone") numbercell = 1000;
555 
556  if( abs(DIETA)<=numbercell && (abs(DIPHI)<=numbercell || ( abs(MaxHit.ietahitm)>=20 && abs(DIPHI)<=numbercell+1)) )
557  {
558  const GlobalPoint pos2 = geo->getPosition(hhit->detid());
559 
560  if(passCuts && hhit->energy()>0)
561  {
562 
563 
564  float factor = 0.0;
565  // factor = CalibFactors[hhit->id()];
566  factor = respRecalib -> getValues(hhit->id())->getValue();
567 
568 
569 
570  //if(i<5){cout<<" calib factors: "<<factor<<" ij "<<100*i+j<<endl;}
571 
572  if (hhit->id().ieta() == MaxHit.ietahitm && hhit->id().iphi() == MaxHit.iphihitm) CentHitFactor=factor;
573 
574  if (hhit->id().ieta() == MaxHit.ietahitm && hhit->id().iphi() == MaxHit.iphihitm) iTime = hhit->time();
575 
576  if(AxB_!="3x3" && AxB_!="5x5" && AxB_!="Cone") LogWarning(" AxB ")<<" Not supported: "<< AxB_;
577 
578 
579  if (abs(DIETA)<=2 && (abs(DIPHI)<=2 || ((abs(MaxHit.ietahitm)>20 && abs(DIPHI)<=4) && !( (abs(MaxHit.ietahitm)==21 || abs(MaxHit.ietahitm)==22) && abs((hhit->id()).ieta())<=20 && abs(DIPHI)>2) )) )
580  {
581 
582  e5x5Before += hhit->energy();
583  e5x5After += hhit->energy()*factor;
584  }
585 
586 
587  if (abs(DIETA)<=1 && (abs(DIPHI)<=1 || ((abs(MaxHit.ietahitm)>20 && abs(DIPHI)<=2) && !(abs(MaxHit.ietahitm)==21 && abs((hhit->id()).ieta())<=20 && abs(DIPHI)>1) )) )
588  {
589 
590  e3x3Before += hhit->energy();
591  e3x3After += hhit->energy()*factor;
592  }
593 
594 
595  if(AxB_=="5x5")
596  {
597 
598  if (abs(DIETA)<=2 && (abs(DIPHI)<=2 || ( abs(MaxHit.ietahitm)>20 && abs(DIPHI)<=4) ) )
599  {
600  if (abs(MaxHit.ietahitm)==21 && abs((hhit->id()).ieta())<=20 && abs(DIPHI)>3) continue;
601 
602  HTime[numHits]= hhit->time();
603  numHits++;
604 
605  eClustBefore += hhit->energy();
606  eClustAfter += hhit->energy()*factor;
607 
608  if ((hhit->id().depth() == 1) && (abs(hhit->id().ieta()) > 17) && (abs(hhit->id().ieta()) < 29))
609  {
610  eBeforeDepth1 += hhit->energy();
611  eAfterDepth1 += hhit->energy()*factor;
612  }
613  else if((hhit->id().depth() == 2) && (abs(hhit->id().ieta()) > 17) && (abs(hhit->id().ieta()) < 29))
614  {
615  eBeforeDepth2 += hhit->energy();
616  eAfterDepth2 += hhit->energy()*factor;
617  }
618  }
619  }//end of 5x5
620 
621 
622 
623  if(AxB_=="3x3")
624  {
625  if (abs(DIETA)<=1 && (abs(DIPHI)<=1 || ( abs(MaxHit.ietahitm)>20 && abs(DIPHI)<=2) ) )
626  {
627  if (abs(MaxHit.ietahitm)==21 && abs((hhit->id()).ieta())<=20 && abs(DIPHI)>2) continue;
628 
629  HTime[numHits]= hhit->time();
630  numHits++;
631 
632  eClustBefore += hhit->energy();
633  eClustAfter += hhit->energy() * factor;
634 
635  if ((hhit->id().depth() == 1) && (abs(hhit->id().ieta()) > 17) && (abs(hhit->id().ieta()) < 29))
636  {
637  eBeforeDepth1 += hhit->energy();
638  eAfterDepth1 += hhit->energy() * factor;
639  }
640  else if((hhit->id().depth() == 2) && (abs(hhit->id().ieta()) > 17) && (abs(hhit->id().ieta()) < 29))
641  {
642  eBeforeDepth2 += hhit->energy();
643  eAfterDepth2 += hhit->energy() * factor;
644  }
645 
646  }
647  }//end of 3x3
648 
649 
650  if (AxB_=="Cone" && getDistInPlaneSimple(gPointHcal, pos2) < calibrationConeSize_) {
651 
652  HTime[numHits]= hhit->time();
653  numHits++;
654 
655  eClustBefore += hhit->energy();
656  eClustAfter += hhit->energy() * factor;
657 
658  if ((hhit->id().depth() == 1) && (abs(hhit->id().ieta()) > 17) && (abs(hhit->id().ieta()) < 29))
659  {
660  eBeforeDepth1 += hhit->energy();
661  eAfterDepth1 += hhit->energy() * factor;
662  }
663  else if((hhit->id().depth() == 2) && (abs(hhit->id().ieta()) > 17) && (abs(hhit->id().ieta()) < 29))
664  {
665  eBeforeDepth2 += hhit->energy();
666  eAfterDepth2 += hhit->energy() * factor;
667  }
668 
669 
670  }//end of Cone
671 
672  }//end of passCuts
673 
674  }//end of DIETA DIPHI
675 
676  } //end of associatedcone HBHE hits cycle
677 
678 
679  int dieta_M_P = 100;
680  int diphi_M_P = 100;
681  if(MaxHit.ietahitm*ietatrue>0) {dieta_M_P = abs (MaxHit.ietahitm-ietatrue);}
682  if(MaxHit.ietahitm*ietatrue<0) {dieta_M_P = abs(MaxHit.ietahitm-ietatrue)-1;}
683  diphi_M_P = abs(MaxHit.iphihitm-iphitrue);
684  diphi_M_P = diphi_M_P>36 ? 72-diphi_M_P : diphi_M_P;
685 
686 
687  if(passCuts)
688 
689  {
690  eventNumber = iEvent.id().event();
691  runNumber = iEvent.id().run();
692 
693  eCentHitBefore = MaxHit.hitenergy;
694  eCentHitAfter = MaxHit.hitenergy*CentHitFactor;
695  eECAL = emEnergy;
696  numValidTrkHits = numVH;
697  numValidTrkStrips = numVS;
698  PtNearBy = ptNear;
699 
700 
701  eTrack = trackE;
702  phiTrack = trackPhi;
703  etaTrack = trackEta;
704 
705  iEta = MaxHit.ietahitm;
706  iPhi = MaxHit.iphihitm;
707 
708  iEtaTr = ietatrue;
709  iPhiTr = iphitrue;
710  iDr = sqrt(diphi_M_P*diphi_M_P+dieta_M_P*dieta_M_P);
711  delR = MaxHit.dr;
712  dietatr = dieta_M_P;
713  diphitr = diphi_M_P;
714 
715  fTree -> Fill();
716  }
717 
718  } //end of isoProdTracks cycle
719 
720 
721 
722 /* ------------------ Some stuff for general tracks ---------- ----*/
723  //cout<<" generalTracks Size: "<< generalTracks->size()<<endl;
724  int n = generalTracks->size();
725  nTracks->Fill(n);
726 
727  if(takeGenTracks_ && iEvent.id().event()%10==1)
728  {
729  gen = generalTracks->size();
730  iso = isoProdTracks->size();
731  pix = isoPixelTracks->size();
732 
733  genPt[0] = -33;
734  genPhi[0] = -33;
735  genEta[0] = -33;
736 
737  isoPt[0] = -33;
738  isoPhi[0] = -33;
739  isoEta[0] = -33;
740 
741  pixPt[0] = -33;
742  pixPhi[0] = -33;
743  pixEta[0] = -33;
744 
745  Int_t gencount=0, isocount=0, pixcount=0;
746  for (reco::TrackCollection::const_iterator gentr=generalTracks->begin(); gentr!=generalTracks->end(); gentr++)
747  {
748  genPt[gencount] = gentr->pt();
749  genPhi[gencount] = gentr->phi();
750  genEta[gencount] = gentr->eta();
751  gencount++;
752  }
753 
754  for (reco::TrackCollection::const_iterator isotr=isoProdTracks->begin(); isotr!=isoProdTracks->end(); isotr++)
755  {
756  isoPt[isocount] = isotr->pt();
757  isoPhi[isocount] = isotr->phi();
758  isoEta[isocount] = isotr->eta();
759  isocount++;
760  }
761 
762  for (reco::IsolatedPixelTrackCandidateCollection::const_iterator pixtr=isoPixelTracks->begin(); pixtr!=isoPixelTracks->end(); pixtr++)
763  {
764  pixPt[pixcount] = pixtr->pt();
765  pixPhi[pixcount] = pixtr->phi();
766  pixEta[pixcount] = pixtr->eta();
767  pixcount++;
768  }
769  }
770 
771  tTree -> Fill();
772 
773 }
774 
775 // ------------ method called once each job just before starting event loop ------------
776 void
778 {
779 
780  // if(!ReadCalibFactors(calibFactorsFileName_.c_str() )) {cout<<"Cant read file with cailib coefficients!! ---"<<endl;}
781 
782 // try{
783 // edm::ESHandle <HcalRespCorrs> recalibCorrs;
784 // iSetup.get<HcalRespCorrsRcd>().get("recalibrate",recalibCorrs);
785 // respRecalib = recalibCorrs.product();
786 //
787 // LogInfo("CalibConstants")<<" Loaded: OK ";
788 //
789 // }catch(const cms::Exception & e) {
790 // LogWarning("CalibConstants")<<" Not Found!! ";
791 // }
792 
793 
794  // rootFile = new TFile(outputFileName_.c_str(),"RECREATE");
795 
796  //@@@@@@@@@@@@@
797  //TFileDirectory ValDir = fs->mkdir("Validation");
798 
799  nTracks = fs->make<TH1F>("nTracks","general;number of general tracks",11,-0.5,10.5);
800 
801 
802  tTree = fs->make<TTree>("tTree", "Tree for gen info");
803 
804  fTree = fs->make<TTree>("fTree", "Tree for IsoTrack Calibration");
805 
806  fTree->Branch("eventNumber", &eventNumber, "eventNumber/I");
807  fTree->Branch("runNumber", &runNumber, "runNumber/I");
808 
809  fTree->Branch("eClustBefore", &eClustBefore,"eClustBefore/F");
810  fTree->Branch("eClustAfter", &eClustAfter,"eClustAfter/F");
811  fTree->Branch("eTrack", &eTrack, "eTrack/F");
812  fTree->Branch("etaTrack", &etaTrack, "etaTrack/F");
813  fTree->Branch("phiTrack", &phiTrack, "phiTrack/F");
814 
815  fTree->Branch("numHits", &numHits, "numHits/I");
816  fTree->Branch("eECAL", &eECAL, "eECAL/F");
817  fTree->Branch("PtNearBy", &PtNearBy, "PtNearBy/F");
818  fTree->Branch("numValidTrkHits", &numValidTrkHits, "numValidTrkHits/F");
819  fTree->Branch("numValidTrkStrips", &numValidTrkStrips, "numValidTrkStrips/F");
820 
821  fTree->Branch("eBeforeDepth1", &eBeforeDepth1, "eBeforeDepth1/F");
822  fTree->Branch("eBeforeDepth2", &eBeforeDepth2, "eBeforeDepth2/F");
823  fTree->Branch("eAfterDepth1", &eAfterDepth1, "eAfterDepth1/F");
824  fTree->Branch("eAfterDepth2", &eAfterDepth2, "eAfterDepth2/F");
825 
826  fTree->Branch("e3x3Before", &e3x3Before, "e3x3Before/F");
827  fTree->Branch("e3x3After", &e3x3After, "e3x3After/F");
828  fTree->Branch("e5x5Before", &e5x5Before, "e5x5Before/F");
829  fTree->Branch("e5x5After", &e5x5After, "e5x5After/F");
830 
831  fTree->Branch("eCentHitAfter", &eCentHitAfter, "eCentHitAfter/F");
832  fTree->Branch("eCentHitBefore", &eCentHitBefore, "eCentHitBefore/F");
833  fTree->Branch("iEta", &iEta, "iEta/I");
834  fTree->Branch("iPhi", &iPhi, "iPhi/I");
835 
836 
837  fTree->Branch("iEtaTr", &iEtaTr, "iEtaTr/I");
838  fTree->Branch("iPhiTr", &iPhiTr, "iPhiTr/I");
839  fTree->Branch("dietatr", &dietatr, "dietatr/I");
840  fTree->Branch("diphitr", &diphitr, "diphitr/I");
841  fTree->Branch("iDr", &iDr, "iDr/F");
842  fTree->Branch("delR", &delR, "delR/F");
843 
844  fTree->Branch("iTime", &iTime, "iTime/F");
845  fTree->Branch("HTime", HTime, "HTime[numHits]/F");
846 
847  fTree->Branch("xTrkEcal", &xTrkEcal, "xTrkEcal/F");
848  fTree->Branch("yTrkEcal", &yTrkEcal, "yTrkEcal/F");
849  fTree->Branch("zTrkEcal", &zTrkEcal, "zTrkEcal/F");
850  fTree->Branch("xTrkHcal", &xTrkHcal, "xTrkHcal/F");
851  fTree->Branch("yTrkHcal", &yTrkHcal, "yTrkHcal/F");
852  fTree->Branch("zTrkHcal", &zTrkHcal, "zTrkHcal/F");
853 
854  if(takeGenTracks_) {
855  tTree->Branch("gen", &gen, "gen/I");
856  tTree->Branch("iso", &iso, "iso/I");
857  tTree->Branch("pix", &pix, "pix/I");
858  tTree->Branch("genPt", genPt, "genPt[gen]/F");
859  tTree->Branch("genPhi", genPhi, "genPhi[gen]/F");
860  tTree->Branch("genEta", genEta, "genEta[gen]/F");
861 
862  tTree->Branch("isoPt", isoPt, "isoPt[iso]/F");
863  tTree->Branch("isoPhi", isoPhi, "isoPhi[iso]/F");
864  tTree->Branch("isoEta", isoEta, "isoEta[iso]/F");
865 
866  tTree->Branch("pixPt", pixPt, "pixPt[pix]/F");
867  tTree->Branch("pixPhi", pixPhi, "pixPhi[pix]/F");
868  tTree->Branch("pixEta", pixEta, "pixEta[pix]/F");
869  }
870 
871 
872 }
873 
874 
875 
876 // ------------ method called once each job just after ending the event loop ------------
877 void
879 {
880  // rootFile->Write();
881  //rootFile->Close();
882 }
883 
884 //define this as a plug-in
RunNumber_t run() const
Definition: EventID.h:39
edm::EDGetTokenT< reco::TrackCollection > tok_genTrack_
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:41
T getUntrackedParameter(std::string const &, T const &) const
const unsigned int nTracks(const reco::Vertex &sv)
static const TGPicture * info(bool iBackgroundIsBlack)
~ValidIsoTrkCalib() override
virtual example_stream void analyze(const edm::Event &, const edm::EventSetup &) override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
void analyze(const edm::Event &, const edm::EventSetup &) override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::vector< HBHERecHit >::const_iterator const_iterator
TrackAssociatorParameters parameters_
math::XYZPoint trkGlobPosAtHcal
void beginJob()
Definition: Breakpoints.cc:15
edm::EDGetTokenT< reco::IsolatedPixelTrackCandidateCollection > tok_track_
int iEvent
Definition: GenABIO.cc:230
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
const CaloGeometry * geo
T sqrt(T t)
Definition: SSEVec.h:18
void beginJob() override
int ieta() const
get the cell ieta
Definition: HcalDetId.h:56
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
TrackDetectorAssociator trackAssociator_
edm::EDGetTokenT< HBHERecHitCollection > tok_hbhe_
const HcalRespCorrs * respRecalib
def gen(fragment, howMuch)
Production test section ####.
ValidIsoTrkCalib(const edm::ParameterSet &)
std::vector< edm::InputTag > genecalLabel_
const_iterator end() const
int iphi() const
get the cell iphi
Definition: HcalDetId.cc:103
T const * product() const
Definition: Handle.h:81
const T & get() const
Definition: EventSetup.h:59
edm::Service< TFileService > fs
edm::EventID id() const
Definition: EventBase.h:60
edm::EDGetTokenT< HORecHitCollection > tok_ho_
fixed size matrix
HLT enums.
math::XYZPoint trkGlobPosAtEcal
Track position at different parts of the calorimeter.
const JetExtendedData & getValue(const Container &, const reco::JetBaseRef &)
get value for the association. Throw exception if no association found
double getDistInPlaneSimple(const GlobalPoint caloPoint, const GlobalPoint rechitPoint)
Definition: ConeDefinition.h:9
edm::EDGetTokenT< reco::TrackCollection > tok_track1_
T const * product() const
Definition: ESHandle.h:86
void endJob() override
DetId getClosestCell(const GlobalPoint &r) const override
const_iterator begin() const