]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG1/TRD/AliTRDcheckDET.cxx
fix error during merging (https://savannah.cern.ch/bugs/index.php?70729)
[u/mrichter/AliRoot.git] / PWG1 / TRD / AliTRDcheckDET.cxx
1 ////////////////////////////////////////////////////////////////////////////
2 //                                                                        //
3 //                                                                        //
4 //  Basic checks for tracking and detector performance                    //
5 //  
6 //     For the moment (15.10.2009) the following checks are implemented    //
7 //       - Number of clusters/track
8 //       - Number of clusters/tracklet
9 //       - Number of tracklets/track from different sources (Barrel, StandAlone)
10 //       - Number of findable tracklets
11 //       - Number of tracks per event and TRD sector
12 //       - <PH>
13 //       - Chi2 distribution for tracks
14 //       - Charge distribution per cluster and tracklet
15 //       - Number of events and tracks per trigger source 
16 //       - Trigger purity
17 //       - Track and Tracklet propagation status
18 //
19 //  Authors:                                                              //
20 //    Anton Andronic <A.Andronic@gsi.de>                                  //
21 //    Alexandru Bercuci <A.Bercuci@gsi.de>                                //
22 //    Markus Fasel <M.Fasel@gsi.de>                                       //
23 //                                                                        //
24 ////////////////////////////////////////////////////////////////////////////
25
26 #include <TAxis.h>
27 #include <TCanvas.h>
28 #include <TFile.h>
29 #include <TH1F.h>
30 #include <TH1I.h>
31 #include <TH2F.h>
32 #include <TF1.h>
33 #include <TGaxis.h>
34 #include <TGraph.h>
35 #include <TGraphErrors.h>
36 #include <TLegend.h>
37 #include <TLinearFitter.h>
38 #include <TMath.h>
39 #include <TMap.h>
40 #include <TProfile2D.h>
41 #include <TObjArray.h>
42 #include <TObject.h>
43 #include <TObjString.h>
44
45 #include <TPad.h>
46 #include <TProfile.h>
47 #include <TProfile2D.h>
48 #include <TROOT.h>
49 #include <TChain.h>
50
51 #include "AliLog.h"
52 #include "AliTRDcluster.h"
53 #include "AliESDHeader.h"
54 #include "AliESDRun.h"
55 #include "AliESDtrack.h"
56 #include "AliExternalTrackParam.h"
57 #include "AliTRDgeometry.h"
58 #include "AliTRDpadPlane.h"
59 #include "AliTRDSimParam.h"
60 #include "AliTRDseedV1.h"
61 #include "AliTRDtrackV1.h"
62 #include "AliTRDtrackerV1.h"
63 #include "AliTRDReconstructor.h"
64 #include "AliTrackReference.h"
65 #include "AliTrackPointArray.h"
66 #include "AliTracker.h"
67 #include "TTreeStream.h"
68
69 #include "info/AliTRDtrackInfo.h"
70 #include "info/AliTRDeventInfo.h"
71 #include "AliTRDcheckDET.h"
72
73 #include <cstdio>
74 #include <iostream>
75
76 ClassImp(AliTRDcheckDET)
77
78 //_______________________________________________________
79 AliTRDcheckDET::AliTRDcheckDET():
80   AliTRDrecoTask()
81   ,fEventInfo(NULL)
82   ,fTriggerNames(NULL)
83   ,fReconstructor(NULL)
84   ,fGeo(NULL)
85   ,fFlags(0)
86 {
87   //
88   // Default constructor
89   //
90   SetNameTitle("checkDET", "Basic TRD data checker");
91 }
92
93 //_______________________________________________________
94 AliTRDcheckDET::AliTRDcheckDET(char* name):
95   AliTRDrecoTask(name, "Basic TRD data checker")
96   ,fEventInfo(NULL)
97   ,fTriggerNames(NULL)
98   ,fReconstructor(NULL)
99   ,fGeo(NULL)
100   ,fFlags(0)
101 {
102   //
103   // Default constructor
104   //
105   DefineInput(2, AliTRDeventInfo::Class());
106
107   fReconstructor = new AliTRDReconstructor;
108   fReconstructor->SetRecoParam(AliTRDrecoParam::GetLowFluxParam());
109   fGeo = new AliTRDgeometry;
110   InitFunctorList();
111 }
112
113
114 //_______________________________________________________
115 AliTRDcheckDET::~AliTRDcheckDET(){
116   //
117   // Destructor
118   // 
119   if(fTriggerNames) delete fTriggerNames;
120   delete fReconstructor;
121   delete fGeo;
122 }
123
124 //_______________________________________________________
125 void AliTRDcheckDET::UserCreateOutputObjects(){
126   //
127   // Create Output Objects
128   //
129   AliTRDrecoTask::UserCreateOutputObjects();
130   if(!fTriggerNames) fTriggerNames = new TMap();
131 }
132
133 //_______________________________________________________
134 void AliTRDcheckDET::UserExec(Option_t *opt){
135   //
136   // Execution function
137   // Filling TRD quality histos
138   //
139
140   fEventInfo = dynamic_cast<AliTRDeventInfo *>(GetInputData(2));
141   AliDebug(2, Form("EventInfo[%p] Header[%p]", (void*)fEventInfo, (void*)(fEventInfo?fEventInfo->GetEventHeader():NULL)));
142
143   AliTRDrecoTask::UserExec(opt);  
144   Int_t nTracks = 0;            // Count the number of tracks per event
145   Int_t triggermask = fEventInfo->GetEventHeader()->GetTriggerMask();
146   TString triggername =  fEventInfo->GetRunInfo()->GetFiredTriggerClasses(triggermask);
147   AliDebug(6, Form("Trigger cluster: %d, Trigger class: %s\n", triggermask, triggername.Data()));
148   dynamic_cast<TH1F *>(fContainer->UncheckedAt(kNeventsTrigger))->Fill(triggermask);
149   for(Int_t iti = 0; iti < fTracks->GetEntriesFast(); iti++){
150     if(!fTracks->UncheckedAt(iti)) continue;
151     AliTRDtrackInfo *fTrackInfo = dynamic_cast<AliTRDtrackInfo *>(fTracks->UncheckedAt(iti));
152     if(!fTrackInfo->GetTrack()) continue;
153     nTracks++;
154   }
155
156   if(nTracks){
157     dynamic_cast<TH1F *>(fContainer->UncheckedAt(kNeventsTriggerTracks))->Fill(triggermask);
158     dynamic_cast<TH1F *>(fContainer->UncheckedAt(kNtracksEvent))->Fill(nTracks);
159   }
160   if(triggermask <= 20 && !fTriggerNames->FindObject(Form("%d", triggermask))){
161     fTriggerNames->Add(new TObjString(Form("%d", triggermask)), new TObjString(triggername));
162     // also set the label for both histograms
163     TH1 *histo = dynamic_cast<TH1F *>(fContainer->UncheckedAt(kNeventsTriggerTracks));
164     histo->GetXaxis()->SetBinLabel(histo->FindBin(triggermask), triggername);
165     histo = dynamic_cast<TH1F *>(fContainer->UncheckedAt(kNeventsTrigger));
166     histo->GetXaxis()->SetBinLabel(histo->FindBin(triggermask), triggername);
167   }
168 }
169
170
171 //_______________________________________________________
172 Bool_t AliTRDcheckDET::PostProcess(){
173   //
174   // Do Postprocessing (for the moment set the number of Reference histograms)
175   //
176   
177   TH1 * h = NULL;
178   
179   // Calculate of the trigger clusters purity
180   h = dynamic_cast<TH1F *>(fContainer->FindObject("hEventsTrigger"));
181   TH1F *h1 = dynamic_cast<TH1F *>(fContainer->FindObject("hEventsTriggerTracks"));
182   h1->Divide(h);
183   Float_t purities[20], val = 0;
184   TString triggernames[20];
185   Int_t nTriggerClasses = 0;
186   for(Int_t ibin = 1; ibin <= h->GetNbinsX(); ibin++){
187     if((val = h1->GetBinContent(ibin))){
188       purities[nTriggerClasses] = val;
189       triggernames[nTriggerClasses] = h1->GetXaxis()->GetBinLabel(ibin);
190       nTriggerClasses++;
191     }
192   }
193   h = dynamic_cast<TH1F *>(fContainer->UncheckedAt(kTriggerPurity));
194   TAxis *ax = h->GetXaxis();
195   for(Int_t itrg = 0; itrg < nTriggerClasses; itrg++){
196     h->Fill(itrg, purities[itrg]);
197     ax->SetBinLabel(itrg+1, triggernames[itrg].Data());
198   }
199   ax->SetRangeUser(-0.5, nTriggerClasses+.5);
200   h->GetYaxis()->SetRangeUser(0,1);
201
202   // track status
203   h=dynamic_cast<TH1F*>(fContainer->At(kTrackStatus));
204   Float_t ok = h->GetBinContent(1);
205   Int_t nerr = h->GetNbinsX();
206   for(Int_t ierr=nerr; ierr--;){
207     h->SetBinContent(ierr+1, ok>0.?1.e2*h->GetBinContent(ierr+1)/ok:0.);
208   }
209   h->SetBinContent(1, 0.);
210
211   // tracklet status
212   
213   TObjArray *arr = dynamic_cast<TObjArray*>(fContainer->UncheckedAt(kTrackletStatus));
214   for(Int_t ily = AliTRDgeometry::kNlayer; ily--;){
215     h=dynamic_cast<TH1F*>(arr->At(ily));
216     Float_t okB = h->GetBinContent(1);
217     Int_t nerrB = h->GetNbinsX();
218     for(Int_t ierr=nerrB; ierr--;){
219       h->SetBinContent(ierr+1, okB>0.?1.e2*h->GetBinContent(ierr+1)/okB:0.);
220     }
221     h->SetBinContent(1, 0.);
222   }
223
224   fNRefFigures = 17;
225
226   return kTRUE;
227 }
228
229 //_______________________________________________________
230 void AliTRDcheckDET::MakeSummary(){
231   //
232   // Create summary plots for TRD check DET
233   // This function creates 2 summary plots:
234   // - General Quantities
235   // - PHS
236   // The function will reuse GetRefFigure
237   //
238   
239   TCanvas *cOut = new TCanvas(Form("summary%s1", GetName()), Form("Summary 1 for task %s", GetName()), 1024, 768);
240   cOut->Divide(3,3);
241
242   // Create figures using GetRefFigure
243   cOut->cd(1); GetRefFigure(kFigNtracksEvent);  
244   cOut->cd(2); GetRefFigure(kFigNtracksSector);
245   cOut->cd(3); GetRefFigure(kFigNclustersTrack);
246   cOut->cd(4); GetRefFigure(kFigNclustersTracklet);
247   cOut->cd(5); GetRefFigure(kFigNtrackletsTrack);
248   cOut->cd(6); GetRefFigure(kFigNTrackletsP);
249   cOut->cd(7); GetRefFigure(kFigChargeCluster);
250   cOut->cd(8); GetRefFigure(kFigChargeTracklet);
251   cOut->SaveAs(Form("TRDsummary%s1.gif", GetName()));
252   delete cOut;
253
254   // Second Plot: PHS
255   cOut = new TCanvas(Form("summary%s2", GetName()), Form("Summary 2 for task %s", GetName()), 1024, 512);
256   cOut->cd(); GetRefFigure(kFigPH);
257   cOut->SaveAs(Form("TRDsummary%s2.gif", GetName())); 
258   delete cOut;
259
260   // Third Plot: Mean Number of clusters as function of eta, phi and layer
261    cOut = new TCanvas(Form("summary%s3", GetName()), Form("Summary 3 for task %s", GetName()), 1024, 768);
262   cOut->cd(); MakePlotMeanClustersLayer();
263   cOut->SaveAs(Form("TRDsummary%s3.gif", GetName())); 
264   delete cOut;
265
266 }
267
268 //_______________________________________________________
269 Bool_t AliTRDcheckDET::GetRefFigure(Int_t ifig){
270   //
271   // Setting Reference Figures
272   //
273   gPad->SetLogy(0);
274   gPad->SetLogx(0);
275   TH1 *h = NULL; TObjArray *arr=NULL;
276   TLegend *leg = NULL;
277   Bool_t kFIRST(1);
278   switch(ifig){
279   case kFigNclustersTrack:
280     (h=(TH1F*)fContainer->FindObject("hNcls"))->Draw("pl");
281     PutTrendValue("NClustersTrack", h->GetMean());
282     PutTrendValue("NClustersTrackRMS", h->GetRMS());
283     return kTRUE;
284   case kFigNclustersTracklet:
285     (h =(TH1F*)fContainer->FindObject("hNclTls"))->Draw("pc");
286     PutTrendValue("NClustersTracklet", h->GetMean());
287     PutTrendValue("NClustersTrackletRMS", h->GetRMS());
288     return kTRUE;
289   case kFigNtrackletsTrack:
290     h=MakePlotNTracklets();
291     PutTrendValue("NTrackletsTrack", h->GetMean());
292     PutTrendValue("NTrackletsTrackRMS", h->GetRMS());
293     return kTRUE;
294   case kFigNTrackletsP:
295     MakePlotnTrackletsVsP();
296     return kTRUE;
297   case kFigNtrackletsCross:
298     h = (TH1F*)fContainer->FindObject("hNtlsCross");
299     if(!MakeBarPlot(h, kRed)) break;
300     PutTrendValue("NTrackletsCross", h->GetMean());
301     PutTrendValue("NTrackletsCrossRMS", h->GetRMS());
302     return kTRUE;
303   case kFigNtrackletsFindable:
304     h = (TH1F*)fContainer->FindObject("hNtlsFindable");
305     if(!MakeBarPlot(h, kGreen)) break;
306     PutTrendValue("NTrackletsFindable", h->GetMean());
307     PutTrendValue("NTrackletsFindableRMS", h->GetRMS());
308     return kTRUE;
309   case kFigNtracksEvent:
310     (h = (TH1F*)fContainer->FindObject("hNtrks"))->Draw("pl");
311     PutTrendValue("NTracksEvent", h->GetMean());
312     PutTrendValue("NTracksEventRMS", h->GetRMS());
313     return kTRUE;
314   case kFigNtracksSector:
315     h = (TH1F*)fContainer->FindObject("hNtrksSector");
316     if(!MakeBarPlot(h, kGreen)) break;
317     PutTrendValue("NTracksSector", h->Integral()/h->GetNbinsX());
318     return kTRUE;
319   case kFigTrackStatus:
320     if(!(h=(TH1F *)fContainer->FindObject("hTrackStatus"))) break;
321     h->GetXaxis()->SetRangeUser(0.5, -1);
322     h->GetYaxis()->CenterTitle();
323     h->Draw("c");
324     PutTrendValue("TrackStatus", h->Integral());
325     gPad->SetLogy(0);
326     return kTRUE;
327   case kFigTrackletStatus:
328     if(!(arr = dynamic_cast<TObjArray*>(fContainer->At(kTrackletStatus)))) break;
329     leg = new TLegend(.68, .7, .97, .97);
330     leg->SetBorderSize(0);leg->SetFillStyle(0);
331     leg->SetHeader("TRD layer");
332     for(Int_t ily=AliTRDgeometry::kNlayer; ily--;){
333       if(!(h=dynamic_cast<TH1F*>(arr->At(ily)))) continue;
334       if(kFIRST){
335         h->Draw("pl");
336         h->GetXaxis()->SetRangeUser(0.5, -1);
337         h->GetYaxis()->CenterTitle();
338         kFIRST = kFALSE;
339       } else h->Draw("samepl");
340       leg->AddEntry(h, Form("ly = %d", ily), "l");
341       PutTrendValue(Form("TrackletStatus%d", ily), h->Integral());
342     }
343     leg->Draw();
344     gPad->SetLogy(0);
345     return kTRUE;
346   case kFigChi2:
347     return kTRUE;
348     MakePlotChi2();
349     return kTRUE;
350   case kFigPH:
351     gPad->SetMargin(0.125, 0.015, 0.1, 0.1);
352     MakePlotPulseHeight();
353     gPad->SetLogy(0);
354     return kTRUE;
355   case kFigChargeCluster:
356     (h = (TH1F*)fContainer->FindObject("hQcl"))->Draw("c");
357     gPad->SetLogy(1);
358     PutTrendValue("ChargeCluster", h->GetMaximumBin());
359     PutTrendValue("ChargeClusterRMS", h->GetRMS());
360     return kTRUE;
361   case kFigChargeTracklet:
362     (h=(TH1F*)fContainer->FindObject("hQtrklt"))->Draw("c");
363     PutTrendValue("ChargeTracklet", h->GetMaximumBin());
364     PutTrendValue("ChargeTrackletRMS", h->GetRMS());
365     return kTRUE;
366   case kFigNeventsTrigger:
367     ((TH1F*)fContainer->FindObject("hEventsTrigger"))->Draw("");
368     return kTRUE;
369   case kFigNeventsTriggerTracks:
370     ((TH1F*)fContainer->FindObject("hEventsTriggerTracks"))->Draw("");
371     return kTRUE;
372   case kFigTriggerPurity: 
373     if(!MakeBarPlot((TH1F*)fContainer->FindObject("hTriggerPurity"), kGreen)) break;
374     break;
375   default:
376     break;
377   }
378   AliInfo(Form("Reference plot [%d] missing result", ifig));
379   return kFALSE;
380 }
381
382 //_______________________________________________________
383 TObjArray *AliTRDcheckDET::Histos(){
384   //
385   // Create QA histograms
386   //
387     
388   if(fContainer) return fContainer;
389   
390   fContainer = new TObjArray(20);
391   //fContainer->SetOwner(kTRUE);
392
393   // Register Histograms
394   TH1 * h = NULL;
395   TAxis *ax = NULL;
396   if(!(h = (TH1F *)gROOT->FindObject("hNcls"))){
397     h = new TH1F("hNcls", "N_{clusters} / track", 181, -0.5, 180.5);
398     h->GetXaxis()->SetTitle("N_{clusters}");
399     h->GetYaxis()->SetTitle("Entries");
400   } else h->Reset();
401   fContainer->AddAt(h, kNclustersTrack);
402
403   TObjArray *arr = new TObjArray(AliTRDgeometry::kNlayer);
404   arr->SetOwner(kTRUE);  arr->SetName("clusters");
405   fContainer->AddAt(arr, kNclustersLayer);
406   for(Int_t ily=AliTRDgeometry::kNlayer; ily--;){
407     if(!(h = (TProfile2D *)gROOT->FindObject(Form("hNcl%d", ily)))){
408       h = new TProfile2D(Form("hNcl%d", ily), Form("Mean Number of clusters in Layer %d", ily), 100, -1.0, 1.0, 50, -1.1*TMath::Pi(), 1.1*TMath::Pi());
409       h->GetXaxis()->SetTitle("#eta");
410       h->GetYaxis()->SetTitle("#phi");
411     } else h->Reset();
412     arr->AddAt(h, ily);
413   }
414
415   if(!(h = (TH1F *)gROOT->FindObject("hNclTls"))){
416     h = new TH1F("hNclTls","N_{clusters} / tracklet", 51, -0.5, 50.5);
417     h->GetXaxis()->SetTitle("N_{clusters}");
418     h->GetYaxis()->SetTitle("Entries");
419   } else h->Reset();
420   fContainer->AddAt(h, kNclustersTracklet);
421
422   if(!(h = (TH1F *)gROOT->FindObject("hNtls"))){
423     h = new TH1F("hNtls", "N_{tracklets} / track", AliTRDgeometry::kNlayer, 0.5, 6.5);
424     h->GetXaxis()->SetTitle("N^{tracklet}");
425     h->GetYaxis()->SetTitle("freq. [%]");
426   } else h->Reset();
427   fContainer->AddAt(h, kNtrackletsTrack);
428
429   if(!(h = (TH1F *)gROOT->FindObject("htlsSTA"))){
430     h = new TH1F("hNtlsSTA", "N_{tracklets} / track (Stand Alone)", AliTRDgeometry::kNlayer, 0.5, 6.5);
431     h->GetXaxis()->SetTitle("N^{tracklet}");
432     h->GetYaxis()->SetTitle("freq. [%]");
433   }
434   fContainer->AddAt(h, kNtrackletsSTA);
435
436   // Binning for momentum dependent tracklet Plots
437   const Int_t kNp(30);
438   Float_t P=0.2;
439   Float_t binsP[kNp+1], binsTrklt[AliTRDgeometry::kNlayer+1];
440   for(Int_t i=0;i<kNp+1; i++,P+=(TMath::Exp(i*i*.001)-1.)) binsP[i]=P;
441   for(Int_t il = 0; il <= AliTRDgeometry::kNlayer; il++) binsTrklt[il] = 0.5 + il;
442   if(!(h = (TH1F *)gROOT->FindObject("htlsBAR"))){
443     // Make tracklets for barrel tracks momentum dependent (if we do not exceed min and max values)
444     h = new TH2F("hNtlsBAR", 
445     "N_{tracklets} / track;p [GeV/c];N^{tracklet};freq. [%]", 
446     kNp, binsP, AliTRDgeometry::kNlayer, binsTrklt);
447   }
448   fContainer->AddAt(h, kNtrackletsBAR);
449
450   // 
451   if(!(h = (TH1F *)gROOT->FindObject("hNtlsCross"))){
452     h = new TH1F("hNtlsCross", "N_{tracklets}^{cross} / track", 7, -0.5, 6.5);
453     h->GetXaxis()->SetTitle("n_{row cross}");
454     h->GetYaxis()->SetTitle("freq. [%]");
455   } else h->Reset();
456   fContainer->AddAt(h, kNtrackletsCross);
457
458   if(!(h = (TH1F *)gROOT->FindObject("hNtlsFindable"))){
459     h = new TH1F("hNtlsFindable", "Found/Findable Tracklets" , 101, -0.005, 1.005);
460     h->GetXaxis()->SetTitle("r [a.u]");
461     h->GetYaxis()->SetTitle("Entries");
462   } else h->Reset();
463   fContainer->AddAt(h, kNtrackletsFindable);
464
465   if(!(h = (TH1F *)gROOT->FindObject("hNtrks"))){
466     h = new TH1F("hNtrks", "N_{tracks} / event", 100, 0, 100);
467     h->GetXaxis()->SetTitle("N_{tracks}");
468     h->GetYaxis()->SetTitle("Entries");
469   } else h->Reset();
470   fContainer->AddAt(h, kNtracksEvent);
471
472   if(!(h = (TH1F *)gROOT->FindObject("hNtrksSector"))){
473     h = new TH1F("hNtrksSector", "N_{tracks} / sector", AliTRDgeometry::kNsector, -0.5, 17.5);
474     h->GetXaxis()->SetTitle("sector");
475     h->GetYaxis()->SetTitle("freq. [%]");
476   } else h->Reset();
477   fContainer->AddAt(h, kNtracksSector);
478
479   if(!(h = (TH1F*)gROOT->FindObject("hTrackStatus"))){
480     const Int_t nerr = 7;
481     h = new TH1F("hTrackStatus", "Track Status", nerr, -0.5, nerr-0.5);
482     const Char_t *label[nerr] = {"OK", "PROL", "PROP", "AJST", "SNP", "TINI", "UPDT"};
483     ax = h->GetXaxis();
484     for(Int_t ierr=nerr; ierr--;) ax->SetBinLabel(ierr+1, label[ierr]);
485     h->SetYTitle("Relative Error to Good [%]");
486   }
487   fContainer->AddAt(h, kTrackStatus);
488
489   arr = new TObjArray(AliTRDgeometry::kNlayer);
490   arr->SetOwner(kTRUE);  arr->SetName("TrackletStatus");
491   fContainer->AddAt(arr, kTrackletStatus);
492   for(Int_t ily=AliTRDgeometry::kNlayer; ily--;){
493     if(!(h = (TH1F *)gROOT->FindObject(Form("hTrackletStatus%d", ily)))){
494       const Int_t nerr = 8;
495       h = new TH1F(Form("hTrackletStatus%d", ily), "Tracklet status", nerr, -0.5, nerr-0.5);
496       h->SetLineColor(ily+1);
497       const Char_t *label[nerr] = {"OK", "Geom", "Bound", "NoCl", "NoAttach", "NoClTr", "NoFit", "Chi2"};
498       ax = h->GetXaxis();
499       for(Int_t ierr=nerr; ierr--;) ax->SetBinLabel(ierr+1, label[ierr]);
500       h->SetYTitle("Relative Error to Good [%]");
501     } else h->Reset();
502     arr->AddAt(h, ily);
503   }
504
505   // <PH> histos
506   arr = new TObjArray(3);
507   arr->SetOwner(kTRUE);  arr->SetName("<PH>");
508   fContainer->AddAt(arr, kPH);
509   if(!(h = (TH1F *)gROOT->FindObject("hPHt"))){
510     h = new TProfile("hPHt", "<PH>", 31, -0.5, 30.5);
511     h->GetXaxis()->SetTitle("Time / 100ns");
512     h->GetYaxis()->SetTitle("<PH> [a.u]");
513   } else h->Reset();
514   arr->AddAt(h, 0);
515   if(!(h = (TH1F *)gROOT->FindObject("hPHx")))
516     h = new TProfile("hPHx", "<PH>", 31, -0.08, 4.88);
517   else h->Reset();
518   arr->AddAt(h, 1);
519   if(!(h = (TH2F *)gROOT->FindObject("hPH2D"))){
520     h = new TH2F("hPH2D", "Charge Distribution / time", 31, -0.5, 30.5, 100, 0, 1024);
521     h->GetXaxis()->SetTitle("Time / 100ns");
522     h->GetYaxis()->SetTitle("Charge / a.u.");
523   } else h->Reset();
524   arr->AddAt(h, 2);
525
526   // Chi2 histos
527   if(!(h = (TH2S*)gROOT->FindObject("hChi2"))){
528     h = new TH2S("hChi2", "#chi^{2} per track", AliTRDgeometry::kNlayer, .5, AliTRDgeometry::kNlayer+.5, 100, 0, 50);
529     h->SetXTitle("ndf");
530     h->SetYTitle("#chi^{2}/ndf");
531     h->SetZTitle("entries");
532   } else h->Reset();
533   fContainer->AddAt(h, kChi2);
534
535   if(!(h = (TH1F *)gROOT->FindObject("hQcl"))){
536     h = new TH1F("hQcl", "Q_{cluster}", 200, 0, 1200);
537     h->GetXaxis()->SetTitle("Q_{cluster} [a.u.]");
538     h->GetYaxis()->SetTitle("Entries");
539   }else h->Reset();
540   fContainer->AddAt(h, kChargeCluster);
541
542   if(!(h = (TH1F *)gROOT->FindObject("hQtrklt"))){
543     h = new TH1F("hQtrklt", "Q_{tracklet}", 6000, 0, 6000);
544     h->GetXaxis()->SetTitle("Q_{tracklet} [a.u.]");
545     h->GetYaxis()->SetTitle("Entries");
546   }else h->Reset();
547   fContainer->AddAt(h, kChargeTracklet);
548
549
550   if(!(h = (TH1F *)gROOT->FindObject("hEventsTrigger")))
551     h = new TH1F("hEventsTrigger", "Trigger Class", 100, 0, 100);
552   else h->Reset();
553   printf("Histos Adding \n");
554   
555   fContainer->AddAt(h, kNeventsTrigger);
556
557   if(!(h = (TH1F *)gROOT->FindObject("hEventsTriggerTracks")))
558     h = new TH1F("hEventsTriggerTracks", "Trigger Class (Tracks)", 100, 0, 100);
559   else h->Reset();
560   fContainer->AddAt(h, kNeventsTriggerTracks);
561
562   if(!(h = (TH1F *)gROOT->FindObject("hTriggerPurity"))){
563     h = new TH1F("hTriggerPurity", "Trigger Purity", 10, -0.5, 9.5);
564     h->GetXaxis()->SetTitle("Trigger Cluster");
565     h->GetYaxis()->SetTitle("freq.");
566   } else h->Reset();
567   fContainer->AddAt(h, kTriggerPurity);
568
569   return fContainer;
570 }
571
572 /*
573 * Plotting Functions
574 */
575
576 //_______________________________________________________
577 TH1 *AliTRDcheckDET::PlotTrackStatus(const AliTRDtrackV1 *track)
578 {
579 //
580 // Plot the track propagation status. The following errors are defined (see AliTRDtrackV1::ETRDtrackError)
581 //   PROL - track prolongation failure
582 //   PROP - track propagation failure
583 //   AJST - crossing sectors failure
584 //   SNP  - too large bending
585 //   TINI - tracklet initialization failure
586 //   UPDT - track position/covariance update failure 
587 //
588 // Performance plot looks as below:
589 //Begin_Html
590 //<img src="TRD/trackStatus.gif">
591 //End_Html
592 //
593   if(track) fkTrack = track;
594   if(!fkTrack){
595     AliDebug(4, "No Track defined.");
596     return NULL;
597   }
598   TH1 *h = NULL;
599   if(!(h = dynamic_cast<TH1F *>(fContainer->At(kTrackStatus)))){
600     AliWarning("No Histogram defined.");
601     return NULL;
602   }
603   h->Fill(fkTrack->GetStatusTRD());
604   return h;
605 }
606
607 //_______________________________________________________
608 TH1 *AliTRDcheckDET::PlotTrackletStatus(const AliTRDtrackV1 *track)
609 {
610 //
611 // Plot the tracklet propagation status. The following errors are defined for tracklet (see AliTRDtrackV1::ETRDlayerError)
612 //   Geom   - 
613 //   Bound  - tracklet too close to chamber walls
614 //   NoCl   - no clusters in the track roads
615 //   NoAttach - fail to attach clusters
616 //   NoClTr - fail to use clusters for fit
617 //   NoFit  - tracklet fit failled
618 //   Chi2   - chi2 tracklet-track over threshold
619 //
620 // Performance plot looks as below:
621 //Begin_Html
622 //<img src="TRD/trackletStatus.gif">
623 //End_Html
624 //
625   if(track) fkTrack = track;
626   if(!fkTrack){
627     AliDebug(4, "No Track defined.");
628     return NULL;
629   }
630   TObjArray *arr =NULL;
631   if(!(arr = dynamic_cast<TObjArray*>(fContainer->At(kTrackletStatus)))){
632     AliWarning("Histograms not defined.");
633     return NULL;
634   }
635
636   TH1 *h = NULL;
637   for(Int_t ily=AliTRDgeometry::kNlayer; ily--;){
638     if(!(h = dynamic_cast<TH1F*>(arr->At(ily)))){
639       AliWarning(Form("Missing histo for layer %d.", ily));
640       continue;
641     }
642     h->Fill(fkTrack->GetStatusTRD(ily));
643   }
644   return h;
645 }
646
647 //_______________________________________________________
648 TH1 *AliTRDcheckDET::PlotNClustersTracklet(const AliTRDtrackV1 *track){
649   //
650   // Plot the mean number of clusters per tracklet
651   //
652   if(track) fkTrack = track;
653   if(!fkTrack){
654     AliDebug(4, "No Track defined.");
655     return NULL;
656   }
657   AliExternalTrackParam *par = fkTrack->GetTrackIn() ? fkTrack->GetTrackIn() : fkTrack->GetTrackOut();
658   TH1 *h = NULL;
659   TProfile2D *hlayer = NULL;
660   Double_t eta = 0., phi = 0.;
661   if(!(h = dynamic_cast<TH1F *>(fContainer->At(kNclustersTracklet)))){
662     AliWarning("No Histogram defined.");
663     return NULL;
664   }
665   AliTRDseedV1 *tracklet = NULL;
666   TObjArray *histosLayer = dynamic_cast<TObjArray *>(fContainer->At(kNclustersLayer));
667   if(!histosLayer){
668     AliWarning("No Histograms for single layer defined");
669   }
670   for(Int_t itl = 0; itl < AliTRDgeometry::kNlayer; itl++){
671     if(!(tracklet = fkTrack->GetTracklet(itl)) || !tracklet->IsOK()) continue;
672     h->Fill(tracklet->GetN2());
673     if(histosLayer && par){
674       if((hlayer = dynamic_cast<TProfile2D *>(histosLayer->At(itl)))){
675         GetEtaPhiAt(par, tracklet->GetX0(), eta, phi);
676         hlayer->Fill(eta, phi, tracklet->GetN2());
677       }
678     }
679   }
680   return h;
681 }
682
683 //_______________________________________________________
684 TH1 *AliTRDcheckDET::PlotNClustersTrack(const AliTRDtrackV1 *track){
685   //
686   // Plot the number of clusters in one track
687   //
688   if(track) fkTrack = track;
689   if(!fkTrack){
690     AliDebug(4, "No Track defined.");
691     return NULL;
692   }
693   TH1 *h = NULL;
694   if(!(h = dynamic_cast<TH1F *>(fContainer->At(kNclustersTrack)))){
695     AliWarning("No Histogram defined.");
696     return NULL;
697   }
698   
699   Int_t nclusters = 0;
700   AliTRDseedV1 *tracklet = NULL;
701   AliExternalTrackParam *par = fkTrack->GetTrackOut() ? fkTrack->GetTrackOut() : fkTrack->GetTrackIn();
702   if(!par) return NULL;
703   Double_t momentumRec = par->P();
704   for(Int_t itl = 0; itl < AliTRDgeometry::kNlayer; itl++){
705     if(!(tracklet = fkTrack->GetTracklet(itl)) || !tracklet->IsOK()) continue;
706     Int_t n(tracklet->GetN());
707     nclusters += n;
708     if(DebugLevel() > 2){
709       Int_t crossing = Int_t(tracklet->IsRowCross());
710       Int_t detector = tracklet->GetDetector();
711       Float_t theta = TMath::ATan(tracklet->GetZref(1));
712       Float_t phi = TMath::ATan(tracklet->GetYref(1));
713       Float_t momentumMC = 0.;
714       Int_t pdg = 0;
715       Int_t kinkIndex = fkESD ? fkESD->GetKinkIndex() : 0;
716       UShort_t nclsTPC = fkESD ? fkESD->GetTPCncls() : 0;
717       if(fkMC){
718         if(fkMC->GetTrackRef()) momentumMC = fkMC->GetTrackRef()->P();
719         pdg = fkMC->GetPDG();
720       }
721       (*DebugStream()) << "NClustersTrack"
722         << "Detector="  << detector
723         << "crossing="  << crossing
724         << "momentumMC="<< momentumMC
725         << "momentumRec=" << momentumRec
726         << "pdg="                               << pdg
727         << "theta="                     << theta
728         << "phi="                               << phi
729         << "kinkIndex=" << kinkIndex
730         << "TPCncls="           << nclsTPC
731         << "TRDncls="   << n
732         << "\n";
733     }
734   }
735   h->Fill(nclusters);
736   return h;
737 }
738
739
740 //_______________________________________________________
741 TH1 *AliTRDcheckDET::PlotNTrackletsTrack(const AliTRDtrackV1 *track){
742   //
743   // Plot the number of tracklets
744   //
745   if(track) fkTrack = track;
746   if(!fkTrack){
747     AliDebug(4, "No Track defined.");
748     return NULL;
749   }
750   TH1 *h = NULL, *hSta = NULL; TH2 *hBarrel = NULL;
751   if(!(h = dynamic_cast<TH1F *>(fContainer->At(kNtrackletsTrack)))){
752     AliWarning("No Histogram defined.");
753     return NULL;
754   }
755   Int_t nTracklets = fkTrack->GetNumberOfTracklets();
756   h->Fill(nTracklets);
757   if(!fkESD) return h;
758   Int_t status = fkESD->GetStatus();
759
760 /*  printf("in/out/refit/pid: TRD[%d|%d|%d|%d]\n", status &AliESDtrack::kTRDin ? 1 : 0, status &AliESDtrack::kTRDout ? 1 : 0, status &AliESDtrack::kTRDrefit ? 1 : 0, status &AliESDtrack::kTRDpid ? 1 : 0);*/
761   Double_t p = 0.;
762   Int_t method = -1;    // to distinguish between stand alone and full barrel tracks in the debugging
763   if((status & AliESDtrack::kTRDin) != 0){
764     method = 1;
765     // Full Barrel Track: Save momentum dependence
766     if(!(hBarrel = dynamic_cast<TH2F *>(fContainer->At(kNtrackletsBAR)))){
767       AliWarning("Method: Barrel.  Histogram not processed!");
768       return NULL;
769     }
770     AliExternalTrackParam *par(fkTrack->GetTrackIn());
771     if(!par){
772       AliError("Input track params missing");
773       return NULL;
774     }
775     p = par->P(); // p needed later in the debug streaming
776     hBarrel->Fill(p, nTracklets);
777   } else {
778     // Stand alone Track: momentum dependence not usefull
779     method = 0;
780     if(!(hSta = dynamic_cast<TH1F *>(fContainer->At(kNtrackletsSTA)))) {
781       AliWarning("Method: StandAlone.  Histogram not processed!");
782       return NULL;
783     }
784     hSta->Fill(nTracklets);
785   }
786
787   if(DebugLevel() > 2){
788     AliTRDseedV1 *tracklet = NULL;
789     Int_t sector = -1, stack = -1, detector;
790     for(Int_t itl = 0; itl < AliTRDgeometry::kNlayer; itl++){
791       if(!(tracklet = fkTrack->GetTracklet(itl)) || !(tracklet->IsOK())) continue;
792       detector = tracklet->GetDetector();
793       sector = fGeo->GetSector(detector);
794       stack = fGeo->GetStack(detector);
795       break;
796     }
797     (*DebugStream()) << "NTrackletsTrack"
798       << "Sector="      << sector
799       << "Stack="        << stack 
800       << "NTracklets="  << nTracklets
801       << "Method="      << method
802       << "p="           << p
803       << "\n";
804   }
805   if(DebugLevel() > 3){
806     AliTRDseedV1 *tracklet = NULL;
807     for(Int_t il = 0; il < AliTRDgeometry::kNlayer; il++){
808       if((tracklet = fkTrack->GetTracklet(il)) && tracklet->IsOK()){
809         (*DebugStream()) << "NTrackletsLayer"
810         << "Layer=" << il
811         << "p=" << p
812         << "\n";
813       }
814     }
815   }
816   return h;
817 }
818
819
820 //_______________________________________________________
821 TH1 *AliTRDcheckDET::PlotNTrackletsRowCross(const AliTRDtrackV1 *track){
822   //
823   // Plot the number of tracklets
824   //
825   if(track) fkTrack = track;
826   if(!fkTrack){
827     AliDebug(4, "No Track defined.");
828     return NULL;
829   }
830   TH1 *h = NULL;
831   if(!(h = dynamic_cast<TH1F *>(fContainer->At(kNtrackletsCross)))){
832     AliWarning("No Histogram defined.");
833     return NULL;
834   }
835
836   Int_t ncross = 0;
837   AliTRDseedV1 *tracklet = NULL;
838   for(Int_t il = 0; il < AliTRDgeometry::kNlayer; il++){
839     if(!(tracklet = fkTrack->GetTracklet(il)) || !tracklet->IsOK()) continue;
840
841     if(tracklet->IsRowCross()) ncross++;
842   }
843   h->Fill(ncross);
844   return h;
845 }
846
847 //_______________________________________________________
848 TH1 *AliTRDcheckDET::PlotFindableTracklets(const AliTRDtrackV1 *track){
849   //
850   // Plots the ratio of number of tracklets vs.
851   // number of findable tracklets
852   //
853   // Findable tracklets are defined as track prolongation
854   // to layer i does not hit the dead area +- epsilon
855   //
856   // In order to check whether tracklet hist active area in Layer i, 
857   // the track is refitted and the fitted position + an uncertainty 
858   // range is compared to the chamber border (also with a different
859   // uncertainty)
860   //
861   // For the track fit two cases are distinguished:
862   // If the track is a stand alone track (defined by the status bit 
863   // encoding, then the track is fitted with the tilted Rieman model
864   // Otherwise the track is fitted with the Kalman fitter in two steps:
865   // Since the track parameters are give at the outer point, we first 
866   // fit in direction inwards. Afterwards we fit again in direction outwards
867   // to extrapolate the track to layers which are not reached by the track
868   // For the Kalman model, the radial track points have to be shifted by
869   // a distance epsilon in the direction that we want to fit
870   //
871   const Float_t epsilon = 0.01;   // dead area tolerance
872   const Float_t epsilonR = 1;    // shift in radial direction of the anode wire position (Kalman filter only)
873   const Float_t deltaY = 0.7;    // Tolerance in the track position in y-direction
874   const Float_t deltaZ = 7.0;    // Tolerance in the track position in z-direction (Padlength)
875   Double_t xAnode[AliTRDgeometry::kNlayer] = {300.2, 312.8, 325.4, 338.0, 350.6, 363.2}; // Take the default X0
876  
877   if(track) fkTrack = track;
878   if(!fkTrack){
879     AliDebug(4, "No Track defined.");
880     return NULL;
881   }
882   TH1 *h = NULL;
883   if(!(h = dynamic_cast<TH1F *>(fContainer->At(kNtrackletsFindable)))){
884     AliWarning("No Histogram defined.");
885     return NULL;
886   }
887   Int_t nFound = 0, nFindable = 0;
888   Int_t stack = -1;
889   Double_t ymin = 0., ymax = 0., zmin = 0., zmax = 0.;
890   Double_t y = 0., z = 0.;
891   AliTRDseedV1 *tracklet = NULL;
892   AliTRDpadPlane *pp;  
893   for(Int_t il = 0; il < AliTRDgeometry::kNlayer; il++){
894     if((tracklet = fkTrack->GetTracklet(il)) && tracklet->IsOK()){
895       tracklet->SetReconstructor(fReconstructor);
896       nFound++;
897     }
898   }
899   // 2 Different cases:
900   // 1st stand alone: here we cannot propagate, but be can do a Tilted Rieman Fit
901   // 2nd barrel track: here we propagate the track to the layers
902   AliTrackPoint points[6];
903   Float_t xyz[3];
904   memset(xyz, 0, sizeof(Float_t) * 3);
905   if(((fkESD->GetStatus() & AliESDtrack::kTRDout) > 0) && !((fkESD->GetStatus() & AliESDtrack::kTRDin) > 0)){
906     // stand alone track
907     for(Int_t il = 0; il < AliTRDgeometry::kNlayer; il++){
908       xyz[0] = xAnode[il];
909       points[il].SetXYZ(xyz);
910     }
911     AliTRDtrackerV1::FitRiemanTilt(const_cast<AliTRDtrackV1 *>(fkTrack), NULL, kTRUE, 6, points);
912   } else {
913     // barrel track
914     //
915     // 2 Steps:
916     // -> Kalman inwards
917     // -> Kalman outwards
918     AliTRDtrackV1 copyTrack(*fkTrack);  // Do Kalman on a (non-constant) copy of the track
919     AliTrackPoint pointsInward[6], pointsOutward[6];
920     for(Int_t il = AliTRDgeometry::kNlayer; il--;){
921       // In order to avoid complications in the Kalman filter if the track points have the same radial
922       // position like the tracklets, we have to shift the radial postion of the anode wire by epsilon
923       // in the direction we want to go
924       // The track points have to be in reverse order for the Kalman Filter inwards
925       xyz[0] = xAnode[AliTRDgeometry::kNlayer - il - 1] - epsilonR;
926       pointsInward[il].SetXYZ(xyz);
927       xyz[0] = xAnode[il] + epsilonR;
928       pointsOutward[il].SetXYZ(xyz);
929     }
930     /*for(Int_t ipt = 0; ipt < AliTRDgeometry::kNlayer; ipt++)
931       printf("%d. X = %f\n", ipt, points[ipt].GetX());*/
932     // Kalman inwards
933     AliTRDtrackerV1::FitKalman(&copyTrack, NULL, kFALSE, 6, pointsInward);
934     memcpy(points, pointsInward, sizeof(AliTrackPoint) * 6); // Preliminary store the inward results in the Array points
935     // Kalman outwards
936     AliTRDtrackerV1::FitKalman(&copyTrack, NULL, kTRUE, 6, pointsInward);
937     memcpy(points, pointsOutward, sizeof(AliTrackPoint) * AliTRDgeometry::kNlayer);
938   }
939   for(Int_t il = 0; il < AliTRDgeometry::kNlayer; il++){
940     y = points[il].GetY();
941     z = points[il].GetZ();
942     if((stack = fGeo->GetStack(z, il)) < 0) continue; // Not findable
943     pp = fGeo->GetPadPlane(il, stack);
944     ymin = pp->GetCol0() + epsilon;
945     ymax = pp->GetColEnd() - epsilon; 
946     zmin = pp->GetRowEnd() + epsilon; 
947     zmax = pp->GetRow0() - epsilon;
948     // ignore y-crossing (material)
949     if((z + deltaZ > zmin && z - deltaZ < zmax) && (y + deltaY > ymin && y - deltaY < ymax)) nFindable++;
950       if(DebugLevel() > 3){
951         Double_t posTracklet[2] = {tracklet ? tracklet->GetYfit(0) : 0, tracklet ? tracklet->GetZfit(0) : 0};
952         Int_t hasTracklet = tracklet ? 1 : 0;
953         (*DebugStream())   << "FindableTracklets"
954           << "layer="     << il
955           << "ytracklet=" << posTracklet[0]
956           << "ytrack="    << y
957           << "ztracklet=" << posTracklet[1]
958           << "ztrack="    << z
959           << "tracklet="  << hasTracklet
960           << "\n";
961       }
962   }
963   
964   h->Fill(nFindable > 0 ? TMath::Min(nFound/static_cast<Double_t>(nFindable), 1.) : 1);
965   AliDebug(2, Form("Findable[Found]: %d[%d|%f]", nFindable, nFound, nFound/static_cast<Float_t>(nFindable > 0 ? nFindable : 1)));
966   return h;
967 }
968
969
970 //_______________________________________________________
971 TH1 *AliTRDcheckDET::PlotChi2(const AliTRDtrackV1 *track){
972   //
973   // Plot the chi2 of the track
974   //
975   if(track) fkTrack = track;
976   if(!fkTrack){
977     AliDebug(4, "No Track defined.");
978     return NULL;
979   }
980   TH1 *h = NULL;
981   if(!(h = dynamic_cast<TH2S*>(fContainer->At(kChi2)))) {
982     AliWarning("No Histogram defined.");
983     return NULL;
984   }
985   Int_t n = fkTrack->GetNumberOfTracklets();
986   if(!n) return NULL;
987
988   h->Fill(n, fkTrack->GetChi2()/n);
989   return h;
990 }
991
992
993 //_______________________________________________________
994 TH1 *AliTRDcheckDET::PlotPHt(const AliTRDtrackV1 *track){
995   //
996   // Plot the average pulse height
997   //
998   if(track) fkTrack = track;
999   if(!fkTrack){
1000     AliDebug(4, "No Track defined.");
1001     return NULL;
1002   }
1003   TProfile *h = NULL; TH2F *phs2D = NULL;
1004   if(!(h = dynamic_cast<TProfile *>(((TObjArray*)(fContainer->At(kPH)))->At(0)))){
1005     AliWarning("No Histogram defined.");
1006     return NULL;
1007   }
1008   if(!(phs2D = dynamic_cast<TH2F *>(((TObjArray*)(fContainer->At(kPH)))->At(2)))){
1009     AliWarning("2D Pulse Height histogram not defined. Histogramm cannot be filled");
1010   }
1011   AliTRDseedV1 *tracklet = NULL;
1012   AliTRDcluster *c = NULL;
1013   for(Int_t itl = 0; itl < AliTRDgeometry::kNlayer; itl++){
1014     if(!(tracklet = fkTrack->GetTracklet(itl)) || !tracklet->IsOK())continue;
1015     Int_t crossing = Int_t(tracklet->IsRowCross());
1016     Int_t detector = tracklet->GetDetector();
1017     tracklet->ResetClusterIter();
1018     while((c = tracklet->NextCluster())){
1019       if(!IsUsingClustersOutsideChamber() && !c->IsInChamber()) continue;
1020       Int_t localtime        = c->GetLocalTimeBin();
1021       Double_t absoluteCharge = TMath::Abs(c->GetQ());
1022       h->Fill(localtime, absoluteCharge);
1023       phs2D->Fill(localtime, absoluteCharge); 
1024       if(DebugLevel() > 3){
1025         Int_t inChamber = c->IsInChamber() ? 1 : 0;
1026         Double_t distance[2];
1027         GetDistanceToTracklet(distance, tracklet, c);
1028         Float_t theta = TMath::ATan(tracklet->GetZref(1));
1029         Float_t phi = TMath::ATan(tracklet->GetYref(1));
1030         AliExternalTrackParam *trdPar = fkTrack->GetTrackIn();
1031         Float_t momentumMC = 0, momentumRec = trdPar ? trdPar->P() : track->P(); // prefer Track Low
1032         Int_t pdg = 0;
1033         Int_t kinkIndex = fkESD ? fkESD->GetKinkIndex() : 0;
1034         UShort_t TPCncls = fkESD ? fkESD->GetTPCncls() : 0;
1035         if(fkMC){
1036           if(fkMC->GetTrackRef()) momentumMC = fkMC->GetTrackRef()->P();
1037           pdg = fkMC->GetPDG();
1038         }
1039         (*DebugStream()) << "PHt"
1040           << "Detector="        << detector
1041           << "crossing="        << crossing
1042           << "inChamber=" << inChamber
1043           << "Timebin="         << localtime
1044           << "Charge="          << absoluteCharge
1045           << "momentumMC="      << momentumMC
1046           << "momentumRec="     << momentumRec
1047           << "pdg="                             << pdg
1048           << "theta="                   << theta
1049           << "phi="                             << phi
1050           << "kinkIndex="       << kinkIndex
1051           << "TPCncls="         << TPCncls
1052           << "dy="        << distance[0]
1053           << "dz="        << distance[1]
1054           << "c.="        << c
1055           << "\n";
1056       }
1057     }
1058   }
1059   return h;
1060 }
1061
1062 //_______________________________________________________
1063 TH1 *AliTRDcheckDET::PlotPHx(const AliTRDtrackV1 *track){
1064   //
1065   // Plots the average pulse height vs the distance from the anode wire
1066   // (plus const anode wire offset)
1067   //
1068   if(track) fkTrack = track;
1069   if(!fkTrack){
1070      AliDebug(4, "No Track defined.");
1071      return NULL;
1072   }
1073   TProfile *h = NULL;
1074   if(!(h = dynamic_cast<TProfile *>(((TObjArray*)(fContainer->At(kPH)))->At(1)))){
1075     AliWarning("No Histogram defined.");
1076     return NULL;
1077   }
1078   AliTRDseedV1 *tracklet(NULL);
1079   AliTRDcluster *c(NULL);
1080   Double_t xd(0.), dqdl(0.);
1081   TVectorD vq(AliTRDseedV1::kNtb), vxd(AliTRDseedV1::kNtb), vdqdl(AliTRDseedV1::kNtb);
1082   for(Int_t itl = 0; itl < AliTRDgeometry::kNlayer; itl++){
1083     if(!(tracklet = fkTrack->GetTracklet(itl)) || !(tracklet->IsOK())) continue;
1084     Int_t det(tracklet->GetDetector());
1085     Bool_t rc(tracklet->IsRowCross());
1086     for(Int_t ic(0); ic<AliTRDseedV1::kNtb; ic++){
1087       Bool_t kFIRST(kFALSE);
1088       if(!(c = tracklet->GetClusters(ic))){
1089          if(!(c = tracklet->GetClusters(AliTRDseedV1::kNtb+ic))) continue;
1090       } else kFIRST=kTRUE;
1091       if(!IsUsingClustersOutsideChamber() && !c->IsInChamber()) continue;
1092       xd = tracklet->GetX0() - c->GetX(); vxd[ic] = xd;
1093       dqdl=tracklet->GetdQdl(ic); vdqdl[ic] = dqdl;
1094       vq[ic]=c->GetQ();
1095       if(kFIRST && (c = tracklet->GetClusters(AliTRDseedV1::kNtb+ic))) vq[ic]+=c->GetQ();
1096       h->Fill(xd, dqdl);
1097     }
1098     if(DebugLevel() > 3){
1099       (*DebugStream()) << "PHx"
1100         << "det="  << det
1101         << "rc="   << rc
1102         << "xd="   << &vxd
1103         << "q="    << &vq
1104         << "dqdl=" << &vdqdl
1105         << "\n";
1106     }
1107   }  
1108   return h;
1109 }
1110
1111 //_______________________________________________________
1112 TH1 *AliTRDcheckDET::PlotChargeCluster(const AliTRDtrackV1 *track){
1113   //
1114   // Plot the cluster charge
1115   //
1116   if(track) fkTrack = track;
1117   if(!fkTrack){
1118     AliDebug(4, "No Track defined.");
1119     return NULL;
1120   }
1121   TH1 *h = NULL;
1122   if(!(h = dynamic_cast<TH1F *>(fContainer->At(kChargeCluster)))){
1123     AliWarning("No Histogram defined.");
1124     return NULL;
1125   }
1126   AliTRDseedV1 *tracklet = NULL;
1127   AliTRDcluster *c = NULL;
1128   for(Int_t itl = 0; itl < AliTRDgeometry::kNlayer; itl++){
1129     if(!(tracklet = fkTrack->GetTracklet(itl)) || !tracklet->IsOK())continue;
1130     for(Int_t ic(0); ic < AliTRDseedV1::kNtb; ic++){
1131       Bool_t kFIRST(kFALSE);
1132       if(!(c = tracklet->GetClusters(ic))) {
1133         if(!(c = tracklet->GetClusters(AliTRDseedV1::kNtb+ic))) continue;
1134       } else kFIRST = kTRUE;
1135       Float_t q(c->GetQ());
1136       if(kFIRST && (c = tracklet->GetClusters(AliTRDseedV1::kNtb+ic))) q+=c->GetQ();
1137       h->Fill(q);
1138     }
1139   }
1140   return h;
1141 }
1142
1143 //_______________________________________________________
1144 TH1 *AliTRDcheckDET::PlotChargeTracklet(const AliTRDtrackV1 *track){
1145   //
1146   // Plot the charge deposit per chamber
1147   //
1148   if(track) fkTrack = track;
1149   if(!fkTrack){
1150     AliDebug(4, "No Track defined.");
1151     return NULL;
1152   }
1153   TH1 *h = NULL;
1154   if(!(h = dynamic_cast<TH1F *>(fContainer->At(kChargeTracklet)))){
1155     AliWarning("No Histogram defined.");
1156     return NULL;
1157   }
1158   AliTRDseedV1 *tracklet = NULL;
1159   AliTRDcluster *c = NULL;
1160   Double_t qTot = 0;
1161   Int_t nTracklets =fkTrack->GetNumberOfTracklets();
1162   for(Int_t itl(0); itl < AliTRDgeometry::kNlayer; itl++){
1163     if(!(tracklet = fkTrack->GetTracklet(itl)) || !tracklet->IsOK()) continue;
1164     qTot = 0.;
1165     for(Int_t ic = AliTRDseedV1::kNclusters; ic--;){
1166       if(!(c = tracklet->GetClusters(ic))) continue;
1167       qTot += TMath::Abs(c->GetQ());
1168     }
1169     h->Fill(qTot);
1170     if(DebugLevel() > 3){
1171       Int_t crossing = (Int_t)tracklet->IsRowCross();
1172       Int_t detector = tracklet->GetDetector();
1173       Float_t theta = TMath::ATan(tracklet->GetZfit(1));
1174       Float_t phi = TMath::ATan(tracklet->GetYfit(1));
1175       Float_t momentum = 0.;
1176       Int_t pdg = 0;
1177       Int_t kinkIndex = fkESD ? fkESD->GetKinkIndex() : 0;
1178       UShort_t nclsTPC = fkESD ? fkESD->GetTPCncls() : 0;
1179       if(fkMC){
1180               if(fkMC->GetTrackRef()) momentum = fkMC->GetTrackRef()->P();
1181         pdg = fkMC->GetPDG();
1182       }
1183       (*DebugStream()) << "ChargeTracklet"
1184         << "Detector="  << detector
1185         << "crossing="  << crossing
1186         << "momentum="  << momentum
1187         << "nTracklets="<< nTracklets
1188         << "pdg="                               << pdg
1189         << "theta="                     << theta
1190         << "phi="                               << phi
1191         << "kinkIndex=" << kinkIndex
1192         << "TPCncls="           << nclsTPC
1193         << "QT="        << qTot
1194         << "\n";
1195     }
1196   }
1197   return h;
1198 }
1199
1200 //_______________________________________________________
1201 TH1 *AliTRDcheckDET::PlotNTracksSector(const AliTRDtrackV1 *track){
1202   //
1203   // Plot the number of tracks per Sector
1204   //
1205   if(track) fkTrack = track;
1206   if(!fkTrack){
1207     AliDebug(4, "No Track defined.");
1208     return NULL;
1209   }
1210   TH1 *h = NULL;
1211   if(!(h = dynamic_cast<TH1F *>(fContainer->At(kNtracksSector)))){
1212     AliWarning("No Histogram defined.");
1213     return NULL;
1214   }
1215
1216   // TODO we should compare with
1217   // sector = Int_t(track->GetAlpha() / AliTRDgeometry::GetAlpha());
1218
1219   AliTRDseedV1 *tracklet = NULL;
1220   Int_t sector = -1;
1221   for(Int_t itl = 0; itl < AliTRDgeometry::kNlayer; itl++){
1222     if(!(tracklet = fkTrack->GetTracklet(itl)) || !tracklet->IsOK()) continue;
1223     sector = static_cast<Int_t>(tracklet->GetDetector()/AliTRDgeometry::kNdets);
1224     break;
1225   }
1226   h->Fill(sector);
1227   return h;
1228 }
1229
1230
1231 //________________________________________________________
1232 void AliTRDcheckDET::SetRecoParam(AliTRDrecoParam *r)
1233 {
1234
1235   fReconstructor->SetRecoParam(r);
1236 }
1237
1238 //________________________________________________________
1239 void AliTRDcheckDET::GetDistanceToTracklet(Double_t *dist, AliTRDseedV1 * const tracklet, AliTRDcluster * const c)
1240 {
1241   Float_t x = c->GetX();
1242   dist[0] = c->GetY() - tracklet->GetYat(x);
1243   dist[1] = c->GetZ() - tracklet->GetZat(x);
1244 }
1245
1246 //________________________________________________________
1247 void AliTRDcheckDET::GetEtaPhiAt(AliExternalTrackParam *track, Double_t x, Double_t &eta, Double_t &phi){
1248   //
1249   // Get phi and eta at a given radial position
1250   // 
1251   AliExternalTrackParam workpar(*track);
1252
1253   Double_t posLocal[3];
1254   Bool_t sucPos = workpar.GetXYZAt(x, fEventInfo->GetRunInfo()->GetMagneticField(), posLocal);
1255   Double_t sagPhi = sucPos ? TMath::ATan2(posLocal[1], posLocal[0]) : 0.;
1256   phi = sagPhi;
1257   eta = workpar.Eta();
1258 }
1259
1260
1261 //_______________________________________________________
1262 TH1* AliTRDcheckDET::MakePlotChi2()
1263 {
1264 // Plot chi2/track normalized to number of degree of freedom 
1265 // (tracklets) and compare with the theoretical distribution.
1266 // 
1267 // Alex Bercuci <A.Bercuci@gsi.de>
1268
1269   return NULL;
1270
1271   TH2S *h2 = (TH2S*)fContainer->At(kChi2);
1272   TF1 f("fChi2", "[0]*pow(x, [1]-1)*exp(-0.5*x)", 0., 50.);
1273   f.SetParLimits(1,1, 1e100);
1274   TLegend *leg = new TLegend(.7,.7,.95,.95);
1275   leg->SetBorderSize(1); leg->SetHeader("Tracklets per Track");
1276   TH1D *h1 = NULL;
1277   Bool_t kFIRST = kTRUE;
1278   for(Int_t il=1; il<=h2->GetNbinsX(); il++){
1279     h1 = h2->ProjectionY(Form("pyChi2%d", il), il, il);
1280     if(h1->Integral()<50) continue;
1281     h1->Scale(1./h1->Integral());
1282     h1->SetMarkerStyle(7);h1->SetMarkerColor(il);
1283     h1->SetLineColor(il);h1->SetLineStyle(2);
1284     f.SetParameter(1, .5*il);f.SetLineColor(il);
1285     h1->Fit(&f, "QW+", kFIRST ? "pc": "pcsame");
1286     leg->AddEntry(h1, Form("%d", il), "l");
1287     if(kFIRST){
1288       h1->GetXaxis()->SetRangeUser(0., 25.);
1289     }
1290     kFIRST = kFALSE;
1291   }
1292   leg->Draw();
1293   gPad->SetLogy();
1294   return h1;
1295 }
1296
1297
1298 //________________________________________________________
1299 TH1* AliTRDcheckDET::MakePlotNTracklets(){
1300   //
1301   // Make nice bar plot of the number of tracklets in each method
1302   //
1303   TH2F *tmp = (TH2F *)fContainer->FindObject("hNtlsBAR");
1304   TH1D *hBAR = tmp->ProjectionY();
1305   TH1F *hSTA = (TH1F *)fContainer->FindObject("hNtlsSTA");
1306   TH1F *hCON = (TH1F *)fContainer->FindObject("hNtls");
1307   TLegend *leg = new TLegend(0.13, 0.75, 0.39, 0.89);
1308   leg->SetBorderSize(1);
1309   leg->SetFillColor(0);
1310
1311   Float_t scale = hCON->Integral();
1312   hCON->Scale(100./scale);
1313   hCON->SetFillColor(kRed);hCON->SetLineColor(kRed);
1314   hCON->SetBarWidth(0.2);
1315   hCON->SetBarOffset(0.6);
1316   hCON->SetStats(kFALSE);
1317   hCON->GetYaxis()->SetRangeUser(0.,40.);
1318   hCON->GetYaxis()->SetTitleOffset(1.2);
1319   hCON->Draw("bar1"); leg->AddEntry(hCON, "Total", "f");
1320   hCON->SetMaximum(55.);
1321
1322   hBAR->Scale(100./scale);
1323   hBAR->SetFillColor(kGreen);hBAR->SetLineColor(kGreen);
1324   hBAR->SetBarWidth(0.2);
1325   hBAR->SetBarOffset(0.2);
1326   hBAR->SetTitle("");
1327   hBAR->SetStats(kFALSE);
1328   hBAR->GetYaxis()->SetRangeUser(0.,40.);
1329   hBAR->GetYaxis()->SetTitleOffset(1.2);
1330   hBAR->Draw("bar1same"); leg->AddEntry(hBAR, "Barrel", "f");
1331
1332   hSTA->Scale(100./scale);
1333   hSTA->SetFillColor(kBlue);hSTA->SetLineColor(kBlue);
1334   hSTA->SetBarWidth(0.2);
1335   hSTA->SetBarOffset(0.4);
1336   hSTA->SetTitle("");
1337   hSTA->SetStats(kFALSE);
1338   hSTA->GetYaxis()->SetRangeUser(0.,40.);
1339   hSTA->GetYaxis()->SetTitleOffset(1.2);
1340   hSTA->Draw("bar1same"); leg->AddEntry(hSTA, "Stand Alone", "f");
1341   leg->Draw();
1342   gPad->Update();
1343   return hCON;
1344 }
1345
1346 //________________________________________________________
1347 void AliTRDcheckDET::MakePlotnTrackletsVsP(){
1348   //
1349   // Plot abundance of tracks with number of tracklets as function of momentum
1350   //
1351
1352
1353
1354
1355   Color_t color[AliTRDgeometry::kNlayer] = {kBlue, kOrange, kBlack, kGreen, kCyan, kRed};
1356   TH1 *h(NULL); TGraphErrors *g[AliTRDgeometry::kNlayer];
1357   for(Int_t itl(0); itl<AliTRDgeometry::kNlayer; itl++){
1358     g[itl] = new TGraphErrors();
1359     g[itl]->SetLineColor(color[itl]);
1360     g[itl]->SetMarkerColor(color[itl]);
1361     g[itl]->SetMarkerStyle(20 + itl);
1362   }
1363
1364   TH2 *hBar = (TH2F *)fContainer->FindObject("hNtlsBAR");
1365   TAxis *ax(hBar->GetXaxis());
1366   Int_t np(ax->GetNbins());
1367   for(Int_t ipBin(1); ipBin<np; ipBin++){
1368     h = hBar->ProjectionY("npBin", ipBin, ipBin);
1369     if(!Int_t(h->Integral())) continue;
1370     h->Scale(100./h->Integral());
1371     Float_t p(ax->GetBinCenter(ipBin)); 
1372     Float_t dp(ax->GetBinWidth(ipBin)); 
1373     Int_t ip(g[0]->GetN());
1374     for(Int_t itl(AliTRDgeometry::kNlayer); itl--;){
1375       g[itl]->SetPoint(ip, p, h->GetBinContent(itl+1));
1376       g[itl]->SetPointError(ip, dp/3.46, h->GetBinError(itl+1));
1377     }
1378   }
1379
1380   TLegend *leg = new TLegend(0.76, 0.6, 1., 0.9);
1381   leg->SetBorderSize(0);
1382   leg->SetHeader("Tracklet/Track");
1383   leg->SetFillStyle(0);
1384   h = hBar->ProjectionX("npxBin"); h->Reset();
1385   h->SetTitle("");
1386   h->GetYaxis()->SetRangeUser(1., 99.);
1387   h->GetYaxis()->SetMoreLogLabels();
1388   h->GetYaxis()->CenterTitle();
1389   h->GetYaxis()->SetTitleOffset(1.2);
1390   h->SetYTitle("Prob. [%]");
1391   h->GetXaxis()->SetRangeUser(0.4, 12.);
1392   h->GetXaxis()->SetMoreLogLabels();
1393   h->GetXaxis()->CenterTitle();
1394   h->Draw("p");
1395   for(Int_t itl(AliTRDgeometry::kNlayer); itl--;){
1396     g[itl]->Draw("pc");
1397     leg->AddEntry(g[itl], Form("n = %d", itl+1),"pl");
1398   }
1399
1400   leg->Draw();
1401   gPad->SetLogx();gPad->SetLogy();
1402 }
1403
1404 //________________________________________________________
1405 Bool_t AliTRDcheckDET::MakePlotPulseHeight(){
1406   //
1407   // Create Plot of the Pluse Height Spectrum
1408   //
1409   TCanvas *output = gPad->GetCanvas();
1410   output->Divide(2);
1411   output->cd(1);
1412   TH1 *h, *h1, *h2;
1413   TObjArray *arr = (TObjArray*)fContainer->FindObject("<PH>");
1414   h = (TH1F*)arr->At(0);
1415   h->SetMarkerStyle(24);
1416   h->SetMarkerColor(kBlack);
1417   h->SetLineColor(kBlack);
1418   h->GetYaxis()->SetTitleOffset(1.5);
1419   h->Draw("e1");
1420   // Trending for the pulse height: plateau value, slope and timebin of the maximum
1421   TLinearFitter fit(1,"pol1");
1422   Double_t time = 0.;
1423   for(Int_t itime = 10; itime <= 20; itime++){
1424     time = static_cast<Double_t>(itime);
1425     fit.AddPoint(&time, h->GetBinContent(itime + 1), h->GetBinError(itime + 1));
1426   }
1427   fit.Eval();
1428   Double_t plateau = fit.GetParameter(0) + 12 * fit.GetParameter(1);
1429   Double_t slope = fit.GetParameter(1);
1430   PutTrendValue("PHplateau", plateau);
1431   PutTrendValue("PHslope", slope);
1432   PutTrendValue("PHamplificationPeak", static_cast<Double_t>(h->GetMaximumBin()-1));
1433   AliDebug(1, Form("plateau %f, slope %f, MaxTime %f", plateau, slope, static_cast<Double_t>(h->GetMaximumBin()-1)));
1434 //   copy the second histogram in a new one with the same x-dimension as the phs with respect to time
1435   h1 = (TH1F *)arr->At(1);
1436   h2 = new TH1F("hphs1","Average PH", 31, -0.5, 30.5);
1437   for(Int_t ibin = h1->GetXaxis()->GetFirst(); ibin < h1->GetNbinsX(); ibin++) 
1438     h2->SetBinContent(ibin, h1->GetBinContent(ibin));
1439   h2->SetMarkerStyle(22);
1440   h2->SetMarkerColor(kBlue);
1441   h2->SetLineColor(kBlue);
1442   h2->Draw("e1same");
1443   gPad->Update();
1444 //   create axis according to the histogram dimensions of the original second histogram
1445   TGaxis *axis = new TGaxis(gPad->GetUxmin(),
1446                     gPad->GetUymax(),
1447                     gPad->GetUxmax(),
1448                     gPad->GetUymax(),
1449                     -0.08, 4.88, 510,"-L");
1450   axis->SetLineColor(kBlue);
1451   axis->SetLabelColor(kBlue);
1452   axis->SetTextColor(kBlue);
1453   axis->SetTitle("x_{0}-x_{c} [cm]");
1454   axis->Draw();
1455
1456   output->cd(2);
1457   TH2 *ph2d = (TH2F *)arr->At(2);
1458   ph2d->GetYaxis()->SetTitleOffset(1.8);
1459   ph2d->SetStats(kFALSE);
1460   ph2d->Draw("colz");
1461   return kTRUE;
1462 }
1463
1464 //________________________________________________________
1465 void AliTRDcheckDET::MakePlotMeanClustersLayer(){
1466   //
1467   // Create Summary plot for the mean number of clusters per layer
1468   //
1469   TCanvas *output = gPad->GetCanvas();
1470   output->Divide(3,2);
1471   TObjArray *histos = (TObjArray *)fContainer->At(kNclustersLayer);
1472   if(!histos){
1473     AliWarning("Histos for each layer not found");
1474     return;
1475   }
1476   TProfile2D *hlayer = NULL;
1477   for(Int_t ily = 0; ily < AliTRDgeometry::kNlayer; ily++){
1478     hlayer = dynamic_cast<TProfile2D *>(histos->At(ily));
1479     output->cd(ily + 1);
1480     gPad->SetGrid(0,0);
1481     hlayer->Draw("colz");
1482   }
1483 }
1484
1485 //________________________________________________________
1486 Bool_t AliTRDcheckDET::MakeBarPlot(TH1 *histo, Int_t color){
1487   //
1488   // Draw nice bar plots
1489   //
1490   if(!histo->GetEntries()) return kFALSE;
1491   histo->Scale(100./histo->Integral());
1492   histo->SetFillColor(color);
1493   histo->SetBarOffset(.2);
1494   histo->SetBarWidth(.6);
1495   histo->Draw("bar1");
1496   return kTRUE;
1497 }