CMS 3D CMS Logo

PixelTrackExtractor.cc
Go to the documentation of this file.
1 #include "PixelTrackExtractor.h"
2 
5 #include "TrackSelector.h"
10 
14 
15 
16 using namespace edm;
17 using namespace std;
18 using namespace reco;
19 using namespace muonisolation;
21 
22 PixelTrackExtractor::PixelTrackExtractor( const ParameterSet& par, edm::ConsumesCollector && iC ) :
23  theTrackCollectionToken(iC.consumes<TrackCollection >(par.getParameter<edm::InputTag>("inputTrackCollection"))),
24  theDepositLabel(par.getUntrackedParameter<string>("DepositLabel")),
25  theDiff_r(par.getParameter<double>("Diff_r")),
26  theDiff_z(par.getParameter<double>("Diff_z")),
27  theDR_Max(par.getParameter<double>("DR_Max")),
28  theDR_Veto(par.getParameter<double>("DR_Veto")),
29  theBeamlineOption(par.getParameter<string>("BeamlineOption")),
30  theBeamSpotToken(iC.mayConsume<BeamSpot>(par.getParameter<edm::InputTag>("BeamSpotLabel"))),
31  theNHits_Min(par.getParameter<unsigned int>("NHits_Min")),
32  theChi2Ndof_Max(par.getParameter<double>("Chi2Ndof_Max")),
33  theChi2Prob_Min(par.getParameter<double>("Chi2Prob_Min")),
34  thePt_Min(par.getParameter<double>("Pt_Min")),
35  thePropagateTracksToRadius(par.getParameter<bool>("PropagateTracksToRadius")),
36  theReferenceRadius(par.getParameter<double>("ReferenceRadius")),
37  theVetoLeadingTrack(par.getParameter<bool>("VetoLeadingTrack")),
38  thePtVeto_Min(par.getParameter<double>("PtVeto_Min")),
39  theDR_VetoPt(par.getParameter<double>("DR_VetoPt"))
40 {
41 }
42 
44  const edm::EventSetup & evSetup, const reco::Track & track) const
45 {
46  Direction dir(track.eta(),track.phi());
47  return reco::IsoDeposit::Vetos(1,veto(dir));
48 }
49 
51 {
53  result.vetoDir = dir;
54  result.dR = theDR_Veto;
55  return result;
56 }
57 
60  return Direction(tk.eta(), tk.phi());
61  }
62 
63  // this should represent a cylinder in global frame at R=refRadius cm, roughly where mid-layer of pixels is
64  double psRadius = theReferenceRadius;
65  double tkDxy = tk.dxy();
66  double s2D = fabs(tk.dxy()) < psRadius ? sqrt(psRadius*psRadius - tkDxy*tkDxy) : 0;
67 
68  // the field we get from the caller is already in units of GeV/cm
69  double dPhi = -s2D*tk.charge()*bz/tk.pt();
70 
71  return Direction(tk.eta(), tk.phi()+dPhi);
72 }
73 
74 IsoDeposit PixelTrackExtractor::deposit(const Event & event, const EventSetup & eventSetup, const Track & muon) const
75 {
76  static const std::string metname = "MuonIsolation|PixelTrackExtractor";
77 
79  eventSetup.get<IdealMagneticFieldRecord>().get(bField);
80  double bz = bField->inInverseGeV(GlobalPoint(0.,0.,0.)).z();
81 
82  Direction muonDir(directionAtPresetRadius(muon, bz));
83  IsoDeposit deposit(muonDir );
85  deposit.setVeto( veto(muonDir) );
86 
87  deposit.addCandEnergy(muon.pt());
88 
90  event.getByToken(theTrackCollectionToken, tracksH);
91  // const TrackCollection tracks = *(tracksH.product());
92  LogTrace(metname)<<"***** TRACK COLLECTION SIZE: "<<tracksH->size();
93 
94  double vtx_z = muon.vz();
95  LogTrace(metname)<<"***** Muon vz: "<<vtx_z;
96  reco::TrackBase::Point beamPoint(0,0, 0);
97 
98  if (theBeamlineOption == "BeamSpotFromEvent"){
99  //pick beamSpot
101  edm::Handle<reco::BeamSpot> beamSpotH;
102 
103  event.getByToken(theBeamSpotToken,beamSpotH);
104 
105  if (beamSpotH.isValid()){
106  beamPoint = beamSpotH->position();
107  LogTrace(metname)<<"Extracted beam point at "<<beamPoint<<std::endl;
108  }
109  }
110 
111  LogTrace(metname)<<"Using beam point at "<<beamPoint<<std::endl;
112 
114  theDiff_r, muonDir, theDR_Max, beamPoint);
115 
116  pars.nHitsMin = theNHits_Min;
117  pars.chi2NdofMax = theChi2Ndof_Max;
118  pars.chi2ProbMin = theChi2Prob_Min;
119  pars.ptMin = thePt_Min;
120 
121  TrackSelector selection(pars);
122  TrackSelector::result_type sel_tracks = selection(*tracksH);
123  LogTrace(metname)<<"all tracks: "<<tracksH->size()<<" selected: "<<sel_tracks.size();
124 
125 
126  double maxPt = -1;
127  Direction maxPtDir;
128  TrackSelector::result_type::const_iterator tkI = sel_tracks.begin();
129  for (; tkI != sel_tracks.end(); ++tkI) {
130  const reco::Track* tk = *tkI;
131  LogTrace(metname) << "This track has: pt= " << tk->pt() << ", eta= "
132  << tk->eta() <<", phi= "<<tk->phi();
133  Direction dirTrk(directionAtPresetRadius(*tk, bz));
134  deposit.addDeposit(dirTrk, tk->pt());
135  double tkDr = (muonDir-dirTrk).deltaR;
136  double tkPt = tk->pt();
137  if (theVetoLeadingTrack && tkPt > thePtVeto_Min
138  && tkDr < theDR_VetoPt
139  && maxPt < tkPt ){
140  maxPt = tkPt;
141  maxPtDir = dirTrk;
142  }
143  }
144  if (maxPt > 0){
145  deposit.setVeto(veto(maxPtDir));
146  LogTrace(metname)<<" Set track veto the leading track with pt "
147  <<maxPt<<" in direction (eta,phi) "
148  <<maxPtDir.eta()<<", "<<maxPtDir.phi();
149  }
150 
151  return deposit;
152 }
double theReferenceRadius
If set to true will compare track eta-phi at ...
void setVeto(const Veto &aVeto)
Set veto.
Definition: IsoDeposit.h:80
double theChi2Ndof_Max
trk.numberOfValidHits >= theNHits_Min
reco::IsoDeposit::Veto veto(const reco::IsoDeposit::Direction &dir) const
const std::string metname
double theDR_Veto
Maximum cone angle for deposits.
std::list< const reco::Track * > result_type
Definition: TrackSelector.h:17
selection
main part
Definition: corrVsCorr.py:98
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:645
bool ev
virtual reco::IsoDeposit::Vetos vetos(const edm::Event &ev, const edm::EventSetup &evSetup, const reco::Track &track) const
void addDeposit(double dr, double deposit)
Add deposit (ie. transverse energy or pT)
Definition: IsoDeposit.cc:23
bool thePropagateTracksToRadius
min track pt to include into iso deposit
GlobalVector inInverseGeV(const GlobalPoint &gp) const
Field value ad specified global point, in 1/Gev.
Definition: MagneticField.h:41
reco::IsoDeposit deposit(const edm::Event &ev, const edm::EventSetup &evSetup, const reco::Track &muon) const override
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:651
vector< ParameterSet > Parameters
double thePtVeto_Min
will veto leading track if
T sqrt(T t)
Definition: SSEVec.h:18
edm::EDGetTokenT< reco::BeamSpot > theBeamSpotToken
"NONE", "BeamSpotFromEvent"
double pt() const
track transverse momentum
Definition: TrackBase.h:621
T z() const
Definition: PV3DBase.h:64
double theDR_Max
z distance to vertex
void addCandEnergy(double et)
Set energy or pT attached to cand trajectory.
Definition: IsoDeposit.h:139
edm::EDGetTokenT< reco::TrackCollection > theTrackCollectionToken
double thePt_Min
ChiSquaredProbability(trk.chi2,trk.ndof) > theChi2Prob_Min.
math::XYZPoint Point
point in the space
Definition: TrackBase.h:83
bool isValid() const
Definition: HandleBase.h:74
#define LogTrace(id)
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
double vz() const
z coordinate of the reference point on track
Definition: TrackBase.h:669
const T & get() const
Definition: EventSetup.h:58
double theDiff_z
transverse distance to vertex
std::string theBeamlineOption
Veto cone angle.
reco::isodeposit::Direction directionAtPresetRadius(const reco::Track &tk, double bz) const
std::vector< Veto > Vetos
Definition: IsoDeposit.h:63
fixed size matrix
HLT enums.
double theDR_VetoPt
.. it is above this threshold
unsigned int theNHits_Min
BeamSpot name.
muonisolation::Range< float > Range
Definition: TrackSelector.h:16
int charge() const
track electric charge
Definition: TrackBase.h:567
const Point & position() const
position
Definition: BeamSpot.h:62
double theChi2Prob_Min
trk.normalizedChi2 < theChi2Ndof_Max
dbl *** dir
Definition: mlp_gen.cc:35
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:591
Definition: event.py:1