Skip to content

Commit 128089b

Browse files
authored
Merge pull request #47370 from AdrianoDee/simple_validation
Adding `SimpleTrackValidation` Analyzer
2 parents 51509a4 + f208846 commit 128089b

File tree

1 file changed

+188
-0
lines changed

1 file changed

+188
-0
lines changed
Lines changed: 188 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,188 @@
1+
// system include files
2+
#include <memory>
3+
4+
#include "TTree.h"
5+
#include "TFile.h"
6+
7+
// user include files
8+
#include "FWCore/Framework/interface/Frameworkfwd.h"
9+
#include "FWCore/Framework/interface/one/EDAnalyzer.h"
10+
11+
#include "FWCore/Framework/interface/Event.h"
12+
#include "FWCore/Framework/interface/MakerMacros.h"
13+
14+
#include "FWCore/ParameterSet/interface/ParameterSet.h"
15+
#include "FWCore/Utilities/interface/InputTag.h"
16+
#include "DataFormats/TrackReco/interface/Track.h"
17+
#include "DataFormats/TrackReco/interface/TrackFwd.h"
18+
19+
#include "SimTracker/TrackerHitAssociation/interface/ClusterTPAssociation.h"
20+
#include "SimDataFormats/PileupSummaryInfo/interface/PileupSummaryInfo.h"
21+
#include "SimDataFormats/Associations/interface/TrackToTrackingParticleAssociator.h"
22+
23+
#include "SimTracker/Common/interface/TrackingParticleSelector.h"
24+
25+
#include "FWCore/ServiceRegistry/interface/Service.h"
26+
#include "CommonTools/UtilAlgos/interface/TFileService.h"
27+
28+
using reco::TrackCollection;
29+
30+
class SimpleTrackValidation : public edm::one::EDAnalyzer<edm::one::SharedResources> {
31+
public:
32+
explicit SimpleTrackValidation(const edm::ParameterSet&);
33+
~SimpleTrackValidation() override = default;
34+
35+
static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
36+
37+
private:
38+
void beginJob() override;
39+
void analyze(const edm::Event&, const edm::EventSetup&) override;
40+
void endJob() override;
41+
42+
// Counters
43+
int global_rt_ = 0;
44+
int global_at_ = 0;
45+
int global_st_ = 0;
46+
int global_dt_ = 0;
47+
int global_ast_ = 0;
48+
49+
TrackingParticleSelector tpSelector;
50+
TTree* output_tree_;
51+
const std::vector<edm::InputTag> trackLabels_;
52+
std::vector<edm::EDGetTokenT<edm::View<reco::Track>>> trackTokens_;
53+
const edm::EDGetTokenT<reco::TrackToTrackingParticleAssociator> trackAssociatorToken_;
54+
const edm::EDGetTokenT<TrackingParticleCollection> trackingParticleToken_;
55+
};
56+
57+
SimpleTrackValidation::SimpleTrackValidation(const edm::ParameterSet& iConfig)
58+
: trackLabels_(iConfig.getParameter<std::vector<edm::InputTag>>("trackLabels")),
59+
trackAssociatorToken_(consumes<reco::TrackToTrackingParticleAssociator>(
60+
iConfig.getUntrackedParameter<edm::InputTag>("trackAssociator"))),
61+
trackingParticleToken_(
62+
consumes<TrackingParticleCollection>(iConfig.getParameter<edm::InputTag>("trackingParticles"))) {
63+
for (auto& itag : trackLabels_) {
64+
trackTokens_.push_back(consumes<edm::View<reco::Track>>(itag));
65+
}
66+
tpSelector = TrackingParticleSelector(iConfig.getParameter<double>("ptMinTP"),
67+
iConfig.getParameter<double>("ptMaxTP"),
68+
iConfig.getParameter<double>("minRapidityTP"),
69+
iConfig.getParameter<double>("maxRapidityTP"),
70+
iConfig.getParameter<double>("tipTP"),
71+
iConfig.getParameter<double>("lipTP"),
72+
iConfig.getParameter<int>("minHitTP"),
73+
iConfig.getParameter<bool>("signalOnlyTP"),
74+
iConfig.getParameter<bool>("intimeOnlyTP"),
75+
iConfig.getParameter<bool>("chargedOnlyTP"),
76+
iConfig.getParameter<bool>("stableOnlyTP"),
77+
iConfig.getParameter<std::vector<int>>("pdgIdTP"),
78+
iConfig.getParameter<bool>("invertRapidityCutTP"),
79+
iConfig.getParameter<double>("minPhiTP"),
80+
iConfig.getParameter<double>("maxPhiTP"));
81+
}
82+
83+
void SimpleTrackValidation::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
84+
using namespace edm;
85+
86+
auto const& associatorByHits = iEvent.get(trackAssociatorToken_);
87+
auto TPCollectionH = iEvent.getHandle(trackingParticleToken_);
88+
TrackingParticleRefVector tpCollection;
89+
90+
for (size_t i = 0, size = TPCollectionH->size(); i < size; ++i) {
91+
auto tp = TrackingParticleRef(TPCollectionH, i);
92+
if (tpSelector(*tp)) {
93+
tpCollection.push_back(tp);
94+
}
95+
}
96+
97+
for (const auto& trackToken : trackTokens_) {
98+
edm::Handle<edm::View<reco::Track>> tracksHandle;
99+
iEvent.getByToken(trackToken, tracksHandle);
100+
const edm::View<reco::Track>& tracks = *tracksHandle;
101+
102+
edm::RefToBaseVector<reco::Track> trackRefs;
103+
for (edm::View<reco::Track>::size_type i = 0; i < tracks.size(); ++i) {
104+
trackRefs.push_back(tracks.refAt(i));
105+
}
106+
107+
reco::RecoToSimCollection recSimColl = associatorByHits.associateRecoToSim(trackRefs, tpCollection);
108+
reco::SimToRecoCollection simRecColl = associatorByHits.associateSimToReco(trackRefs, tpCollection);
109+
110+
int rt = 0; // number of reconstructed tracks;
111+
int at = 0; // number of reconstructed tracks associated to a tracking particle
112+
int ast = 0; // number of tracking particles associated to at least a reconstructed track.
113+
int dt = 0; // number of duplicates (number of sim associated to more than one reco);
114+
int st = tpCollection.size(); // number of tracking particles;
115+
116+
for (const auto& track : trackRefs) {
117+
rt++;
118+
auto foundTP = recSimColl.find(track);
119+
if (foundTP != recSimColl.end()) {
120+
const auto& tp = foundTP->val;
121+
if (!tp.empty()) {
122+
at++;
123+
}
124+
if (simRecColl.find(tp[0].first) != simRecColl.end()) {
125+
if (simRecColl[tp[0].first].size() > 1) {
126+
dt++;
127+
}
128+
}
129+
}
130+
}
131+
for (const TrackingParticleRef& tpr : tpCollection) {
132+
auto foundTrack = simRecColl.find(tpr);
133+
if (foundTrack != simRecColl.end()) {
134+
ast++;
135+
}
136+
}
137+
138+
LogPrint("TrackValidator") << "Tag " << trackLabels_[0].label() << " Total simulated " << st
139+
<< " Associated tracks " << at << " Total reconstructed " << rt;
140+
global_rt_ += rt;
141+
global_st_ += st;
142+
global_at_ += at;
143+
global_dt_ += dt;
144+
global_ast_ += ast;
145+
}
146+
}
147+
148+
void SimpleTrackValidation::beginJob() {
149+
edm::Service<TFileService> fs;
150+
output_tree_ = fs->make<TTree>("output", "Simple Track Validation TTree");
151+
152+
output_tree_->Branch("rt", &global_rt_);
153+
output_tree_->Branch("at", &global_at_);
154+
output_tree_->Branch("st", &global_st_);
155+
output_tree_->Branch("dt", &global_dt_);
156+
output_tree_->Branch("ast", &global_ast_);
157+
}
158+
159+
void SimpleTrackValidation::endJob() { output_tree_->Fill(); }
160+
161+
void SimpleTrackValidation::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
162+
edm::ParameterSetDescription desc;
163+
164+
desc.add<std::vector<edm::InputTag>>("trackLabels", {edm::InputTag("generalTracks")});
165+
desc.add<edm::InputTag>("trackingParticles", edm::InputTag("mix", "MergedTrackTruth"));
166+
desc.add<edm::InputTag>("trackAssociator", edm::InputTag("trackingParticleRecoTrackAsssociation"));
167+
168+
// TP Selector parameters
169+
desc.add<double>("ptMinTP", 0.9);
170+
desc.add<double>("ptMaxTP", 1e100);
171+
desc.add<double>("minRapidityTP", -2.4);
172+
desc.add<double>("maxRapidityTP", 2.4);
173+
desc.add<double>("tipTP", 3.5);
174+
desc.add<double>("lipTP", 30.0);
175+
desc.add<int>("minHitTP", 0);
176+
desc.add<bool>("signalOnlyTP", true);
177+
desc.add<bool>("intimeOnlyTP", false);
178+
desc.add<bool>("chargedOnlyTP", true);
179+
desc.add<bool>("stableOnlyTP", false);
180+
desc.add<std::vector<int>>("pdgIdTP", {});
181+
desc.add<bool>("invertRapidityCutTP", false);
182+
desc.add<double>("minPhiTP", -3.2);
183+
desc.add<double>("maxPhiTP", 3.2);
184+
185+
descriptions.addWithDefaultLabel(desc);
186+
}
187+
188+
DEFINE_FWK_MODULE(SimpleTrackValidation);

0 commit comments

Comments
 (0)