Skip to content

Commit

Permalink
Merge pull request cms-sw#5665 from jbrinson/mydev
Browse files Browse the repository at this point in the history
RecoTracker/DeDx -- Added configurable parameters for nValidHits, nMissIn, nMissMid, maxRelTrkIsoDeltaRp3, relTrkIsoDeltaRSize
  • Loading branch information
nclopezo committed Oct 8, 2014
2 parents 5ed993e + 3455dc5 commit a4a5e7d
Show file tree
Hide file tree
Showing 2 changed files with 45 additions and 4 deletions.
44 changes: 40 additions & 4 deletions RecoTracker/DeDx/plugins/HLTDeDxFilter.cc
Original file line number Diff line number Diff line change
Expand Up @@ -21,6 +21,7 @@
#include "FWCore/Framework/interface/EventSetup.h"
#include "DataFormats/TrackReco/interface/DeDxData.h"
//#include "DataFormats/Math/interface/deltaPhi.h"
#include "DataFormats/Math/interface/deltaR.h"
#include "DataFormats/Common/interface/ValueMap.h"
#include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
#include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
Expand All @@ -38,6 +39,11 @@ HLTDeDxFilter::HLTDeDxFilter(const edm::ParameterSet& iConfig) : HLTFilter(iConf
minPT_ = iConfig.getParameter<double> ("minPT");
minNOM_ = iConfig.getParameter<double> ("minNOM");
maxETA_ = iConfig.getParameter<double> ("maxETA");
minNumValidHits_ = iConfig.getParameter<double> ("minNumValidHits");
maxNHitMissIn_ = iConfig.getParameter<double> ("maxNHitMissIn");
maxNHitMissMid_ = iConfig.getParameter<double> ("maxNHitMissMid");
maxRelTrkIsoDeltaRp3_ = iConfig.getParameter<double> ("maxRelTrkIsoDeltaRp3");
relTrkIsoDeltaRSize_ = iConfig.getParameter<double> ("relTrkIsoDeltaRSize");
inputTracksTag_ = iConfig.getParameter< edm::InputTag > ("inputTracksTag");
inputdedxTag_ = iConfig.getParameter< edm::InputTag > ("inputDeDxTag");
inputTracksToken_ = consumes<reco::TrackCollection>(iConfig.getParameter< edm::InputTag > ("inputTracksTag"));
Expand All @@ -58,6 +64,11 @@ void HLTDeDxFilter::fillDescriptions(edm::ConfigurationDescriptions& description
desc.add<double>("minPT",0.0);
desc.add<double>("minNOM",0.0);
desc.add<double>("maxETA",5.5);
desc.add<double>("minNumValidHits",0);
desc.add<double>("maxNHitMissIn",99);
desc.add<double>("maxNHitMissMid",99);
desc.add<double>("maxRelTrkIsoDeltaRp3", -1);
desc.add<double>("relTrkIsoDeltaRSize", 0.3);
desc.add<edm::InputTag>("inputTracksTag",edm::InputTag("hltL3Mouns"));
desc.add<edm::InputTag>("inputDeDxTag",edm::InputTag("HLTdedxHarm2"));
descriptions.add("hltDeDxFilter",desc);
Expand Down Expand Up @@ -94,12 +105,22 @@ bool

bool accept=false;
int NTracks = 0;
//fill local arrays for eta, phi, and pt
float eta[trackCollection.size()], phi[trackCollection.size()], pt[trackCollection.size()];
for(unsigned int i=0; i<trackCollection.size(); i++){
reco::TrackRef track = reco::TrackRef( trackCollectionHandle, i );
if(track->pt()>minPT_ && fabs(track->eta())<maxETA_ && dEdxTrack[track].numberOfMeasurements()>minNOM_ && dEdxTrack[track].dEdx()>minDEDx_){
NTracks++;
eta[i] = trackCollection[i].eta();
phi[i] = trackCollection[i].phi();
pt[i] = trackCollection[i].pt();
}
for(unsigned int i=0; i<trackCollection.size(); i++){
reco::TrackRef track = reco::TrackRef( trackCollectionHandle, i );
if(pt[i]>minPT_ && fabs(eta[i])<maxETA_ && dEdxTrack[track].numberOfMeasurements()>minNOM_ && dEdxTrack[track].dEdx()>minDEDx_){
NTracks++;
if(track->numberOfValidHits() < minNumValidHits_) continue;
if(track->hitPattern().trackerLayersWithoutMeasurement( reco::HitPattern::MISSING_INNER_HITS) > maxNHitMissIn_) continue;
if(track->hitPattern().trackerLayersWithoutMeasurement( reco::HitPattern::TRACK_HITS) > maxNHitMissMid_) continue;
if (saveTags()){
Particle::Charge q = track->charge();
Particle::Charge q = track->charge();
//SAVE DEDX INFORMATION AS IF IT WAS THE MASS OF THE PARTICLE
Particle::LorentzVector p4(track->px(), track->py(), track->pz(), sqrt(pow(track->p(),2) + pow(dEdxTrack[track].dEdx(),2)));
Particle::Point vtx(track->vx(),track->vy(), track->vz());
Expand All @@ -109,6 +130,21 @@ bool
cand.setTrack(track);
chargedCandidates->push_back(cand);
}

//calculate relative trk isolation only if parameter maxRelTrkIsoDeltaRp3 is less than 0
if(maxRelTrkIsoDeltaRp3_ >= 0){
auto ptCone = trackCollection[i].pt();
for(unsigned int j=0; j<trackCollection.size(); j++){
reco::TrackRef track2 = reco::TrackRef( trackCollectionHandle, j );
if (i==j) continue; // do not compare track to itself
auto trkDeltaR2 = deltaR2(eta[i], phi[i], eta[j], phi[j]);
if (trkDeltaR2 < relTrkIsoDeltaRSize_ * relTrkIsoDeltaRSize_){
ptCone+=pt[j];
}
}
double relTrkIso = (ptCone - pt[i])/(pt[i]);
if (relTrkIso > maxRelTrkIsoDeltaRp3_) continue;
}
accept=true;
}
}
Expand Down
5 changes: 5 additions & 0 deletions RecoTracker/DeDx/plugins/HLTDeDxFilter.h
Original file line number Diff line number Diff line change
Expand Up @@ -33,6 +33,11 @@ class HLTDeDxFilter : public HLTFilter {
double minPT_;
double minNOM_;
double maxETA_;
double minNumValidHits_;
double maxNHitMissIn_;
double maxNHitMissMid_;
double maxRelTrkIsoDeltaRp3_;
double relTrkIsoDeltaRSize_;
edm::EDGetToken inputTracksToken_;
edm::EDGetToken inputdedxToken_;
edm::InputTag inputTracksTag_;
Expand Down

0 comments on commit a4a5e7d

Please sign in to comment.