Project
Loading...
Searching...
No Matches
checkStack.cxx
Go to the documentation of this file.
1// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
2// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3// All rights not expressly granted are reserved.
4//
5// This software is distributed under the terms of the GNU General Public
6// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7//
8// In applying this license CERN does not waive the privileges and immunities
9// granted to it by virtue of its status as an Intergovernmental Organization
10// or submit itself to any jurisdiction.
11
12// Executable to check functioning of stack
13// Analyses kinematics and track references of a kinematics file
14
18#include "DetectorsBase/Stack.h"
21#include "TFile.h"
22#include "TTree.h"
23#ifdef NDEBUG
24#undef NDEBUG
25#endif
26#include <cassert>
27#include <fairlogger/Logger.h>
30#include <unordered_map>
31
32int main(int argc, char** argv)
33{
34 const char* nameprefix = argv[1];
35
36 fair::Logger::SetConsoleSeverity("DEBUG");
37 TFile f(o2::base::NameConf::getMCKinematicsFileName(nameprefix).c_str());
38
39 LOG(debug) << "Checking input file :" << f.GetPath();
40
41 std::vector<o2::MCTrack>* mctracks = nullptr;
42 auto tr = (TTree*)f.Get("o2sim");
43 assert(tr);
44
45 auto mcbr = tr->GetBranch("MCTrack");
46 assert(mcbr);
47 mcbr->SetAddress(&mctracks);
48
49 std::vector<o2::TrackReference>* trackrefs = nullptr;
50 auto refbr = tr->GetBranch("TrackRefs");
51 assert(refbr);
52 refbr->SetAddress(&trackrefs);
53
55
56 // when present we also read some hits for ITS to test consistency of trackID assignments
58 auto hittr = (TTree*)hitf.Get("o2sim");
59 auto hitbr = hittr ? hittr->GetBranch("ITSHit") : nullptr;
60 std::vector<o2::itsmft::Hit>* hits = nullptr;
61 if (hitbr) {
62 hitbr->SetAddress(&hits);
63 }
64
65 for (int eventID = 0; eventID < mcbr->GetEntries(); ++eventID) {
66 mcbr->GetEntry(eventID);
67 refbr->GetEntry(eventID);
68 LOG(debug) << "-- Entry --" << eventID;
69 LOG(debug) << "Have " << mctracks->size() << " tracks";
70
71 std::unordered_map<int, bool> trackidsinITS_fromhits;
72 if (hitbr) {
73 hitbr->GetEntry(eventID);
74 LOG(debug) << "Have " << hits->size() << " hits";
75
76 // check that trackIDs from the hits are within range
77 int maxid = 0;
78 for (auto& h : *hits) {
79 maxid = std::max(maxid, h.GetTrackID());
80 trackidsinITS_fromhits[h.GetTrackID()] = true;
81 assert(maxid < mctracks->size());
82 }
83 }
84
85 int ti = 0;
86
87 // record tracks that left a hit in TPC
88 // (we know that these tracks should then have a TrackRef)
89 std::vector<int> trackidsinTPC;
90 std::vector<int> trackidsinITS;
91
92 // fetch the encoding of DetIDs to bits for the hit properties (can be fetched from any MCEventHeader)
93 auto& mcEventHeader = mcreader.getMCEventHeader(0, 0);
94
95 int primaries = 0;
96 int physicalprimaries = 0;
97 int secondaries = 0;
98 for (auto& t : *mctracks) {
99 // perform checks on the mass
100 if (t.GetMass() < 0) {
101 LOG(info) << "Mass not found for PDG " << t.GetPdgCode();
102 }
103
104 if (t.isSecondary()) {
105 // check that mother indices are monotonic
106 // for primaries, this may be different (for instance with Pythia8)
107 assert(ti > t.getMotherTrackId());
108 }
109
110 if (t.leftTrace(o2::detectors::DetID::TPC, mcEventHeader.getDetId2HitBitLUT())) {
111 trackidsinTPC.emplace_back(ti);
112 }
113 if (t.leftTrace(o2::detectors::DetID::ITS, mcEventHeader.getDetId2HitBitLUT())) {
114 trackidsinITS.emplace_back(ti);
115 }
116
117 bool physicalPrim = o2::mcutils::MCTrackNavigator::isPhysicalPrimary(t, *mctracks);
118 LOG(debug) << " track " << ti << "\t" << t.getMotherTrackId() << " hits " << t.hasHits() << " isPhysicalPrimary " << physicalPrim;
119 if (t.isPrimary()) {
120 primaries++;
121 } else {
122 secondaries++;
123 }
124 if (physicalPrim) {
125 physicalprimaries++;
126 }
127 ti++;
128 }
129
130 if (hitbr) {
131 assert(trackidsinITS.size() == trackidsinITS_fromhits.size());
132 for (auto id : trackidsinITS) {
133 assert(trackidsinITS_fromhits[id] == true);
134 }
135 }
136
137 LOG(debug) << "Have " << trackidsinTPC.size() << " tracks with hits in TPC";
138 LOG(debug) << "Have " << trackrefs->size() << " track refs";
139 LOG(info) << "Have " << primaries << " primaries and " << physicalprimaries << " physical primaries";
140
141 // check correct working of MCKinematicsReader
142 bool havereferences = trackrefs->size();
143 if (havereferences) {
144 for (auto& trackID : trackidsinTPC) {
145 auto tpc_trackrefs = mcreader.getTrackRefs(eventID, trackID);
146 LOG(debug) << " Track " << trackID << " has " << tpc_trackrefs.size() << " TrackRefs";
147 assert(tpc_trackrefs.size() > 0);
148 for (auto& ref : tpc_trackrefs) {
149 assert(ref.getTrackID() == trackID);
150 }
151 }
152 }
153 }
154 LOG(info) << "STACK TEST SUCCESSFULL\n";
155 return 0;
156}
Definition of the Names Generator class.
Definition of the Stack class.
Definition of the ITSMFT Hit class.
Definition of the MCTrack class.
Definition of a container to keep Monte Carlo truth external to simulation objects.
Utility functions for MC particles.
std::ostringstream debug
Class for time synchronization of RawReader instances.
static std::string getHitsFileName(DId d, const std::string_view prefix=STANDARDSIMPREFIX)
static std::string getMCKinematicsFileName(const std::string_view prefix=STANDARDSIMPREFIX)
Definition NameConf.h:46
static constexpr ID ITS
Definition DetID.h:63
static constexpr ID TPC
Definition DetID.h:64
static bool isPhysicalPrimary(o2::MCTrack const &p, std::vector< o2::MCTrack > const &pcontainer)
Definition MCUtils.cxx:73
o2::dataformats::MCEventHeader const & getMCEventHeader(int source, int event) const
retrieves the MCEventHeader for a given eventID and sourceID
gsl::span< o2::TrackReference > getTrackRefs(int source, int event, int track) const
get all primaries for a certain event
GLsizeiptr size
Definition glcorearb.h:659
GLdouble f
Definition glcorearb.h:310
#define main
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"