1+ // Copyright 2019-2026 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+ /// \file CheckClusterSizeVsEta.C
13+ /// \brief A post-processing macro to draw average cluster size vs eta
14+
15+ #include <TCanvas.h>
16+ #include <TFile.h>
17+ #include <TH1F.h>
18+ #include <TH2F.h>
19+ #include <TNtuple.h>
20+ #include <TString.h>
21+ #include <TTree.h>
22+ #include <TROOT.h>
23+ #include <TStyle.h>
24+
25+ // ### required input file: CheckClusters.root, which is the output of CheckClusters.C macro
26+ void CheckClusterSizeVsEta (const std ::string & strFileInput = "CheckClusters.root" )
27+ {
28+ gStyle -> SetOptStat (0 );
29+
30+ TFile * fileInput = new TFile (strFileInput .c_str ());
31+ TTree * tree = (TTree * )fileInput -> Get ("ntc" );
32+
33+ std ::cout << "Opened tree: " << tree -> GetName ()
34+ << ", entries = " << tree -> GetEntries () << std ::endl ;
35+
36+ // set branch addresses
37+ Float_t event ;
38+ Float_t mcTrackID ;
39+ Float_t hitLocX , hitLocZ ;
40+ Float_t hitGlobX , hitGlobY , hitGlobZ ;
41+ Float_t clusGlobX , clusGlobY , clusGlobZ ;
42+ Float_t clusLocX , clusLocZ ;
43+ Float_t rofFrame ;
44+ Float_t clusSize ;
45+ Float_t chipID ;
46+ Float_t layer ;
47+ Float_t disk ;
48+ Float_t subdet ;
49+ Float_t row , col ;
50+ Float_t pt ;
51+
52+ // set branch addresses
53+ tree -> SetBranchAddress ("event" , & event );
54+ tree -> SetBranchAddress ("mcTrackID" , & mcTrackID );
55+ tree -> SetBranchAddress ("hitLocX" , & hitLocX );
56+ tree -> SetBranchAddress ("hitLocZ" , & hitLocZ );
57+ tree -> SetBranchAddress ("hitGlobX" , & hitGlobX );
58+ tree -> SetBranchAddress ("hitGlobY" , & hitGlobY );
59+ tree -> SetBranchAddress ("hitGlobZ" , & hitGlobZ );
60+ tree -> SetBranchAddress ("clusGlobX" , & clusGlobX );
61+ tree -> SetBranchAddress ("clusGlobY" , & clusGlobY );
62+ tree -> SetBranchAddress ("clusGlobZ" , & clusGlobZ );
63+ tree -> SetBranchAddress ("clusLocX" , & clusLocX );
64+ tree -> SetBranchAddress ("clusLocZ" , & clusLocZ );
65+ tree -> SetBranchAddress ("rofFrame" , & rofFrame );
66+ tree -> SetBranchAddress ("clusSize" , & clusSize );
67+ tree -> SetBranchAddress ("chipID" , & chipID );
68+ tree -> SetBranchAddress ("layer" , & layer );
69+ tree -> SetBranchAddress ("disk" , & disk );
70+ tree -> SetBranchAddress ("subdet" , & subdet );
71+ tree -> SetBranchAddress ("row" , & row );
72+ tree -> SetBranchAddress ("col" , & col );
73+ tree -> SetBranchAddress ("pt" , & pt );
74+
75+ // Some QA histograms
76+ TH1F * hPt = new TH1F ("hPt" , "p_{T};p_{T};Entries" , 100 , 0. , 10. );
77+ TH1F * hClusSize = new TH1F ("hClusSize" , "Cluster size;clusSize;Entries" , 20 , 0. , 20. );
78+ TH1F * hLayer = new TH1F ("hLayer" , "Layer;layer;Entries" , 20 , -0.5 , 19.5 );
79+ TH1F * hDxGlob = new TH1F ("hDxGlob" , "clusGlobX - hitGlobX;#DeltaX [global];Entries" , 200 , -1. , 1. );
80+ TH1F * hDzGlob = new TH1F ("hDzGlob" , "clusGlobZ - hitGlobZ;#DeltaZ [global];Entries" , 200 , -1. , 1. );
81+ TH2F * hHitXY = new TH2F ("hHitXY" , "Hit global XY;hitGlobX;hitGlobY" , 200 , -20. , 20. , 200 , -20. , 20. );
82+ TH2F * hClusVsHitX = new TH2F ("hClusVsHitX" , "clusGlobX vs hitGlobX;hitGlobX;clusGlobX" , 200 , -20. , 20. , 200 , -20. , 20. );
83+
84+ // histograms for cluster size vs eta for each layer (ML and OT barrel layers):
85+ const int nLayers = 8 ; // ML and OT barrel layers only
86+ TH2F * hClustSizePerLayerVsEta [nLayers ];
87+ for (int i = 0 ; i < nLayers ; i ++ )
88+ {
89+ hClustSizePerLayerVsEta [i ] = new TH2F (Form ("hClustSizePerLayerVsEta_Lay%d" , i ), Form ("Cluster size vs eta for layer %d;#eta;Cluster size" , i ), 160 , -4 , 4 , 101 , -0.5 , 100.5 );
90+ }
91+
92+ // Loop over entries
93+ const Long64_t nEntries = tree -> GetEntries ();
94+ for (Long64_t i = 0 ; i < nEntries ; ++ i )
95+ {
96+ tree -> GetEntry (i );
97+
98+ // Fill QA histograms
99+ float dXGlob = clusGlobX - hitGlobX ;
100+ float dZGlob = clusGlobZ - hitGlobZ ;
101+ hPt -> Fill (pt );
102+ hClusSize -> Fill (clusSize );
103+ hLayer -> Fill (layer );
104+ hDxGlob -> Fill (dXGlob );
105+ hDzGlob -> Fill (dZGlob );
106+ hHitXY -> Fill (hitGlobX , hitGlobY );
107+ hClusVsHitX -> Fill (hitGlobX , clusGlobX );
108+
109+ // cls size vs eta:
110+ float clustR = sqrt (clusGlobX * clusGlobX + clusGlobY * clusGlobY );
111+ float clustPhi = atan2 (clusGlobY , clusGlobX );
112+ float clustTheta = atan2 (clustR , clusGlobZ );
113+ float clustEta = - log (tan (clustTheta / 2 ));
114+
115+ // !!! important: to avoid VD layers (numeration for ML starts from 0, while VD layers are also numbered as 0,1,2)
116+ if (clustR > 5 ) // cm
117+ hClustSizePerLayerVsEta [(int )layer ]-> Fill (clustEta , clusSize );
118+
119+ // progress print
120+ if ((i + 1 ) % 200000 == 0 )
121+ {
122+ std ::cout << "Processed " << (i + 1 ) << " / " << nEntries << " entries" << std ::endl ;
123+ }
124+ }
125+
126+ // Save histograms to file
127+ TFile * fout = TFile ::Open ("clusterSizes_vs_eta.root" , "RECREATE" );
128+ hPt -> Write ();
129+ hClusSize -> Write ();
130+ hLayer -> Write ();
131+ hDxGlob -> Write ();
132+ hDzGlob -> Write ();
133+ hHitXY -> Write ();
134+ hClusVsHitX -> Write ();
135+
136+ // draw some QA histograms
137+ TCanvas * c1 = new TCanvas ("canv_clusters_QA" , "Clusters QA" , 1200 , 800 );
138+ c1 -> Divide (2 , 2 );
139+ c1 -> cd (1 );
140+ hPt -> Draw ();
141+ c1 -> cd (2 );
142+ hClusSize -> Draw ();
143+ c1 -> cd (3 );
144+ hDxGlob -> Draw ();
145+ c1 -> cd (4 );
146+ hHitXY -> Draw ("COLZ" );
147+
148+ TCanvas * canv_clsSize_vs_eta [nLayers ];
149+ TProfile * profPerLayerVsEta [nLayers ];
150+ for (int i = 0 ; i < nLayers ; i ++ )
151+ {
152+ canv_clsSize_vs_eta [i ] = new TCanvas (Form ("canv_clsSize_vs_eta_Lay%d" , i ), Form ("Cluster size vs eta for layer %d" , i ), 800 , 600 );
153+ hClustSizePerLayerVsEta [i ]-> Draw ("COLZ" );
154+ gPad -> SetLogz ();
155+ profPerLayerVsEta [i ] = hClustSizePerLayerVsEta [i ]-> ProfileX ();
156+ profPerLayerVsEta [i ]-> DrawCopy ("same" );
157+
158+ hClustSizePerLayerVsEta [i ]-> Write ();
159+ profPerLayerVsEta [i ]-> Write ();
160+ }
161+
162+ // canvas with profiles for all layers together
163+ TCanvas * canv_av_clsSize_vs_eta_all_layers = new TCanvas ("canv_clsSize_vs_eta_all_layers" , "Cluster size vs eta for all layers" , 800 , 600 );
164+ TLegend * legLayers = new TLegend (0.3 , 0.52 , 0.65 , 0.89 );
165+ int colors [] = {kRed , kBlue , kMagenta + 1 , kCyan + 1 , kGray + 2 , kRed , kBlue , kMagenta + 1 , kCyan , kAzure + 1 , kOrange - 9 , kRed + 2 , kBlue + 2 , kMagenta + 2 };
166+ for (int i = 0 ; i < nLayers ; i ++ )
167+ {
168+ profPerLayerVsEta [i ]-> SetLineColor (colors [i ]);
169+ profPerLayerVsEta [i ]-> SetMarkerColor (colors [i ]);
170+ profPerLayerVsEta [i ]-> SetMarkerStyle (i < 5 ? 20 : 24 );
171+ profPerLayerVsEta [i ]-> SetTitle (";#eta;#LTcluster size#GT" );
172+ profPerLayerVsEta [i ]-> GetYaxis ()-> SetRangeUser (0. , 12.5 );
173+ profPerLayerVsEta [i ]-> DrawCopy (i == 0 ? "P" : "P same" );
174+ legLayers -> AddEntry (profPerLayerVsEta [i ], Form ("Layer %d" , i ), "P" );
175+ }
176+
177+ legLayers -> Draw ();
178+ gPad -> SetGrid ();
179+ canv_av_clsSize_vs_eta_all_layers -> SaveAs ("clsSize_vs_eta_all_layers.png" );
180+ canv_av_clsSize_vs_eta_all_layers -> Write ();
181+
182+ // fout->Close();
183+ }
0 commit comments