04fa961a |
1 | /************************************************************************** |
2 | * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. * |
3 | * * |
4 | * Author: The ALICE Off-line Project. * |
5 | * Contributors are mentioned in the code where appropriate. * |
6 | * * |
7 | * Permission to use, copy, modify and distribute this software and its * |
8 | * documentation strictly for non-commercial purposes is hereby granted * |
9 | * without fee, provided that the above copyright notice appears in all * |
10 | * copies and that both the copyright notice and this permission notice * |
11 | * appear in the supporting documentation. The authors make no claims * |
12 | * about the suitability of this software for any purpose. It is * |
13 | * provided "as is" without express or implied warranty. * |
14 | **************************************************************************/ |
15 | |
16 | /* $Id$ */ |
17 | |
18 | /////////////////////////////////////////////////////////////////////////////// |
19 | // // |
20 | // This class creates and fills the monitor histograms for the TPC // |
21 | // // |
22 | /////////////////////////////////////////////////////////////////////////////// |
23 | |
24 | |
25 | #include "AliMonitorTPC.h" |
26 | #include "AliTPCRawStream.h" |
27 | #include "AliTPCClustersRow.h" |
28 | #include "AliTPCclusterMI.h" |
29 | #include "AliTPCtrack.h" |
30 | |
31 | |
32 | ClassImp(AliMonitorDataTPC) |
33 | |
34 | |
35 | //_____________________________________________________________________________ |
36 | AliMonitorDataTPC::AliMonitorDataTPC() |
37 | { |
38 | fPt = fEta = fPhi = NULL; |
39 | fSize = 0; |
40 | } |
41 | |
42 | //_____________________________________________________________________________ |
43 | AliMonitorDataTPC::AliMonitorDataTPC(Int_t size) |
44 | { |
45 | fPt = new Float_t[size]; |
46 | fEta = new Float_t[size]; |
47 | fPhi = new Float_t[size]; |
48 | fSize = size; |
49 | } |
50 | |
51 | //_____________________________________________________________________________ |
52 | AliMonitorDataTPC::~AliMonitorDataTPC() |
53 | { |
54 | delete[] fPt; |
55 | delete[] fEta; |
56 | delete[] fPhi; |
57 | } |
58 | |
59 | //_____________________________________________________________________________ |
60 | void AliMonitorDataTPC::SetSize(Int_t size) |
61 | { |
62 | if (size > fSize) { |
63 | delete[] fPt; |
64 | delete[] fEta; |
65 | delete[] fPhi; |
66 | fPt = new Float_t[size]; |
67 | fEta = new Float_t[size]; |
68 | fPhi = new Float_t[size]; |
69 | fSize = size; |
70 | } |
71 | } |
72 | |
73 | |
74 | |
75 | ClassImp(AliMonitorTPC) |
76 | |
77 | |
78 | //_____________________________________________________________________________ |
79 | AliMonitorTPC::AliMonitorTPC(AliTPCParam* param) |
80 | { |
81 | // create a TPC monitor object with the given parameters |
82 | |
83 | fParam = param; |
84 | fData = new AliMonitorDataTPC(10000); |
85 | } |
86 | |
87 | //_____________________________________________________________________________ |
88 | AliMonitorTPC::~AliMonitorTPC() |
89 | { |
90 | delete fData; |
91 | } |
92 | |
93 | |
94 | //_____________________________________________________________________________ |
95 | void AliMonitorTPC::CreateHistos(TFolder* folder) |
96 | { |
97 | // create the TPC monitor histograms |
98 | |
99 | fFolder = folder->AddFolder("TPC", "TPC"); |
100 | |
101 | fPadsCharge = CreateHisto1("PadsCharge", "charge distribution of pads", |
102 | 100, 0, 200, "charge", "#Delta N/N", |
103 | AliMonitorHisto::kNormEvents); |
104 | |
105 | fClustersCharge = CreateHisto1("ClustersCharge", |
106 | "charge distribution of clusters", |
107 | 100, 0, 1000, "charge", "#Delta N/N", |
108 | AliMonitorHisto::kNormEvents); |
109 | |
110 | Int_t nRows = fParam->GetNRowLow() + fParam->GetNRowUp(); |
111 | fNClustersVsRow = CreateHisto1("NClustersVsRow", |
112 | "mean number of clusters per pad row", |
113 | nRows, -0.5, nRows-0.5, |
114 | "pad row", "<N_{clusters}>", |
115 | AliMonitorHisto::kNormEvents); |
116 | |
117 | Int_t nSector = fParam->GetNInnerSector(); |
118 | fNClustersVsSector = CreateHisto1("NClustersVsSector", |
119 | "mean number of clusters per sector", |
120 | nSector, -0.5, nSector-0.5, |
121 | "sector", "<N_{clusters}>", |
122 | AliMonitorHisto::kNormEvents); |
123 | |
124 | fNTracks = CreateTrend("NTracks", "number of tracks per event", |
125 | "N_{tracks}"); |
126 | |
127 | fTrackPt = CreateHisto1("TrackPt", "pt distribution of tracks", |
128 | 90, 0, 3, "p_{t} [GeV/c]", "#Delta N/N", |
129 | AliMonitorHisto::kNormEntries); |
130 | |
131 | fTrackEta = CreateHisto1("TrackEta", "eta distribution of tracks", |
132 | 100, -2, 2, "#eta", "#Delta N/N", |
133 | AliMonitorHisto::kNormEntries); |
134 | |
135 | fTrackPhi = CreateHisto1("TrackPhi", "phi distribution of tracks", |
136 | 120, 0, 360, "#phi [#circ]", "#Delta N/N", |
137 | AliMonitorHisto::kNormEntries); |
138 | fTrackPhi->SetDescription("The phi distribution should be flat on average.\nIf it is not flat check for dead TPC sectors."); |
139 | |
140 | fTrackDEdxVsP = CreateHisto2("TrackDEdxVsP", "dE/dx of tracks", |
141 | 100, 0, 3, 100, 0, 200, |
142 | "p [GeV/c]", "dE/dx", "#Delta N/N", |
143 | AliMonitorHisto::kNormEntries); |
144 | } |
145 | |
146 | |
147 | //_____________________________________________________________________________ |
148 | void AliMonitorTPC::CreateBranches(TTree* tree) |
149 | { |
150 | // create a branch with TPC variables |
151 | |
152 | tree->Branch("TPC", "AliMonitorDataTPC", &fData); |
153 | } |
154 | |
155 | |
156 | //_____________________________________________________________________________ |
157 | void AliMonitorTPC::FillHistos(AliRunLoader* runLoader, |
158 | AliRawReader* rawReader) |
159 | { |
160 | // fill the TPC monitor histogrms |
161 | |
162 | rawReader->Reset(); |
163 | AliTPCRawStream input(rawReader); |
164 | while (input.Next()) { |
165 | fPadsCharge->Fill(input.GetSignal()); |
166 | } |
167 | |
168 | AliLoader* tpcLoader = runLoader->GetLoader("TPCLoader"); |
169 | if (!tpcLoader) return; |
170 | |
171 | |
172 | tpcLoader->LoadRecPoints(); |
173 | TTree* clusters = tpcLoader->TreeR(); |
174 | if (!clusters) return; |
175 | AliTPCClustersRow* clustersRow = new AliTPCClustersRow; |
176 | clustersRow->SetClass("AliTPCclusterMI"); |
177 | clusters->SetBranchAddress("Segment", &clustersRow); |
178 | |
179 | for (Int_t i = 0; i < clusters->GetEntries(); i++) { |
180 | clusters->GetEntry(i); |
181 | Int_t iSector, iRow; |
182 | fParam->AdjustSectorRow(clustersRow->GetID(), iSector, iRow); |
183 | if (iSector >= fParam->GetNInnerSector()) { |
184 | iRow += fParam->GetNRowLow(); |
185 | iSector -= fParam->GetNInnerSector(); |
186 | } |
187 | |
188 | TObjArray* array = clustersRow->GetArray(); |
189 | for (Int_t j = 0; j < array->GetEntriesFast(); j++) { |
190 | AliTPCclusterMI* cluster = (AliTPCclusterMI*) array->At(j); |
191 | fClustersCharge->Fill(cluster->GetQ()); |
192 | fNClustersVsRow->Fill(iRow); |
193 | fNClustersVsSector->Fill(iSector); |
194 | } |
195 | } |
196 | fNClustersVsSector->ScaleErrorBy(100.); |
197 | |
198 | delete clustersRow; |
199 | tpcLoader->UnloadRecPoints(); |
200 | |
201 | |
202 | tpcLoader->LoadTracks(); |
203 | TTree* tracks = tpcLoader->TreeT(); |
204 | if (!tracks) return; |
205 | AliTPCtrack* track = new AliTPCtrack; |
206 | tracks->SetBranchAddress("tracks", &track); |
207 | |
208 | fNTracks->Fill(tracks->GetEntries()); |
209 | fData->fNTracks = (Int_t) tracks->GetEntries(); |
210 | fData->SetSize(fData->fNTracks); |
211 | for (Int_t i = 0; i < tracks->GetEntries(); i++) { |
212 | tracks->GetEntry(i); |
213 | fTrackPt->Fill(track->Pt()); |
214 | fTrackEta->Fill(track->Eta()); |
215 | fTrackPhi->Fill(track->Phi() * TMath::RadToDeg()); |
216 | fTrackDEdxVsP->Fill(track->P(), track->GetdEdx()); |
217 | |
218 | fData->fPt[i] = track->Pt(); |
219 | fData->fEta[i] = track->Eta(); |
220 | fData->fPhi[i] = track->Phi() * TMath::RadToDeg(); |
221 | } |
222 | |
223 | delete track; |
224 | tpcLoader->UnloadTracks(); |
225 | } |