File indexing completed on 2024-04-06 12:31:02
0001 #ifndef SimTracker_TrackAssociatorProducers_TrackAssociatorByHitsImpl_h
0002 #define SimTracker_TrackAssociatorProducers_TrackAssociatorByHitsImpl_h
0003
0004 #include "DataFormats/Common/interface/Handle.h"
0005 #include "DataFormats/Common/interface/Ref.h"
0006 #include "SimTracker/TrackerHitAssociation/interface/TrackerHitAssociator.h"
0007
0008
0009 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0010 #include "DataFormats/TrackReco/interface/Track.h"
0011 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
0012
0013 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h"
0014 #include "SimDataFormats/EncodedEventId/interface/EncodedEventId.h"
0015 #include "SimDataFormats/Associations/interface/TrackToTrackingParticleAssociatorBaseImpl.h"
0016
0017 class TrackerTopology;
0018
0019 namespace edm {
0020 class EDProductGetter;
0021 }
0022
0023 class TrackAssociatorByHitsImpl : public reco::TrackToTrackingParticleAssociatorBaseImpl {
0024 public:
0025 enum SimToRecoDenomType { denomnone, denomsim, denomreco };
0026
0027 typedef std::pair<TrackingParticleRef, TrackPSimHitRef> SimHitTPPair;
0028 typedef std::vector<SimHitTPPair> SimHitTPAssociationList;
0029
0030 TrackAssociatorByHitsImpl(edm::EDProductGetter const& productGetter,
0031 std::unique_ptr<TrackerHitAssociator> iAssociate,
0032 TrackerTopology const* iTopo,
0033 SimHitTPAssociationList const* iSimHitsTPAssoc,
0034 SimToRecoDenomType iSimToRecoDenominator,
0035 double iQuality_SimToReco,
0036 double iPurity_SimToReco,
0037 double iCut_RecoToSim,
0038 bool iUsePixels,
0039 bool iUseGrouped,
0040 bool iUseSplitting,
0041 bool ThreeHitTracksAreSpecial,
0042 bool AbsoluteNumberOfHits);
0043
0044
0045
0046
0047
0048 reco::RecoToSimCollection associateRecoToSim(const edm::RefToBaseVector<reco::Track>&,
0049 const edm::RefVector<TrackingParticleCollection>&) const override;
0050
0051
0052 reco::SimToRecoCollection associateSimToReco(const edm::RefToBaseVector<reco::Track>&,
0053 const edm::RefVector<TrackingParticleCollection>&) const override;
0054
0055
0056
0057 reco::RecoToSimCollection associateRecoToSim(const edm::Handle<edm::View<reco::Track> >& tCH,
0058 const edm::Handle<TrackingParticleCollection>& tPCH) const override {
0059 return TrackToTrackingParticleAssociatorBaseImpl::associateRecoToSim(tCH, tPCH);
0060 }
0061
0062
0063
0064 reco::SimToRecoCollection associateSimToReco(const edm::Handle<edm::View<reco::Track> >& tCH,
0065 const edm::Handle<TrackingParticleCollection>& tPCH) const override {
0066 return TrackToTrackingParticleAssociatorBaseImpl::associateSimToReco(tCH, tPCH);
0067 }
0068
0069
0070
0071 reco::RecoToSimCollectionSeed associateRecoToSim(const edm::Handle<edm::View<TrajectorySeed> >&,
0072 const edm::Handle<TrackingParticleCollection>&) const override;
0073
0074 reco::SimToRecoCollectionSeed associateSimToReco(const edm::Handle<edm::View<TrajectorySeed> >&,
0075 const edm::Handle<TrackingParticleCollection>&) const override;
0076
0077 private:
0078 template <typename iter>
0079 void getMatchedIds(std::vector<SimHitIdpr>&, std::vector<SimHitIdpr>&, int&, iter, iter, TrackerHitAssociator*) const;
0080
0081 int getShared(std::vector<SimHitIdpr>&, std::vector<SimHitIdpr>&, TrackingParticle const&) const;
0082
0083 template <typename iter>
0084 int getDoubleCount(iter, iter, TrackerHitAssociator*, TrackingParticle const&) const;
0085
0086
0087 edm::EDProductGetter const* productGetter_;
0088 std::unique_ptr<TrackerHitAssociator> associate;
0089 TrackerTopology const* tTopo;
0090 SimHitTPAssociationList const* simHitsTPAssoc;
0091
0092 SimToRecoDenomType SimToRecoDenominator;
0093 const double quality_SimToReco;
0094 const double purity_SimToReco;
0095 const double cut_RecoToSim;
0096 const bool UsePixels;
0097 const bool UseGrouped;
0098 const bool UseSplitting;
0099 const bool ThreeHitTracksAreSpecial;
0100 const bool AbsoluteNumberOfHits;
0101
0102 const TrackingRecHit* getHitPtr(edm::OwnVector<TrackingRecHit>::const_iterator iter) const { return &*iter; }
0103 const TrackingRecHit* getHitPtr(trackingRecHit_iterator iter) const { return &**iter; }
0104
0105
0106 };
0107
0108 #endif