CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
AlCaHOCalibProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 // Dec 2015 : Added bool m_cosmic to choose cosmic or collision run through python file
3 
4 //integrate the code with some 8_0_X or 7_6_X recent IB and run the
5 // following tests: 4.22, 8.0, 25.0, 140.53. You can always activate them using
6 //runTheMatrix.py -l 4.22
7 
8 // 7th Nov 2015 : tmpHOCalib.ecal03 = iso05.sumPt; // iso03.emEt+muonenr.em;
9 // tmpHOCalib.inslumi=lumiScale->begin()->pileup();
10 //
11 // April 2015 : Remove all digi part
12 // Also look for HO geometry in CMSSW in parallel with stanalone one.
13 // Official one has problem in reco geometry, particularly tiles at the edge of wheel
14 // Remove all histogrammes except occupancy one
15 // Remove Trigger bits
16 // But addition of these variables, ilumi (analyser), inslumi (analyser), nprim
17 
18 // Feb09 2009
19 // Move the initialisation of SteppingHelixPropagator from ::beginJob() to ::produce()
20 //
21 // Oct3 2008
22 // Difference in tag V00-02-45 with previous code
23 
24 // 1. One new object on data format, which was realised in
25 // CRUZET data analysis.
26 //2. Remove all histogram and cout in the code
27 //3. An upgrade in code, which increases the acceptance of
28 // muon near the edge (this also realised in CRUZET data).
29 // Difference in wrt V00-02-45
30 // 1. initialisation tmpHOCalib.htime = -1000;
31 // 2. By mistake HLT was commented out
32 
33 // Package: AlCaHOCalibProducer
34 // Class: AlCaHOCalibProducer
35 //
68 //
69 // Original Author: Gobinda Majumder
70 // Created: Fri Jul 6 17:17:21 CEST 2007
71 //
72 //
73 
74 
75 // system include files
76 #include <memory>
77 
78 // user include files
81 
89 
93 
98 
102 
104 
107 
110 
113 
116 #include "CLHEP/Vector/LorentzVector.h"
117 
120 
125 
132 
133 // Necessary includes for identify severity of flagged problems in HO rechits
134 //#include "RecoLocalCalo/HcalRecAlgos/interface/HcalCaloFlagLabels.h"
139 #include "CLHEP/Units/GlobalPhysicalConstants.h"
140 #include "CLHEP/Units/GlobalSystemOfUnits.h"
141 
142 #include "TH2F.h"
143 
144 /* C++ Headers */
145 #include <string>
146 
147 #include <iostream>
148 #include <fstream>
149 //
150 // class decleration
151 //
152 
153 
155 public:
156  explicit AlCaHOCalibProducer(const edm::ParameterSet&);
158 
162 
163 
164 private:
165  virtual void produce(edm::Event&, const edm::EventSetup&) override;
166  virtual void beginJob() override ;
167  virtual void endJob() override ;
168  virtual void beginRun(edm::Run const &, edm::EventSetup const &) override;
169  void fillHOStore(const reco::TrackRef& ncosm,
170  HOCalibVariables& tmpHOCalib,
171  std::auto_ptr<HOCalibVariableCollection> &hostore,
172  int Noccu_old, int indx,
175  const edm::Event& iEvent, const edm::EventSetup& iSetup);
176  void findHOEtaPhi(int iphsect, int& ietaho, int& iphiho);
177  // virtual void endRun(edm::Run const &, edm::EventSetup const &) override;
178  // ----------member data ---------------------------
179 
180  float xhor0; //x-position in ring 0
181  float yhor0; //y-position in ring 0
182  float xhor1; //x-position in ring 1
183  float yhor1; //y-position in ring 1
184  int iring; //Ring number -2,-1,0,1,2
185 
186  float localxhor0; //local x-distance from edege in ring 0
187  float localyhor0; //local y-distance from edege in ring 0
188  float localxhor1; //local x-distance from edege in ring 1
189  float localyhor1; //local y-distance from edege in ring 1
190 
191  TH2F* ho_occupency[5];
193  bool m_cosmic;
194 
195  const int netabin= 16;
196  const int nphimx = 72;
197  const int netamx = 32;
198  const int ncidmx = 5;
199  const double rHOL0 = 382.0;
200  const double rHOL1 = 407.0;
201 
202  edm::InputTag muonTags_; // cosmicMuons (for cosmic run) or muons (for collision run)
203 
207  // edm::EDGetTokenT<LumiDetails> tok_lumi_;
209 
213 
214  bool m_hbinfo;
216  int m_endTS;
217  double m_sigma;
218 
220  int Noccu;
221  int nRuns;
222 
223  // SteppingHelixPropagator* stepProp;
224  FreeTrajectoryState getFreeTrajectoryState( const reco::Track& tk, const MagneticField* field, int itag, bool dir);
225 
226  unsigned int Ntp; // # of HLT trigger paths (should be the same for all events!)
227  std::map<std::string, bool> fired;
228 
229  //hcal severity ES
231  // edm::ESHandle<HcalChannelQuality> theHcalChStatus;
233  int Nevents;
234 };
235 
236 //
237 // constants, enums and typedefs
238 //
239 
240 //
241 // static data member definitions
242 //
243 
244 //
245 // constructors and destructor
246 //
248  //register your products
249 
250  m_hbinfo = iConfig.getUntrackedParameter<bool>("hbinfo", false);
251  m_sigma = iConfig.getUntrackedParameter<double>("sigma", 0.05);
252  m_occupancy = iConfig.getUntrackedParameter<bool>("plotOccupancy", false);
253  m_cosmic = iConfig.getUntrackedParameter<bool>("CosmicData", false);
254 
255  // keep InputTag muonTags_ since it is used below. - cowden
256  muonTags_ = iConfig.getUntrackedParameter<edm::InputTag>("muons");
257  tok_muonsCosmic_ = consumes<reco::TrackCollection>(muonTags_);
258  tok_muons_ = consumes<edm::View<reco::Muon> >(muonTags_);
259  tok_vertex_ = consumes<reco::VertexCollection >(iConfig.getParameter<edm::InputTag>("vertexTags"));
260  // tok_lumi_ = consumes<LumiDetails ,edm::InLumi>(iConfig.getParameter<edm::InputTag>("lumiTags"));
261  tok_lumi_ = consumes<LumiScalersCollection>(iConfig.getParameter<edm::InputTag>("lumiTags"));
262  tok_ho_ = consumes<HORecHitCollection>(iConfig.getParameter<edm::InputTag>("hoInput"));
263  tok_hbhe_ = consumes<HBHERecHitCollection>(iConfig.getParameter<edm::InputTag>("hbheInput"));
264  tok_tower_ = consumes<CaloTowerCollection>(iConfig.getParameter<edm::InputTag>("towerInput"));
265 
266  produces<HOCalibVariableCollection>("HOCalibVariableCollection").setBranchAlias("HOCalibVariableCollection");
267 
268  if (m_occupancy) {
270 
271  char title[200];
272 
273  for (int ij=0; ij<5; ij++) {
274  sprintf(title, "ho_occupency (>%i #sigma)", ij+2);
275  ho_occupency[ij] = fs->make<TH2F>(title, title, netamx+1, -netamx-0.5, netamx/2+0.5, nphimx, 0.5, nphimx+0.5);
276  }
277  }
278 }
279 
280 
282 {
283 
284  // do anything here that needs to be done at desctruction time
285  // (e.g. close files, deallocate resources etc.)
286 }
287 
288 
289 //
290 // member functions
291 //
292 
293 // ------------ method called to produce the data ------------
294 void
296 {
297 
298  int irun = iEvent.id().run();
299  // int ilumi = iEvent.luminosityBlock();
300 
301  Nevents++;
302 
303  if (Nevents%5000==1) edm::LogInfo("HOCalib") <<"AlCaHOCalibProducer Processing event # "<<Nevents<<" "<<Noccu<<" "<<irun<<" "<<iEvent.id().event();
304 
305  std::auto_ptr<HOCalibVariableCollection> hostore (new HOCalibVariableCollection);
306 
308  edm::Handle<edm::View<reco::Muon> > collisionmuon;
309 
310  bool muonOK(true);
311  HOCalibVariables tmpHOCalib;
312  tmpHOCalib.nprim = -1;
313  tmpHOCalib.inslumi=-1.;
314 
315  if (m_cosmic) {
316  iEvent.getByToken(tok_muonsCosmic_, cosmicmuon);
317  muonOK = (cosmicmuon.isValid() && cosmicmuon->size()>0);
318  } else {
319  iEvent.getByToken(tok_muons_,collisionmuon);
320  muonOK = (collisionmuon.isValid() && collisionmuon->size()>0);
321 
322  if (iEvent.isRealData()) {
324  iEvent.getByToken(tok_vertex_, primaryVertices);
325  if (primaryVertices.isValid()) { tmpHOCalib.nprim = primaryVertices->size();}
326 
327  tmpHOCalib.inslumi=0.;
328 
330  iEvent.getByToken(tok_lumi_, lumiScale);
331 
332  if (lumiScale.isValid()) {
333  tmpHOCalib.inslumi=lumiScale->begin()->pileup();
334  }
335  }
336  }
337 
338  if (muonOK) {
339 
340  int Noccu_old = Noccu;
342  if (m_cosmic) {
343  int indx(0);
344  for(reco::TrackCollection::const_iterator ncosm = cosmicmuon->begin();
345  ncosm != cosmicmuon->end(); ++ncosm,++indx) {
346  if ((*ncosm).ndof() < 15) continue;
347  if ((*ncosm).normalizedChi2() >30.0) continue;
348  reco::TrackRef tRef = reco::TrackRef(cosmicmuon,indx);
349  fillHOStore(tRef,tmpHOCalib,hostore,Noccu_old,indx,cosmicmuon,muon1,
350  iEvent, iSetup);
351  }
352  } else {
353  for( muon1 = collisionmuon->begin(); muon1 < collisionmuon->end(); muon1++ ) {
354  if ((!muon1->isGlobalMuon()) || (!muon1->isTrackerMuon())) continue;
355  reco::TrackRef ncosm = muon1->innerTrack();
356  fillHOStore(ncosm,tmpHOCalib,hostore,Noccu_old,0,cosmicmuon,muon1,
357  iEvent, iSetup);
358  }
359  }
360  }
361 
362  iEvent.put(hostore, "HOCalibVariableCollection");
363 
364 }
365 
366 // ------------ method called once each job just before starting event loop ------------
367 void
369 {
370  Nevents = 0;
371  nRuns = 0;
372  Noccu = 0;
373 }
374 
375 // ------------ method called once each job just after ending the event loop ------------
377  if (m_occupancy) {
378  for (int ij=0; ij<5; ij++) {
379  ho_occupency[ij]->Scale(1./std::max(1,Noccu));
380  }
381  }
382  edm::LogInfo("HOCalib") <<" AlCaHOCalibProducer processed event "<< Nevents;
383 }
384 
385 // ------------ method called once each job just after ending the event loop ------------
386 //void
387 // AlCaHOCalibProducer::endRun(edm::Run const & run, edm::EventSetup & es) {}
388 // ------------ method called once each job just before starting event loop ------------
389 void
391  const edm::EventSetup & es) {
392 
393  // HCAL channel status map ****************************************
394  edm::ESHandle<HcalChannelQuality> hcalChStatus;
395  es.get<HcalChannelQualityRcd>().get("withTopo", hcalChStatus );
396  theHcalChStatus = hcalChStatus.product();
397 }
398 
399 
401  HOCalibVariables& tmpHOCalib,
402  std::auto_ptr<HOCalibVariableCollection> &hostore,
403  int Noccu_old, int indx,
406  const edm::Event& iEvent,
407  const edm::EventSetup& iSetup) {
408 
410  iSetup.get<CaloGeometryRecord>().get(pG);
411  const CaloGeometry* geo = pG.product();
412  const CaloSubdetectorGeometry* gHO =
414 
415  // Get Hcal Severity Level Computer, so that the severity of each rechit flag/status may be determined
416  // edm::ESHandle<HcalSeverityLevelComputer> hcalSevLvlComputerHndl;
418  const HcalSeverityLevelComputer* hcalSevLvlComputer = hcalSevLvlComputerHndl.product();
419 
420  int charge = ncosm->charge();
421 
422  double innerr = (*ncosm).innerPosition().Perp2();
423  double outerr = (*ncosm).outerPosition().Perp2();
424  int iiner = (innerr <outerr) ? 1 : 0;
425 
426  //---------------------------------------------------
427  // in_to_out Dir in_to_out Dir
428  // StandAlone ^ ^ Cosmic ^ |
429  // | | | v
430  //---------------------------------------------------Y=0
431  // StandAlone | | Cosmic ^ |
432  // v v | v
433  //----------------------------------------------------
434 
435  double posx, posy, posz;
436  double momx, momy, momz;
437 
438  if (iiner==1) {
439  posx = (*ncosm).innerPosition().X();
440  posy = (*ncosm).innerPosition().Y();
441  posz = (*ncosm).innerPosition().Z();
442 
443  momx = (*ncosm).innerMomentum().X();
444  momy = (*ncosm).innerMomentum().Y();
445  momz = (*ncosm).innerMomentum().Z();
446 
447  } else {
448  posx = (*ncosm).outerPosition().X();
449  posy = (*ncosm).outerPosition().Y();
450  posz = (*ncosm).outerPosition().Z();
451 
452  momx = (*ncosm).outerMomentum().X();
453  momy = (*ncosm).outerMomentum().Y();
454  momz = (*ncosm).outerMomentum().Z();
455  }
456 
457 
458  PositionType trkpos(posx, posy, posz);
459 
460  CLHEP::Hep3Vector tmpmuon3v(posx, posy, posz);
461  CLHEP::Hep3Vector tmpmuondir(momx, momy, momz);
462 
463  bool samedir = (tmpmuon3v.dot(tmpmuondir) >0) ? true : false;
464  for (int ij=0; ij<3; ij++) {tmpHOCalib.caloen[ij] = 0.0;}
465  int inearbymuon = 0;
466  localxhor0 = localyhor0 = 20000; //GM for 22OCT07 data
467 
468  if (m_cosmic) {
469  int ind(0);
470  for(reco::TrackCollection::const_iterator ncosmcor=cosmicmuon->begin();
471  ncosmcor != cosmicmuon->end(); ++ncosmcor,++ind) {
472  if (indx==ind) continue;
473  CLHEP::Hep3Vector tmpmuon3vcor;
474  CLHEP::Hep3Vector tmpmom3v;
475  if (iiner==1) {
476  tmpmuon3vcor = CLHEP::Hep3Vector((*ncosmcor).innerPosition().X(),(*ncosmcor).innerPosition().Y(),(*ncosmcor).innerPosition().Z());
477  tmpmom3v = CLHEP::Hep3Vector((*ncosmcor).innerMomentum().X(),(*ncosmcor).innerMomentum().Y(),(*ncosmcor).innerMomentum().Z());
478  } else {
479  tmpmuon3vcor = CLHEP::Hep3Vector((*ncosmcor).outerPosition().X(),(*ncosmcor).outerPosition().Y(),(*ncosmcor).outerPosition().Z());
480  tmpmom3v = CLHEP::Hep3Vector((*ncosmcor).outerMomentum().X(),(*ncosmcor).outerMomentum().Y(),(*ncosmcor).outerMomentum().Z());
481 
482  }
483 
484  if (tmpmom3v.mag()<0.2 || (*ncosmcor).ndof()<5) continue;
485 
486  double angle = tmpmuon3v.angle(tmpmuon3vcor);
487  if (angle < 7.5*CLHEP::deg) {inearbymuon=1;} // break;}
488 
489  // if (muonTagsi_.label() =="cosmicMuons") {
490  if (angle <7.5*CLHEP::deg) { tmpHOCalib.caloen[0] +=1.;}
491  if (angle <15.0*CLHEP::deg) { tmpHOCalib.caloen[1] +=1.;}
492  if (angle <35.0*CLHEP::deg) { tmpHOCalib.caloen[2] +=1.;}
493  }
494  } else {
495  // if (muonTags_.label() =="muons") {
497  iEvent.getByToken(tok_tower_, calotower);
498 
499  for (CaloTowerCollection::const_iterator calt = calotower->begin();
500  calt !=calotower->end(); calt++) {
501  //CMSSW_2_1_x const math::XYZVector towermom = (*calt).momentum();
502  double ith = (*calt).momentum().theta();
503  double iph = (*calt).momentum().phi();
504 
505  CLHEP::Hep3Vector calo3v(sin(ith)*cos(iph), sin(ith)*sin(iph), cos(ith));
506 
507  double angle = tmpmuon3v.angle(calo3v);
508 
509  if (angle < 7.5*CLHEP::deg) {tmpHOCalib.caloen[0] += calt->emEnergy()+calt->hadEnergy();}
510  if (angle < 15*CLHEP::deg) {tmpHOCalib.caloen[1] += calt->emEnergy()+calt->hadEnergy();}
511  if (angle < 35*CLHEP::deg) {tmpHOCalib.caloen[2] += calt->emEnergy()+calt->hadEnergy();}
512  }
513  }
514  if ((m_cosmic) || (tmpHOCalib.caloen[0] <=10.0)) {
515 
516  GlobalPoint glbpt(posx, posy, posz);
517 
518  double mom = sqrt(momx*momx + momy*momy +momz*momz);
519 
520  momx /= mom;
521  momy /= mom;
522  momz /= mom;
523 
524  DirectionType trkdir(momx, momy, momz);
525 
526  tmpHOCalib.trkdr = (*ncosm).d0();
527  tmpHOCalib.trkdz = (*ncosm).dz();
528  tmpHOCalib.nmuon = (m_cosmic) ? cosmicmuon->size() : 1;
529  tmpHOCalib.trkvx = glbpt.x();
530  tmpHOCalib.trkvy = glbpt.y();
531  tmpHOCalib.trkvz = glbpt.z();
532  tmpHOCalib.trkmm = mom*charge;
533  tmpHOCalib.trkth = trkdir.theta();
534  tmpHOCalib.trkph = trkdir.phi();
535  tmpHOCalib.isect2 = -2;
536  tmpHOCalib.isect = -2;
537  tmpHOCalib.hodx = -100;
538  tmpHOCalib.hody = -100;
539  tmpHOCalib.hoang = -2.0;
540  tmpHOCalib.momatho = -2;
541  tmpHOCalib.ndof = (inearbymuon ==0) ? (int)(*ncosm).ndof() : -(int)(*ncosm).ndof();
542  tmpHOCalib.chisq = (*ncosm).normalizedChi2(); // max(1.,tmpHOCalib.ndof);
543  if (!m_cosmic) {
544  reco::MuonEnergy muonenr = muon1->calEnergy();
545  reco::MuonIsolation iso03 = muon1->isolationR03();
546  reco::MuonIsolation iso05 = muon1->isolationR05();
547 
548  tmpHOCalib.tkpt03 = iso03.sumPt;
549  tmpHOCalib.ecal03 = iso05.sumPt; // iso03.emEt+muonenr.em;
550  tmpHOCalib.hcal03 = iso03.hadEt+muonenr.had;
551  }
552  tmpHOCalib.therr = 0.;
553  tmpHOCalib.pherr = 0.;
554  if (iiner==1) {
555  reco::TrackBase::CovarianceMatrix innercov = (*ncosm).innerStateCovariance();
556  tmpHOCalib.therr = innercov(1,1); //thetaError();
557  tmpHOCalib.pherr = innercov(2,2); //phi0Error();
558  } else {
559  reco::TrackBase::CovarianceMatrix outercov = (*ncosm).outerStateCovariance();
560  tmpHOCalib.therr = outercov(1,1); //thetaError();
561  tmpHOCalib.pherr = outercov(2,2); //phi0Error();
562  }
563  edm::ESHandle<MagneticField> theMagField;
564  iSetup.get<IdealMagneticFieldRecord>().get(theMagField );
565 
566  SteppingHelixPropagator myHelix(&*theMagField,anyDirection);
567  myHelix.setMaterialMode(false);
568  myHelix.applyRadX0Correction(true);
569  double phiho = trkpos.phi();
570  if (phiho<0) phiho +=CLHEP::twopi;
571 
572  int iphisect_dt=int(6*(phiho+10.0*CLHEP::deg)/CLHEP::pi); //for u 18/12/06
573  if (iphisect_dt>=12) iphisect_dt=0;
574 
575  int iphisect = -1;
576  bool ipath = false;
577  for (int kl = 0; kl<=2; kl++) {
578 
579  int iphisecttmp = (kl<2) ? iphisect_dt + kl : iphisect_dt - 1;
580  if (iphisecttmp <0) iphisecttmp = 11;
581  if (iphisecttmp >=12) iphisecttmp = 0;
582 
583  double phipos = iphisecttmp*CLHEP::pi/6.;
584  double phirot = phipos;
585 
586  GlobalVector xLocal(-sin(phirot), cos(phirot), 0.);
587  GlobalVector yLocal(0., 0., 1.);
588  GlobalVector zLocal = xLocal.cross(yLocal).unit();
589  // GlobalVector zLocal(cos(phirot), sin(phirot), 0.0);
590 
591 
592  FreeTrajectoryState freetrajectorystate_ = getFreeTrajectoryState(*ncosm,&(*theMagField), iiner, samedir);
593 
594  Surface::RotationType rot(xLocal, yLocal, zLocal);
595 
596  for (int ik=1; ik>=0; ik--) { //propagate track in two HO layers
597 
598  double radial = rHOL1;
599  if (ik==0) radial = rHOL0;
600 
601  Surface::PositionType pos(radial*cos(phipos), radial*sin(phipos), 0.);
602  PlaneBuilder::ReturnType aPlane = PlaneBuilder().plane(pos,rot);
603 
604  auto aPlane2 = new Plane(pos,rot);
605 
606  SteppingHelixStateInfo steppingHelixstateinfo_;
607  myHelix.propagate(SteppingHelixStateInfo(freetrajectorystate_), (*aPlane2), steppingHelixstateinfo_);
608 
609  if (steppingHelixstateinfo_.isValid()) {
610 
611  GlobalPoint hotrkpos2xx(steppingHelixstateinfo_.position().x(), steppingHelixstateinfo_.position().y(), steppingHelixstateinfo_.position().z());
612 
613  if (ik==1) {
614  HcalDetId ClosestCell = (HcalDetId) gHO->getClosestCell(hotrkpos2xx);
615  int ixeta = ClosestCell.ieta();
616  int ixphi = ClosestCell.iphi();
617  tmpHOCalib.isect2 = 100*std::abs(ixeta+50)+std::abs(ixphi);
618  }
619 
620 
621  GlobalVector hotrkpos2(steppingHelixstateinfo_.position().x(), steppingHelixstateinfo_.position().y(), steppingHelixstateinfo_.position().z());
622  CLHEP::Hep3Vector hotrkdir2(steppingHelixstateinfo_.momentum().x(), steppingHelixstateinfo_.momentum().y(),steppingHelixstateinfo_.momentum().z());
623 
624  LocalVector lclvt0 = (*aPlane).toLocal(hotrkpos2);
625 
626  double xx = lclvt0.x();
627  double yy = lclvt0.y();
628 
629  if (ik ==1) {
630  if ((std::abs(yy) < 130 && xx >-64.7 && xx <138.2) //Ring-0
631  ||(std::abs(yy) > 130 && std::abs(yy) <700 && xx >-76.3 && xx <140.5)) { //Ring +-1,2
632  ipath = true; //Only look for tracks which as hits in layer 1
633  iphisect = iphisecttmp;
634  }
635  }
636 
637  if (iphisect != iphisecttmp) continue; //Look for ring-0 only when ring1 is accepted for that sector
638 
639  switch (ik)
640  {
641  case 0 :
642  xhor0 = xx; //lclvt0.x();
643  yhor0 = yy; //lclvt0.y();
644  break;
645  case 1 :
646  xhor1 = xx; //lclvt0.x();
647  yhor1 = yy; //lclvt0.y();
648  tmpHOCalib.momatho = hotrkdir2.mag();
649  tmpHOCalib.hoang = CLHEP::Hep3Vector(zLocal.x(),zLocal.y(),zLocal.z()).dot(hotrkdir2.unit());
650  break;
651  default : break;
652  }
653  } else {
654  break;
655  }
656  }
657  if (ipath) break;
658  }
659  if (ipath) { //If muon crossed HO laeyrs
660 
661  int ietaho = 50;
662  int iphiho = -1;
663 
664  for (int ij=0; ij<9; ij++) {tmpHOCalib.hosig[ij]=-100.0;}
665  for (int ij=0; ij<18; ij++) {tmpHOCalib.hocorsig[ij]=-100.0;}
666  for (int ij=0; ij<9; ij++) {tmpHOCalib.hbhesig[ij]=-100.0;}
667  tmpHOCalib.hocro = -100;
668  tmpHOCalib.htime = -1000;
669 
670  int isect = 0;
671 
672  findHOEtaPhi(iphisect, ietaho, iphiho);
673 
674  if (ietaho !=0 && iphiho !=0 && std::abs(iring)<=2) { //Muon passed through a tower
675  isect = 100*std::abs(ietaho+50)+std::abs(iphiho);
676  if (std::abs(ietaho) >=netabin || iphiho<0) isect *=-1; //Not extrapolated to any tower
677  if (std::abs(ietaho) >=netabin) isect -=1000000; //not matched with eta
678  if (iphiho<0) isect -=2000000; //not matched with phi
679  tmpHOCalib.isect = isect;
680 
681  tmpHOCalib.hodx = localxhor1;
682  tmpHOCalib.hody = localyhor1;
683 
684  if (iring==0) {
685  tmpHOCalib.hocorsig[8] = localxhor0;
686  tmpHOCalib.hocorsig[9] = localyhor0;
687  }
688 
689  int etamn=-4;
690  int etamx=4;
691  if (iring==1) {etamn=5; etamx = 10;}
692  if (iring==2) {etamn=11; etamx = 16;}
693  if (iring==-1){etamn=-10; etamx = -5;}
694  if (iring==-2){etamn=-16; etamx = -11;}
695 
696  int phimn = 1;
697  int phimx = 2;
698  if (iring ==0) {
699  phimx =2*int((iphiho+1)/2.);
700  phimn = phimx - 1;
701  } else {
702  phimn = 3*int((iphiho+1)/3.) - 1;
703  phimx = phimn + 2;
704  }
705 
706  if (phimn <1) phimn += nphimx;
707  if (phimx >72) phimx -= nphimx;
708 
709  if (m_hbinfo) {
710  for (int ij=0; ij<9; ij++) {tmpHOCalib.hbhesig[ij]=-100.0;}
711 
712  edm::Handle<HBHERecHitCollection> hbheht;// iEvent.getByType(hbheht);
713  iEvent.getByToken(tok_hbhe_,hbheht);
714 
715  if ((*hbheht).size()>0) {
716  if(!(*hbheht).size()) throw (int)(*hbheht).size();
717 
718  for (HBHERecHitCollection::const_iterator jk=(*hbheht).begin(); jk!=(*hbheht).end(); jk++){
719  HcalDetId id =(*jk).id();
720  int tmpeta= id.ieta();
721  int tmpphi= id.iphi();
722 
723  int deta = tmpeta-ietaho;
724  if (tmpeta<0 && ietaho>0) deta += 1;
725  if (tmpeta>0 && ietaho<0) deta -= 1;
726 
727  // if (tmpeta==-1 && ietaho== 1) deta = -1;
728  // if (tmpeta== 1 && ietaho==-1) deta = 1;
729 
730  int dphi = tmpphi-iphiho;
731  if (dphi>nphimx/2) { dphi -=nphimx;}
732  if (dphi<-nphimx/2) { dphi +=nphimx;}
733 
734  // if (phimn >phimx) {
735  // if (dphi==71) dphi=-1;
736  // if (dphi==-71) dphi=1;
737  // }
738 
739  if (m_occupancy) {
740  float signal = (*jk).energy();
741  // int tmpeta1 = (tmpeta>0) ? tmpeta -1 : -tmpeta +14;
742  if (signal >-100 && Noccu == Noccu_old) {
743  for (int ij=0; ij<5; ij++) {
744  if (signal >(ij+2)*m_sigma) {
745  ho_occupency[ij]->Fill(tmpeta, tmpphi);
746  }
747  }
748  }
749  }
750 
751  int ipass2 = (std::abs(deta) <=1 && std::abs(dphi)<=1) ? 1 : 0; //NEED correction in full CMS detector
752  if ( ipass2 ==0 ) continue;
753 
754  float signal = (*jk).energy();
755 
756  if (3*(deta+1)+dphi+1<9) tmpHOCalib.hbhesig[3*(deta+1)+dphi+1] = signal;
757  }
758  }
759  } //m_hbinfo #endif
760 
762  iEvent.getByToken(tok_ho_,hoht);
763 
764  if ((*hoht).size()>0) {
765  for (HORecHitCollection::const_iterator jk=(*hoht).begin(); jk!=(*hoht).end(); jk++){
766  HcalDetId id =(*jk).id();
767  int tmpeta= id.ieta();
768  int tmpphi= id.iphi();
769 
770  int ipass1 =0;
771  if (tmpeta >=etamn && tmpeta <=etamx) {
772  if (phimn < phimx) {
773  ipass1 = (tmpphi >=phimn && tmpphi <=phimx ) ? 1 : 0;
774  } else {
775  ipass1 = (tmpphi==71 || tmpphi ==72 || tmpphi==1) ? 1 : 0;
776  }
777  }
778 
779  int deta = tmpeta-ietaho;
780  int dphi = tmpphi -iphiho;
781 
782  if (tmpeta<0 && ietaho>0) deta += 1;
783  if (tmpeta>0 && ietaho<0) deta -= 1;
784  // if (tmpeta==-1 && ietaho== 1) deta = -1;
785  // if (tmpeta== 1 && ietaho==-1) deta = 1;
786 
787  if (dphi>nphimx/2) { dphi -=nphimx;}
788  if (dphi<-nphimx/2) { dphi +=nphimx;}
789  // if (phimn>phimx) {
790  // if (dphi==71) dphi=-1;
791  // if (dphi==-71) dphi=1;
792  // }
793 
794  float signal = (*jk).energy();
795 
796  int ipass2 = (std::abs(deta) <=1 && std::abs(dphi)<=1) ? 1 : 0;
797 
798  if (ipass1 ==0 && ipass2 ==0 ) continue;
799 
800  if (ipass1 ==1) {
801  int tmpdph = tmpphi-phimn;
802  if (tmpdph<0) tmpdph = 2; //only case of iphi==1, where phimn=71
803 
804  int ilog = 2*(tmpeta-etamn)+tmpdph;
805  if (iring !=0) {
806  if (iring >0) {
807  ilog = 3*(tmpeta-etamn)+tmpdph; //Again CMS correction
808  } else {
809  ilog = 3*(etamx-tmpeta)+tmpdph; //Again CMS correction
810  }
811  }
812  if (ilog>-1 && ilog<18) {
813  tmpHOCalib.hocorsig[ilog] = signal;
814  }
815  }
816 
817  if (ipass2 ==1) {
818 
819  if (3*(deta+1)+dphi+1<9) {
820  tmpHOCalib.hosig[3*(deta+1)+dphi+1] = signal; //Again CMS azimuthal near phi 1&72
821  }
822  }
823 
824  if (deta==0 && dphi ==0) {
825  tmpHOCalib.htime = (*jk).time();
826  tmpHOCalib.hoflag = (*jk).flags();
827 
828  // Get Channel Quality information for the given detID
829  unsigned theStatusValue = theHcalChStatus->getValues(id)->getValue();
830  // Now get severity of problems for the given detID, based on the rechit flag word and the channel quality status value
831  int hitSeverity=hcalSevLvlComputer->getSeverityLevel(id, (*jk).flags(),theStatusValue);
832  tmpHOCalib.hoflag = hitSeverity;
833  int crphi = tmpphi + 6;
834  if (crphi >72) crphi -=72;
835 
836  for (HORecHitCollection::const_iterator jcr=(*hoht).begin(); jcr!=(*hoht).end(); jcr++){
837  const HORecHit reccr = (const HORecHit)(*jcr);
838  HcalDetId idcr =reccr.id();
839  int etacr= idcr.ieta();
840  int phicr= idcr.iphi();
841  if (tmpeta==etacr && crphi ==phicr) {
842 
843  tmpHOCalib.hocro = reccr.energy();
844 
845  }
846  }
847  }
848  }
849  }
850  }
851 
852  //GMA Npass++;
853  if (Noccu == Noccu_old) Noccu++;
854  hostore->push_back(tmpHOCalib);
855  } // if (ipath)
856  } // Cut on calo energy
857 }
858 
859 void AlCaHOCalibProducer::findHOEtaPhi(int iphisect, int& ietaho, int& iphiho) {
860 
861  //18/12/06 : use only position, not angle phi
862 
863  const double etalow[16]={ 0.025, 35.195, 70.625, 106.595, 141.565, 180.765, 220.235, 261.385, 304.525, 349.975, 410.025, 452.085, 506.645, 565.025, 627.725, 660.25};
864  const double etahgh[16]={ 35.145, 70.575, 106.545, 125.505, 180.715, 220.185, 261.335, 304.475, 349.925, 392.575, 452.035, 506.595, 564.975, 627.675, 661.075, 700.25};
865 
866  const double philow[6]={-76.27, -35.11, 0.35, 35.81, 71.77, 108.93}; //Ring+/-1 & 2
867  const double phihgh[6]={-35.81, -0.35, 35.11, 71.07, 108.23, 140.49};
868 
869  const double philow00[6]={-60.27, -32.91, 0.35, 33.61, 67.37, 102.23}; //Ring0 L0
870  const double phihgh00[6]={-33.61, -0.35, 32.91, 66.67, 101.53, 129.49};
871 
872  const double philow01[6]={-64.67, -34.91, 0.35, 35.61, 71.37, 108.33}; //Ring0 L1
873  const double phihgh01[6]={-35.61, -0.35, 34.91, 70.67, 107.63, 138.19};
874 
875  iring = -10;
876 
877  double tmpdy = std::abs(yhor1);
878  for (int ij=0; ij<netabin; ij++) {
879  if (tmpdy >etalow[ij] && tmpdy <etahgh[ij]) {
880  ietaho = ij+1;
881  float tmp1 = fabs(tmpdy-etalow[ij]);
882  float tmp2 = fabs(tmpdy-etahgh[ij]);
883 
884  localyhor1 = (tmp1 < tmp2) ? -tmp1 : tmp2;
885  if (yhor1 <0) localyhor1 *=-1.;
886 
887  if (ij<4) iring =0;
888  if (ij>=4 && ij<10) iring=1;
889  if (ij>=10 && ij<netabin) iring=2;
890  break;
891  }
892  }
893 
894  int tmpphi = 0;
895  int tmpphi0 = 0;
896 
897  if (ietaho >4) { //Ring 1 and 2
898  for (int ij=0; ij<6; ij++) {
899  if (xhor1 >philow[ij] && xhor1 <phihgh[ij]) {
900  tmpphi=ij+1;
901  float tmp1 = fabs(xhor1-philow[ij]);
902  float tmp2 = fabs(xhor1-phihgh[ij]);
903  localxhor1 = (tmp1 < tmp2) ? -tmp1 : tmp2;
904  break;
905  }
906  }
907  } else { //Ring 0
908  for (int ij=0; ij<6; ij++) {
909  if (xhor1 >philow01[ij] && xhor1 <phihgh01[ij]) {
910  tmpphi=ij+1;
911  float tmp1 = fabs(xhor1-philow01[ij]);
912  float tmp2 = fabs(xhor1-phihgh01[ij]);
913  localxhor1 = (tmp1 < tmp2) ? -tmp1 : tmp2;
914  break;
915  }
916  }
917 
918  for (int ij=0; ij<6; ij++) {
919  if (xhor0 >philow00[ij] && xhor0 <phihgh00[ij]) {
920  tmpphi0=ij+1;
921  float tmp1 = fabs(xhor0-philow00[ij]);
922  float tmp2 = fabs(xhor0-phihgh00[ij]);
923  localxhor0 = (tmp1 < tmp2) ? -tmp1 : tmp2;
924  if (tmpphi !=tmpphi0) localxhor0 +=10000.;
925  break;
926  }
927  }
928 
929  double tmpdy = std::abs(yhor0);
930  for (int ij=0; ij<4; ij++) {
931  if (tmpdy >etalow[ij] && tmpdy <etahgh[ij]) {
932  float tmp1 = fabs(tmpdy-etalow[ij]);
933  float tmp2 = fabs(tmpdy-etahgh[ij]);
934  localyhor0 = (tmp1 < tmp2) ? -tmp1 : tmp2;
935  if (yhor0 <0) localyhor0 *=-1.;
936  if (ij+1 != ietaho) localyhor0 +=10000.;
937  break;
938  }
939  }
940  }
941 
942  if (tmpphi!=0) {
943  iphiho = 6*iphisect -2 + tmpphi;
944  if (iphiho <=0) iphiho +=nphimx;
945  if (iphiho >nphimx) iphiho -=nphimx;
946  }
947 
948  // isect2 = 15*iring+iphisect+1;
949 
950  if (yhor1 <0) {
951  if (std::abs(ietaho) >netabin) { //Initialised with 50
952  ietaho +=1;
953  } else {
954  ietaho *=-1;
955  }
956  // isect2 *=-1;
957  iring *=-1;
958  }
959 }
960 
962 {
963 
964  if (iiner ==0) {
965  GlobalPoint gpos( tk.outerX(), tk.outerY(), tk.outerZ());
966  GlobalVector gmom( tk.outerPx(), tk.outerPy(), tk.outerPz());
967  if (dir) gmom *=-1.;
968  GlobalTrajectoryParameters par( gpos, gmom, tk.charge(), field);
969  CurvilinearTrajectoryError err( tk.extra()->outerStateCovariance());
970  return FreeTrajectoryState( par, err);
971  } else {
972  GlobalPoint gpos( tk.innerPosition().X(), tk.innerPosition().Y(), tk.innerPosition().Z());
973  GlobalVector gmom( tk.innerMomentum().X(), tk.innerMomentum().Y(), tk.innerMomentum().Z());
974  if (dir) gmom *=-1.;
975  GlobalTrajectoryParameters par( gpos, -gmom, tk.charge(), field);
976  CurvilinearTrajectoryError err( tk.extra()->innerStateCovariance());
977  return FreeTrajectoryState( par, err);
978  }
979 
980 }
981 
983 
984 //define this as a plug-in
986 
987 
RunNumber_t run() const
Definition: EventID.h:39
float hadEt
hcal sum-Et
Definition: MuonIsolation.h:9
T getParameter(std::string const &) const
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:43
EventNumber_t event() const
Definition: EventID.h:41
Basic3DVector< float > DirectionType
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< CaloTowerCollection > tok_tower_
FreeTrajectoryState getFreeTrajectoryState(const reco::Track &tk, const MagneticField *field, int itag, bool dir)
float sumPt
sum-pt of tracks
Definition: MuonIsolation.h:7
double outerPy() const
y coordinate of momentum vector at the outermost hit position
Definition: Track.h:129
const HcalChannelQuality * theHcalChStatus
const TrackExtraRef & extra() const
reference to &quot;extra&quot; object
Definition: Track.h:189
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
tuple primaryVertices
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Theta< T > theta() const
ReturnType plane(const PositionType &pos, const RotationType &rot) const
Definition: PlaneBuilder.h:22
std::vector< CaloTower >::const_iterator const_iterator
T y() const
Definition: PV3DBase.h:63
virtual void produce(edm::Event &, const edm::EventSetup &) override
std::map< std::string, bool > fired
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
ErrorD< N >::type type
Definition: Error.h:33
edm::EDGetTokenT< reco::TrackCollection > tok_muonsCosmic_
const Item * getValues(DetId fId, bool throwOnFail=true) const
Geom::Phi< T > phi() const
virtual void endJob() override
GlobalVector momentum() const
void findHOEtaPhi(int iphsect, int &ietaho, int &iphiho)
bool isRealData() const
Definition: EventBase.h:63
Definition: Plane.h:17
unsigned int hoflag
double outerZ() const
z coordinate of the outermost hit position
Definition: Track.h:149
const Double_t pi
const math::XYZPoint & innerPosition() const
position of the innermost hit
Definition: Track.h:55
int iEvent
Definition: GenABIO.cc:230
void fillHOStore(const reco::TrackRef &ncosm, HOCalibVariables &tmpHOCalib, std::auto_ptr< HOCalibVariableCollection > &hostore, int Noccu_old, int indx, edm::Handle< reco::TrackCollection > cosmicmuon, edm::View< reco::Muon >::const_iterator muon1, const edm::Event &iEvent, const edm::EventSetup &iSetup)
void applyRadX0Correction(bool applyRadX0Correction)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
float energy() const
Definition: CaloRecHit.h:17
GlobalPoint position() const
T sqrt(T t)
Definition: SSEVec.h:18
Basic3DVector< float > RotationType
Vector3DBase< typename PreciseFloatType< T, U >::Type, FrameTag > cross(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:119
T z() const
Definition: PV3DBase.h:64
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
int ieta() const
get the cell ieta
Definition: HcalDetId.h:56
math::Error< 5 >::type CovarianceMatrix
virtual void beginRun(edm::Run const &, edm::EventSetup const &) override
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void propagate(const SteppingHelixStateInfo &ftsStart, const Surface &sDest, SteppingHelixStateInfo &out) const
Propagate to Plane given a starting point.
virtual void beginJob() override
Basic3DVector< float > PositionType
double outerX() const
x coordinate of the outermost hit position
Definition: Track.h:139
bool isValid() const
Definition: HandleBase.h:75
double outerPz() const
z coordinate of momentum vector at the outermost hit position
Definition: Track.h:134
virtual DetId getClosestCell(const GlobalPoint &r) const
edm::EDGetTokenT< LumiScalersCollection > tok_lumi_
int iphi() const
get the cell iphi
Definition: HcalDetId.cc:101
edm::EDGetTokenT< reco::VertexCollection > tok_vertex_
edm::EDGetTokenT< edm::View< reco::Muon > > tok_muons_
const T & get() const
Definition: EventSetup.h:56
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:20
T const * product() const
Definition: ESHandle.h:86
std::vector< HOCalibVariables > HOCalibVariableCollection
collection of HOcalibration variabale
int getSeverityLevel(const DetId &myid, const uint32_t &myflag, const uint32_t &mystatus) const
edm::EDGetTokenT< HORecHitCollection > tok_ho_
edm::EventID id() const
Definition: EventBase.h:59
AlCaHOCalibProducer(const edm::ParameterSet &)
T dot(const Basic3DVector &v) const
Scalar product, or &quot;dot&quot; product, with a vector of same type.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
const math::XYZVector & innerMomentum() const
momentum vector at the innermost hit position
Definition: Track.h:60
double outerY() const
y coordinate of the outermost hit position
Definition: Track.h:144
int charge() const
track electric charge
Definition: TrackBase.h:562
dbl *** dir
Definition: mlp_gen.cc:35
uint32_t getValue() const
T x() const
Definition: PV3DBase.h:62
edm::ESHandle< HcalSeverityLevelComputer > hcalSevLvlComputerHndl
void setMaterialMode(bool noMaterial)
Switch for material effects mode: no material effects if true.
HcalDetId id() const
get the id
Definition: HORecHit.h:19
double outerPx() const
x coordinate of momentum vector at the outermost hit position
Definition: Track.h:124
edm::EDGetTokenT< HBHERecHitCollection > tok_hbhe_
Definition: Run.h:43
math::Error< dimension >::type CovarianceMatrix
5 parameter covariance matrix
Definition: TrackBase.h:77
T angle(T x1, T y1, T z1, T x2, T y2, T z2)
Definition: angle.h:11