Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Added configurable parameters for nValidHits, nMissIn, nMissMid, maxRelTrkIsoDeltaRp3, relTrkIsoDeltaRSize #5665

Merged
merged 5 commits into from Oct 8, 2014
Merged
Show file tree
Hide file tree
Changes from 3 commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Jump to
Jump to file
Failed to load files.
Diff view
Diff view
27 changes: 27 additions & 0 deletions RecoTracker/DeDx/plugins/HLTDeDxFilter.cc
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,10 @@ 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");
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 +63,10 @@ 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<edm::InputTag>("inputTracksTag",edm::InputTag("hltL3Mouns"));
desc.add<edm::InputTag>("inputDeDxTag",edm::InputTag("HLTdedxHarm2"));
descriptions.add("hltDeDxFilter",desc);
Expand Down Expand Up @@ -98,6 +107,9 @@ bool
reco::TrackRef track = reco::TrackRef( trackCollectionHandle, i );
if(track->pt()>minPT_ && fabs(track->eta())<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();
//SAVE DEDX INFORMATION AS IF IT WAS THE MASS OF THE PARTICLE
Expand All @@ -109,6 +121,21 @@ bool
cand.setTrack(track);
chargedCandidates->push_back(cand);
}

//calculate relative trk isolation only if parameter maxRelTrkIsoDeltaRp3 is different from default value of 99
if(maxRelTrkIsoDeltaRp3_ >= 0){
double ptCone = track->pt();
for(unsigned int j=0; j<trackCollection.size(); j++){
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Add here if(i==j) continue;
Then you don't need eta and phi comparisons below (except if there are special cases of two tracks same parameters)

reco::TrackRef track2 = reco::TrackRef( trackCollectionHandle, j );
if (track->eta() == track2->eta() && track->phi() == track2->phi()) continue; // do not compare track to itself
double trkDeltaR = deltaR(track->eta(), track->phi(), track2->eta(), track2->phi());
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Use deltaR2() and not deltaR().

Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

tracks do not store eta phi pt (candidate do)
so here you must first create local arrays
float eta[trackCollection.size()], phi[trackCollection.size()], pt[trackCollection.size()];
fill in one loop
and then proceed with the combinatorial.

(there is selection on i so ptcone optimization is irrelevant)

Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Hi Vincenzo,

Could you clarify for me what your concerns are about the code as it
written? Do you think that it will not work properly in this release or
that it will take too much time?

I am not very familiar with the syntax you suggested so if you could point
me to a piece of code where it is done this way, I am happy to implement
these changes. It's just not immediately clear to me why we should do it
this way.

Thanks for your help!

Cheers,
Jessica

On Fri, Oct 3, 2014 at 2:14 PM, Vincenzo Innocente <notifications@github.com

wrote:

In RecoTracker/DeDx/plugins/HLTDeDxFilter.cc:

@@ -109,6 +121,21 @@ bool
cand.setTrack(track);
chargedCandidates->push_back(cand);
}
+

  •   //calculate relative trk isolation only if parameter maxRelTrkIsoDeltaRp3 is different from default value of 99
    
  •   if(maxRelTrkIsoDeltaRp3_ >= 0){
    
  •   double ptCone = track->pt();
    
  •   for(unsigned int j=0; j<trackCollection.size(); j++){
    
  • reco::TrackRef track2 = reco::TrackRef( trackCollectionHandle, j );
  •     if (track->eta() == track2->eta() && track->phi() == track2->phi()) continue; // do not compare track to itself
    
  •     double trkDeltaR = deltaR(track->eta(), track->phi(), track2->eta(), track2->phi());
    

tracks do not store eta phi pt (candidate do)
so here you must first create local arrays
float eta[trackCollection.size()], phi[trackCollection.size()],
pt[trackCollection.size()];
float ptcone[trackCollection.size()]
fill in one loop
and then proceed with the combinatorial.
you can loop from i+1 if you fill both ptcone[i] and ptcone[j]
when you arrive at j (unless selection occurred) the sum over [0,j-1] have
been done in the previous iteration over i...


Reply to this email directly or view it on GitHub
https://github.com/cms-sw/cmssw/pull/5665/files#r18391493.

Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

On 3 Oct, 2014, at 3:39 PM, Jessica Brinson jbrins1@gmail.com wrote:

Hi Vincenzo,

Could you clarify for me what your concerns are about the code as it written? Do you think that it will not work properly in this release or that it will take too much time?
it will take time.

I am not very familiar with the syntax you suggested so if you could point me to a piece of code where it is done this way, I am happy to implement these changes. It's just not immediately clear to me why we should do it this way.

I am sure I did it somewhere else
for instance here
https://cmssdt.cern.ch/SDT/lxr/source/RecoEcal/EgammaClusterAlgos/src/Multi5x5BremRecoveryClusterAlgo.cc#0050
used in the double loop at
https://cmssdt.cern.ch/SDT/lxr/source/RecoEcal/EgammaClusterAlgos/src/Multi5x5BremRecoveryClusterAlgo.cc#0100
(this is a square window, for you is deltaR2)

v.=

if (trkDeltaR < 0.3){
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Better not have such a parameter (0.3) in the middle of the code, better name it properly.

ptCone+=track2->pt();
}
}
double relTrkIso = (ptCone - track->pt())/(track->pt());
if (relTrkIso > maxRelTrkIsoDeltaRp3_) continue;
}
accept=true;
}
}
Expand Down
4 changes: 4 additions & 0 deletions RecoTracker/DeDx/plugins/HLTDeDxFilter.h
Expand Up @@ -33,6 +33,10 @@ class HLTDeDxFilter : public HLTFilter {
double minPT_;
double minNOM_;
double maxETA_;
double minNumValidHits_;
double maxNHitMissIn_;
double maxNHitMissMid_;
double maxRelTrkIsoDeltaRp3_;
edm::EDGetToken inputTracksToken_;
edm::EDGetToken inputdedxToken_;
edm::InputTag inputTracksTag_;
Expand Down