CMS 3D CMS Logo

EgammaTrackExtractor.cc
Go to the documentation of this file.
2 
21 
22 using namespace edm;
23 using namespace std;
24 using namespace reco;
25 using namespace egammaisolation;
27 
28 EgammaTrackExtractor::EgammaTrackExtractor( const ParameterSet& par, edm::ConsumesCollector & iC ) :
29  theTrackCollectionToken(iC.consumes<View<Track> >(par.getParameter<edm::InputTag>("inputTrackCollection"))),
30  theDepositLabel(par.getUntrackedParameter<std::string>("DepositLabel")),
31  theDiff_r(par.getParameter<double>("Diff_r")),
32  theDiff_z(par.getParameter<double>("Diff_z")),
33  theDR_Max(par.getParameter<double>("DR_Max")),
34  theDR_Veto(par.getParameter<double>("DR_Veto")),
35  theBeamlineOption(par.getParameter<std::string>("BeamlineOption")),
36  theBeamSpotToken(iC.mayConsume<reco::BeamSpot>(par.getParameter<edm::InputTag>("BeamSpotLabel"))),
37  theNHits_Min(par.getParameter<unsigned int>("NHits_Min")),
38  theChi2Ndof_Max(par.getParameter<double>("Chi2Ndof_Max")),
39  theChi2Prob_Min(par.getParameter<double>("Chi2Prob_Min")),
40  thePt_Min(par.getParameter<double>("Pt_Min")),
41  dzOptionString(par.getParameter<std::string>("dzOption"))
42 {
43  if( ! dzOptionString.compare("dz") ) dzOption = EgammaTrackSelector::dz;
44  else if( ! dzOptionString.compare("vz") ) dzOption = EgammaTrackSelector::vz;
45  else if( ! dzOptionString.compare("bs") ) dzOption = EgammaTrackSelector::bs;
46  else if( ! dzOptionString.compare("vtx") )dzOption = EgammaTrackSelector::vtx;
48 
49 }
50 
52  const edm::EventSetup & evSetup, const reco::Track & track) const
53 {
54  reco::isodeposit::Direction dir(track.eta(),track.phi());
55  return reco::IsoDeposit::Vetos(1,veto(dir));
56 }
57 
59 {
61  result.vetoDir = dir;
62  result.dR = theDR_Veto;
63  return result;
64 }
65 
66 IsoDeposit EgammaTrackExtractor::deposit(const Event & event, const EventSetup & eventSetup, const Candidate & candTk) const
67 {
68  static std::string metname = "EgammaIsolationAlgos|EgammaTrackExtractor";
69 
71  double dzCut=0;
72 
73  reco::TrackBase::Point beamPoint(0,0, 0);
74  if (theBeamlineOption == "BeamSpotFromEvent"){
75  //pick beamSpot
78 
79  event.getByToken(theBeamSpotToken,beamSpotH);
80 
81  if (beamSpotH.isValid()){
82  beamPoint = beamSpotH->position();
83  }
84  }
85 
86  Handle<View<Track> > tracksH;
87  event.getByToken(theTrackCollectionToken, tracksH);
88 
89  if( candTk.isElectron() ){
90  const reco::GsfElectron* elec = dynamic_cast<const reco::GsfElectron*>(&candTk);
91  candDir = reco::isodeposit::Direction(elec->gsfTrack()->eta(), elec->gsfTrack()->phi());
92  } else {
93  candDir = reco::isodeposit::Direction(candTk.eta(), candTk.phi());
94  }
95 
96  IsoDeposit deposit(candDir );
97  deposit.setVeto( veto(candDir) );
98  deposit.addCandEnergy(candTk.et());
99 
100  View<Track>::const_iterator itrTr = tracksH->begin();
101  View<Track>::const_iterator trEnd = tracksH->end();
102  for (itrTr = tracksH->begin();itrTr != trEnd; ++itrTr) {
103 
104  if(candDir.deltaR( reco::isodeposit::Direction(itrTr->eta(),itrTr->phi()) ) > theDR_Max )
105  continue;
106 
107  if(itrTr->normalizedChi2() > theChi2Ndof_Max)
108  continue;
109 
110  if(itrTr->pt() < thePt_Min)
111  continue;
112 
113  if(theChi2Prob_Min > 0 && ChiSquaredProbability(itrTr->chi2(), itrTr->ndof()) < theChi2Prob_Min )
114  continue;
115 
116  if(theNHits_Min > 0 && itrTr->numberOfValidHits() < theNHits_Min)
117  continue;
118 
119  if( candTk.isElectron() ){
120  const reco::GsfElectron* elec = dynamic_cast<const reco::GsfElectron*>(&candTk);
121  switch(dzOption) {
122  case EgammaTrackSelector::dz : dzCut = elec->gsfTrack()->dz() - itrTr->dz() ; break;
123  case EgammaTrackSelector::vz : dzCut = elec->gsfTrack()->vz() - itrTr->vz() ; break;
124  case EgammaTrackSelector::bs : dzCut = elec->gsfTrack()->dz(beamPoint) - itrTr->dz(beamPoint) ; break;
125  case EgammaTrackSelector::vtx: dzCut = itrTr->dz(elec->gsfTrack()->vertex()) ; break;
126  default: dzCut = elec->gsfTrack()->vz() - itrTr->vz() ; break;
127  }
128  } else {
129  switch(dzOption) {
130  case EgammaTrackSelector::dz : dzCut = (*itrTr).dz() - candTk.vertex().z() ; break;
131  case EgammaTrackSelector::vz : dzCut = (*itrTr).vz() - candTk.vertex().z() ; break;
132  case EgammaTrackSelector::bs : dzCut = (*itrTr).dz(beamPoint) - candTk.vertex().z() ; break;
133  case EgammaTrackSelector::vtx: dzCut = (*itrTr).dz(candTk.vertex()); break;
134  default : dzCut = (*itrTr).vz() - candTk.vertex().z(); break;
135  }
136  }
137 
138  if(fabs(dzCut) > theDiff_z)
139  continue;
140 
141  if(fabs(itrTr->dxy(beamPoint) ) > theDiff_r)
142  continue;
143 
144  deposit.addDeposit(reco::isodeposit::Direction(itrTr->eta(), itrTr->phi()), itrTr->pt());
145 
146  }
147 
148  return deposit;
149 }
double theDR_Veto
Maximum cone angle for deposits.
edm::EDGetTokenT< edm::View< reco::Track > > theTrackCollectionToken
virtual reco::IsoDeposit deposit(const edm::Event &ev, const edm::EventSetup &evSetup, const reco::Track &muon) const
std::string dzOptionString
Endcap requirements to determine if isolated for selective filling.
const std::string metname
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:645
reco::IsoDeposit::Veto veto(const reco::IsoDeposit::Direction &dir) const
bool ev
double theChi2Prob_Min
trk.normalizedChi2 < theChi2Ndof_Max
double theDiff_r
minimum candidate et
virtual double et() const =0
transverse energy
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:651
double thePt_Min
ChiSquaredProbability(trk.chi2,trk.ndof) > theChi2Prob_Min.
float ChiSquaredProbability(double chiSquared, double nrDOF)
virtual bool isElectron() const =0
math::XYZPoint Point
point in the space
Definition: TrackBase.h:83
bool isValid() const
Definition: HandleBase.h:74
virtual double eta() const =0
momentum pseudorapidity
double theDiff_z
transverse distance to vertex
unsigned int theNHits_Min
BeamSpot name.
virtual reco::IsoDeposit::Vetos vetos(const edm::Event &ev, const edm::EventSetup &evSetup, const reco::Track &track) const
std::vector< Veto > Vetos
Definition: IsoDeposit.h:63
fixed size matrix
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
edm::EDGetTokenT< reco::BeamSpot > theBeamSpotToken
double theDR_Max
z distance to vertex
virtual const Point & vertex() const =0
vertex position
const Point & position() const
position
Definition: BeamSpot.h:62
dbl *** dir
Definition: mlp_gen.cc:35
double theChi2Ndof_Max
trk.numberOfValidHits >= theNHits_Min
virtual double phi() const =0
momentum azimuthal angle
std::string theBeamlineOption
Veto cone angle.
Definition: event.py:1
double deltaR(const Direction &dir2) const
virtual GsfTrackRef gsfTrack() const
reference to a GsfTrack
Definition: GsfElectron.h:185