Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:05:53

0001 #include <iostream>
0002 
0003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0004 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0005 #include "DataFormats/Candidate/interface/CandidateFwd.h"
0006 #include "DataFormats/RecoCandidate/interface/RecoCandidate.h"
0007 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
0008 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidateFwd.h"
0009 #include "DataFormats/Candidate/interface/Candidate.h"
0010 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0011 #include "DataFormats/HLTReco/interface/TriggerRefsCollections.h"
0012 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0013 #include "DataFormats/TrackReco/interface/Track.h"
0014 #include "DataFormats/VertexReco/interface/Vertex.h"
0015 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0016 #include "DataFormats/Common/interface/RefToBase.h"
0017 #include "RecoVertex/KalmanVertexFit/interface/KalmanVertexFitter.h"
0018 #include "RecoVertex/VertexPrimitives/interface/TransientVertex.h"
0019 #include "HLTDisplacedmumumuVtxProducer.h"
0020 
0021 using namespace edm;
0022 using namespace reco;
0023 using namespace std;
0024 using namespace trigger;
0025 //
0026 // constructors and destructor
0027 //
0028 HLTDisplacedmumumuVtxProducer::HLTDisplacedmumumuVtxProducer(const edm::ParameterSet& iConfig)
0029     : transientTrackRecordToken_(esConsumes(edm::ESInputTag("", "TransientTrackBuilder"))),
0030       srcTag_(iConfig.getParameter<edm::InputTag>("Src")),
0031       srcToken_(consumes<reco::RecoChargedCandidateCollection>(srcTag_)),
0032       previousCandTag_(iConfig.getParameter<edm::InputTag>("PreviousCandTag")),
0033       previousCandToken_(consumes<trigger::TriggerFilterObjectWithRefs>(previousCandTag_)),
0034       maxEta_(iConfig.getParameter<double>("MaxEta")),
0035       minPt_(iConfig.getParameter<double>("MinPt")),
0036       minPtTriplet_(iConfig.getParameter<double>("MinPtTriplet")),
0037       minInvMass_(iConfig.getParameter<double>("MinInvMass")),
0038       maxInvMass_(iConfig.getParameter<double>("MaxInvMass")),
0039       chargeOpt_(iConfig.getParameter<int>("ChargeOpt")) {
0040   produces<VertexCollection>();
0041 }
0042 
0043 HLTDisplacedmumumuVtxProducer::~HLTDisplacedmumumuVtxProducer() = default;
0044 
0045 void HLTDisplacedmumumuVtxProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0046   edm::ParameterSetDescription desc;
0047   desc.add<edm::InputTag>("Src", edm::InputTag("hltL3MuonCandidates"));
0048   desc.add<edm::InputTag>("PreviousCandTag", edm::InputTag(""));
0049   desc.add<double>("MaxEta", 2.5);
0050   desc.add<double>("MinPt", 0.0);
0051   desc.add<double>("MinPtTriplet", 0.0);
0052   desc.add<double>("MinInvMass", 1.0);
0053   desc.add<double>("MaxInvMass", 20.0);
0054   desc.add<int>("ChargeOpt", -1);
0055   descriptions.add("hltDisplacedmumumuVtxProducer", desc);
0056 }
0057 
0058 // ------------ method called on each new Event  ------------
0059 void HLTDisplacedmumumuVtxProducer::produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0060   double const MuMass = 0.106;
0061   double const MuMass2 = MuMass * MuMass;
0062 
0063   // get hold of muon trks
0064   Handle<RecoChargedCandidateCollection> mucands;
0065   iEvent.getByToken(srcToken_, mucands);
0066 
0067   // get the transient track builder
0068   auto const& theB = iSetup.getHandle(transientTrackRecordToken_);
0069 
0070   std::unique_ptr<VertexCollection> vertexCollection(new VertexCollection());
0071 
0072   // look at all mucands,  check cuts and make vertices
0073   double e1, e2, e3;
0074   Particle::LorentzVector p, p1, p2, p3;
0075 
0076   RecoChargedCandidateCollection::const_iterator cand1;
0077   RecoChargedCandidateCollection::const_iterator cand2;
0078   RecoChargedCandidateCollection::const_iterator cand3;
0079 
0080   // get the objects passing the previous filter
0081   Handle<TriggerFilterObjectWithRefs> previousCands;
0082   iEvent.getByToken(previousCandToken_, previousCands);
0083 
0084   vector<RecoChargedCandidateRef> vPrevCands;
0085   previousCands->getObjects(TriggerMuon, vPrevCands);
0086 
0087   for (cand1 = mucands->begin(); cand1 != mucands->end(); cand1++) {
0088     TrackRef tk1 = cand1->get<TrackRef>();
0089     LogDebug("HLTDisplacedMumumuFilter") << " 1st muon in loop: q*pt= " << cand1->charge() * cand1->pt()
0090                                          << ", eta= " << cand1->eta() << ", hits= " << tk1->numberOfValidHits();
0091 
0092     //first check if this muon passed the previous filter
0093     if (!checkPreviousCand(tk1, vPrevCands))
0094       continue;
0095 
0096     // cuts
0097     if (fabs(cand1->eta()) > maxEta_)
0098       continue;
0099     if (cand1->pt() < minPt_)
0100       continue;
0101 
0102     cand2 = cand1;
0103     cand2++;
0104     for (; cand2 != mucands->end(); cand2++) {
0105       TrackRef tk2 = cand2->get<TrackRef>();
0106 
0107       // eta cut
0108       LogDebug("HLTMuonDimuonFilter") << " 2nd muon in loop: q*pt= " << cand2->charge() * cand2->pt()
0109                                       << ", eta= " << cand2->eta() << ", hits= " << tk2->numberOfValidHits()
0110                                       << ", d0= " << tk2->d0();
0111       //first check if this muon passed the previous filter
0112       if (!checkPreviousCand(tk2, vPrevCands))
0113         continue;
0114 
0115       // cuts
0116       if (fabs(cand2->eta()) > maxEta_)
0117         continue;
0118       if (cand2->pt() < minPt_)
0119         continue;
0120 
0121       cand3 = cand2;
0122       cand3++;
0123       for (; cand3 != mucands->end(); cand3++) {
0124         TrackRef tk3 = cand3->get<TrackRef>();
0125 
0126         // eta cut
0127         LogDebug("HLTMuonDimuonFilter") << " 3rd muon in loop: q*pt= " << cand3->charge() * cand3->pt()
0128                                         << ", eta= " << cand3->eta() << ", hits= " << tk3->numberOfValidHits()
0129                                         << ", d0= " << tk3->d0();
0130         //first check if this muon passed the previous filter
0131         if (!checkPreviousCand(tk3, vPrevCands))
0132           continue;
0133 
0134         // cuts
0135         if (fabs(cand3->eta()) > maxEta_)
0136           continue;
0137         if (cand3->pt() < minPt_)
0138           continue;
0139 
0140         // opposite sign or same sign
0141         if (chargeOpt_ > 0) {
0142           if (fabs(cand1->charge() + cand2->charge() + cand3->charge()) != chargeOpt_)
0143             continue;
0144         }
0145 
0146         // Combined dimuon system
0147         e1 = sqrt(cand1->momentum().Mag2() + MuMass2);
0148         e2 = sqrt(cand2->momentum().Mag2() + MuMass2);
0149         e3 = sqrt(cand3->momentum().Mag2() + MuMass2);
0150         p1 = Particle::LorentzVector(cand1->px(), cand1->py(), cand1->pz(), e1);
0151         p2 = Particle::LorentzVector(cand2->px(), cand2->py(), cand2->pz(), e2);
0152         p3 = Particle::LorentzVector(cand3->px(), cand3->py(), cand3->pz(), e3);
0153         p = p1 + p2 + p3;
0154 
0155         if (p.pt() < minPtTriplet_)
0156           continue;
0157 
0158         double invmass = abs(p.mass());
0159         LogDebug("HLTDisplacedMumumuFilter") << " ... 1-2 invmass= " << invmass;
0160 
0161         if (invmass < minInvMass_)
0162           continue;
0163         if (invmass > maxInvMass_)
0164           continue;
0165 
0166         // do the vertex fit
0167         vector<TransientTrack> t_tks;
0168         TransientTrack ttkp1 = (*theB).build(&tk1);
0169         TransientTrack ttkp2 = (*theB).build(&tk2);
0170         TransientTrack ttkp3 = (*theB).build(&tk3);
0171         t_tks.push_back(ttkp1);
0172         t_tks.push_back(ttkp2);
0173         t_tks.push_back(ttkp3);
0174 
0175         if (t_tks.size() != 3)
0176           continue;
0177 
0178         KalmanVertexFitter kvf;
0179         TransientVertex tv = kvf.vertex(t_tks);
0180 
0181         if (!tv.isValid())
0182           continue;
0183 
0184         Vertex vertex = tv;
0185 
0186         // put vertex in the event
0187         vertexCollection->push_back(vertex);
0188       }
0189     }
0190   }
0191   iEvent.put(std::move(vertexCollection));
0192 }
0193 
0194 bool HLTDisplacedmumumuVtxProducer::checkPreviousCand(const TrackRef& trackref,
0195                                                       const vector<RecoChargedCandidateRef>& refVect) const {
0196   bool ok = false;
0197   for (auto& i : refVect) {
0198     if (i->get<TrackRef>() == trackref) {
0199       ok = true;
0200       break;
0201     }
0202   }
0203   return ok;
0204 }