]>
Commit | Line | Data |
---|---|---|
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 | ||
c650b50c | 140 | fTrackNCl = CreateHisto1("TrackNCl", "Number of clusters per track", |
141 | 200, 0, 200, "N_{clusters}", "#Delta N/N", | |
142 | AliMonitorHisto::kNormEntries); | |
143 | ||
04fa961a | 144 | fTrackDEdxVsP = CreateHisto2("TrackDEdxVsP", "dE/dx of tracks", |
145 | 100, 0, 3, 100, 0, 200, | |
146 | "p [GeV/c]", "dE/dx", "#Delta N/N", | |
147 | AliMonitorHisto::kNormEntries); | |
148 | } | |
149 | ||
150 | ||
151 | //_____________________________________________________________________________ | |
152 | void AliMonitorTPC::CreateBranches(TTree* tree) | |
153 | { | |
154 | // create a branch with TPC variables | |
155 | ||
156 | tree->Branch("TPC", "AliMonitorDataTPC", &fData); | |
157 | } | |
158 | ||
159 | ||
160 | //_____________________________________________________________________________ | |
161 | void AliMonitorTPC::FillHistos(AliRunLoader* runLoader, | |
162 | AliRawReader* rawReader) | |
163 | { | |
164 | // fill the TPC monitor histogrms | |
165 | ||
166 | rawReader->Reset(); | |
167 | AliTPCRawStream input(rawReader); | |
168 | while (input.Next()) { | |
169 | fPadsCharge->Fill(input.GetSignal()); | |
170 | } | |
171 | ||
172 | AliLoader* tpcLoader = runLoader->GetLoader("TPCLoader"); | |
173 | if (!tpcLoader) return; | |
174 | ||
175 | ||
176 | tpcLoader->LoadRecPoints(); | |
177 | TTree* clusters = tpcLoader->TreeR(); | |
178 | if (!clusters) return; | |
179 | AliTPCClustersRow* clustersRow = new AliTPCClustersRow; | |
180 | clustersRow->SetClass("AliTPCclusterMI"); | |
181 | clusters->SetBranchAddress("Segment", &clustersRow); | |
182 | ||
183 | for (Int_t i = 0; i < clusters->GetEntries(); i++) { | |
184 | clusters->GetEntry(i); | |
185 | Int_t iSector, iRow; | |
186 | fParam->AdjustSectorRow(clustersRow->GetID(), iSector, iRow); | |
187 | if (iSector >= fParam->GetNInnerSector()) { | |
188 | iRow += fParam->GetNRowLow(); | |
189 | iSector -= fParam->GetNInnerSector(); | |
190 | } | |
191 | ||
192 | TObjArray* array = clustersRow->GetArray(); | |
193 | for (Int_t j = 0; j < array->GetEntriesFast(); j++) { | |
194 | AliTPCclusterMI* cluster = (AliTPCclusterMI*) array->At(j); | |
195 | fClustersCharge->Fill(cluster->GetQ()); | |
196 | fNClustersVsRow->Fill(iRow); | |
197 | fNClustersVsSector->Fill(iSector); | |
198 | } | |
199 | } | |
97d6eb66 | 200 | fNClustersVsSector->ScaleErrorBy(10.); |
04fa961a | 201 | |
202 | delete clustersRow; | |
203 | tpcLoader->UnloadRecPoints(); | |
204 | ||
205 | ||
206 | tpcLoader->LoadTracks(); | |
207 | TTree* tracks = tpcLoader->TreeT(); | |
208 | if (!tracks) return; | |
209 | AliTPCtrack* track = new AliTPCtrack; | |
210 | tracks->SetBranchAddress("tracks", &track); | |
211 | ||
212 | fNTracks->Fill(tracks->GetEntries()); | |
213 | fData->fNTracks = (Int_t) tracks->GetEntries(); | |
214 | fData->SetSize(fData->fNTracks); | |
215 | for (Int_t i = 0; i < tracks->GetEntries(); i++) { | |
216 | tracks->GetEntry(i); | |
217 | fTrackPt->Fill(track->Pt()); | |
218 | fTrackEta->Fill(track->Eta()); | |
219 | fTrackPhi->Fill(track->Phi() * TMath::RadToDeg()); | |
c650b50c | 220 | fTrackNCl->Fill(track->GetNumberOfClusters()); |
04fa961a | 221 | fTrackDEdxVsP->Fill(track->P(), track->GetdEdx()); |
222 | ||
223 | fData->fPt[i] = track->Pt(); | |
224 | fData->fEta[i] = track->Eta(); | |
225 | fData->fPhi[i] = track->Phi() * TMath::RadToDeg(); | |
226 | } | |
227 | ||
228 | delete track; | |
229 | tpcLoader->UnloadTracks(); | |
230 | } |