File indexing completed on 2024-04-06 12:31:01
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020 #include <memory>
0021
0022
0023 #include "FWCore/Framework/interface/Frameworkfwd.h"
0024 #include "FWCore/Framework/interface/global/EDProducer.h"
0025
0026 #include "FWCore/Framework/interface/Event.h"
0027 #include "FWCore/Framework/interface/MakerMacros.h"
0028
0029 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0030 #include "FWCore/Utilities/interface/EDGetToken.h"
0031
0032 #include "SimDataFormats/Associations/interface/TrackToTrackingParticleAssociator.h"
0033 #include "QuickTrackAssociatorByHitsImpl.h"
0034 #include "SimTracker/TrackerHitAssociation/interface/TrackerHitAssociator.h"
0035
0036
0037
0038
0039 namespace {}
0040
0041 class QuickTrackAssociatorByHitsProducer : public edm::global::EDProducer<> {
0042 public:
0043 explicit QuickTrackAssociatorByHitsProducer(const edm::ParameterSet&);
0044 ~QuickTrackAssociatorByHitsProducer() override;
0045
0046 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0047
0048 private:
0049 void beginJob() override;
0050 void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
0051 void endJob() override;
0052 edm::ParameterSet makeHitAssociatorParameters(const edm::ParameterSet&);
0053
0054
0055 TrackerHitAssociator::Config trackerHitAssociatorConfig_;
0056 edm::EDGetTokenT<ClusterTPAssociation> cluster2TPToken_;
0057 double qualitySimToReco_;
0058 double puritySimToReco_;
0059 double pixelHitWeight_;
0060 double cutRecoToSim_;
0061 QuickTrackAssociatorByHitsImpl::SimToRecoDenomType simToRecoDenominator_;
0062 bool threeHitTracksAreSpecial_;
0063 bool useClusterTPAssociation_;
0064 bool absoluteNumberOfHits_;
0065 };
0066
0067
0068
0069
0070
0071
0072
0073
0074
0075
0076
0077
0078 QuickTrackAssociatorByHitsProducer::QuickTrackAssociatorByHitsProducer(const edm::ParameterSet& iConfig)
0079 : qualitySimToReco_(iConfig.getParameter<double>("Quality_SimToReco")),
0080 puritySimToReco_(iConfig.getParameter<double>("Purity_SimToReco")),
0081 pixelHitWeight_(iConfig.getParameter<double>("PixelHitWeight")),
0082 cutRecoToSim_(iConfig.getParameter<double>("Cut_RecoToSim")),
0083 threeHitTracksAreSpecial_(iConfig.getParameter<bool>("ThreeHitTracksAreSpecial")),
0084 useClusterTPAssociation_(iConfig.getParameter<bool>("useClusterTPAssociation")),
0085 absoluteNumberOfHits_(iConfig.getParameter<bool>("AbsoluteNumberOfHits")) {
0086
0087
0088
0089
0090 std::string denominatorString = iConfig.getParameter<std::string>("SimToRecoDenominator");
0091 if (denominatorString == "sim")
0092 simToRecoDenominator_ = QuickTrackAssociatorByHitsImpl::denomsim;
0093 else if (denominatorString == "reco")
0094 simToRecoDenominator_ = QuickTrackAssociatorByHitsImpl::denomreco;
0095 else
0096 throw cms::Exception("QuickTrackAssociatorByHitsImpl") << "SimToRecoDenominator not specified as sim or reco";
0097
0098
0099
0100
0101
0102 bool useGrouped, useSplitting;
0103 if (iConfig.exists("UseGrouped"))
0104 useGrouped = iConfig.getParameter<bool>("UseGrouped");
0105 else
0106 useGrouped = true;
0107
0108 if (iConfig.exists("UseSplitting"))
0109 useSplitting = iConfig.getParameter<bool>("UseSplitting");
0110 else
0111 useSplitting = true;
0112
0113
0114
0115 if (!(useGrouped && useSplitting)) {
0116 LogDebug("QuickTrackAssociatorByHitsImpl")
0117 << "UseGrouped and/or UseSplitting has been set to false, but this associator ignores that setting.";
0118 }
0119
0120
0121 produces<reco::TrackToTrackingParticleAssociator>();
0122
0123 if (useClusterTPAssociation_) {
0124 cluster2TPToken_ = consumes<ClusterTPAssociation>(iConfig.getParameter<edm::InputTag>("cluster2TPSrc"));
0125 } else {
0126 trackerHitAssociatorConfig_ =
0127 TrackerHitAssociator::Config(makeHitAssociatorParameters(iConfig), consumesCollector());
0128 }
0129 }
0130
0131 QuickTrackAssociatorByHitsProducer::~QuickTrackAssociatorByHitsProducer() {
0132
0133
0134 }
0135
0136
0137
0138
0139
0140
0141 edm::ParameterSet QuickTrackAssociatorByHitsProducer::makeHitAssociatorParameters(const edm::ParameterSet& iConfig) {
0142 edm::ParameterSet hitAssociatorParameters;
0143 hitAssociatorParameters.addParameter<bool>("associatePixel", iConfig.getParameter<bool>("associatePixel"));
0144 hitAssociatorParameters.addParameter<bool>("associateStrip", iConfig.getParameter<bool>("associateStrip"));
0145
0146
0147
0148 hitAssociatorParameters.addParameter<bool>("associateRecoTracks", true);
0149
0150 hitAssociatorParameters.addParameter<edm::InputTag>("pixelSimLinkSrc",
0151 iConfig.getParameter<edm::InputTag>("pixelSimLinkSrc"));
0152 hitAssociatorParameters.addParameter<edm::InputTag>("stripSimLinkSrc",
0153 iConfig.getParameter<edm::InputTag>("stripSimLinkSrc"));
0154
0155 return hitAssociatorParameters;
0156 }
0157
0158
0159 void QuickTrackAssociatorByHitsProducer::produce(edm::StreamID,
0160 edm::Event& iEvent,
0161 const edm::EventSetup& iSetup) const {
0162 using namespace edm;
0163
0164 const ClusterTPAssociation* clusterAssoc = nullptr;
0165 std::unique_ptr<TrackerHitAssociator> trackAssoc;
0166 if (useClusterTPAssociation_) {
0167 edm::Handle<ClusterTPAssociation> clusterAssocHandle;
0168 iEvent.getByToken(cluster2TPToken_, clusterAssocHandle);
0169 clusterAssoc = clusterAssocHandle.product();
0170 } else {
0171
0172
0173 trackAssoc = std::make_unique<TrackerHitAssociator>(iEvent, trackerHitAssociatorConfig_);
0174 }
0175
0176 auto impl = std::make_unique<QuickTrackAssociatorByHitsImpl>(iEvent.productGetter(),
0177 std::move(trackAssoc),
0178 clusterAssoc,
0179 absoluteNumberOfHits_,
0180 qualitySimToReco_,
0181 puritySimToReco_,
0182 pixelHitWeight_,
0183 cutRecoToSim_,
0184 threeHitTracksAreSpecial_,
0185 simToRecoDenominator_);
0186
0187 auto toPut = std::make_unique<reco::TrackToTrackingParticleAssociator>(std::move(impl));
0188 iEvent.put(std::move(toPut));
0189 }
0190
0191
0192 void QuickTrackAssociatorByHitsProducer::beginJob() {}
0193
0194
0195 void QuickTrackAssociatorByHitsProducer::endJob() {}
0196
0197
0198 void QuickTrackAssociatorByHitsProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0199
0200
0201 edm::ParameterSetDescription desc;
0202 desc.setUnknown();
0203 descriptions.addDefault(desc);
0204 }
0205
0206
0207 DEFINE_FWK_MODULE(QuickTrackAssociatorByHitsProducer);