File indexing completed on 2024-04-06 12:31:06
0001
0002
0003
0004
0005
0006
0007
0008 #include "TH1F.h"
0009
0010 #include "FWCore/Framework/interface/Frameworkfwd.h"
0011 #include "FWCore/Framework/interface/MakerMacros.h"
0012 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0013 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0014 #include "FWCore/ServiceRegistry/interface/Service.h"
0015
0016 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0017
0018 #include "SimTracker/TrackHistory/interface/TrackClassifier.h"
0019
0020
0021
0022
0023
0024 class TrackingParticleCategoriesAnalyzer : public edm::one::EDAnalyzer<edm::one::SharedResources> {
0025 public:
0026 explicit TrackingParticleCategoriesAnalyzer(const edm::ParameterSet &);
0027 ~TrackingParticleCategoriesAnalyzer() override;
0028
0029 private:
0030 void analyze(const edm::Event &, const edm::EventSetup &) override;
0031
0032
0033
0034 edm::EDGetTokenT<TrackingParticleCollection> trackingTruth_;
0035
0036 std::size_t totalTrakingParticles_;
0037
0038 TrackClassifier classifier_;
0039
0040 TH1F *trackingParticleCategories_;
0041
0042 Int_t numberTrackingParticleCategories_;
0043 };
0044
0045 TrackingParticleCategoriesAnalyzer::TrackingParticleCategoriesAnalyzer(const edm::ParameterSet &config)
0046 : classifier_(config, consumesCollector()) {
0047
0048 trackingTruth_ = consumes<TrackingParticleCollection>(config.getUntrackedParameter<edm::InputTag>("trackingTruth"));
0049
0050
0051 usesResource("TFileService");
0052 edm::Service<TFileService> fs;
0053
0054
0055 TFileDirectory directory = fs->mkdir("TrackingParticleCategoriesAnalyzer");
0056
0057
0058 numberTrackingParticleCategories_ = TrackCategories::Unknown + 1;
0059
0060
0061 trackingParticleCategories_ = fs->make<TH1F>("Frequency",
0062 "Frequency for the different track categories",
0063 numberTrackingParticleCategories_,
0064 -0.5,
0065 numberTrackingParticleCategories_ - 0.5);
0066
0067
0068 for (Int_t i = 0; i < numberTrackingParticleCategories_; ++i)
0069 trackingParticleCategories_->GetXaxis()->SetBinLabel(i + 1, TrackCategories::Names[i]);
0070 }
0071
0072 TrackingParticleCategoriesAnalyzer::~TrackingParticleCategoriesAnalyzer() {}
0073
0074 void TrackingParticleCategoriesAnalyzer::analyze(const edm::Event &event, const edm::EventSetup &setup) {
0075
0076 edm::Handle<TrackingParticleCollection> TPCollection;
0077 event.getByToken(trackingTruth_, TPCollection);
0078
0079
0080 classifier_.newEvent(event, setup);
0081
0082
0083 for (std::size_t index = 0; index < TPCollection->size(); index++) {
0084 TrackingParticleRef trackingParticle(TPCollection, index);
0085
0086
0087 classifier_.evaluate(trackingParticle);
0088
0089
0090 for (Int_t i = 0; i != numberTrackingParticleCategories_; ++i)
0091 if (classifier_.is((TrackCategories::Category)i))
0092 trackingParticleCategories_->Fill(i);
0093 }
0094 }
0095
0096 DEFINE_FWK_MODULE(TrackingParticleCategoriesAnalyzer);