]>
Commit | Line | Data |
---|---|---|
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 | // An analysis task to check the TRD data in simulated data | |
20 | // | |
21 | //*-- Sylwester Radomski | |
22 | ////////////////////////////////////////////////////////////////////////////// | |
23 | // track type codding | |
24 | // | |
25 | // tpci = kTPCin | |
26 | // tpco = kTPCout | |
27 | // tpcz = kTPCout && !kTRDout | |
28 | // trdo = kTRDout | |
29 | // trdr = kTRDref | |
30 | // trdz = kTRDout && !kTRDref | |
31 | // | |
32 | ||
33 | #include <TCanvas.h> | |
34 | #include <TChain.h> | |
35 | #include <TFile.h> | |
36 | #include <TGaxis.h> | |
37 | #include <TH1D.h> | |
38 | #include <TH2D.h> | |
39 | #include <TROOT.h> | |
40 | #include <TStyle.h> | |
41 | ||
42 | #include "AliTRDQATask.h" | |
43 | #include "AliESD.h" | |
44 | #include "AliLog.h" | |
45 | ||
46 | //______________________________________________________________________________ | |
47 | AliTRDQATask::AliTRDQATask(const char *name) : | |
48 | AliAnalysisTask(name,""), | |
49 | fChain(0), | |
50 | fESD(0) | |
51 | { | |
52 | // Constructor. | |
53 | // Input slot #0 works with an Ntuple | |
54 | DefineInput(0, TChain::Class()); | |
55 | // Output slot #0 writes into a TH1 container | |
56 | DefineOutput(0, TObjArray::Class()) ; | |
57 | } | |
58 | ||
59 | //______________________________________________________________________________ | |
60 | void AliTRDQATask::ConnectInputData(const Option_t *) | |
61 | { | |
62 | // Initialisation of branch container and histograms | |
63 | ||
64 | AliInfo(Form("*** Initialization of %s", GetName())) ; | |
65 | ||
66 | // Get input data | |
67 | fChain = dynamic_cast<TChain *>(GetInputData(0)) ; | |
68 | if (!fChain) { | |
69 | AliError(Form("Input 0 for %s not found\n", GetName())); | |
70 | return ; | |
71 | } | |
72 | ||
73 | // One should first check if the branch address was taken by some other task | |
74 | char ** address = (char **)GetBranchAddress(0, "ESD"); | |
75 | if (address) { | |
76 | fESD = (AliESD*)(*address); | |
77 | } else { | |
78 | fESD = new AliESD(); | |
79 | SetBranchAddress(0, "ESD", &fESD); | |
80 | } | |
81 | } | |
82 | ||
83 | //________________________________________________________________________ | |
84 | void AliTRDQATask::CreateOutputObjects() | |
85 | { | |
86 | // create histograms | |
87 | ||
88 | OpenFile(0) ; | |
89 | ||
90 | fNTracks = new TH1D("ntracks", ";number of all tracks", 500, -0.5, 499.5); | |
91 | fEventSize = new TH1D("evSize", ";event size (MB)", 100, 0, 5); | |
92 | ||
93 | fTrackStatus = new TH1D("trackStatus", ";status bit", 32, -0.5, 31.5); | |
94 | fKinkIndex = new TH1D("kinkIndex", ";kink index", 41, -20.5, 20.5); | |
95 | ||
96 | fParIn = new TH1D("parIn", "Inner Plane", 2, -0.5, 1.5); | |
97 | fParOut = new TH1D("parOut", "Outer Plane", 2, -0.5, 1.5); | |
98 | ||
99 | fXIn = new TH1D("xIn", ";X at the inner plane (cm)", 200, 50, 250); | |
100 | fXOut = new TH1D("xOut", ";X at the outer plane (cm)", 300, 50, 400); | |
101 | ||
102 | const int knNameAlpha = 4; | |
103 | const char *namesAlpha[knNameAlpha] = {"alphaTPCi", "alphaTPCo", "alphaTRDo", "alphaTRDr"}; | |
104 | //TH1D *fAlpha[4]; | |
105 | for(int i=0; i<knNameAlpha; i++) { | |
106 | fAlpha[i] = new TH1D(namesAlpha[i], "alpha", 100, -4, 4); | |
107 | } | |
108 | fSectorTRD = new TH1D ("sectorTRD", ";sector TRD", 20, -0.5, 19.5); | |
109 | ||
110 | ||
111 | // track parameters | |
112 | const int knbits = 6; | |
113 | const char *suf[knbits] = {"TPCi", "TPCo", "TPCz", "TRDo", "TRDr", "TRDz"}; | |
114 | for(int i=0; i<knbits; i++) { | |
115 | fPt[i] = new TH1D(Form("pt%s",suf[i]), ";p_{T} (GeV/c);entries TPC", 50, 0, 10); | |
116 | fTheta[i] = new TH1D(Form("theta%s", suf[i]), "theta (rad)", 100, -4, 4); | |
117 | fSigmaY[i] = new TH1D(Form("sigmaY%s",suf[i]), ";sigma Y (cm)", 200, 0, 1); | |
118 | fChi2[i] = new TH1D(Form("Chi2%s", suf[i]), ";#chi2 / ndf", 100, 0, 10); | |
119 | fPlaneYZ[i] = new TH2D(Form("planeYZ%s", suf[i]), Form("%sy (cm);z (cm)", suf[i]), | |
120 | 100, -60, 60, 500, -500, 500); | |
121 | } | |
122 | ||
123 | // efficiency | |
124 | fEffPt[0] = (TH1D*) fPt[0]->Clone(Form("eff_%s_%s", suf[0], suf[1])); | |
125 | fEffPt[1] = (TH1D*) fPt[0]->Clone(Form("eff_%s_%s", suf[1], suf[3])); | |
126 | fEffPt[2] = (TH1D*) fPt[0]->Clone(Form("eff_%s_%s", suf[3], suf[4])); | |
127 | fEffPt[3] = (TH1D*) fPt[0]->Clone(Form("eff_%s_%s", suf[1], suf[4])); | |
128 | ||
129 | for(int i=0; i<4; i++) { | |
130 | fEffPt[i]->Sumw2(); | |
131 | fEffPt[i]->SetMarkerStyle(20); | |
132 | fEffPt[i]->SetMinimum(0); | |
133 | fEffPt[i]->SetMaximum(1.1); | |
134 | } | |
135 | ||
136 | // track features | |
137 | fClustersTRD[0] = new TH1D("clsTRDo", "TRDo;number of clusters", 130, -0.5, 129.5);; | |
138 | fClustersTRD[1] = new TH1D("clsTRDr", "TRDr;number of clusters", 130, -0.5, 129.5);; | |
139 | fClustersTRD[2] = new TH1D("clsTRDz", "TRDz;number of clusters", 130, -0.5, 129.5);; | |
140 | ||
141 | // for good refitted tracks only | |
142 | fTime = new TH1D("time", ";time bin", 25, -0.5, 24.5); | |
143 | fBudget = new TH1D("budget", ";material budget", 100, 0, 100); | |
144 | fQuality = new TH1D("quality", ";track quality", 100, 0, 1.1); | |
145 | fSignal = new TH1D("signal", ";signal", 100, 0, 1e3); | |
146 | ||
147 | // dEdX and PID | |
148 | fTrdSigMom = new TH2D("trdSigMom", ";momentum (GeV/c);signal", 100, 0, 3, 100, 0, 1e3); | |
149 | fTpcSigMom = new TH2D("tpcSigMom", ";momentum (GeV/c);signal", 100, 0, 3, 100, 0, 200); | |
150 | ||
151 | const char *pidName[6] = {"El", "Mu", "Pi", "K", "P", "Ch"}; | |
152 | for(int i=0; i<6; i++) { | |
153 | ||
154 | // TPC | |
155 | fTpcPID[i] = new TH1D(Form("tpcPid%s",pidName[i]), pidName[i], 100, 0, 1.5); | |
156 | fTpcPID[i]->GetXaxis()->SetTitle("probability"); | |
157 | ||
158 | fTpcSigMomPID[i] = new TH2D(Form("tpcSigMom%s",pidName[i]), "", 100, 0, 3, 100, 0, 200); | |
159 | fTpcSigMomPID[i]->SetTitle(Form("%s;momentum (GeV/c);signal",pidName[i])); | |
160 | ||
161 | // TRD | |
162 | fTrdPID[i] = new TH1D(Form("trdPid%s",pidName[i]), pidName[i], 100, 0, 1.5); | |
163 | fTrdPID[i]->GetXaxis()->SetTitle("probability"); | |
164 | ||
165 | fTrdSigMomPID[i] = new TH2D(Form("trdSigMom%s",pidName[i]), "", 100, 0, 3, 100, 0, 1e3); | |
166 | fTrdSigMomPID[i]->SetTitle(Form("%s;momentum (GeV/c);signal",pidName[i])); | |
167 | } | |
168 | ||
169 | ||
170 | // create output container | |
171 | fOutputContainer = new TObjArray(150); | |
172 | ||
173 | // register histograms to the container | |
174 | int counter = 0; | |
175 | ||
176 | fOutputContainer->AddAt(fNTracks, counter++); | |
177 | fOutputContainer->AddAt(fEventSize, counter++); | |
178 | fOutputContainer->AddAt(fTrackStatus, counter++); | |
179 | fOutputContainer->AddAt(fKinkIndex, counter++); | |
180 | fOutputContainer->AddAt(fParIn, counter++); | |
181 | fOutputContainer->AddAt(fParOut, counter++); | |
182 | fOutputContainer->AddAt(fXIn, counter++); | |
183 | fOutputContainer->AddAt(fXOut, counter++); | |
184 | fOutputContainer->AddAt(fAlpha[0], counter++); | |
185 | fOutputContainer->AddAt(fAlpha[1], counter++); | |
186 | fOutputContainer->AddAt(fAlpha[2], counter++); | |
187 | fOutputContainer->AddAt(fAlpha[3], counter++); | |
188 | ||
189 | fOutputContainer->AddAt(fSectorTRD, counter++); | |
190 | for(int i=0; i<knbits; i++) { | |
191 | fOutputContainer->AddAt(fPt[i], counter++); | |
192 | fOutputContainer->AddAt(fTheta[i], counter++); | |
193 | fOutputContainer->AddAt(fSigmaY[i], counter++); | |
194 | fOutputContainer->AddAt(fChi2[i], counter++); | |
195 | fOutputContainer->AddAt(fPlaneYZ[i], counter++); | |
196 | } | |
197 | fOutputContainer->AddAt(fEffPt[0], counter++); | |
198 | fOutputContainer->AddAt(fEffPt[1], counter++); | |
199 | fOutputContainer->AddAt(fEffPt[2], counter++); | |
200 | fOutputContainer->AddAt(fEffPt[3], counter++); | |
201 | ||
202 | fOutputContainer->AddAt(fClustersTRD[0], counter++); | |
203 | fOutputContainer->AddAt(fClustersTRD[1], counter++); | |
204 | fOutputContainer->AddAt(fClustersTRD[2], counter++); | |
205 | fOutputContainer->AddAt(fTime, counter++); | |
206 | fOutputContainer->AddAt(fBudget, counter++); | |
207 | fOutputContainer->AddAt(fQuality, counter++); | |
208 | fOutputContainer->AddAt(fSignal, counter++); | |
209 | fOutputContainer->AddAt(fTrdSigMom, counter++); | |
210 | fOutputContainer->AddAt(fTpcSigMom, counter++); | |
211 | for(int i=0; i<6; i++) { | |
212 | fOutputContainer->AddAt(fTpcPID[i], counter++); | |
213 | fOutputContainer->AddAt(fTpcSigMomPID[i], counter++); | |
214 | fOutputContainer->AddAt(fTrdPID[i], counter++); | |
215 | fOutputContainer->AddAt(fTrdSigMomPID[i], counter++); | |
216 | } | |
217 | ||
218 | //AliInfo(Form("Number of histograms = %d", counter)); | |
219 | ||
220 | } | |
221 | ||
222 | //______________________________________________________________________________ | |
223 | void AliTRDQATask::Exec(Option_t *) | |
224 | { | |
225 | // Process one event | |
226 | ||
227 | Long64_t entry = fChain->GetReadEntry() ; | |
228 | ||
229 | // Processing of one event | |
230 | ||
231 | if (!fESD) { | |
232 | AliError("fESD is not connected to the input!") ; | |
233 | return ; | |
234 | } | |
235 | ||
236 | if ( !((entry-1)%100) ) | |
237 | AliInfo(Form("%s ----> Processing event # %lld", (dynamic_cast<TChain *>(fChain))->GetFile()->GetName(), entry)) ; | |
238 | ||
239 | int nTracks = fESD->GetNumberOfTracks(); | |
240 | fNTracks->Fill(nTracks); | |
241 | ||
242 | // track loop | |
243 | for(int i=0; i<nTracks; i++) { | |
244 | ||
245 | AliESDtrack *track = fESD->GetTrack(i); | |
246 | const AliExternalTrackParam *paramOut = track->GetOuterParam(); | |
247 | const AliExternalTrackParam *paramIn = track->GetInnerParam(); | |
248 | ||
249 | fParIn->Fill(!!paramIn); | |
250 | if (!paramIn) continue; | |
251 | fXIn->Fill(paramIn->GetX()); | |
252 | ||
253 | fParOut->Fill(!!paramOut); | |
254 | if (!paramOut) continue; | |
255 | fXOut->Fill(paramOut->GetX()); | |
256 | ||
257 | int sector = GetSector(paramOut->GetAlpha()); | |
258 | if (!CheckSector(sector)) continue; | |
259 | fSectorTRD->Fill(sector); | |
260 | ||
261 | fKinkIndex->Fill(track->GetKinkIndex(0)); | |
262 | if (track->GetKinkIndex(0)) continue; | |
263 | ||
264 | UInt_t u = 1; | |
265 | UInt_t status = track->GetStatus(); | |
266 | for(int bit=0; bit<32; bit++) | |
267 | if (u<<bit & status) fTrackStatus->Fill(bit); | |
268 | ||
269 | const int knbits = 6; | |
270 | int bit[6] = {0,0,0,0,0,0}; | |
271 | bit[0] = status & AliESDtrack::kTPCin; | |
272 | bit[1] = status & AliESDtrack::kTPCout; | |
273 | bit[2] = (status & AliESDtrack::kTPCout) && !(status & AliESDtrack::kTRDout); | |
274 | bit[3] = status & AliESDtrack::kTRDout; | |
275 | bit[4] = status & AliESDtrack::kTRDrefit; | |
276 | bit[5] = (status & AliESDtrack::kTRDout) && !(status & AliESDtrack::kTRDrefit); | |
277 | ||
278 | ||
279 | // transverse momentum | |
280 | const double *val = track->GetParameter(); // parameters at the vertex | |
281 | double pt = 1./TMath::Abs(val[4]); | |
282 | ||
283 | for(int b=0; b<knbits; b++) { | |
284 | if (bit[b]) { | |
285 | fPt[b]->Fill(pt); | |
286 | fTheta[b]->Fill(val[3]); | |
287 | fSigmaY[b]->Fill(TMath::Sqrt(paramOut->GetSigmaY2())); | |
288 | fChi2[b]->Fill(track->GetTRDchi2()/track->GetTRDncls()); | |
289 | fPlaneYZ[b]->Fill(paramOut->GetY(), paramOut->GetZ()); | |
290 | } | |
291 | } | |
292 | ||
293 | // sectors | |
294 | if (bit[1]) { | |
295 | fAlpha[0]->Fill(paramIn->GetAlpha()); | |
296 | fAlpha[1]->Fill(paramOut->GetAlpha()); | |
297 | } | |
298 | ||
299 | if (bit[3]) fAlpha[2]->Fill(paramOut->GetAlpha()); | |
300 | if (bit[4]) fAlpha[3]->Fill(paramOut->GetAlpha()); | |
301 | ||
302 | // clusters | |
303 | for(int b=0; b<3; b++) | |
304 | if (bit[3+b]) fClustersTRD[b]->Fill(track->GetTRDncls()); | |
305 | ||
306 | // refitted only | |
307 | if (!bit[4]) continue; | |
308 | ||
309 | fQuality->Fill(track->GetTRDQuality()); | |
310 | fBudget->Fill(track->GetTRDBudget()); | |
311 | fSignal->Fill(track->GetTRDsignal()); | |
312 | ||
313 | fTrdSigMom->Fill(track->GetP(), track->GetTRDsignal()); | |
314 | fTpcSigMom->Fill(track->GetP(), track->GetTPCsignal()); | |
315 | ||
316 | // PID only | |
317 | if (status & AliESDtrack::kTRDpid) { | |
318 | ||
319 | for(int l=0; l<6; l++) fTime->Fill(track->GetTRDTimBin(l)); | |
320 | ||
321 | // fill pid histograms | |
322 | double trdr0 = 0, tpcr0 = 0; | |
323 | int trdBestPid = 5, tpcBestPid = 5; // charged | |
324 | const double kminPidValue = 0.9; | |
325 | ||
326 | double pp[5]; | |
327 | track->GetTPCpid(pp); // ESD inconsequence | |
328 | ||
329 | for(int pid=0; pid<5; pid++) { | |
330 | ||
331 | trdr0 += track->GetTRDpid(pid); | |
332 | tpcr0 += pp[pid]; | |
333 | ||
334 | fTrdPID[pid]->Fill(track->GetTRDpid(pid)); | |
335 | fTpcPID[pid]->Fill(pp[pid]); | |
336 | ||
337 | if (track->GetTRDpid(pid) > kminPidValue) trdBestPid = pid; | |
338 | if (pp[pid] > kminPidValue) tpcBestPid = pid; | |
339 | } | |
340 | ||
341 | fTrdPID[5]->Fill(trdr0); // check unitarity | |
342 | fTrdSigMomPID[trdBestPid]->Fill(track->GetP(), track->GetTRDsignal()); | |
343 | ||
344 | fTpcPID[5]->Fill(tpcr0); // check unitarity | |
345 | fTpcSigMomPID[tpcBestPid]->Fill(track->GetP(), track->GetTPCsignal()); | |
346 | } | |
347 | ||
348 | } | |
349 | ||
350 | CalculateEff(); | |
351 | PostData(0, fOutputContainer); | |
352 | } | |
353 | ||
354 | //______________________________________________________________________________ | |
355 | void AliTRDQATask::Terminate(Option_t *) | |
356 | { | |
357 | // Processing when the event loop is ended | |
358 | fOutputContainer = (TObjArray*)GetOutputData(0); | |
359 | int counter = 0; | |
360 | fNTracks = (TH1D*)fOutputContainer->At(counter++); | |
361 | fEventSize = (TH1D*)fOutputContainer->At(counter++); | |
362 | fTrackStatus = (TH1D*)fOutputContainer->At(counter++); | |
363 | fKinkIndex = (TH1D*)fOutputContainer->At(counter++); | |
364 | fParIn = (TH1D*)fOutputContainer->At(counter++); | |
365 | fParOut = (TH1D*)fOutputContainer->At(counter++); | |
366 | fXIn = (TH1D*)fOutputContainer->At(counter++); | |
367 | fXOut = (TH1D*)fOutputContainer->At(counter++); | |
368 | fAlpha[0] = (TH1D*)fOutputContainer->At(counter++); | |
369 | fAlpha[1] = (TH1D*)fOutputContainer->At(counter++); | |
370 | fAlpha[2] = (TH1D*)fOutputContainer->At(counter++); | |
371 | fAlpha[3] = (TH1D*)fOutputContainer->At(counter++); | |
372 | ||
373 | fSectorTRD = (TH1D*)fOutputContainer->At(counter++); | |
374 | const int knbits = 6; | |
375 | for(int i=0; i<knbits; i++) { | |
376 | fPt[i] = (TH1D*)fOutputContainer->At(counter++); | |
377 | fTheta[i] = (TH1D*)fOutputContainer->At(counter++); | |
378 | fSigmaY[i] = (TH1D*)fOutputContainer->At(counter++); | |
379 | fChi2[i] = (TH1D*)fOutputContainer->At(counter++); | |
380 | fPlaneYZ[i] = (TH2D*)fOutputContainer->At(counter++); | |
381 | } | |
382 | fEffPt[0] = (TH1D*)fOutputContainer->At(counter++); | |
383 | fEffPt[1] = (TH1D*)fOutputContainer->At(counter++); | |
384 | fEffPt[2] = (TH1D*)fOutputContainer->At(counter++); | |
385 | fEffPt[3] = (TH1D*)fOutputContainer->At(counter++); | |
386 | ||
387 | fClustersTRD[0] = (TH1D*)fOutputContainer->At(counter++); | |
388 | fClustersTRD[1] = (TH1D*)fOutputContainer->At(counter++); | |
389 | fClustersTRD[2] = (TH1D*)fOutputContainer->At(counter++); | |
390 | fTime = (TH1D*)fOutputContainer->At(counter++); | |
391 | fBudget = (TH1D*)fOutputContainer->At(counter++); | |
392 | fQuality = (TH1D*)fOutputContainer->At(counter++); | |
393 | fSignal = (TH1D*)fOutputContainer->At(counter++); | |
394 | fTrdSigMom = (TH2D*)fOutputContainer->At(counter++); | |
395 | fTpcSigMom = (TH2D*)fOutputContainer->At(counter++); | |
396 | for(int i=0; i<6; i++) { | |
397 | fTpcPID[i] = (TH1D*)fOutputContainer->At(counter++); | |
398 | fTpcSigMomPID[i] = (TH2D*)fOutputContainer->At(counter++); | |
399 | fTrdPID[i] = (TH1D*)fOutputContainer->At(counter++); | |
400 | fTrdSigMomPID[i] = (TH2D*)fOutputContainer->At(counter++); | |
401 | } | |
402 | ||
403 | // create efficiency histograms | |
404 | Bool_t problem = kFALSE ; | |
405 | AliInfo(Form(" *** %s Report:", GetName())) ; | |
406 | ||
407 | CalculateEff(); | |
408 | ||
409 | DrawESD() ; | |
410 | DrawGeoESD() ; | |
411 | //DrawConvESD() ; | |
412 | DrawPidESD() ; | |
413 | ||
414 | char line[1024] ; | |
415 | sprintf(line, ".!tar -zcf %s.tar.gz *.eps", GetName()) ; | |
416 | gROOT->ProcessLine(line); | |
417 | ||
418 | AliInfo(Form("!!! All the eps files are in %s.tar.gz !!!", GetName())) ; | |
419 | ||
420 | char * report ; | |
421 | if(problem) | |
422 | report="Problems found, please check!!!"; | |
423 | else | |
424 | report="OK"; | |
425 | ||
426 | AliInfo(Form("*** %s Summary Report: %s\n",GetName(), report)) ; | |
427 | ||
428 | } | |
429 | ||
430 | //______________________________________________________________________________ | |
431 | const int AliTRDQATask::GetSector(const double alpha) const | |
432 | { | |
433 | // Gets the sector number | |
434 | ||
435 | double size = TMath::DegToRad() * 20.; | |
436 | int sector = (int)((alpha + TMath::Pi())/size); | |
437 | return sector; | |
438 | } | |
439 | ||
440 | //______________________________________________________________________________ | |
441 | const int AliTRDQATask::CheckSector(const int sector) const | |
442 | { | |
443 | // Checks the sector number | |
444 | const int knSec = 8; | |
445 | int sec[] = {2,3,5,6,11,12,13,15}; | |
446 | ||
447 | for(int i=0; i<knSec; i++) | |
448 | if (sector == sec[i]) return 1; | |
449 | ||
450 | return 0; | |
451 | } | |
452 | ||
453 | //______________________________________________________________________________ | |
454 | void AliTRDQATask::CalculateEff() | |
455 | { | |
456 | // calculates the efficiency | |
457 | ||
458 | for(int i=0; i<4; i++) fEffPt[i]->Reset(); | |
459 | ||
460 | fEffPt[0]->Add(fPt[1]); | |
461 | fEffPt[0]->Divide(fPt[0]); | |
462 | ||
463 | fEffPt[1]->Add(fPt[3]); | |
464 | fEffPt[1]->Divide(fPt[1]); | |
465 | ||
466 | fEffPt[2]->Add(fPt[4]); | |
467 | fEffPt[2]->Divide(fPt[3]); | |
468 | ||
469 | fEffPt[3]->Add(fPt[4]); | |
470 | fEffPt[3]->Divide(fPt[1]); | |
471 | } | |
472 | ||
473 | //______________________________________________________________________________ | |
474 | void AliTRDQATask::DrawESD() const | |
475 | { | |
476 | // Makes a few plots | |
477 | ||
478 | TCanvas * cTRD = new TCanvas("cTRD", "TRD ESD Test", 400, 10, 600, 700) ; | |
479 | cTRD->Divide(6,3) ; | |
480 | ||
481 | gROOT->SetStyle("Plain"); | |
482 | gStyle->SetPalette(1); | |
483 | gStyle->SetOptStat(0); | |
484 | ||
485 | TGaxis::SetMaxDigits(3); | |
486 | ||
487 | gStyle->SetLabelFont(52, "XYZ"); | |
488 | gStyle->SetTitleFont(62, "XYZ"); | |
489 | gStyle->SetPadRightMargin(0.02); | |
490 | ||
491 | // draw all | |
492 | ||
493 | const int knplots = 18; | |
494 | const int knover[knplots] = {1,1,1,4,1,1,1,1,1,1,2,1,1,3,1,1,1,1}; | |
495 | const int knnames = 24; | |
496 | const char *names[knnames] = { | |
497 | "ntracks", "kinkIndex", "trackStatus", | |
498 | "ptTPCi", "ptTPCo", "ptTRDo", "ptTRDr", "ptTPCz", "ptTRDz", | |
499 | "eff_TPCi_TPCo", "eff_TPCo_TRDo", "eff_TRDo_TRDr", "eff_TPCo_TRDr", | |
500 | "clsTRDo", "clsTRDr", "clsTRDz", | |
501 | "alphaTPCi", "alphaTPCo", "alphaTRDo", "alphaTRDr", "sectorTRD", | |
502 | "time", "budget", "signal" | |
503 | }; | |
504 | ||
505 | const int klogy[knnames] = { | |
506 | 1,1,1, | |
507 | 1,1,1, | |
508 | 0,0,0,0, | |
509 | 1,1, | |
510 | 0,0,0,0,0, | |
511 | 0,1,1 | |
512 | }; | |
513 | ||
514 | int nhist=0; | |
515 | for(int i=0; i<knplots; i++) { | |
516 | cTRD->cd(i+1) ; | |
517 | ||
518 | // new TCanvas(names[i], names[nhist], 500, 300); | |
519 | ||
520 | for(int j=0; j<knover[i]; j++) { | |
521 | TH1D *hist = dynamic_cast<TH1D*>(gDirectory->FindObject(names[nhist++])); | |
522 | if (!hist) continue; | |
523 | if (hist->GetMaximum() > 0 ) | |
524 | gPad->SetLogy(klogy[i]); | |
525 | if (strstr(hist->GetName(), "eff")) { | |
526 | hist->SetMarkerStyle(20); | |
527 | hist->SetMinimum(0); | |
528 | hist->SetMaximum(1.2); | |
529 | } | |
530 | ||
531 | if (!j) hist->Draw(); | |
532 | else hist->Draw("SAME"); | |
533 | } | |
534 | } | |
535 | cTRD->Print("TRD_ESD.eps"); | |
536 | } | |
537 | ||
538 | //______________________________________________________________________________ | |
539 | void AliTRDQATask::DrawGeoESD() const | |
540 | { | |
541 | // Makes a few plots | |
542 | ||
543 | TCanvas * cTRDGeo = new TCanvas("cTRDGeo", "TRD ESDGeo Test", 400, 10, 600, 700) ; | |
544 | cTRDGeo->Divide(4,2) ; | |
545 | ||
546 | gStyle->SetOptStat(0); | |
547 | TGaxis::SetMaxDigits(3); | |
548 | ||
549 | gStyle->SetLabelFont(52, "XYZ"); | |
550 | gStyle->SetTitleFont(62, "XYZ"); | |
551 | ||
552 | gStyle->SetPadTopMargin(0.06); | |
553 | gStyle->SetTitleBorderSize(0); | |
554 | ||
555 | // draw all | |
556 | const int knnames = 7; | |
557 | const char *names[knnames] = { | |
558 | "xIn", "xOut", | |
559 | "planeYZTPCo", "planeYZTPCz", "planeYZTRDo", "planeYZTRDr", "planeYZTRDz", | |
560 | }; | |
561 | ||
562 | const char *opt[knnames] = { | |
563 | "", "", | |
564 | "colz","colz", "colz", "colz", "colz" | |
565 | }; | |
566 | ||
567 | const int klogy[knnames] = { | |
568 | 1,1, | |
569 | 0,0,0,0,0 | |
570 | }; | |
571 | ||
572 | for(int i=0; i<knnames; i++) { | |
573 | cTRDGeo->cd(i+1) ; | |
574 | TH1D *hist = dynamic_cast<TH1D*>(gDirectory->FindObject(names[i])); | |
575 | if (!hist) continue; | |
576 | ||
577 | //if (i<2) new TCanvas(names[i], names[i], 500, 300); | |
578 | //else new TCanvas(names[i], names[i], 300, 900); | |
579 | ||
580 | if (hist->GetMaximum() > 0 ) | |
581 | gPad->SetLogy(klogy[i]); | |
582 | if (strstr(opt[i],"colz")) gPad->SetRightMargin(0.1); | |
583 | ||
584 | hist->Draw(opt[i]); | |
585 | AliInfo(Form("%s\t%d", names[i], hist->GetEntries())); | |
586 | } | |
587 | ||
588 | cTRDGeo->Print("TRD_Geo.eps"); | |
589 | } | |
590 | ||
591 | //______________________________________________________________________________ | |
592 | void AliTRDQATask::DrawConvESD() const | |
593 | { | |
594 | // Makes a few plots | |
595 | ||
596 | AliInfo("Plotting....") ; | |
597 | TCanvas * cTRDConv = new TCanvas("cTRDConv", "TRD ESDConv Test", 400, 10, 600, 700) ; | |
598 | cTRDConv->Divide(3,2) ; | |
599 | ||
600 | gROOT->SetStyle("Plain"); | |
601 | gROOT->ForceStyle(); | |
602 | gStyle->SetPalette(1); | |
603 | ||
604 | TGaxis::SetMaxDigits(3); | |
605 | ||
606 | gStyle->SetLabelFont(52, "XYZ"); | |
607 | gStyle->SetTitleFont(62, "XYZ"); | |
608 | gStyle->SetPadRightMargin(0.02); | |
609 | ||
610 | const int knnames = 9; | |
611 | const int knplots = 5; | |
612 | const int knover[knplots] = {3,1,1,3,1}; | |
613 | ||
614 | const char *names[knnames] = { | |
615 | "sigmaYTPCo","sigmaYTRDo", "sigmaYTRDr", "sigmaYTPCz", "sigmaYTRDz", | |
616 | "Chi2TPCo", "Chi2TRDo", "Chi2TRDr", "Chi2TRDz" | |
617 | }; | |
618 | ||
619 | const char *opt[knplots] = { | |
620 | "", "", "","","", | |
621 | }; | |
622 | ||
623 | const int klogy[knplots] = { | |
624 | 0,0,0,1,1 | |
625 | }; | |
626 | ||
627 | int nhist = 0; | |
628 | for(int i=0; i<knplots; i++) { | |
629 | cTRDConv->cd(i+1) ; | |
630 | //new TCanvas(names[i], names[i], 500, 300); | |
631 | if (strstr(opt[i],"colz")) gPad->SetRightMargin(0.1); | |
632 | ||
633 | for(int j=0; j<knover[i]; j++) { | |
634 | TH1D *hist = dynamic_cast<TH1D*>(gDirectory->FindObject(names[nhist++])); | |
635 | if ( hist->GetMaximum() > 0 ) | |
636 | gPad->SetLogy(klogy[i]); | |
637 | if (!j) hist->Draw(opt[i]); | |
638 | else hist->Draw("same"); | |
639 | } | |
640 | ||
641 | } | |
642 | cTRDConv->Print("TRD_Conv.eps"); | |
643 | } | |
644 | ||
645 | //______________________________________________________________________________ | |
646 | void AliTRDQATask::DrawPidESD() const | |
647 | { | |
648 | // Makes a few plots | |
649 | ||
650 | TCanvas * cTRDPid = new TCanvas("cTRDPid", "TRD ESDPid Test", 400, 10, 600, 700) ; | |
651 | cTRDPid->Divide(9,3) ; | |
652 | ||
653 | gROOT->SetStyle("Plain"); | |
654 | gROOT->ForceStyle(); | |
655 | gStyle->SetPalette(1); | |
656 | gStyle->SetOptStat(0); | |
657 | ||
658 | TGaxis::SetMaxDigits(3); | |
659 | ||
660 | gStyle->SetLabelFont(52, "XYZ"); | |
661 | gStyle->SetTitleFont(62, "XYZ"); | |
662 | ||
663 | gStyle->SetPadTopMargin(0.05); | |
664 | gStyle->SetPadRightMargin(0.02); | |
665 | ||
666 | // draw all | |
667 | ||
668 | const int knnames = 27; | |
669 | const char *names[knnames] = { | |
670 | ||
671 | "signal", "trdSigMom", "tpcSigMom", | |
672 | ||
673 | "trdPidEl", "trdPidMu", "trdPidPi", "trdPidK", "trdPidP", "trdPidCh", | |
674 | "trdSigMomEl", "trdSigMomMu", "trdSigMomPi", "trdSigMomK", "trdSigMomP", "trdSigMomCh", | |
675 | ||
676 | "tpcPidEl", "tpcPidMu", "tpcPidPi", "tpcPidK", "tpcPidP", "tpcPidCh", | |
677 | "tpcSigMomEl", "tpcSigMomMu", "tpcSigMomPi", "tpcSigMomK", "tpcSigMomP", "tpcSigMomCh" | |
678 | ||
679 | }; | |
680 | ||
681 | const char *opt[knnames] = { | |
682 | ||
683 | "", "colz", "colz", | |
684 | ||
685 | "", "", "", "", "", "" , | |
686 | "colz", "colz", "colz", "colz", "colz", "colz", | |
687 | ||
688 | "", "", "", "", "", "" , | |
689 | "colz", "colz", "colz", "colz", "colz", "colz" | |
690 | }; | |
691 | ||
692 | const int klogy[knnames] = { | |
693 | ||
694 | 0,0,0, | |
695 | ||
696 | 1,1,1,1,1,1, | |
697 | 0,0,0,0,0,0, | |
698 | ||
699 | 1,1,1,1,1,1, | |
700 | 0,0,0,0,0,0 | |
701 | }; | |
702 | ||
703 | for(int i=0; i<knnames; i++) { | |
704 | cTRDPid->cd(i+1) ; | |
705 | ||
706 | TH1D *hist = dynamic_cast<TH1D*>(gDirectory->FindObject(names[i])); | |
707 | if (!hist) continue; | |
708 | ||
709 | //new TCanvas(names[i], names[i], 500, 300); | |
710 | if ( hist->GetMaximum() > 0 ) | |
711 | gPad->SetLogy(klogy[i]); | |
712 | if (strstr(opt[i],"colz")) gPad->SetRightMargin(0.1); | |
713 | ||
714 | if (strstr(names[i],"sigMom")) gPad->SetLogz(1); | |
715 | if (strstr(names[i],"SigMom")) gPad->SetLogz(1); | |
716 | ||
717 | hist->Draw(opt[i]); | |
718 | AliInfo(Form("%s\t%d", names[i], hist->GetEntries())); | |
719 | } | |
720 | cTRDPid->Print("TRD_Pid.eps"); | |
721 | } |