Go to the documentation of this file.00001 #include <iostream>
00002
00003 #include "FWCore/Framework/interface/ESHandle.h"
00004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00005
00006 #include "DataFormats/Candidate/interface/CandidateFwd.h"
00007 #include "DataFormats/RecoCandidate/interface/RecoCandidate.h"
00008 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
00009 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidateFwd.h"
00010 #include "DataFormats/Candidate/interface/Candidate.h"
00011 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
00012 #include "DataFormats/HLTReco/interface/TriggerRefsCollections.h"
00013
00014 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00015 #include "DataFormats/TrackReco/interface/Track.h"
00016 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
00017 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
00018
00019 #include "RecoVertex/KalmanVertexFit/interface/KalmanVertexFitter.h"
00020 #include "RecoVertex/VertexPrimitives/interface/TransientVertex.h"
00021 #include "DataFormats/VertexReco/interface/Vertex.h"
00022 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00023 #include "DataFormats/Common/interface/RefToBase.h"
00024
00025
00026 #include "HLTDisplacedmumumuVtxProducer.h"
00027
00028 using namespace edm;
00029 using namespace reco;
00030 using namespace std;
00031 using namespace trigger;
00032
00033
00034
00035 HLTDisplacedmumumuVtxProducer::HLTDisplacedmumumuVtxProducer(const edm::ParameterSet& iConfig):
00036 src_ (iConfig.getParameter<edm::InputTag>("Src")),
00037 previousCandTag_(iConfig.getParameter<edm::InputTag>("PreviousCandTag")),
00038 maxEta_ (iConfig.getParameter<double>("MaxEta")),
00039 minPt_ (iConfig.getParameter<double>("MinPt")),
00040 minPtTriplet_ (iConfig.getParameter<double>("MinPtTriplet")),
00041 minInvMass_ (iConfig.getParameter<double>("MinInvMass")),
00042 maxInvMass_ (iConfig.getParameter<double>("MaxInvMass")),
00043 chargeOpt_ (iConfig.getParameter<int>("ChargeOpt"))
00044 {
00045 produces<VertexCollection>();
00046 }
00047
00048
00049 HLTDisplacedmumumuVtxProducer::~HLTDisplacedmumumuVtxProducer()
00050 {
00051
00052 }
00053
00054
00055
00056 void HLTDisplacedmumumuVtxProducer::beginJob()
00057 {
00058
00059 }
00060
00061
00062 void HLTDisplacedmumumuVtxProducer::endJob()
00063 {
00064
00065 }
00066
00067
00068 void HLTDisplacedmumumuVtxProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00069 {
00070 double const MuMass = 0.106;
00071 double const MuMass2 = MuMass*MuMass;
00072
00073
00074
00075 Handle<RecoChargedCandidateCollection> mucands;
00076 iEvent.getByLabel (src_,mucands);
00077
00078
00079 edm::ESHandle<TransientTrackBuilder> theB;
00080 iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",theB);
00081
00082 std::auto_ptr<VertexCollection> vertexCollection(new VertexCollection());
00083
00084
00085 double e1,e2,e3;
00086 Particle::LorentzVector p,p1,p2,p3;
00087
00088 RecoChargedCandidateCollection::const_iterator cand1;
00089 RecoChargedCandidateCollection::const_iterator cand2;
00090 RecoChargedCandidateCollection::const_iterator cand3;
00091
00092
00093 Handle<TriggerFilterObjectWithRefs> previousCands;
00094 iEvent.getByLabel (previousCandTag_,previousCands);
00095
00096 vector<RecoChargedCandidateRef> vPrevCands;
00097 previousCands->getObjects(TriggerMuon,vPrevCands);
00098
00099 for (cand1=mucands->begin(); cand1!=mucands->end(); cand1++) {
00100 TrackRef tk1 = cand1->get<TrackRef>();
00101 LogDebug("HLTDisplacedMumumuFilter") << " 1st muon in loop: q*pt= " << cand1->charge()*cand1->pt() << ", eta= " << cand1->eta() << ", hits= " << tk1->numberOfValidHits();
00102
00103
00104 if( ! checkPreviousCand( tk1, vPrevCands) ) continue;
00105
00106
00107 if (fabs(cand1->eta())>maxEta_) continue;
00108 if (cand1->pt() < minPt_) continue;
00109
00110 cand2 = cand1; cand2++;
00111 for (; cand2!=mucands->end(); cand2++) {
00112 TrackRef tk2 = cand2->get<TrackRef>();
00113
00114
00115 LogDebug("HLTMuonDimuonFilter") << " 2nd muon in loop: q*pt= " << cand2->charge()*cand2->pt() << ", eta= " << cand2->eta() << ", hits= " << tk2->numberOfValidHits() << ", d0= " << tk2->d0();
00116
00117 if( ! checkPreviousCand( tk2, vPrevCands) ) continue;
00118
00119
00120 if (fabs(cand2->eta())>maxEta_) continue;
00121 if (cand2->pt() < minPt_) continue;
00122
00123 cand3 = cand2; cand3++;
00124 for (; cand3!=mucands->end(); cand3++) {
00125 TrackRef tk3 = cand3->get<TrackRef>();
00126
00127
00128 LogDebug("HLTMuonDimuonFilter") << " 3rd muon in loop: q*pt= " << cand3->charge()*cand3->pt() << ", eta= " << cand3->eta() << ", hits= " << tk3->numberOfValidHits() << ", d0= " << tk3->d0();
00129
00130 if( ! checkPreviousCand( tk3, vPrevCands) ) continue;
00131
00132
00133 if (fabs(cand3->eta())>maxEta_) continue;
00134 if (cand3->pt() < minPt_) continue;
00135
00136
00137 if (chargeOpt_>0) {
00138 if (fabs (cand1->charge() + cand2->charge() + cand3->charge()) != chargeOpt_) continue;
00139 }
00140
00141
00142 e1 = sqrt(cand1->momentum().Mag2()+MuMass2);
00143 e2 = sqrt(cand2->momentum().Mag2()+MuMass2);
00144 e3 = sqrt(cand3->momentum().Mag2()+MuMass2);
00145 p1 = Particle::LorentzVector(cand1->px(),cand1->py(),cand1->pz(),e1);
00146 p2 = Particle::LorentzVector(cand2->px(),cand2->py(),cand2->pz(),e2);
00147 p3 = Particle::LorentzVector(cand3->px(),cand3->py(),cand3->pz(),e3);
00148 p = p1+p2+p3;
00149
00150
00151 if (p.pt()<minPtTriplet_) continue;
00152
00153 double invmass = abs(p.mass());
00154 LogDebug("HLTDisplacedMumumuFilter") << " ... 1-2 invmass= " << invmass;
00155
00156 if (invmass<minInvMass_) continue;
00157 if (invmass>maxInvMass_) continue;
00158
00159
00160 vector<TransientTrack> t_tks;
00161 TransientTrack ttkp1 = (*theB).build(&tk1);
00162 TransientTrack ttkp2 = (*theB).build(&tk2);
00163 TransientTrack ttkp3 = (*theB).build(&tk3);
00164 t_tks.push_back(ttkp1);
00165 t_tks.push_back(ttkp2);
00166 t_tks.push_back(ttkp3);
00167
00168
00169 if (t_tks.size()!=3) continue;
00170
00171 KalmanVertexFitter kvf;
00172 TransientVertex tv = kvf.vertex(t_tks);
00173
00174 if (!tv.isValid()) continue;
00175
00176 Vertex vertex = tv;
00177
00178
00179 vertexCollection->push_back(vertex);
00180 }
00181 }
00182 }
00183 iEvent.put(vertexCollection);
00184 }
00185
00186
00187
00188 bool HLTDisplacedmumumuVtxProducer::checkPreviousCand(const TrackRef& trackref, vector<RecoChargedCandidateRef> & refVect){
00189 bool ok=false;
00190 for (unsigned int i=0; i<refVect.size(); i++) {
00191 if ( refVect[i]->get<TrackRef>() == trackref ) {
00192 ok=true;
00193 break;
00194 }
00195 }
00196 return ok;
00197 }