Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include <iostream>
0002 
0003 #include "FWCore/Framework/interface/ESHandle.h"
0004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0005 
0006 #include "DataFormats/Candidate/interface/CandidateFwd.h"
0007 #include "DataFormats/RecoCandidate/interface/RecoCandidate.h"
0008 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
0009 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidateFwd.h"
0010 #include "DataFormats/Candidate/interface/Candidate.h"
0011 
0012 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0013 #include "DataFormats/TrackReco/interface/Track.h"
0014 
0015 #include "DataFormats/VertexReco/interface/Vertex.h"
0016 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0017 #include "DataFormats/Common/interface/RefToBase.h"
0018 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0019 #include "DataFormats/HLTReco/interface/TriggerRefsCollections.h"
0020 
0021 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0022 
0023 #include "DataFormats/Math/interface/Error.h"
0024 #include "DataFormats/Math/interface/Point3D.h"
0025 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
0026 #include "DataFormats/GeometryCommonDetAlgo/interface/GlobalError.h"
0027 
0028 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0029 
0030 #include "HLTDisplacedtktktkFilter.h"
0031 #include "TMath.h"
0032 
0033 //
0034 // constructors and destructor
0035 //
0036 HLTDisplacedtktktkFilter::HLTDisplacedtktktkFilter(const edm::ParameterSet& iConfig)
0037     : HLTFilter(iConfig),
0038       fastAccept_(iConfig.getParameter<bool>("FastAccept")),
0039       minLxySignificance_(iConfig.getParameter<double>("MinLxySignificance")),
0040       maxLxySignificance_(iConfig.getParameter<double>("MaxLxySignificance")),
0041       maxNormalisedChi2_(iConfig.getParameter<double>("MaxNormalisedChi2")),
0042       minVtxProbability_(iConfig.getParameter<double>("MinVtxProbability")),
0043       minCosinePointingAngle_(iConfig.getParameter<double>("MinCosinePointingAngle")),
0044       triggerTypeDaughters_(iConfig.getParameter<int>("triggerTypeDaughters")),
0045       DisplacedVertexTag_(iConfig.getParameter<edm::InputTag>("DisplacedVertexTag")),
0046       DisplacedVertexToken_(consumes<reco::VertexCollection>(DisplacedVertexTag_)),
0047       beamSpotTag_(iConfig.getParameter<edm::InputTag>("BeamSpotTag")),
0048       beamSpotToken_(consumes<reco::BeamSpot>(beamSpotTag_)),
0049       TrackTag_(iConfig.getParameter<edm::InputTag>("TrackTag")),
0050       TrackToken_(consumes<reco::RecoChargedCandidateCollection>(TrackTag_)) {}
0051 
0052 HLTDisplacedtktktkFilter::~HLTDisplacedtktktkFilter() = default;
0053 
0054 void HLTDisplacedtktktkFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0055   edm::ParameterSetDescription desc;
0056   makeHLTFilterDescription(desc);
0057   desc.add<bool>("FastAccept", false);
0058   desc.add<double>("MinLxySignificance", 0.0);
0059   desc.add<double>("MaxLxySignificance", 0.0);
0060   desc.add<double>("MaxNormalisedChi2", 10.0);
0061   desc.add<double>("MinVtxProbability", 0.0);
0062   desc.add<double>("MinCosinePointingAngle", -2.0);
0063   desc.add<int>("triggerTypeDaughters", 0);
0064   desc.add<edm::InputTag>("DisplacedVertexTag", edm::InputTag("hltDisplacedtktktkVtx"));
0065   desc.add<edm::InputTag>("BeamSpotTag", edm::InputTag("hltOnlineBeamSpot"));
0066   desc.add<edm::InputTag>("TrackTag", edm::InputTag("hltL3MuonCandidates"));
0067   descriptions.add("hltDisplacedtktktkFilter", desc);
0068 }
0069 
0070 // ------------ method called once each job just before starting event loop  ------------
0071 void HLTDisplacedtktktkFilter::beginJob() {}
0072 
0073 // ------------ method called once each job just after ending the event loop  ------------
0074 void HLTDisplacedtktktkFilter::endJob() {}
0075 
0076 // ------------ method called on each new Event  ------------
0077 bool HLTDisplacedtktktkFilter::hltFilter(edm::Event& iEvent,
0078                                          const edm::EventSetup& iSetup,
0079                                          trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0080   // get beam spot
0081   reco::BeamSpot vertexBeamSpot;
0082   edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
0083   iEvent.getByToken(beamSpotToken_, recoBeamSpotHandle);
0084   vertexBeamSpot = *recoBeamSpotHandle;
0085 
0086   // get displaced vertices
0087   reco::VertexCollection displacedVertexColl;
0088   edm::Handle<reco::VertexCollection> displacedVertexCollHandle;
0089   bool foundVertexColl = iEvent.getByToken(DisplacedVertexToken_, displacedVertexCollHandle);
0090   if (foundVertexColl)
0091     displacedVertexColl = *displacedVertexCollHandle;
0092 
0093   // get track collection
0094   edm::Handle<reco::RecoChargedCandidateCollection> trackcands;
0095   iEvent.getByToken(TrackToken_, trackcands);
0096 
0097   // All HLT filters must create and fill an HLT filter object,
0098   // recording any reconstructed physics objects satisfying (or not)
0099   // this HLT filter, and place it in the Event.
0100 
0101   if (saveTags())
0102     filterproduct.addCollectionTag(TrackTag_);
0103 
0104   bool triggered = false;
0105 
0106   // loop over vertex collection
0107   for (auto const& displacedVertex : displacedVertexColl) {
0108     // check if the vertex actually consists of exactly three track tracks, reject the event if not
0109     if (displacedVertex.tracksSize() != 3) {
0110       edm::LogError("HLTDisplacedtktktkFilter") << "HLTDisplacedtktktkFilter: ERROR: the tktktk vertex must have "
0111                                                    "exactly three tracks by definition. It now has n tracks = "
0112                                                 << displacedVertex.tracksSize() << std::endl;
0113       return false;
0114     }
0115 
0116     float normChi2 = displacedVertex.normalizedChi2();
0117     if (normChi2 > maxNormalisedChi2_)
0118       continue;
0119 
0120     double vtxProb = 0.0;
0121     if ((displacedVertex.chi2() >= 0.0) && (displacedVertex.ndof() > 0))
0122       vtxProb = TMath::Prob(displacedVertex.chi2(), displacedVertex.ndof());
0123     if (vtxProb < minVtxProbability_)
0124       continue;
0125 
0126     // get the three tracks from the vertex
0127     auto trackIt = displacedVertex.tracks_begin();
0128     reco::TrackRef vertextkRef1 = (*trackIt).castTo<reco::TrackRef>();
0129     // the second one
0130     trackIt++;
0131     reco::TrackRef vertextkRef2 = (*trackIt).castTo<reco::TrackRef>();
0132     // the third one
0133     trackIt++;
0134     reco::TrackRef vertextkRef3 = (*trackIt).castTo<reco::TrackRef>();
0135 
0136     // first find these three tracks in the track collection
0137     reco::RecoChargedCandidateCollection::const_iterator cand1;
0138     reco::RecoChargedCandidateCollection::const_iterator cand2;
0139     reco::RecoChargedCandidateCollection::const_iterator cand3;
0140 
0141     int iFoundRefs = 0;
0142     for (auto cand = trackcands->begin(); cand != trackcands->end(); cand++) {
0143       reco::TrackRef tkRef = cand->get<reco::TrackRef>();
0144       if (tkRef == vertextkRef1) {
0145         cand1 = cand;
0146         iFoundRefs++;
0147       }
0148       if (tkRef == vertextkRef2) {
0149         cand2 = cand;
0150         iFoundRefs++;
0151       }
0152       if (tkRef == vertextkRef3) {
0153         cand3 = cand;
0154         iFoundRefs++;
0155       }
0156     }
0157 
0158     if (iFoundRefs < 3) {
0159       edm::LogError("HLTDisplacedtktktkFilter") << "HLTDisplacedtktktkFilter: ERROR: the tracks matched with the "
0160                                                    "tktktk vertex tracks should be at least three by definition."
0161                                                 << std::endl;
0162       return false;
0163     }
0164     // calculate three-track transverse momentum
0165     math::XYZVector pperp(cand1->px() + cand2->px() + cand3->px(), cand1->py() + cand2->py() + cand3->py(), 0.);
0166 
0167     const reco::Vertex::Point& vpoint = displacedVertex.position();
0168     reco::Vertex::Error verr = displacedVertex.error();
0169     // translate to global error, should be improved
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(
0173         -1 * ((vertexBeamSpot.x0() - vpoint.x()) + (vpoint.z() - vertexBeamSpot.z0()) * vertexBeamSpot.dxdz()),
0174         -1 * ((vertexBeamSpot.y0() - vpoint.y()) + (vpoint.z() - vertexBeamSpot.z0()) * vertexBeamSpot.dydz()),
0175         0);
0176 
0177     float lxy = displacementFromBeamspot.perp();
0178     float lxyerr = sqrt(err.rerr(displacementFromBeamspot));
0179 
0180     //calculate the angle between the decay length and the tktk momentum
0181     reco::Vertex::Point vperp(displacementFromBeamspot.x(), displacementFromBeamspot.y(), 0.);
0182 
0183     float cosAlpha = vperp.Dot(pperp) / (vperp.R() * pperp.R());
0184 
0185     // check thresholds
0186     if (cosAlpha < minCosinePointingAngle_)
0187       continue;
0188     if (minLxySignificance_ > 0. && lxy / lxyerr < minLxySignificance_)
0189       continue;
0190     if (maxLxySignificance_ > 0. && lxy / lxyerr > maxLxySignificance_)
0191       continue;
0192     triggered = true;
0193 
0194     // now add the tracks that passed to the filter object
0195 
0196     // Ref to Candidate object to be recorded in filter object
0197     reco::RecoChargedCandidateRef ref1;
0198     reco::RecoChargedCandidateRef ref2;
0199     reco::RecoChargedCandidateRef ref3;
0200 
0201     ref1 = reco::RecoChargedCandidateRef(
0202         edm::Ref<reco::RecoChargedCandidateCollection>(trackcands, distance(trackcands->begin(), cand1)));
0203     filterproduct.addObject(triggerTypeDaughters_, ref1);
0204 
0205     ref2 = reco::RecoChargedCandidateRef(
0206         edm::Ref<reco::RecoChargedCandidateCollection>(trackcands, distance(trackcands->begin(), cand2)));
0207     filterproduct.addObject(triggerTypeDaughters_, ref2);
0208 
0209     ref3 = reco::RecoChargedCandidateRef(
0210         edm::Ref<reco::RecoChargedCandidateCollection>(trackcands, distance(trackcands->begin(), cand3)));
0211     filterproduct.addObject(triggerTypeDaughters_, ref3);
0212   }
0213 
0214   LogDebug("HLTDisplacedtktktkFilter") << " >>>>> Result of HLTDisplacedtktktkFilter is " << triggered << std::endl;
0215 
0216   return triggered;
0217 }