CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_13_patch3/src/HLTrigger/btau/src/HLTDisplacedmumuVtxProducer.cc

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 "HLTDisplacedmumuVtxProducer.h"
00027 
00028 using namespace edm;
00029 using namespace reco;
00030 using namespace std; 
00031 using namespace trigger;
00032 //
00033 // constructors and destructor
00034 //
00035 HLTDisplacedmumuVtxProducer::HLTDisplacedmumuVtxProducer(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         minPtPair_ (iConfig.getParameter<double>("MinPtPair")),
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 HLTDisplacedmumuVtxProducer::~HLTDisplacedmumuVtxProducer()
00050 {
00051 
00052 }
00053 
00054 
00055 // ------------ method called once each job just before starting event loop  ------------
00056 void HLTDisplacedmumuVtxProducer::beginJob()
00057 {
00058 
00059 }
00060 
00061 // ------------ method called once each job just after ending the event loop  ------------
00062 void HLTDisplacedmumuVtxProducer::endJob() 
00063 {
00064         
00065 }
00066 
00067 // ------------ method called on each new Event  ------------
00068 void HLTDisplacedmumuVtxProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00069 {
00070         double const MuMass = 0.106;
00071         double const MuMass2 = MuMass*MuMass;
00072         
00073         
00074         // get hold of muon trks
00075         Handle<RecoChargedCandidateCollection> mucands;
00076         iEvent.getByLabel (src_,mucands);
00077         
00078         //get the transient track builder:
00079         edm::ESHandle<TransientTrackBuilder> theB;
00080         iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",theB);
00081 
00082         std::auto_ptr<VertexCollection> vertexCollection(new VertexCollection());
00083 
00084         // look at all mucands,  check cuts and make vertices
00085         double e1,e2;
00086         Particle::LorentzVector p,p1,p2;
00087         
00088         RecoChargedCandidateCollection::const_iterator cand1;
00089         RecoChargedCandidateCollection::const_iterator cand2;
00090         
00091         // get the objects passing the previous filter
00092         Handle<TriggerFilterObjectWithRefs> previousCands;
00093         iEvent.getByLabel (previousCandTag_,previousCands);
00094 
00095         vector<RecoChargedCandidateRef> vPrevCands;
00096         previousCands->getObjects(TriggerMuon,vPrevCands);
00097 
00098         for (cand1=mucands->begin(); cand1!=mucands->end(); cand1++) {
00099                TrackRef tk1 = cand1->get<TrackRef>();
00100                LogDebug("HLTDisplacedMumuFilter") << " 1st muon in loop: q*pt= " << cand1->charge()*cand1->pt() << ", eta= " << cand1->eta() << ", hits= " << tk1->numberOfValidHits();
00101              
00102                //first check if this muon passed the previous filter
00103                if( ! checkPreviousCand( tk1, vPrevCands) ) continue;
00104         
00105                // cuts
00106                if (fabs(cand1->eta())>maxEta_) continue;
00107                if (cand1->pt() < minPt_) continue;
00108               
00109                cand2 = cand1; cand2++;
00110                for (; cand2!=mucands->end(); cand2++) {
00111                          TrackRef tk2 = cand2->get<TrackRef>();
00112                  
00113                          // eta cut
00114                          LogDebug("HLTDisplacedmumuVtxProducer") << " 2nd muon in loop: q*pt= " << cand2->charge()*cand2->pt() << ", eta= " << cand2->eta() << ", hits= " << tk2->numberOfValidHits() << ", d0= " << tk2->d0();
00115                          //first check if this muon passed the previous filter
00116                          if( ! checkPreviousCand( tk2, vPrevCands) ) continue;
00117                          
00118                          // cuts
00119                          if (fabs(cand2->eta())>maxEta_) continue;
00120                          if (cand2->pt() < minPt_) continue;
00121                          
00122                          // opposite sign or same sign
00123                          if (chargeOpt_<0) {
00124                            if (cand1->charge()*cand2->charge()>0) continue;
00125                          } else if (chargeOpt_>0) {
00126                            if (cand1->charge()*cand2->charge()<0) continue;
00127                          }
00128                          
00129                          // Combined dimuon system
00130                          e1 = sqrt(cand1->momentum().Mag2()+MuMass2);
00131                          e2 = sqrt(cand2->momentum().Mag2()+MuMass2);
00132                          p1 = Particle::LorentzVector(cand1->px(),cand1->py(),cand1->pz(),e1);
00133                          p2 = Particle::LorentzVector(cand2->px(),cand2->py(),cand2->pz(),e2);
00134                          p = p1+p2;
00135                          
00136                          
00137                          if (p.pt()<minPtPair_) continue;
00138                          
00139                          double invmass = abs(p.mass());
00140                          LogDebug("HLTDisplacedMumuFilter") << " ... 1-2 invmass= " << invmass;
00141                          
00142                          if (invmass<minInvMass_) continue;
00143                          if (invmass>maxInvMass_) continue;
00144                          
00145                          // do the vertex fit
00146                          vector<TransientTrack> t_tks;
00147                          TransientTrack ttkp1 = (*theB).build(&tk1);
00148                          TransientTrack ttkp2 = (*theB).build(&tk2);
00149                          t_tks.push_back(ttkp1);
00150                          t_tks.push_back(ttkp2);
00151                          
00152                 
00153                          if (t_tks.size()!=2) continue;
00154                         
00155                          KalmanVertexFitter kvf;
00156                          TransientVertex tv = kvf.vertex(t_tks);
00157 
00158                          if (!tv.isValid()) continue;
00159                          
00160                          Vertex vertex = tv;
00161 
00162                          // put vertex in the event
00163                          vertexCollection->push_back(vertex);
00164                }
00165         }
00166         iEvent.put(vertexCollection);
00167 }
00168 
00169 
00170 
00171 bool HLTDisplacedmumuVtxProducer::checkPreviousCand(const TrackRef& trackref, vector<RecoChargedCandidateRef> & refVect){
00172   bool ok=false;
00173   for (unsigned int i=0; i<refVect.size(); i++) {
00174     if ( refVect[i]->get<TrackRef>() == trackref ) {
00175       ok=true;
00176       break;
00177     }
00178   }
00179   return ok;
00180 }