Project
Loading...
Searching...
No Matches
TimeFrame.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.
15
17#include "TRKSimulation/Hit.h"
19#include "Framework/Logger.h"
21#include <TTree.h>
22#include <TRandom3.h>
23#include <vector>
24#include <array>
25
26namespace o2::trk
27{
28
29template <int nLayers>
30int TimeFrame<nLayers>::loadROFsFromHitTree(TTree* hitsTree, GeometryTGeo* gman, const nlohmann::json& config)
31{
32 constexpr std::array<int, 2> startLayer{0, 3};
33 const Long64_t nEvents = hitsTree->GetEntries();
34
35 gman->fillMatrixCache(o2::math_utils::bit2Mask(o2::math_utils::TransformType::T2L) | o2::math_utils::bit2Mask(o2::math_utils::TransformType::L2G));
36
37 std::vector<o2::trk::Hit>* trkHit = nullptr;
38 hitsTree->SetBranchAddress("TRKHit", &trkHit);
39
40 const int inROFpileup{config.contains("inROFpileup") ? config["inROFpileup"].get<int>() : 1};
41
42 // Calculate number of ROFs and initialize data structures
43 this->mNrof = (nEvents + inROFpileup - 1) / inROFpileup;
44
45 // Reset and prepare ROF data structures
46 for (int iLayer{0}; iLayer < nLayers; ++iLayer) {
47 this->mMinR[iLayer] = std::numeric_limits<float>::max();
48 this->mMaxR[iLayer] = std::numeric_limits<float>::lowest();
49 this->mROFramesClusters[iLayer].clear();
50 this->mROFramesClusters[iLayer].resize(this->mNrof + 1, 0);
51 this->mUnsortedClusters[iLayer].clear();
52 this->mTrackingFrameInfo[iLayer].clear();
53 this->mClusterExternalIndices[iLayer].clear();
54 }
55
56 // Pre-count hits to reserve memory efficiently
57 int totalNHits{0};
58 std::array<int, nLayers> clusterCountPerLayer{};
59 for (Long64_t iEvent = 0; iEvent < nEvents; ++iEvent) {
60 hitsTree->GetEntry(iEvent);
61 for (const auto& hit : *trkHit) {
62 if (gman->getDisk(hit.GetDetectorID()) != -1) {
63 continue; // skip non-barrel hits
64 }
65 int subDetID = gman->getSubDetID(hit.GetDetectorID());
66 const int layer = startLayer[subDetID] + gman->getLayer(hit.GetDetectorID());
67 if (layer >= nLayers) {
68 continue;
69 }
70 ++clusterCountPerLayer[layer];
71 totalNHits++;
72 }
73 }
74
75 // Reserve memory for all layers
76 for (int iLayer{0}; iLayer < nLayers; ++iLayer) {
77 this->mUnsortedClusters[iLayer].reserve(clusterCountPerLayer[iLayer]);
78 this->mTrackingFrameInfo[iLayer].reserve(clusterCountPerLayer[iLayer]);
79 this->mClusterExternalIndices[iLayer].reserve(clusterCountPerLayer[iLayer]);
80 }
81 clearResizeBoundedVector(this->mClusterSize, totalNHits, this->mMemoryPool.get());
82
83 std::array<float, 11> resolution{0.001, 0.001, 0.001, 0.001, 0.004, 0.004, 0.004, 0.004, 0.004, 0.004, 0.004};
84 if (config["geometry"]["pitch"].size() == nLayers) {
85 for (int iLayer{0}; iLayer < config["geometry"]["pitch"].size(); ++iLayer) {
86 LOGP(info, "Setting resolution for layer {} from config", iLayer);
87 LOGP(info, "Layer {} pitch {} cm", iLayer, config["geometry"]["pitch"][iLayer].get<float>());
88 resolution[iLayer] = config["geometry"]["pitch"][iLayer].get<float>() / std::sqrt(12.f);
89 }
90 }
91 LOGP(info, "Number of active parts in VD: {}", gman->getNumberOfActivePartsVD());
92
93 int hitCounter{0};
95
96 int iRof{0}; // Current ROF index
97 for (Long64_t iEvent = 0; iEvent < nEvents; ++iEvent) {
98 hitsTree->GetEntry(iEvent);
99
100 for (auto& hit : *trkHit) {
101 if (gman->getDisk(hit.GetDetectorID()) != -1) {
102 continue; // skip non-barrel hits for this test
103 }
104 int subDetID = gman->getSubDetID(hit.GetDetectorID());
105 const int layer = startLayer[subDetID] + gman->getLayer(hit.GetDetectorID());
106
107 float alpha{0.f};
110 float r{0.f};
111 if (layer >= nLayers) {
112 continue;
113 }
114 if (layer >= 3) {
115 int chipID = hit.GetDetectorID();
116 alpha = gman->getSensorRefAlphaMLOT(chipID);
117 const o2::math_utils::Transform3D& l2g = gman->getMatrixL2G(chipID);
118 auto locXYZ = l2g ^ (hit.GetPos());
119 locXYZ.SetX(locXYZ.X() + gRandom->Gaus(0.0, resolution[layer]));
120 locXYZ.SetZ(locXYZ.Z() + gRandom->Gaus(0.0, resolution[layer]));
121 gloXYZ = gman->getMatrixL2G(chipID) * locXYZ;
122 trkXYZ = gman->getMatrixT2L(chipID - gman->getNumberOfActivePartsVD()) ^ locXYZ;
123 r = std::hypot(gloXYZ.X(), gloXYZ.Y());
124 } else {
125 const auto& hitPos = hit.GetPos();
126 r = std::hypot(hitPos.X(), hitPos.Y());
127 alpha = std::atan2(hitPos.Y(), hitPos.X()) + gRandom->Gaus(0.0, resolution[layer] / r);
128 o2::math_utils::bringTo02Pi(alpha);
129 gloXYZ.SetX(r * std::cos(alpha));
130 gloXYZ.SetY(r * std::sin(alpha));
131 gloXYZ.SetZ(hitPos.Z() + gRandom->Gaus(0.0, resolution[layer]));
132 trkXYZ.SetX(r);
133 trkXYZ.SetY(0.f);
134 trkXYZ.SetZ(gloXYZ.Z());
135 }
136 this->mMinR[layer] = std::min(this->mMinR[layer], r);
137 this->mMaxR[layer] = std::max(this->mMaxR[layer], r);
138 this->addTrackingFrameInfoToLayer(layer, gloXYZ.x(), gloXYZ.y(), gloXYZ.z(), trkXYZ.x(), alpha,
139 std::array<float, 2>{trkXYZ.y(), trkXYZ.z()},
140 std::array<float, 3>{resolution[layer] * resolution[layer], 0., resolution[layer] * resolution[layer]});
142 this->addClusterToLayer(layer, gloXYZ.x(), gloXYZ.y(), gloXYZ.z(), this->mUnsortedClusters[layer].size());
143 this->addClusterExternalIndexToLayer(layer, hitCounter);
144 MCCompLabel label{hit.GetTrackID(), static_cast<int>(iEvent), 0};
145 labels->addElement(hitCounter, label);
146 this->mClusterSize[hitCounter] = 1; // For compatibility with cluster-based tracking, set cluster size to 1 for hits
147 hitCounter++;
148 }
149 trkHit->clear();
150
151 // Update ROF structure when we complete an ROF or reach the last event
152 if ((iEvent + 1) % inROFpileup == 0 || iEvent == nEvents - 1) {
153 iRof++;
154 for (unsigned int iLayer{0}; iLayer < this->mUnsortedClusters.size(); ++iLayer) {
155 this->mROFramesClusters[iLayer][iRof] = this->mUnsortedClusters[iLayer].size(); // effectively calculating an exclusive sum
156 }
157 // Update primary vertices ROF structure
158 }
159 this->mClusterLabels = labels;
160 }
161 return this->mNrof;
162}
163
164template <int nLayers>
165void TimeFrame<nLayers>::getPrimaryVerticesFromMC(TTree* mcHeaderTree, int nRofs, Long64_t nEvents, int inROFpileup)
166{
167 auto mcheader = new o2::dataformats::MCEventHeader;
168 mcHeaderTree->SetBranchAddress("MCEventHeader.", &mcheader);
169
170 this->mROFramesPV.clear();
171 this->mROFramesPV.resize(nRofs + 1, 0);
172 this->mPrimaryVertices.clear();
173
174 int iRof{0};
175 for (Long64_t iEvent = 0; iEvent < nEvents; ++iEvent) {
176 mcHeaderTree->GetEntry(iEvent);
178 vertex.setXYZ(mcheader->GetX(), mcheader->GetY(), mcheader->GetZ());
179 vertex.setNContributors(30);
180 vertex.setChi2(0.f);
181 LOGP(debug, "ROF {}: Added primary vertex at ({}, {}, {})", iRof, mcheader->GetX(), mcheader->GetY(), mcheader->GetZ());
182 this->mPrimaryVertices.push_back(vertex);
183 if ((iEvent + 1) % inROFpileup == 0 || iEvent == nEvents - 1) {
184 iRof++;
185 this->mROFramesPV[iRof] = this->mPrimaryVertices.size(); // effectively calculating an exclusive sum
186 }
187 }
188 this->mMultiplicityCutMask.resize(nRofs, true);
189}
190
191// Explicit template instantiation for TRK with 11 layers
192template class TimeFrame<11>;
193
194} // namespace o2::trk
std::vector< std::string > labels
TRK TimeFrame class derived from ITS TimeFrame.
Definition of the TRK Hit class.
std::ostringstream debug
uint64_t vertex
Definition RawEventData.h:9
A container to hold and manage MC truth information/labels.
const Mat3D & getMatrixT2L(int sensID) const
const Mat3D & getMatrixL2G(int sensID) const
int getSubDetID(int index) const
float getSensorRefAlphaMLOT(int index) const
int getLayer(int index) const
void fillMatrixCache(int mask)
int getNumberOfActivePartsVD() const
int getDisk(int index) const
int loadROFsFromHitTree(TTree *hitsTree, GeometryTGeo *gman, const nlohmann::json &config)
Definition TimeFrame.cxx:30
void getPrimaryVerticesFromMC(TTree *mcHeaderTree, int nRofs, Long64_t nEvents, int inROFpileup)
GLfloat GLfloat GLfloat alpha
Definition glcorearb.h:279
GLsizeiptr size
Definition glcorearb.h:659
GLuint GLsizei const GLchar * label
Definition glcorearb.h:2519
GLenum GLuint GLint GLint layer
Definition glcorearb.h:1310
GLboolean r
Definition glcorearb.h:1233
static constexpr int L2G
Definition Cartesian.h:54
static constexpr int T2L
Definition Cartesian.h:55
const int nEvents
Definition test_Fifo.cxx:27