Project
Loading...
Searching...
No Matches
TrackClusters.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#define _USE_MATH_DEFINES
13
14#include <cmath>
15#include <memory>
16
17// root includes
18#include "TFile.h"
19#include "TObjArray.h"
20
21// o2 includes
23#include "TPCQC/TrackClusters.h"
24#include "TPCQC/Tracks.h"
25#include "TPCQC/Helpers.h"
26#include "GPUO2InterfaceRefit.h"
28
30using namespace o2::tpc::qc;
31
32struct binning {
33 int bins;
34 double min;
35 double max;
36};
37
38const binning binsClusters{160, 0., 160.};
39const binning binsSharedClusters{160, 0., 160.};
40const binning binsFoundClusters{160, 0., 160.};
41const binning binsCrossedRows{160, 0., 160.};
42const binning binsRatio{150, 0., 1.5};
43
44//______________________________________________________________________________
46{
47 TH1::AddDirectory(false);
48 mMapHist["clusters"].emplace_back(std::make_unique<TH1F>("clusters", "Clusters;NClusters;Entries", binsClusters.bins, binsClusters.min, binsClusters.max));
49 mMapHist["sharedClusters"].emplace_back(std::make_unique<TH1F>("sharedClusters", "sharedClusters;NSharedClusters;Entries", binsSharedClusters.bins, binsSharedClusters.min, binsSharedClusters.max));
50 mMapHist["crossedRows"].emplace_back(std::make_unique<TH1F>("crossedRows", "crossedRows;NCrossedRows;Entries", binsCrossedRows.bins, binsCrossedRows.min, binsCrossedRows.max));
51 mMapHist["sharedClustersOverClusters"].emplace_back(std::make_unique<TH1F>("sharedClustersOverClusters", "sharedClustersOverClusters;NSharedClusters/NClusters;Entries", binsRatio.bins, binsRatio.min, binsRatio.max));
52 mMapHist["clustersOverCrossedRow"].emplace_back(std::make_unique<TH1F>("clustersOverCrossedRow", "clustersOverCrossedRow;NClusters/NCrossedRows;Entries", binsRatio.bins, binsRatio.min, binsRatio.max));
53}
54
55//______________________________________________________________________________
57{
58 for (const auto& pair : mMapHist) {
59 for (auto& hist : pair.second) {
60 hist->Reset();
61 }
62 }
63}
64
65//______________________________________________________________________________
66bool TrackClusters::processTrackAndClusters(const std::vector<o2::tpc::TrackTPC>* tracks, const o2::tpc::ClusterNativeAccess* clusterIndex, std::vector<o2::tpc::TPCClRefElem>* clusRefs)
67{
68
69 std::vector<unsigned char> mBufVec;
70 mBufVec.resize(clusterIndex->nClustersTotal);
71
72 o2::gpu::GPUO2InterfaceRefit::fillSharedClustersAndOccupancyMap(clusterIndex, *tracks, clusRefs->data(), mBufVec.data());
73
74 for (auto const& track : (*tracks)) {
75 const auto dEdxTot = track.getdEdx().dEdxTotTPC;
76 const auto nCls = uint8_t(track.getNClusters());
77 const auto eta = track.getEta();
78
79 if (nCls < mCutMinNCls || dEdxTot < mCutMindEdxTot || std::fabs(eta) > mCutAbsEta) {
80 continue;
81 }
82
83 uint8_t shared = 200, found = 0, crossed = 0;
84
85 o2::TrackMethods::countTPCClusters(track, *clusRefs, mBufVec, *clusterIndex, shared, found, crossed);
86
87 mMapHist["clusters"][0]->Fill(found);
88 mMapHist["sharedClusters"][0]->Fill(shared);
89 mMapHist["crossedRows"][0]->Fill(crossed);
90 mMapHist["sharedClustersOverClusters"][0]->Fill(static_cast<float>(shared) / static_cast<float>(found));
91 mMapHist["clustersOverCrossedRow"][0]->Fill(static_cast<float>(found) / static_cast<float>(crossed));
92 }
93
94 return true;
95}
96
97//______________________________________________________________________________
98void TrackClusters::dumpToFile(const std::string filename)
99{
100 auto f = std::unique_ptr<TFile>(TFile::Open(filename.c_str(), "recreate"));
101 for (const auto& [name, histos] : mMapHist) {
102 TObjArray arr;
103 arr.SetName(name.data());
104 for (auto& hist : histos) {
105 arr.Add(hist.get());
106 }
107 arr.Write(arr.GetName(), TObject::kSingleKey);
108 }
109 f->Close();
110}
const binning binsCrossedRows
const binning binsClusters
const binning binsFoundClusters
const binning binsSharedClusters
ClassImp(o2::tpc::qc::TrackClusters)
const binning binsRatio
static void countTPCClusters(const o2::tpc::TrackTPC &track, const gsl::span< const o2::tpc::TPCClRefElem > &tpcClusRefs, const gsl::span< const unsigned char > &tpcClusShMap, const o2::tpc::ClusterNativeAccess &tpcClusAcc, uint8_t &shared, uint8_t &found, uint8_t &crossed)
static void fillSharedClustersAndOccupancyMap(const o2::tpc::ClusterNativeAccess *cl, const gsl::span< const o2::tpc::TrackTPC > trks, const o2::tpc::TPCClRefElem *trackRef, uint8_t *shmap, uint32_t *ocmap=nullptr, uint32_t nHbfPerTf=0, const GPUParam *param=nullptr)
Shared cluster and crossed rows TPC quality control task.
void dumpToFile(std::string filename)
Dump results to a file.
void initializeHistograms()
Initialize all histograms.
void resetHistograms()
Reset all histograms.
bool processTrackAndClusters(const std::vector< o2::tpc::TrackTPC > *tracks, const o2::tpc::ClusterNativeAccess *clusterIndex, std::vector< o2::tpc::TPCClRefElem > *clusRefs)
GLuint const GLchar * name
Definition glcorearb.h:781
GLdouble f
Definition glcorearb.h:310
std::string filename()
int bins
Definition PID.cxx:34
double max
Definition PID.cxx:36
double min
Definition PID.cxx:35