Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:29:02

0001 #include <memory>
0002 #include <set>
0003 
0004 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
0005 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
0006 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
0007 #include "TrackingTools/IPTools/interface/IPTools.h"
0008 
0009 #include "FWCore/Framework/interface/stream/EDProducer.h"
0010 #include "FWCore/Framework/interface/Event.h"
0011 #include "FWCore/Framework/interface/MakerMacros.h"
0012 #include "FWCore/Utilities/interface/InputTag.h"
0013 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0014 #include "FWCore/Utilities/interface/ESGetToken.h"
0015 
0016 #include "DataFormats/Common/interface/Handle.h"
0017 #include "DataFormats/TrackReco/interface/Track.h"
0018 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0019 #include "DataFormats/VertexReco/interface/Vertex.h"
0020 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0021 #include "RecoVertex/VertexTools/interface/VertexDistance3D.h"
0022 
0023 #include "DataFormats/TrackReco/interface/Track.h"
0024 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0025 #include "DataFormats/VertexReco/interface/Vertex.h"
0026 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0027 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0028 
0029 #include "RecoVertex/AdaptiveVertexFit/interface/AdaptiveVertexFitter.h"
0030 #include "RecoVertex/KalmanVertexFit/interface/KalmanVertexUpdator.h"
0031 #include "RecoVertex/KalmanVertexFit/interface/KalmanVertexTrackCompatibilityEstimator.h"
0032 #include "RecoVertex/KalmanVertexFit/interface/KalmanVertexSmoother.h"
0033 #include "DataFormats/Math/interface/deltaR.h"
0034 
0035 #include "RecoVertex/ConfigurableVertexReco/interface/ConfigurableVertexReconstructor.h"
0036 #include "RecoVertex/AdaptiveVertexFinder/interface/TrackVertexArbitration.h"
0037 #include "DataFormats/Candidate/interface/VertexCompositePtrCandidate.h"
0038 
0039 #include "RecoVertex/AdaptiveVertexFinder/interface/TTHelpers.h"
0040 
0041 //#define VTXDEBUG
0042 
0043 inline const unsigned int nTracks(const reco::Vertex &sv) { return sv.nTracks(); }
0044 inline const unsigned int nTracks(const reco::VertexCompositePtrCandidate &sv) {
0045   return sv.numberOfSourceCandidatePtrs();
0046 }
0047 
0048 template <class InputContainer, class VTX>
0049 class TemplatedVertexArbitrator : public edm::stream::EDProducer<> {
0050 public:
0051   typedef std::vector<VTX> Product;
0052   TemplatedVertexArbitrator(const edm::ParameterSet &params);
0053 
0054   static void fillDescriptions(edm::ConfigurationDescriptions &cdesc) {
0055     edm::ParameterSetDescription pdesc;
0056     pdesc.add<edm::InputTag>("beamSpot", edm::InputTag("offlineBeamSpot"));
0057     pdesc.add<edm::InputTag>("primaryVertices", edm::InputTag("offlinePrimaryVertices"));
0058     if (std::is_same<VTX, reco::Vertex>::value) {
0059       pdesc.add<edm::InputTag>("tracks", edm::InputTag("generalTracks"));
0060       pdesc.add<edm::InputTag>("secondaryVertices", edm::InputTag("vertexMerger"));
0061     } else if (std::is_same<VTX, reco::VertexCompositePtrCandidate>::value) {
0062       pdesc.add<edm::InputTag>("tracks", edm::InputTag("particleFlow"));
0063       pdesc.add<edm::InputTag>("secondaryVertices", edm::InputTag("candidateVertexMerger"));
0064     } else {
0065       pdesc.add<edm::InputTag>("tracks", edm::InputTag("generalTracks"));
0066       pdesc.add<edm::InputTag>("secondaryVertices", edm::InputTag("vertexMerger"));
0067     }
0068     pdesc.add<double>("dLenFraction", 0.333);
0069     pdesc.add<double>("dRCut", 0.4);
0070     pdesc.add<double>("distCut", 0.04);
0071     pdesc.add<double>("sigCut", 5.0);
0072     pdesc.add<double>("fitterSigmacut", 3.0);
0073     pdesc.add<double>("fitterTini", 256);
0074     pdesc.add<double>("fitterRatio", 0.25);
0075     pdesc.add<int>("trackMinLayers", 4);
0076     pdesc.add<double>("trackMinPt", 0.4);
0077     pdesc.add<int>("trackMinPixels", 1);
0078     pdesc.add<double>("maxTimeSignificance", 3.5);
0079     if (std::is_same<VTX, reco::Vertex>::value) {
0080       cdesc.add("trackVertexArbitratorDefault", pdesc);
0081     } else if (std::is_same<VTX, reco::VertexCompositePtrCandidate>::value) {
0082       cdesc.add("candidateVertexArbitratorDefault", pdesc);
0083     } else {
0084       cdesc.addDefault(pdesc);
0085     }
0086   }
0087 
0088   void produce(edm::Event &event, const edm::EventSetup &es) override;
0089 
0090 private:
0091   bool trackFilter(const reco::TrackRef &track) const;
0092 
0093   edm::EDGetTokenT<reco::VertexCollection> token_primaryVertex;
0094   edm::EDGetTokenT<Product> token_secondaryVertex;
0095   edm::EDGetTokenT<InputContainer> token_tracks;
0096   edm::EDGetTokenT<reco::BeamSpot> token_beamSpot;
0097   edm::ESGetToken<TransientTrackBuilder, TransientTrackRecord> token_trackBuilder;
0098   std::unique_ptr<TrackVertexArbitration<VTX> > theArbitrator;
0099 };
0100 
0101 template <class InputContainer, class VTX>
0102 TemplatedVertexArbitrator<InputContainer, VTX>::TemplatedVertexArbitrator(const edm::ParameterSet &params) {
0103   token_primaryVertex = consumes<reco::VertexCollection>(params.getParameter<edm::InputTag>("primaryVertices"));
0104   token_secondaryVertex = consumes<Product>(params.getParameter<edm::InputTag>("secondaryVertices"));
0105   token_beamSpot = consumes<reco::BeamSpot>(params.getParameter<edm::InputTag>("beamSpot"));
0106   token_tracks = consumes<InputContainer>(params.getParameter<edm::InputTag>("tracks"));
0107   token_trackBuilder =
0108       esConsumes<TransientTrackBuilder, TransientTrackRecord>(edm::ESInputTag("", "TransientTrackBuilder"));
0109   produces<Product>();
0110   theArbitrator.reset(new TrackVertexArbitration<VTX>(params));
0111 }
0112 
0113 template <class InputContainer, class VTX>
0114 void TemplatedVertexArbitrator<InputContainer, VTX>::produce(edm::Event &event, const edm::EventSetup &es) {
0115   using namespace reco;
0116 
0117   edm::Handle<Product> secondaryVertices;
0118   event.getByToken(token_secondaryVertex, secondaryVertices);
0119   Product theSecVertexColl = *(secondaryVertices.product());
0120 
0121   edm::Handle<VertexCollection> primaryVertices;
0122   event.getByToken(token_primaryVertex, primaryVertices);
0123 
0124   auto recoVertices = std::make_unique<Product>();
0125   if (!primaryVertices->empty()) {
0126     const reco::Vertex &pv = (*primaryVertices)[0];
0127 
0128     edm::Handle<InputContainer> tracks;
0129     event.getByToken(token_tracks, tracks);
0130 
0131     edm::ESHandle<TransientTrackBuilder> trackBuilder = es.getHandle(token_trackBuilder);
0132 
0133     edm::Handle<BeamSpot> beamSpot;
0134     event.getByToken(token_beamSpot, beamSpot);
0135 
0136     std::vector<TransientTrack> selectedTracks;
0137     for (typename InputContainer::const_iterator track = tracks->begin(); track != tracks->end(); ++track) {
0138       selectedTracks.push_back(tthelpers::buildTT(tracks, trackBuilder, track - tracks->begin()));
0139     }
0140 
0141     //        const edm::RefVector< TrackCollection > tracksForArbitration= selectedTracks;
0142     Product theRecoVertices = theArbitrator->trackVertexArbitrator(beamSpot, pv, selectedTracks, theSecVertexColl);
0143 
0144     for (unsigned int ivtx = 0; ivtx < theRecoVertices.size(); ivtx++) {
0145       if (!(nTracks(theRecoVertices[ivtx]) > 1))
0146         continue;
0147       recoVertices->push_back(theRecoVertices[ivtx]);
0148     }
0149   }
0150   event.put(std::move(recoVertices));
0151 }
0152 
0153 typedef TemplatedVertexArbitrator<reco::TrackCollection, reco::Vertex> TrackVertexArbitrator;
0154 typedef TemplatedVertexArbitrator<edm::View<reco::Candidate>, reco::VertexCompositePtrCandidate>
0155     CandidateVertexArbitrator;
0156 
0157 DEFINE_FWK_MODULE(TrackVertexArbitrator);
0158 DEFINE_FWK_MODULE(CandidateVertexArbitrator);