00001 #ifndef MuonAssociatorByHits_h
00002 #define MuonAssociatorByHits_h
00003
00004 #include "FWCore/Framework/interface/Frameworkfwd.h"
00005 #include "FWCore/Framework/interface/Event.h"
00006 #include "FWCore/Framework/interface/EventSetup.h"
00007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00008 #include "FWCore/Utilities/interface/InputTag.h"
00009 #include "DataFormats/Common/interface/Handle.h"
00010 #include "DataFormats/Common/interface/Ref.h"
00011 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00012 #include "DataFormats/TrackReco/interface/Track.h"
00013 #include "DataFormats/MuonReco/interface/Muon.h"
00014 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
00015 #include "DataFormats/RecoCandidate/interface/TrackAssociation.h"
00016 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h"
00017 #include "SimMuon/MCTruth/interface/DTHitAssociator.h"
00018 #include "SimMuon/MCTruth/interface/MuonTruth.h"
00019 #include "SimMuon/MCTruth/interface/RPCHitAssociator.h"
00020 #include "SimTracker/TrackerHitAssociation/interface/TrackerHitAssociator.h"
00021 #include "SimTracker/TrackAssociation/interface/TrackAssociatorBase.h"
00022
00023
00024 #include <boost/ptr_container/ptr_vector.hpp>
00025
00026 class MuonAssociatorByHits : public TrackAssociatorBase {
00027
00028 public:
00029 typedef std::pair <uint32_t, EncodedEventId> SimHitIdpr;
00030
00031 typedef std::pair<unsigned int,std::vector<SimHitIdpr> > uint_SimHitIdpr_pair;
00032 typedef boost::ptr_vector<uint_SimHitIdpr_pair> MapOfMatchedIds;
00033
00034 MuonAssociatorByHits( const edm::ParameterSet& );
00035 ~MuonAssociatorByHits();
00036
00037
00038 using TrackAssociatorBase::associateRecoToSim;
00039 using TrackAssociatorBase::associateSimToReco;
00040
00041
00043 reco::RecoToSimCollection associateRecoToSim(const edm::RefToBaseVector<reco::Track>&,
00044 const edm::RefVector<TrackingParticleCollection>&,
00045 const edm::Event * event = 0, const edm::EventSetup * setup = 0) const;
00046
00048 reco::SimToRecoCollection associateSimToReco(const edm::RefToBaseVector<reco::Track>&,
00049 const edm::RefVector<TrackingParticleCollection>&,
00050 const edm::Event * event = 0, const edm::EventSetup * setup = 0) const ;
00051
00052
00053 void getMatchedIds
00054 (MapOfMatchedIds & tracker_matchedIds_valid, MapOfMatchedIds & muon_matchedIds_valid,
00055 MapOfMatchedIds & tracker_matchedIds_INVALID, MapOfMatchedIds & muon_matchedIds_INVALID,
00056 int& n_tracker_valid, int& n_dt_valid, int& n_csc_valid, int& n_rpc_valid,
00057 int& n_tracker_matched_valid, int& n_dt_matched_valid, int& n_csc_matched_valid, int& n_rpc_matched_valid,
00058 int& n_tracker_INVALID, int& n_dt_INVALID, int& n_csc_INVALID, int& n_rpc_INVALID,
00059 int& n_tracker_matched_INVALID, int& n_dt_matched_INVALID, int& n_csc_matched_INVALID, int& n_rpc_matched_INVALID,
00060 trackingRecHit_iterator begin, trackingRecHit_iterator end,
00061 TrackerHitAssociator* trackertruth, DTHitAssociator& dttruth, MuonTruth& csctruth, RPCHitAssociator& rpctruth,
00062 bool printRts) const;
00063
00064 int getShared(MapOfMatchedIds & matchedIds, TrackingParticleCollection::const_iterator trpart) const;
00065
00066
00067 enum MuonTrackType { InnerTk, OuterTk, GlobalTk, Segments };
00068 struct RefToBaseSort {
00069 template<typename T> bool operator()(const edm::RefToBase<T> &r1, const edm::RefToBase<T> &r2) const {
00070 return (r1.id() == r2.id() ? r1.key() < r2.key() : r1.id() < r2.id());
00071 }
00072 };
00073 typedef std::map<edm::RefToBase<reco::Muon>, std::vector<std::pair<TrackingParticleRef, double> >, RefToBaseSort> MuonToSimCollection;
00074 typedef std::map<TrackingParticleRef, std::vector<std::pair<edm::RefToBase<reco::Muon>, double> > > SimToMuonCollection;
00075
00076
00077 void associateMuons(MuonToSimCollection & recoToSim, SimToMuonCollection & simToReco,
00078 const edm::RefToBaseVector<reco::Muon> &, MuonTrackType ,
00079 const edm::RefVector<TrackingParticleCollection>&,
00080 const edm::Event * event = 0, const edm::EventSetup * setup = 0) const ;
00081
00082 void associateMuons(MuonToSimCollection & recoToSim, SimToMuonCollection & simToReco,
00083 const edm::Handle<edm::View<reco::Muon> > &, MuonTrackType ,
00084 const edm::Handle<TrackingParticleCollection>&,
00085 const edm::Event * event = 0, const edm::EventSetup * setup = 0) const ;
00086
00087 private:
00088 const bool includeZeroHitMuons;
00089 const bool acceptOneStubMatchings;
00090 bool UseTracker;
00091 bool UseMuon;
00092 const bool AbsoluteNumberOfHits_track;
00093 unsigned int NHitCut_track;
00094 double EfficiencyCut_track;
00095 double PurityCut_track;
00096 const bool AbsoluteNumberOfHits_muon;
00097 unsigned int NHitCut_muon;
00098 double EfficiencyCut_muon;
00099 double PurityCut_muon;
00100 const bool UsePixels;
00101 const bool UseGrouped;
00102 const bool UseSplitting;
00103 const bool ThreeHitTracksAreSpecial;
00104 const bool dumpDT;
00105 const bool dumpInputCollections;
00106 const bool crossingframe;
00107 edm::InputTag simtracksTag;
00108 edm::InputTag simtracksXFTag;
00109 const edm::ParameterSet& conf_;
00110
00111 int LayerFromDetid(const DetId&) const;
00112 const TrackingRecHit* getHitPtr(edm::OwnVector<TrackingRecHit>::const_iterator iter) const {return &*iter;}
00113 const TrackingRecHit* getHitPtr(trackingRecHit_iterator iter) const {return &**iter;}
00114
00115 std::string write_matched_simtracks(const std::vector<SimHitIdpr>&) const;
00116
00117
00118 typedef std::vector<std::pair<trackingRecHit_iterator, trackingRecHit_iterator> > TrackHitsCollection;
00119 struct IndexMatch {
00120 IndexMatch(size_t index, double global_quality) : idx(index), quality(global_quality) {}
00121 size_t idx; double quality;
00122 bool operator<(const IndexMatch &other) const { return other.quality < quality; }
00123 };
00124 typedef std::map<size_t, std::vector<IndexMatch> > IndexAssociation;
00125
00126 IndexAssociation associateSimToRecoIndices(const TrackHitsCollection &,
00127 const edm::RefVector<TrackingParticleCollection>&,
00128 const edm::Event * event = 0, const edm::EventSetup * setup = 0) const ;
00129 IndexAssociation associateRecoToSimIndices(const TrackHitsCollection &,
00130 const edm::RefVector<TrackingParticleCollection>&,
00131 const edm::Event * event = 0, const edm::EventSetup * setup = 0) const ;
00132
00133
00134
00135 };
00136
00137 #endif