Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:15:50

0001 #include <algorithm>
0002 #include <cmath>
0003 #include <vector>
0004 
0005 #include "FWCore/Framework/interface/ESHandle.h"
0006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0007 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0008 
0009 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0010 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0011 #include "DataFormats/RecoCandidate/interface/RecoCandidate.h"
0012 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
0013 #include "DataFormats/VertexReco/interface/Vertex.h"
0014 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0015 
0016 #include "RecoVertex/KalmanVertexFit/interface/KalmanVertexFitter.h"
0017 #include "RecoVertex/VertexPrimitives/interface/TransientVertex.h"
0018 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
0019 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
0020 
0021 #include "HLTmumutktkFilter.h"
0022 #include "TMath.h"
0023 
0024 using namespace edm;
0025 using namespace reco;
0026 using namespace std;
0027 using namespace trigger;
0028 
0029 // ----------------------------------------------------------------------
0030 HLTmumutktkFilter::HLTmumutktkFilter(const edm::ParameterSet& iConfig)
0031     : HLTFilter(iConfig),
0032       muCandTag_(iConfig.getParameter<edm::InputTag>("MuonTag")),
0033       muCandToken_(consumes<reco::RecoChargedCandidateCollection>(muCandTag_)),
0034       trkCandTag_(iConfig.getParameter<edm::InputTag>("TrackTag")),
0035       trkCandToken_(consumes<reco::RecoChargedCandidateCollection>(trkCandTag_)),
0036       MuMuTkVertexTag_(iConfig.getParameter<edm::InputTag>("MuMuTkVertexTag")),
0037       MuMuTkVertexToken_(consumes<reco::VertexCollection>(MuMuTkVertexTag_)),
0038       beamSpotTag_(iConfig.getParameter<edm::InputTag>("BeamSpotTag")),
0039       beamSpotToken_(consumes<reco::BeamSpot>(beamSpotTag_)),
0040       maxEta_(iConfig.getParameter<double>("MaxEta")),
0041       minPt_(iConfig.getParameter<double>("MinPt")),
0042       maxNormalisedChi2_(iConfig.getParameter<double>("MaxNormalisedChi2")),
0043       minVtxProbability_(iConfig.getParameter<double>("MinVtxProbability")),
0044       minLxySignificance_(iConfig.getParameter<double>("MinLxySignificance")),
0045       minCosinePointingAngle_(iConfig.getParameter<double>("MinCosinePointingAngle")) {}
0046 
0047 // ----------------------------------------------------------------------
0048 HLTmumutktkFilter::~HLTmumutktkFilter() = default;
0049 
0050 void HLTmumutktkFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0051   edm::ParameterSetDescription desc;
0052   makeHLTFilterDescription(desc);
0053   desc.add<double>("MaxEta", 2.5);
0054   desc.add<double>("MinPt", 0.0);
0055   desc.add<double>("MaxNormalisedChi2", 10.0);
0056   desc.add<double>("MinVtxProbability", 0.0);
0057   desc.add<double>("MinLxySignificance", 3.0);
0058   desc.add<double>("MinCosinePointingAngle", 0.9);
0059   desc.add<edm::InputTag>("MuonTag", edm::InputTag("hltL3MuonCandidates"));
0060   desc.add<edm::InputTag>("TrackTag", edm::InputTag("hltMumukAllConeTracks"));
0061   desc.add<edm::InputTag>("MuMuTkVertexTag", edm::InputTag("hltDisplacedmumuVtxProducerDoubleMu4Jpsi"));
0062   desc.add<edm::InputTag>("BeamSpotTag", edm::InputTag("hltOnineBeamSpot"));
0063   descriptions.add("HLTmumutktkFilter", desc);
0064 }
0065 
0066 // ----------------------------------------------------------------------
0067 bool HLTmumutktkFilter::hltFilter(edm::Event& iEvent,
0068                                   const edm::EventSetup& iSetup,
0069                                   trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0070   //get the beamspot position
0071   edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
0072   iEvent.getByToken(beamSpotToken_, recoBeamSpotHandle);
0073   const reco::BeamSpot& vertexBeamSpot = *recoBeamSpotHandle;
0074 
0075   // get vertices
0076   reco::VertexCollection displacedVertexColl;
0077   edm::Handle<reco::VertexCollection> displacedVertexCollHandle;
0078   bool foundVertexColl = iEvent.getByToken(MuMuTkVertexToken_, displacedVertexCollHandle);
0079   if (foundVertexColl)
0080     displacedVertexColl = *displacedVertexCollHandle;
0081 
0082   // get muon collection
0083   Handle<RecoChargedCandidateCollection> mucands;
0084   iEvent.getByToken(muCandToken_, mucands);
0085 
0086   // get track candidates around displaced muons
0087   Handle<RecoChargedCandidateCollection> trkcands;
0088   iEvent.getByToken(trkCandToken_, trkcands);
0089 
0090   // Ref to Candidate object to be recorded in filter object
0091   RecoChargedCandidateRef refMu1;
0092   RecoChargedCandidateRef refMu2;
0093   RecoChargedCandidateRef refTrk1;
0094   RecoChargedCandidateRef refTrk2;
0095 
0096   if (saveTags()) {
0097     filterproduct.addCollectionTag(muCandTag_);
0098     filterproduct.addCollectionTag(trkCandTag_);
0099   }
0100 
0101   bool triggered = false;
0102 
0103   // loop over vertex collection
0104   reco::VertexCollection::iterator it;
0105   for (it = displacedVertexColl.begin(); it != displacedVertexColl.end(); it++) {
0106     reco::Vertex displacedVertex = *it;
0107 
0108     // check if the vertex actually consists of exactly two muon + 1 track, throw exception if not
0109     if (displacedVertex.tracksSize() != 4)
0110       throw cms::Exception("BadLogic") << "HLTmumutktkFilter: ERROR: the Jpsi+trk vertex must have "
0111                                        << "exactly two muons + 2 trk by definition. It now has n trakcs = "
0112                                        << displacedVertex.tracksSize() << std::endl;
0113 
0114     float normChi2 = displacedVertex.normalizedChi2();
0115     if (normChi2 > maxNormalisedChi2_)
0116       continue;
0117 
0118     double vtxProb = 0.0;
0119     if ((displacedVertex.chi2() >= 0.0) && (displacedVertex.ndof() > 0))
0120       vtxProb = TMath::Prob(displacedVertex.chi2(), displacedVertex.ndof());
0121     if (vtxProb < minVtxProbability_)
0122       continue;
0123 
0124     // get the four tracks from the vertex
0125     std::vector<reco::TrackRef> vertexTrksRefVec;
0126     auto trackIt = displacedVertex.tracks_begin();
0127     for (trackIt = displacedVertex.tracks_begin(); trackIt != displacedVertex.tracks_end(); trackIt++) {
0128       vertexTrksRefVec.push_back((*trackIt).castTo<reco::TrackRef>());
0129     }
0130 
0131     // first find the tracks in the input muon/track collection
0132     std::vector<reco::RecoChargedCandidateCollection::const_iterator> mucandVec;
0133     std::vector<reco::RecoChargedCandidateCollection::const_iterator> trkcandVec;
0134     for (auto cand = mucands->begin(); cand != mucands->end(); cand++) {
0135       reco::TrackRef tkRef = cand->get<reco::TrackRef>();
0136       for (auto& iVec : vertexTrksRefVec) {
0137         if (tkRef == iVec) {
0138           mucandVec.push_back(cand);
0139           break;
0140         }
0141       }
0142     }
0143     if (mucandVec.size() < 2)
0144       throw cms::Exception("BadLogic") << "HLTmumutktkFilterr: ERROR: the vertex must have "
0145                                        << " at least two muons by definition." << std::endl;
0146 
0147     for (auto cand = trkcands->begin(); cand != trkcands->end(); cand++) {
0148       reco::TrackRef tkRef = cand->get<reco::TrackRef>();
0149       for (auto& iVec : vertexTrksRefVec) {
0150         if (tkRef == iVec) {
0151           trkcandVec.push_back(cand);
0152           break;
0153         }
0154       }
0155     }
0156     if (trkcandVec.size() < 2)
0157       throw cms::Exception("BadLogic") << "HLTmumutktkFilterr: ERROR: the vertex must have "
0158                                        << " at least two tracks by definition." << std::endl;
0159 
0160     // calculate four-track transverse momentum
0161     math::XYZVector pperp(
0162         mucandVec.at(0)->px() + mucandVec.at(1)->px() + trkcandVec.at(0)->px() + trkcandVec.at(1)->px(),
0163         mucandVec.at(0)->py() + mucandVec.at(1)->py() + trkcandVec.at(0)->px() + trkcandVec.at(1)->py(),
0164         0.);
0165 
0166     // get vertex position and error to calculate the decay length significance
0167     const reco::Vertex::Point& vpoint = displacedVertex.position();
0168     reco::Vertex::Error verr = displacedVertex.error();
0169     GlobalPoint secondaryVertex(vpoint.x(), vpoint.y(), vpoint.z());
0170     GlobalError err(verr.At(0, 0), verr.At(1, 0), verr.At(1, 1), verr.At(2, 0), verr.At(2, 1), verr.At(2, 2));
0171 
0172     GlobalPoint displacementFromBeamspot(-1 * ((vertexBeamSpot.x0() - secondaryVertex.x()) +
0173                                                (secondaryVertex.z() - vertexBeamSpot.z0()) * vertexBeamSpot.dxdz()),
0174                                          -1 * ((vertexBeamSpot.y0() - secondaryVertex.y()) +
0175                                                (secondaryVertex.z() - vertexBeamSpot.z0()) * vertexBeamSpot.dydz()),
0176                                          0);
0177     float lxy = displacementFromBeamspot.perp();
0178     float lxyerr = sqrt(err.rerr(displacementFromBeamspot));
0179     float lxysign = 0;
0180     if (lxyerr != 0)
0181       lxysign = lxy / lxyerr;
0182 
0183     //calculate the angle between the decay length and the four-track momentum
0184     Vertex::Point vperp(displacementFromBeamspot.x(), displacementFromBeamspot.y(), 0.);
0185     float cosAlpha = vperp.Dot(pperp) / (vperp.Rho() * pperp.Rho());
0186 
0187     if (pperp.Rho() < minPt_)
0188       continue;
0189     if (lxysign < minLxySignificance_)
0190       continue;
0191     if (cosAlpha < minCosinePointingAngle_)
0192       continue;
0193     triggered = true;
0194 
0195     refMu1 = RecoChargedCandidateRef(
0196         Ref<RecoChargedCandidateCollection>(mucands, distance(mucands->begin(), mucandVec.at(0))));
0197     filterproduct.addObject(TriggerMuon, refMu1);
0198     refMu2 = RecoChargedCandidateRef(
0199         Ref<RecoChargedCandidateCollection>(mucands, distance(mucands->begin(), mucandVec.at(1))));
0200     filterproduct.addObject(TriggerMuon, refMu2);
0201     refTrk1 = RecoChargedCandidateRef(
0202         Ref<RecoChargedCandidateCollection>(trkcands, distance(trkcands->begin(), trkcandVec.at(0))));
0203     filterproduct.addObject(TriggerTrack, refTrk1);
0204     refTrk2 = RecoChargedCandidateRef(
0205         Ref<RecoChargedCandidateCollection>(trkcands, distance(trkcands->begin(), trkcandVec.at(1))));
0206     filterproduct.addObject(TriggerTrack, refTrk2);
0207 
0208   }  //end loop vertices
0209 
0210   LogDebug("HLTDisplacedMumuTrkFilter") << " >>>>> Result of HLTDisplacedMuMuTrkFilter is " << triggered;
0211   return triggered;
0212 }
0213 
0214 bool HLTmumutktkFilter::triggerdByPreviousLevel(const reco::RecoChargedCandidateRef& candref,
0215                                                 const std::vector<reco::RecoChargedCandidateRef>& vcands) {
0216   unsigned int i = 0;
0217   unsigned int i_max = vcands.size();
0218   for (; i != i_max; ++i) {
0219     if (candref == vcands[i])
0220       return true;
0221   }
0222   return false;
0223 }