CMS 3D CMS Logo

HLTPixlMBFilt Class Reference

See header file for documentation. More...

#include <HLTrigger/special/interface/HLTPixlMBFilt.h>

Inheritance diagram for HLTPixlMBFilt:

HLTFilter edm::EDFilter edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

virtual bool filter (edm::Event &, const edm::EventSetup &)
 HLTPixlMBFilt (const edm::ParameterSet &)
 ~HLTPixlMBFilt ()

Private Attributes

double min_Pt_
float min_sep_
unsigned int min_trks_
edm::InputTag pixlTag_


Detailed Description

See header file for documentation.

Date
2008/01/10 07:52:24
Revision
1.3

Author:
Mika Huhtinen

Definition at line 22 of file HLTPixlMBFilt.h.


Constructor & Destructor Documentation

HLTPixlMBFilt::HLTPixlMBFilt ( const edm::ParameterSet iConfig  )  [explicit]

Definition at line 34 of file HLTPixlMBFilt.cc.

References edm::InputTag::encode(), LogDebug, min_Pt_, min_sep_, min_trks_, and pixlTag_.

00034                                                            :
00035     pixlTag_ (iConfig.getParameter<edm::InputTag>("pixlTag")),
00036     min_Pt_  (iConfig.getParameter<double>("MinPt")),
00037     min_trks_  (iConfig.getParameter<unsigned int>("MinTrks")),
00038     min_sep_  (iConfig.getParameter<double>("MinSep"))
00039 
00040 {
00041   LogDebug("") << "MinPt cut " << min_Pt_   << "pixl: " << pixlTag_.encode();
00042   LogDebug("") << "Requesting : " << min_trks_ << " tracks from same vertex ";
00043   LogDebug("") << "Requesting tracks from same vertex eta-phi separation by " << min_sep_;
00044 
00045    //register your products
00046   produces<trigger::TriggerFilterObjectWithRefs>();
00047 }

HLTPixlMBFilt::~HLTPixlMBFilt (  ) 

Definition at line 49 of file HLTPixlMBFilt.cc.

00050 {
00051 }


Member Function Documentation

bool HLTPixlMBFilt::filter ( edm::Event iEvent,
const edm::EventSetup iSetup 
) [virtual]

Implements HLTFilter.

Definition at line 58 of file HLTPixlMBFilt.cc.

References edm::Event::getByLabel(), k, LogDebug, min_sep_, min_trks_, module(), path(), pixlTag_, edm::Event::put(), HcalSimpleRecAlgoImpl::reco(), funct::sqrt(), std, tracks, and trigger::TriggerTrack.

00059 {
00060    using namespace std;
00061    using namespace edm;
00062    using namespace reco;
00063    using namespace trigger;
00064 
00065    // All HLT filters must create and fill an HLT filter object,
00066    // recording any reconstructed physics objects satisfying (or not)
00067    // this HLT filter, and place it in the Event.
00068 
00069 
00070 
00071    // Specific filter code
00072 
00073    // get hold of products from Event
00074 
00075    Handle<RecoChargedCandidateCollection> tracks;
00076 
00077    iEvent.getByLabel(pixlTag_,tracks);
00078 
00079    // pixel tracks
00080    int npixl_tot = 0;
00081    vector<double> etastore;
00082    vector<double> phistore;
00083    vector<int> itstore;
00084    bool accept;
00085    RecoChargedCandidateCollection::const_iterator apixl(tracks->begin());
00086    RecoChargedCandidateCollection::const_iterator epixl(tracks->end());
00087    RecoChargedCandidateCollection::const_iterator ipixl, jpixl;
00088    unsigned int nsame_vtx=0;
00089    int itrk = -1;
00090    if (tracks->size() >= min_trks_) {
00091      for (ipixl=apixl; ipixl!=epixl; ipixl++){ 
00092        itrk++;
00093        const double& ztrk1 = ipixl->vz();                   
00094        const double& etatrk1 = ipixl->momentum().eta();
00095        const double& phitrk1 = ipixl->momentum().phi();
00096        nsame_vtx=1;
00097        etastore.clear();
00098        phistore.clear();
00099        itstore.clear();
00100        etastore.push_back(etatrk1);
00101        phistore.push_back(phitrk1);
00102        itstore.push_back(itrk);
00103        if (fabs(ztrk1) < 15.0) {
00104          //  check this track against all others to see if others start from same point
00105          int jtrk=-1;
00106          for (jpixl=apixl; jpixl!=epixl; jpixl++) {
00107            jtrk++;
00108            if (jpixl==ipixl) continue;
00109            const double& ztrk2 = jpixl->vz();               
00110            const double& etatrk2 = jpixl->momentum().eta();
00111            const double& phitrk2 = jpixl->momentum().phi();
00112            double eta_dist=etatrk2-etatrk1;
00113            double phi_dist=phitrk2-phitrk1;
00114            double etaphi_dist=sqrt(eta_dist*eta_dist + phi_dist*phi_dist);
00115            if (fabs(ztrk2-ztrk1) < 1.0 && etaphi_dist > min_sep_) {
00116               if (min_trks_ <= 2 || itstore.size() <= 1) {
00117                 etastore.push_back(etatrk2);
00118                 phistore.push_back(phitrk2);
00119                 itstore.push_back(jtrk);
00120                 nsame_vtx++;
00121               } else {
00122                 // check also separation to already found 'second' tracks
00123                 LogDebug("") << "HLTPixlMBFilt: with mintrks=2 we should not be here...";
00124                 bool isok = true;
00125                 for (unsigned int k=1; k < itstore.size(); k++) {
00126                   eta_dist=etatrk2-etastore.at(k);
00127                   phi_dist=phitrk2-phistore.at(k);
00128                   etaphi_dist=sqrt(eta_dist*eta_dist + phi_dist*phi_dist);
00129                   if (etaphi_dist < min_sep_) {
00130                     isok=false;
00131                     break;
00132                   }
00133                 }
00134                 if (isok) {
00135                   etastore.push_back(etatrk2);
00136                   phistore.push_back(phitrk2);
00137                   itstore.push_back(jtrk);
00138                   nsame_vtx++;
00139                 }
00140               }
00141            }
00142            if (nsame_vtx >= min_trks_) break;
00143          }
00144        }
00145        npixl_tot++;
00146 
00147        if (nsame_vtx >= min_trks_) break;
00148      }
00149 
00150      //   final filter decision:
00151      //   request at least min_trks_ tracks compatible with vertex-region
00152      accept = (nsame_vtx >= min_trks_ ) ;
00153 
00154    } else {
00155      accept = false;
00156    }
00157 
00158    // At this point we have the indices of the accepted tracks stored in itstore
00159    // we now move them to the filterobject
00160 
00161    // The filter object
00162    auto_ptr<TriggerFilterObjectWithRefs> filterobject (new TriggerFilterObjectWithRefs(path(),module()));
00163 
00164    if (accept) {
00165      for (unsigned int ipos=0; ipos < itstore.size(); ipos++) {
00166        int iaddr=itstore.at(ipos);
00167        filterobject->addObject(TriggerTrack,RecoChargedCandidateRef(tracks,iaddr));
00168      }
00169    }
00170    // put filter object into the Event
00171    iEvent.put(filterobject);
00172 
00173 
00174   LogDebug("") << "Number of pixel-track objects accepted:"
00175                << " " << npixl_tot;
00176 
00177   // return with final filter decision
00178   return accept;
00179 
00180 }


Member Data Documentation

double HLTPixlMBFilt::min_Pt_ [private]

Definition at line 32 of file HLTPixlMBFilt.h.

Referenced by HLTPixlMBFilt().

float HLTPixlMBFilt::min_sep_ [private]

Definition at line 34 of file HLTPixlMBFilt.h.

Referenced by filter(), and HLTPixlMBFilt().

unsigned int HLTPixlMBFilt::min_trks_ [private]

Definition at line 33 of file HLTPixlMBFilt.h.

Referenced by filter(), and HLTPixlMBFilt().

edm::InputTag HLTPixlMBFilt::pixlTag_ [private]

Definition at line 30 of file HLTPixlMBFilt.h.

Referenced by filter(), and HLTPixlMBFilt().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:24:40 2009 for CMSSW by  doxygen 1.5.4