CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
HcalCorrPFCalculation Class Reference
Inheritance diagram for HcalCorrPFCalculation:
edm::EDAnalyzer

Public Member Functions

virtual void analyze (edm::Event const &ev, edm::EventSetup const &c)
 
virtual void beginJob ()
 
virtual void endJob ()
 
 HcalCorrPFCalculation (edm::ParameterSet const &conf)
 
 ~HcalCorrPFCalculation ()
 
- Public Member Functions inherited from edm::EDAnalyzer
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 

Private Member Functions

double getDistInPlaneSimple (const GlobalPoint caloPoint, const GlobalPoint rechitPoint)
 
double RecalibFactor (HcalDetId id)
 

Private Attributes

Bool_t AddRecalib
 
bool Conecorr_
 
Bool_t doHF
 
double eEcalCone
 
double eHcalCone
 
double eHcalConeNoise
 
TH1F * enEcalB
 
TH1F * enEcalE
 
double energyECALmip
 
TProfile * enHcal
 
TProfile * enHcalNoise
 
edm::Service< TFileServicefs
 
const CaloGeometrygeo
 
Int_t iEta
 
Int_t iPhi
 
TProfile * nCells
 
TProfile * nCellsNoise
 
int nevtot
 
TrackAssociatorParameters parameters_
 
bool PFcorr_
 
const HcalPFCorrspfRecalib
 
TTree * pfTree
 
double radius_
 
bool Respcorr_
 
const HcalRespCorrsrespRecalib
 
TFile * rootFile
 
SteppingHelixPropagatorstepPropF
 
double taECALCone_
 
double taHCALCone_
 
MagneticFieldtheMagField
 
TrackDetectorAssociator trackAssociator_
 
int UsedCells
 
int UsedCellsNoise
 
Float_t xTrkEcal
 
Float_t xTrkHcal
 
Float_t yTrkEcal
 
Float_t yTrkHcal
 
Float_t zTrkEcal
 
Float_t zTrkHcal
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
typedef WorkerT< EDAnalyzerWorkerType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDAnalyzer
CurrentProcessingContext const * currentContext () const
 

Detailed Description

Definition at line 50 of file HcalCorrPFCalculation.cc.

Constructor & Destructor Documentation

HcalCorrPFCalculation::HcalCorrPFCalculation ( edm::ParameterSet const &  conf)

Definition at line 112 of file HcalCorrPFCalculation.cc.

References edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), and ExpressReco_HICollisions_FallBack::parameters.

112  {
113 
114  // outputFile_ = conf.getUntrackedParameter<std::string>("outputFile", "myfile.root");
115 
116  Respcorr_ = conf.getUntrackedParameter<bool>("RespcorrAdd", false);
117  PFcorr_ = conf.getUntrackedParameter<bool>("PFcorrAdd", false);
118  Conecorr_ = conf.getUntrackedParameter<bool>("ConeCorrAdd", true);
119  radius_ = conf.getUntrackedParameter<double>("ConeRadiusCm", 40.);
120  energyECALmip = conf.getParameter<double>("energyECALmip");
121 
122  edm::ParameterSet parameters = conf.getParameter<edm::ParameterSet>("TrackAssociatorParameters");
125 
126  taECALCone_=conf.getUntrackedParameter<double>("TrackAssociatorECALCone",0.5);
127  taHCALCone_=conf.getUntrackedParameter<double>("TrackAssociatorHCALCone",0.6);
128 
129 }
void useDefaultPropagator()
use the default propagator
TrackAssociatorParameters parameters_
tuple conf
Definition: dbtoconf.py:185
TrackDetectorAssociator trackAssociator_
void loadParameters(const edm::ParameterSet &)
HcalCorrPFCalculation::~HcalCorrPFCalculation ( )

Definition at line 183 of file HcalCorrPFCalculation.cc.

183  {
184 
185 
186 }

Member Function Documentation

void HcalCorrPFCalculation::analyze ( edm::Event const &  ev,
edm::EventSetup const &  c 
)
virtual

Implements edm::EDAnalyzer.

Definition at line 188 of file HcalCorrPFCalculation.cc.

References abs, edm::SortedCollection< T, SORT >::begin(), DeDxDiscriminatorTools::charge(), gather_cfg::cout, MaxHit_struct::depthhit, MaxHit_struct::dr, ExpressReco_HICollisions_FallBack::e, edm::SortedCollection< T, SORT >::end(), edm::EventSetup::get(), edm::Event::getByLabel(), edm::Event::getByType(), CaloSubdetectorGeometry::getClosestCell(), getDistInPlaneSimple(), DetId::Hcal, HcalBarrel, HcalEndcap, HcalForward, MaxHit_struct::hitenergy, HcalDetId::ieta(), MaxHit_struct::ietahitm, info, HcalDetId::iphi(), MaxHit_struct::iphihitm, edm::HandleBase::isValid(), npart, L1TEmulatorMonitor_cff::p, PV3DBase< T, PVType, FrameType >::phi(), pos, edm::Handle< T >::product(), edm::ESHandle< class >::product(), ExpressReco_HICollisions_FallBack::pt, mathSSE::sqrt(), TrackDetMatchInfo::trkGlobPosAtEcal, and TrackDetMatchInfo::trkGlobPosAtHcal.

188  {
189 
190  AddRecalib=kFALSE;
191 
192  try{
193 
194  edm::ESHandle <HcalRespCorrs> recalibCorrs;
195  c.get<HcalRespCorrsRcd>().get("recalibrate",recalibCorrs);
196  respRecalib = recalibCorrs.product();
197 
199  c.get<HcalPFCorrsRcd>().get("recalibrate",pfCorrs);
200  pfRecalib = pfCorrs.product();
201 
202  AddRecalib = kTRUE;;
203  // LogMessage("CalibConstants")<<" OK ";
204 
205  }catch(const cms::Exception & e) {
206  LogWarning("CalibConstants")<<" Not Found!! ";
207  }
208 
209 
210 
212  c.get<CaloGeometryRecord>().get(pG);
213  geo = pG.product();
214 
215  parameters_.useEcal = true;
216  parameters_.useHcal = true;
217  parameters_.useCalo = false;
218  parameters_.useMuon = false;
221 
222 
223  // double eta_bin[42]={0.,.087,.174,.261,.348,.435,.522,.609,.696,.783,
224  //.870,.957,1.044,1.131,1.218,1.305,1.392,1.479,1.566,1.653,1.740,1.830,1.930,2.043,2.172,
225  //2.322,2.500,2.650,2.853,3.000,3.139,3.314,3.489,3.664,3.839,4.013,4.191,4.363,4.538,4.716,4.889,5.191};
226 
227  // MC info
228  double phi_MC = -999999.; // phi of initial particle from HepMC
229  double eta_MC = -999999.; // eta of initial particle from HepMC
230  double mom_MC = 50.; // P of initial particle from HepMC
231  bool MC = false;
232 
233  // MC information
234 
235 
237  // ev.getByLabel("VtxSmeared",evtMC);
238  ev.getByLabel("generator",evtMC);
239  if (!evtMC.isValid())
240  {
241  std::cout << "no HepMCProduct found" << std::endl;
242  }
243  else
244  {
245  MC=true;
246  // std::cout << "*** source HepMCProduct found"<< std::endl;
247  }
248 
249  // MC particle with highest pt is taken as a direction reference
250  double maxPt = -99999.;
251  int npart = 0;
252 
253  GlobalPoint pos (0,0,0);
254 
255  HepMC::GenEvent * myGenEvent = new HepMC::GenEvent(*(evtMC->GetEvent()));
256  for ( HepMC::GenEvent::particle_iterator p = myGenEvent->particles_begin();
257  p != myGenEvent->particles_end(); ++p )
258  {
259  double phip = (*p)->momentum().phi();
260  double etap = (*p)->momentum().eta();
261  double pt = (*p)->momentum().perp();
262  mom_MC = (*p)->momentum().rho();
263  if(pt > maxPt) { npart++; maxPt = pt; phi_MC = phip; eta_MC = etap; }
264  GlobalVector mom ((*p)->momentum().x(),(*p)->momentum().y(),(*p)->momentum().z());
265  int charge = -1;
266 
267  if(abs((*p)->pdg_id())==211) charge = (*p)->pdg_id()/abs((*p)->pdg_id()); // pions only !!!
268  else continue;
269 
270  const FreeTrajectoryState *freetrajectorystate_ =
271  new FreeTrajectoryState(pos, mom ,charge , &(*theMagField));
272 
273  TrackDetMatchInfo info = trackAssociator_.associate(ev, c, *freetrajectorystate_ , parameters_);
274  // TrackDetMatchInfo info = trackAssociator_.associate(iEvent, iSetup,trackAssociator_.getFreeTrajectoryState(iSetup, *trit), parameters_);
275 
276  float etahcal=info.trkGlobPosAtHcal.eta();
277  float phihcal=info.trkGlobPosAtHcal.phi();
278 
279  float etaecal=info.trkGlobPosAtEcal.eta();
280  // float phiecal=info.trkGlobPosAtEcal.phi();
281 
282 
283  xTrkEcal=info.trkGlobPosAtEcal.x();
284  yTrkEcal=info.trkGlobPosAtEcal.y();
285  zTrkEcal=info.trkGlobPosAtEcal.z();
286 
287  xTrkHcal=info.trkGlobPosAtHcal.x();
288  yTrkHcal=info.trkGlobPosAtHcal.y();
289  zTrkHcal=info.trkGlobPosAtHcal.z();
290 
292 
294 
295  if (etahcal>2.6) doHF = kTRUE;
296 
297 
298 
300  ev.getByType(hbhe);
301  const HBHERecHitCollection Hithbhe = *(hbhe.product());
302 
304  ev.getByType(hfcoll);
305  const HFRecHitCollection Hithf = *(hfcoll.product());
306 
307 
309  ev.getByType(hocoll);
310  const HORecHitCollection Hitho = *(hocoll.product());
311 
312 
314  ev.getByLabel("ecalRecHit","EcalRecHitsEE",ecalEE);
315  const EERecHitCollection HitecalEE = *(ecalEE.product());
316 
318  ev.getByLabel("ecalRecHit","EcalRecHitsEB",ecalEB);
319  const EBRecHitCollection HitecalEB = *(ecalEB.product());
320 
321 
322 
323  // energy in ECAL
324  eEcalCone = 0.;
325  // int numrechitsEcal = 0;
326 
327  //Hcal:
328  eHcalCone = 0.;
329  eHcalConeNoise = 0.;
330  UsedCells = 0;
331  UsedCellsNoise = 0;
332 
333 
334  Int_t iphitrue = -10;
335  Int_t ietatrue = 100;
336 
337  if (etahcal<1.392)
338  {
340  // const GlobalPoint tempPoint(newx, newy, newz);
341  //const DetId tempId = gHB->getClosestCell(tempPoint);
342  const HcalDetId tempId = gHB->getClosestCell(gPointHcal);
343  ietatrue = tempId.ieta();
344  iphitrue = tempId.iphi();
345  }
346 
347  if (etahcal>1.392 && etahcal<3.0)
348  {
350  const HcalDetId tempId = gHE->getClosestCell(gPointHcal);
351  ietatrue = tempId.ieta();
352  iphitrue = tempId.iphi();
353  }
354 
355  if (etahcal>3.0 && etahcal<5.0)
356  {
358  const HcalDetId tempId = gHF->getClosestCell(gPointHcal);
359  ietatrue = tempId.ieta();
360  iphitrue = tempId.iphi();
361  }
362 
363  //Calculate Ecal energy:
364  for (EBRecHitCollection::const_iterator ehit=HitecalEB.begin(); ehit!=HitecalEB.end(); ehit++)
365  {
366 
367  GlobalPoint pos = geo->getPosition(ehit->detid());
368  float dr = getDistInPlaneSimple(gPointEcal,pos);
369  if (dr < 10.) eEcalCone += ehit->energy();
370  }
371 
372  for (EERecHitCollection::const_iterator ehit=HitecalEE.begin(); ehit!=HitecalEE.end(); ehit++)
373  {
374 
375  GlobalPoint pos = geo->getPosition(ehit->detid());
376  float dr = getDistInPlaneSimple(gPointEcal,pos);
377  if (dr < 10.) eEcalCone += ehit->energy();
378  }
379  if(abs(etaecal)<1.5) enEcalB -> Fill(eEcalCone);
380  if(abs(etaecal)>1.5 && abs(etaecal)<3.1) enEcalE -> Fill(eEcalCone);
381 
382 
383  MaxHit_struct MaxHit;
384 
385  MaxHit.hitenergy=-100.;
386 
387 
388  Float_t recal = 1.0;
389 
390 
391  for (HBHERecHitCollection::const_iterator hhit=Hithbhe.begin(); hhit!=Hithbhe.end(); hhit++)
392  {
393 
394  recal = RecalibFactor(hhit->detid());
395  //cout<<"recal: "<<recal<<endl;
396 
397  GlobalPoint pos = geo->getPosition(hhit->detid());
398  float phihit = pos.phi();
399  float etahit = pos.eta();
400 
401  int iphihit = (hhit->id()).iphi();
402  int ietahit = (hhit->id()).ieta();
403  int depthhit = (hhit->id()).depth();
404  float enehit = hhit->energy()* recal;
405 
406  //Set noise RecHit opposite to track hits
407  int iphihitNoise = iphihit >36 ? iphihit-36 : iphihit+36;
408  int ietahitNoise = ietahit;
409  int depthhitNoise = depthhit;
410 
411  double dphi = fabs(phihcal - phihit);
412  if(dphi > 4.*atan(1.)) dphi = 8.*atan(1.) - dphi;
413  double deta = fabs(etahcal - etahit);
414  double dr = sqrt(dphi*dphi + deta*deta);
415 
416  //dr = getDistInPlaneSimple(gPointHcal,pos);
417 
418  if(dr<0.5)
419  {
420 
421  for (HBHERecHitCollection::const_iterator hhit2=Hithbhe.begin(); hhit2!=Hithbhe.end(); hhit2++)
422  {
423  int iphihit2 = (hhit2->id()).iphi();
424  int ietahit2 = (hhit2->id()).ieta();
425  int depthhit2 = (hhit2->id()).depth();
426  float enehit2 = hhit2->energy() * recal;
427 
428  if (iphihit==iphihit2 && ietahit==ietahit2 && depthhit!=depthhit2) enehit = enehit+enehit2;
429 
430  }
431 
432  //Find a Hit with Maximum Energy
433 
434  if(enehit > MaxHit.hitenergy)
435  {
436  MaxHit.hitenergy = enehit;
437  MaxHit.ietahitm = (hhit->id()).ieta();
438  MaxHit.iphihitm = (hhit->id()).iphi();
439  MaxHit.dr = dr;
440  //MaxHit.depthhit = (hhit->id()).depth();
441  MaxHit.depthhit = 1;
442  }
443  }
444 
445  if(dr<radius_ && enehit>0.01)
446  {
447  eHcalCone += enehit;
448  UsedCells++;
449 
450  // cout<<"track: ieta "<<ietahit<<" iphi: "<<iphihit<<" depth: "<<depthhit<<" energydepos: "<<enehit<<endl;
451 
452  for (HBHERecHitCollection::const_iterator hhit2=Hithbhe.begin(); hhit2!=Hithbhe.end(); hhit2++)
453  {
454  recal = RecalibFactor(hhit2->detid());
455  int iphihit2 = (hhit2->id()).iphi();
456  int ietahit2 = (hhit2->id()).ieta();
457  int depthhit2 = (hhit2->id()).depth();
458  float enehit2 = hhit2->energy()* recal;
459 
460  if (iphihitNoise == iphihit2 && ietahitNoise == ietahit2 && depthhitNoise == depthhit2 && enehit2>0.01)
461  {
462 
463  eHcalConeNoise += hhit2->energy()*recal;
464  UsedCellsNoise++;
465  //cout<<"Noise: ieta "<<ietahit2<<" iphi: "<<iphihit2<<" depth: "<<depthhit2<<" energydepos: "<<enehit2<<endl;
466  }
467  }
468  }
469  } //end of all HBHE hits cycle
470 
471  if(doHF){
472  for (HFRecHitCollection::const_iterator hhit=Hithf.begin(); hhit!=Hithf.end(); hhit++)
473  {
474 
475  recal = RecalibFactor(hhit->detid());
476 
477  GlobalPoint pos = geo->getPosition(hhit->detid());
478  float phihit = pos.phi();
479  float etahit = pos.eta();
480 
481  int iphihit = (hhit->id()).iphi();
482  int ietahit = (hhit->id()).ieta();
483  int depthhit = (hhit->id()).depth();
484  float enehit = hhit->energy()* recal;
485 
486  //Set noise RecHit opposite to track hits
487  int iphihitNoise = iphihit >36 ? iphihit-36 : iphihit+36;
488  int ietahitNoise = ietahit;
489  int depthhitNoise = depthhit;
490 
491 
492  double dphi = fabs(phihcal - phihit);
493  if(dphi > 4.*atan(1.)) dphi = 8.*atan(1.) - dphi;
494  double deta = fabs(etahcal - etahit);
495  double dr = sqrt(dphi*dphi + deta*deta);
496 
497  dr = getDistInPlaneSimple(gPointHcal,pos);
498 
499 
500  if(dr<60.)
501  {
502  //Find a Hit with Maximum Energy
503 
504  if(enehit > MaxHit.hitenergy)
505  {
506  MaxHit.hitenergy = enehit;
507  MaxHit.ietahitm = (hhit->id()).ieta();
508  MaxHit.iphihitm = (hhit->id()).iphi();
509  MaxHit.dr = dr;
510  MaxHit.depthhit = 1;
511  }
512  }
513 
514  if(dr<radius_ && enehit>0.01)
515  {
516 
517  eHcalCone += enehit;
518  UsedCells++;
519 
520  for (HFRecHitCollection::const_iterator hhit2=Hithf.begin(); hhit2!=Hithf.end(); hhit2++)
521  {
522  recal = RecalibFactor(hhit2->detid());
523 
524  int iphihit2 = (hhit2->id()).iphi();
525  int ietahit2 = (hhit2->id()).ieta();
526  int depthhit2 = (hhit2->id()).depth();
527  float enehit2 = hhit2->energy()* recal;
528 
529  if (iphihitNoise == iphihit2 && ietahitNoise == ietahit2 && depthhitNoise == depthhit2 && enehit2>0.01)
530  {
531  eHcalConeNoise += hhit2->energy()*recal;
532  UsedCellsNoise++;
533  }
534  }
535 
536  }
537  } //end of all HF hits cycle
538  } //end of doHF
539 
540  int dieta_M_P = 100;
541  int diphi_M_P = 100;
542  if(MaxHit.ietahitm*ietatrue>0) {dieta_M_P = abs (MaxHit.ietahitm-ietatrue);}
543  if(MaxHit.ietahitm*ietatrue<0) {dieta_M_P = abs(MaxHit.ietahitm-ietatrue)-1;}
544  diphi_M_P = abs(MaxHit.iphihitm-iphitrue);
545  diphi_M_P = diphi_M_P>36 ? 72-diphi_M_P : diphi_M_P;
546 
547  float iDr = sqrt(diphi_M_P*diphi_M_P+dieta_M_P*dieta_M_P);
548 
549 
550  Bool_t passCuts = kFALSE;
551  //passCuts=kTRUE;
552  if(eEcalCone < energyECALmip && iDr<2.) passCuts = kTRUE;
553 
554  if(passCuts)
555  {
556  enHcal -> Fill(ietatrue, eHcalCone);
557  nCells -> Fill(ietatrue, UsedCells);
558  enHcalNoise -> Fill(ietatrue, eHcalConeNoise);
559  nCellsNoise -> Fill(ietatrue, UsedCellsNoise);
560 
561  iEta = ietatrue;
562  iPhi = iphitrue;
563 
564  pfTree->Fill();
565  }
566  }
567 }
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:43
const HcalPFCorrs * pfRecalib
Geom::Phi< T > phi() const
Definition: PV3DBase.h:63
TrackAssociatorParameters parameters_
std::vector< T >::const_iterator const_iterator
double getDistInPlaneSimple(const GlobalPoint caloPoint, const GlobalPoint rechitPoint)
#define abs(x)
Definition: mlp_lapack.h:159
double npart
Definition: HydjetWrapper.h:45
double charge(const std::vector< uint8_t > &Ampls)
math::XYZPoint trkGlobPosAtHcal
T sqrt(T t)
Definition: SSEVec.h:28
int ieta() const
get the cell ieta
Definition: HcalDetId.h:38
const GlobalPoint & getPosition(const DetId &id) const
Get the position of a given detector id.
Definition: CaloGeometry.cc:68
const CaloGeometry * geo
bool isValid() const
Definition: HandleBase.h:76
const HcalRespCorrs * respRecalib
const_iterator end() const
virtual DetId getClosestCell(const GlobalPoint &r) const
int iphi() const
get the cell iphi
Definition: HcalDetId.h:40
T const * product() const
Definition: ESHandle.h:62
T const * product() const
Definition: Handle.h:74
TrackDetectorAssociator trackAssociator_
TrackDetMatchInfo associate(const edm::Event &, const edm::EventSetup &, const FreeTrajectoryState &, const AssociatorParameters &)
double RecalibFactor(HcalDetId id)
math::XYZPoint trkGlobPosAtEcal
Track position at different parts of the calorimeter.
tuple cout
Definition: gather_cfg.py:41
const_iterator begin() const
void HcalCorrPFCalculation::beginJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 569 of file HcalCorrPFCalculation.cc.

569  {
570 
571  // TProfile *nCells, *nCellsNoise, *en, *enNoise;
572  //TFile *rootFile;
573 
574  //rootFile = new TFile(outputFile_.c_str(),"RECREATE");
575 
576 
577  nCells = fs->make<TProfile>("nCells", "nCells", 83, -41.5, 41.5);
578  nCellsNoise = fs->make<TProfile>("nCellsNoise", "nCellsNoise", 83, -41.5, 41.5);
579 
580  enHcal = fs->make<TProfile>("enHcal", "enHcal", 83, -41.5, 41.5);
581  enHcalNoise = fs->make<TProfile>("enHcalNoise", "enHcalNoise", 83, -41.5, 41.5);
582 
583  enEcalB = fs->make<TH1F>("enEcalB", "enEcalB", 500, -5,50);
584  enEcalE = fs->make<TH1F>("enEcalE", "enEcalE", 500, -5,50);
585 
586  pfTree = new TTree("pfTree", "Tree for pf info");
587 
588  pfTree->Branch("eEcalCone", &eEcalCone, "eEcalCone/F");
589  pfTree->Branch("eHcalCone", &eHcalCone, "eHcalCone/F");
590  pfTree->Branch("eHcalConeNoise", &eHcalConeNoise, "eHcalConeNoise/F");
591 
592  pfTree->Branch("UsedCellsNoise", &UsedCellsNoise, "UsedCellsNoise/I");
593  pfTree->Branch("UsedCells", &UsedCells, "UsedCells/I");
594 
595 
596  // pfTree->Branch("etaTrack", &etaTrack, "etaTrack/F");
597  //pfTree->Branch("phiTrack", &phiTrack, "phiTrack/F");
598 
599  pfTree->Branch("iEta", &iEta, "iEta/I");
600  pfTree->Branch("iPhi", &iPhi, "iPhi/I");
601 
602 
603 }
T * make() const
make new ROOT object
edm::Service< TFileService > fs
void HcalCorrPFCalculation::endJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 604 of file HcalCorrPFCalculation.cc.

605 {
606 
607  /*
608  nCells -> Write();
609  nCellsNoise -> Write();
610  enHcal -> Write();
611  enHcalNoise -> Write();
612 
613  enEcalB -> Write();
614  enEcalE -> Write();
615 
616  rootFile->Close();
617  */
618 }
double HcalCorrPFCalculation::getDistInPlaneSimple ( const GlobalPoint  caloPoint,
const GlobalPoint  rechitPoint 
)
private

Definition at line 131 of file HcalCorrPFCalculation.cc.

References Vector3DBase< T, FrameTag >::dot(), PV3DBase< T, PVType, FrameType >::mag(), Vector3DBase< T, FrameTag >::unit(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

133 {
134 
135  // Simplified version of getDistInPlane
136  // Assume track direction is origin -> point of hcal intersection
137 
138  const GlobalVector caloIntersectVector(caloPoint.x(),
139  caloPoint.y(),
140  caloPoint.z());
141 
142  const GlobalVector caloIntersectUnitVector = caloIntersectVector.unit();
143 
144  const GlobalVector rechitVector(rechitPoint.x(),
145  rechitPoint.y(),
146  rechitPoint.z());
147 
148  const GlobalVector rechitUnitVector = rechitVector.unit();
149  double dotprod = caloIntersectUnitVector.dot(rechitUnitVector);
150  double rechitdist = caloIntersectVector.mag()/dotprod;
151 
152 
153  const GlobalVector effectiveRechitVector = rechitdist*rechitUnitVector;
154  const GlobalPoint effectiveRechitPoint(effectiveRechitVector.x(),
155  effectiveRechitVector.y(),
156  effectiveRechitVector.z());
157 
158 
159  GlobalVector distance_vector = effectiveRechitPoint-caloPoint;
160 
161  if (dotprod > 0.)
162  {
163  return distance_vector.mag();
164  }
165  else
166  {
167  return 999999.;
168  }
169 }
T y() const
Definition: PV3DBase.h:57
PreciseFloatType< T, U >::Type dot(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:107
T mag() const
Definition: PV3DBase.h:61
T z() const
Definition: PV3DBase.h:58
Vector3DBase unit() const
Definition: Vector3DBase.h:57
T x() const
Definition: PV3DBase.h:56
double HcalCorrPFCalculation::RecalibFactor ( HcalDetId  id)
private

Definition at line 171 of file HcalCorrPFCalculation.cc.

172 {
173  Float_t resprecal = 1.;
174  Float_t pfrecal = 1.;
175  if(AddRecalib) {
176  if(Respcorr_) resprecal = respRecalib -> getValues(id)->getValue();
177  if(PFcorr_) pfrecal = pfRecalib -> getValues(id)->getValue();
178  }
179  Float_t factor = resprecal*pfrecal;
180  return factor;
181 }
const HcalPFCorrs * pfRecalib
const HcalRespCorrs * respRecalib

Member Data Documentation

Bool_t HcalCorrPFCalculation::AddRecalib
private

Definition at line 70 of file HcalCorrPFCalculation.cc.

bool HcalCorrPFCalculation::Conecorr_
private

Definition at line 64 of file HcalCorrPFCalculation.cc.

Bool_t HcalCorrPFCalculation::doHF
private

Definition at line 69 of file HcalCorrPFCalculation.cc.

double HcalCorrPFCalculation::eEcalCone
private

Definition at line 101 of file HcalCorrPFCalculation.cc.

double HcalCorrPFCalculation::eHcalCone
private

Definition at line 101 of file HcalCorrPFCalculation.cc.

double HcalCorrPFCalculation::eHcalConeNoise
private

Definition at line 101 of file HcalCorrPFCalculation.cc.

TH1F* HcalCorrPFCalculation::enEcalB
private

Definition at line 82 of file HcalCorrPFCalculation.cc.

TH1F * HcalCorrPFCalculation::enEcalE
private

Definition at line 82 of file HcalCorrPFCalculation.cc.

double HcalCorrPFCalculation::energyECALmip
private

Definition at line 67 of file HcalCorrPFCalculation.cc.

TProfile * HcalCorrPFCalculation::enHcal
private

Definition at line 81 of file HcalCorrPFCalculation.cc.

TProfile * HcalCorrPFCalculation::enHcalNoise
private

Definition at line 81 of file HcalCorrPFCalculation.cc.

edm::Service<TFileService> HcalCorrPFCalculation::fs
private

Definition at line 79 of file HcalCorrPFCalculation.cc.

const CaloGeometry* HcalCorrPFCalculation::geo
private

Definition at line 91 of file HcalCorrPFCalculation.cc.

Int_t HcalCorrPFCalculation::iEta
private

Definition at line 107 of file HcalCorrPFCalculation.cc.

Int_t HcalCorrPFCalculation::iPhi
private

Definition at line 107 of file HcalCorrPFCalculation.cc.

TProfile* HcalCorrPFCalculation::nCells
private

Definition at line 81 of file HcalCorrPFCalculation.cc.

TProfile * HcalCorrPFCalculation::nCellsNoise
private

Definition at line 81 of file HcalCorrPFCalculation.cc.

int HcalCorrPFCalculation::nevtot
private

Definition at line 71 of file HcalCorrPFCalculation.cc.

TrackAssociatorParameters HcalCorrPFCalculation::parameters_
private
bool HcalCorrPFCalculation::PFcorr_
private

Definition at line 63 of file HcalCorrPFCalculation.cc.

const HcalPFCorrs* HcalCorrPFCalculation::pfRecalib
private

Definition at line 74 of file HcalCorrPFCalculation.cc.

TTree* HcalCorrPFCalculation::pfTree
private

Definition at line 83 of file HcalCorrPFCalculation.cc.

double HcalCorrPFCalculation::radius_
private

Definition at line 65 of file HcalCorrPFCalculation.cc.

bool HcalCorrPFCalculation::Respcorr_
private

Definition at line 62 of file HcalCorrPFCalculation.cc.

const HcalRespCorrs* HcalCorrPFCalculation::respRecalib
private

Definition at line 73 of file HcalCorrPFCalculation.cc.

TFile* HcalCorrPFCalculation::rootFile
private

Definition at line 84 of file HcalCorrPFCalculation.cc.

SteppingHelixPropagator* HcalCorrPFCalculation::stepPropF
private

Definition at line 76 of file HcalCorrPFCalculation.cc.

double HcalCorrPFCalculation::taECALCone_
private

Definition at line 88 of file HcalCorrPFCalculation.cc.

double HcalCorrPFCalculation::taHCALCone_
private

Definition at line 89 of file HcalCorrPFCalculation.cc.

MagneticField* HcalCorrPFCalculation::theMagField
private

Definition at line 77 of file HcalCorrPFCalculation.cc.

TrackDetectorAssociator HcalCorrPFCalculation::trackAssociator_
private

Definition at line 86 of file HcalCorrPFCalculation.cc.

int HcalCorrPFCalculation::UsedCells
private

Definition at line 104 of file HcalCorrPFCalculation.cc.

int HcalCorrPFCalculation::UsedCellsNoise
private

Definition at line 104 of file HcalCorrPFCalculation.cc.

Float_t HcalCorrPFCalculation::xTrkEcal
private

Definition at line 93 of file HcalCorrPFCalculation.cc.

Float_t HcalCorrPFCalculation::xTrkHcal
private

Definition at line 97 of file HcalCorrPFCalculation.cc.

Float_t HcalCorrPFCalculation::yTrkEcal
private

Definition at line 94 of file HcalCorrPFCalculation.cc.

Float_t HcalCorrPFCalculation::yTrkHcal
private

Definition at line 98 of file HcalCorrPFCalculation.cc.

Float_t HcalCorrPFCalculation::zTrkEcal
private

Definition at line 95 of file HcalCorrPFCalculation.cc.

Float_t HcalCorrPFCalculation::zTrkHcal
private

Definition at line 99 of file HcalCorrPFCalculation.cc.