00001
00012 #include "HLTrigger/special/interface/HLTPixlMBFilt.h"
00013
00014 #include "FWCore/Framework/interface/Event.h"
00015 #include "FWCore/Framework/interface/EventSetup.h"
00016 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00017
00018 #include "DataFormats/Common/interface/Handle.h"
00019
00020 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
00021 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidateFwd.h"
00022
00023 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
00024
00025 #include "DataFormats/TrackReco/interface/Track.h"
00026 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00027
00028 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00029
00030
00031
00032
00033
00034 HLTPixlMBFilt::HLTPixlMBFilt(const edm::ParameterSet& iConfig) :
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
00046 produces<trigger::TriggerFilterObjectWithRefs>();
00047 }
00048
00049 HLTPixlMBFilt::~HLTPixlMBFilt()
00050 {
00051 }
00052
00053
00054
00055
00056
00057
00058 bool HLTPixlMBFilt::filter(edm::Event& iEvent, const edm::EventSetup& iSetup)
00059 {
00060 using namespace std;
00061 using namespace edm;
00062 using namespace reco;
00063 using namespace trigger;
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075 Handle<RecoChargedCandidateCollection> tracks;
00076
00077 iEvent.getByLabel(pixlTag_,tracks);
00078
00079
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
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
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
00151
00152 accept = (nsame_vtx >= min_trks_ ) ;
00153
00154 } else {
00155 accept = false;
00156 }
00157
00158
00159
00160
00161
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
00171 iEvent.put(filterobject);
00172
00173
00174 LogDebug("") << "Number of pixel-track objects accepted:"
00175 << " " << npixl_tot;
00176
00177
00178 return accept;
00179
00180 }