]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG1/TRD/AliTRDcheckDET.cxx
fix coverity
[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 "AliTRDinfoGen.h"
72 #include "AliTRDcheckDET.h"
73
74 #include <cstdio>
75 #include <iostream>
76
77 ClassImp(AliTRDcheckDET)
78
79 //_______________________________________________________
80 AliTRDcheckDET::AliTRDcheckDET():
81   AliTRDrecoTask()
82   ,fEventInfo(NULL)
83   ,fTriggerNames(NULL)
84   ,fFlags(0)
85 {
86   //
87   // Default constructor
88   //
89   SetNameTitle("TRDcheckDET", "Basic TRD data checker");
90 }
91
92 //_______________________________________________________
93 AliTRDcheckDET::AliTRDcheckDET(char* name):
94   AliTRDrecoTask(name, "Basic TRD data checker")
95   ,fEventInfo(NULL)
96   ,fTriggerNames(NULL)
97   ,fFlags(0)
98 {
99   //
100   // Default constructor
101   //
102   DefineInput(2, AliTRDeventInfo::Class());
103
104   InitFunctorList();
105 }
106
107
108 //_______________________________________________________
109 AliTRDcheckDET::~AliTRDcheckDET(){
110   //
111   // Destructor
112   // 
113   if(fTriggerNames) delete fTriggerNames;
114 }
115
116 //_______________________________________________________
117 void AliTRDcheckDET::UserCreateOutputObjects(){
118   //
119   // Create Output Objects
120   //
121   AliTRDrecoTask::UserCreateOutputObjects();
122   if(!fTriggerNames) fTriggerNames = new TMap();
123 }
124
125 //_______________________________________________________
126 void AliTRDcheckDET::UserExec(Option_t *opt){
127   //
128   // Execution function
129   // Filling TRD quality histos
130   //
131   fEventInfo = dynamic_cast<AliTRDeventInfo *>(GetInputData(2));
132   AliDebug(2, Form("EventInfo[%p] Header[%p]", (void*)fEventInfo, (void*)(fEventInfo?fEventInfo->GetEventHeader():NULL)));
133
134   AliTRDrecoTask::UserExec(opt);  
135
136   TH1F *histo(NULL); AliTRDtrackInfo *fTrackInfo(NULL); Int_t nTracks(0);               // Count the number of tracks per event
137   for(Int_t iti = 0; iti < fTracks->GetEntriesFast(); iti++){
138     if(!fTracks->UncheckedAt(iti)) continue;
139     if(!(fTrackInfo = dynamic_cast<AliTRDtrackInfo *>(fTracks->UncheckedAt(iti)))) continue;
140     if(!fTrackInfo->GetTrack()) continue;
141     nTracks++;
142   }
143   if(nTracks)
144     dynamic_cast<TH1F *>(fContainer->UncheckedAt(kNtracksEvent))->Fill(nTracks);
145
146   if(!fEventInfo->GetEventHeader()) return; // For trigger statistics event header is essential
147   Int_t triggermask = fEventInfo->GetEventHeader()->GetTriggerMask();
148   TString triggername =  fEventInfo->GetRunInfo()->GetFiredTriggerClasses(triggermask);
149   AliDebug(6, Form("Trigger cluster: %d, Trigger class: %s\n", triggermask, triggername.Data()));
150   if((histo = dynamic_cast<TH1F *>(fContainer->UncheckedAt(kNeventsTrigger)))) histo->Fill(triggermask);
151
152   if(nTracks){
153     if((histo = dynamic_cast<TH1F *>(fContainer->UncheckedAt(kNeventsTriggerTracks)))) histo->Fill(triggermask);
154   }
155   if(triggermask <= 20 && !fTriggerNames->FindObject(Form("%d", triggermask))){
156     fTriggerNames->Add(new TObjString(Form("%d", triggermask)), new TObjString(triggername));
157     // also set the label for both histograms
158     if((histo = dynamic_cast<TH1F *>(fContainer->UncheckedAt(kNeventsTriggerTracks))))
159       histo->GetXaxis()->SetBinLabel(histo->FindBin(triggermask), triggername);
160     if((histo = dynamic_cast<TH1F *>(fContainer->UncheckedAt(kNeventsTrigger))))
161       histo->GetXaxis()->SetBinLabel(histo->FindBin(triggermask), triggername);
162   }
163 }
164
165
166 //_______________________________________________________
167 Bool_t AliTRDcheckDET::PostProcess(){
168   //
169   // Do Postprocessing (for the moment set the number of Reference histograms)
170   //
171   
172   TH1F *h(NULL), *h1(NULL);
173
174   // Calculate of the trigger clusters purity
175   if((h  = dynamic_cast<TH1F *>(fContainer->FindObject("hEventsTrigger"))) &&
176      (h1 = dynamic_cast<TH1F *>(fContainer->FindObject("hEventsTriggerTracks")))) {
177     h1->Divide(h);
178     Float_t purities[20], val = 0; memset(purities, 0, 20*sizeof(Float_t));
179     TString triggernames[20];
180     Int_t nTriggerClasses = 0;
181     for(Int_t ibin = 1; ibin <= h->GetNbinsX(); ibin++){
182       if((val = h1->GetBinContent(ibin))){
183         purities[nTriggerClasses] = val;
184         triggernames[nTriggerClasses] = h1->GetXaxis()->GetBinLabel(ibin);
185         nTriggerClasses++;
186       }
187     }
188
189     if((h = dynamic_cast<TH1F *>(fContainer->UncheckedAt(kTriggerPurity)))){
190       TAxis *ax = h->GetXaxis();
191       for(Int_t itrg = 0; itrg < nTriggerClasses; itrg++){
192         h->Fill(itrg, purities[itrg]);
193         ax->SetBinLabel(itrg+1, triggernames[itrg].Data());
194       }
195       ax->SetRangeUser(-0.5, nTriggerClasses+.5);
196       h->GetYaxis()->SetRangeUser(0,1);
197     }
198   }
199
200   // track status
201   if((h=dynamic_cast<TH1F*>(fContainer->At(kTrackStatus)))){
202     Float_t ok = h->GetBinContent(1);
203     Int_t nerr = h->GetNbinsX();
204     for(Int_t ierr=nerr; ierr--;){
205       h->SetBinContent(ierr+1, ok>0.?1.e2*h->GetBinContent(ierr+1)/ok:0.);
206     }
207     h->SetBinContent(1, 0.);
208   }
209   // tracklet status
210   
211   TObjArray *arr(NULL);
212   if(( arr = dynamic_cast<TObjArray*>(fContainer->UncheckedAt(kTrackletStatus)))){
213     for(Int_t ily = AliTRDgeometry::kNlayer; ily--;){
214       if(!(h=dynamic_cast<TH1F*>(arr->At(ily)))) continue;
215       Float_t okB = h->GetBinContent(1);
216       Int_t nerrB = h->GetNbinsX();
217       for(Int_t ierr=nerrB; ierr--;){
218         h->SetBinContent(ierr+1, okB>0.?1.e2*h->GetBinContent(ierr+1)/okB:0.);
219       }
220       h->SetBinContent(1, 0.);
221     }
222   }
223   fNRefFigures = 17;
224
225   return kTRUE;
226 }
227
228 //_______________________________________________________
229 void AliTRDcheckDET::MakeSummary(){
230   //
231   // Create summary plots for TRD check DET
232   // This function creates 2 summary plots:
233   // - General Quantities
234   // - PHS
235   // The function will reuse GetRefFigure
236   //
237   
238   TCanvas *cOut = new TCanvas(Form("summary%s1", GetName()), Form("Summary 1 for task %s", GetName()), 1024, 768);
239   cOut->Divide(3,3);
240
241   // Create figures using GetRefFigure
242   cOut->cd(1); GetRefFigure(kFigNtracksEvent);  
243   cOut->cd(2); GetRefFigure(kFigNtracksSector);
244   cOut->cd(3); GetRefFigure(kFigNclustersTrack);
245   cOut->cd(4); GetRefFigure(kFigNclustersTracklet);
246   cOut->cd(5); GetRefFigure(kFigNtrackletsTrack);
247   cOut->cd(6); GetRefFigure(kFigNTrackletsP);
248   cOut->cd(7); GetRefFigure(kFigChargeCluster);
249   cOut->cd(8); GetRefFigure(kFigChargeTracklet);
250   cOut->SaveAs("TRD_TrackQuality.gif");
251   delete cOut;
252
253   // Second Plot: PHS
254   cOut = new TCanvas(Form("summary%s2", GetName()), Form("Summary 2 for task %s", GetName()), 1024, 512);
255   cOut->cd(); GetRefFigure(kFigPH);
256   cOut->SaveAs("TRD_PH.gif"); 
257   delete cOut;
258
259   // Third Plot: Mean Number of clusters as function of eta, phi and layer
260    cOut = new TCanvas(Form("summary%s3", GetName()), Form("Summary 3 for task %s", GetName()), 1024, 768);
261   cOut->cd(); MakePlotMeanClustersLayer();
262   cOut->SaveAs("TRD_MeanNclusters.gif"); 
263   delete cOut;
264
265 }
266
267 //_______________________________________________________
268 Bool_t AliTRDcheckDET::GetRefFigure(Int_t ifig){
269   //
270   // Setting Reference Figures
271   //
272   gPad->SetLogy(0);
273   gPad->SetLogx(0);
274   TH1 *h = NULL; TObjArray *arr=NULL;
275   TLegend *leg = NULL;
276   Bool_t kFIRST(1);
277   switch(ifig){
278   case kFigNclustersTrack:
279     (h=(TH1F*)fContainer->FindObject("hNcls"))->Draw("pl");
280     PutTrendValue("NClustersTrack", h->GetMean());
281     PutTrendValue("NClustersTrackRMS", h->GetRMS());
282     return kTRUE;
283   case kFigNclustersTracklet:
284     (h =(TH1F*)fContainer->FindObject("hNclTls"))->Draw("pc");
285     PutTrendValue("NClustersTracklet", h->GetMean());
286     PutTrendValue("NClustersTrackletRMS", h->GetRMS());
287     return kTRUE;
288   case kFigNtrackletsTrack:
289     h=MakePlotNTracklets();
290     if(h){
291       PutTrendValue("NTrackletsTrack", h->GetMean());
292       PutTrendValue("NTrackletsTrackRMS", h->GetRMS());
293     }
294     return kTRUE;
295   case kFigNTrackletsP:
296     MakePlotnTrackletsVsP();
297     return kTRUE;
298   case kFigNtrackletsCross:
299     h = (TH1F*)fContainer->FindObject("hNtlsCross");
300     if(!MakeBarPlot(h, kRed)) break;
301     PutTrendValue("NTrackletsCross", h->GetMean());
302     PutTrendValue("NTrackletsCrossRMS", h->GetRMS());
303     return kTRUE;
304   case kFigNtrackletsFindable:
305     h = (TH1F*)fContainer->FindObject("hNtlsFindable");
306     if(!MakeBarPlot(h, kGreen)) break;
307     PutTrendValue("NTrackletsFindable", h->GetMean());
308     PutTrendValue("NTrackletsFindableRMS", h->GetRMS());
309     return kTRUE;
310   case kFigNtracksEvent:
311     (h = (TH1F*)fContainer->FindObject("hNtrks"))->Draw("pl");
312     PutTrendValue("NTracksEvent", h->GetMean());
313     PutTrendValue("NTracksEventRMS", h->GetRMS());
314     return kTRUE;
315   case kFigNtracksSector:
316     h = (TH1F*)fContainer->FindObject("hNtrksSector");
317     if(!MakeBarPlot(h, kGreen)) break;
318     PutTrendValue("NTracksSector", h->Integral()/h->GetNbinsX());
319     return kTRUE;
320   case kFigTrackStatus:
321     if(!(h=(TH1F *)fContainer->FindObject("hTrackStatus"))) break;
322     h->GetXaxis()->SetRangeUser(0.5, -1);
323     h->GetYaxis()->CenterTitle();
324     h->Draw("c");
325     PutTrendValue("TrackStatus", h->Integral());
326     gPad->SetLogy(0);
327     return kTRUE;
328   case kFigTrackletStatus:
329     if(!(arr = dynamic_cast<TObjArray*>(fContainer->At(kTrackletStatus)))) break;
330     leg = new TLegend(.68, .7, .97, .97);
331     leg->SetBorderSize(0);leg->SetFillStyle(0);
332     leg->SetHeader("TRD layer");
333     for(Int_t ily=AliTRDgeometry::kNlayer; ily--;){
334       if(!(h=dynamic_cast<TH1F*>(arr->At(ily)))) continue;
335       if(kFIRST){
336         h->Draw("pl");
337         h->GetXaxis()->SetRangeUser(0.5, -1);
338         h->GetYaxis()->CenterTitle();
339         kFIRST = kFALSE;
340       } else h->Draw("samepl");
341       leg->AddEntry(h, Form("ly = %d", ily), "l");
342       PutTrendValue(Form("TrackletStatus%d", ily), h->Integral());
343     }
344     leg->Draw();
345     gPad->SetLogy(0);
346     return kTRUE;
347   case kFigChi2:
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   
554   fContainer->AddAt(h, kNeventsTrigger);
555
556   if(!(h = (TH1F *)gROOT->FindObject("hEventsTriggerTracks")))
557     h = new TH1F("hEventsTriggerTracks", "Trigger Class (Tracks)", 100, 0, 100);
558   else h->Reset();
559   fContainer->AddAt(h, kNeventsTriggerTracks);
560
561   if(!(h = (TH1F *)gROOT->FindObject("hTriggerPurity"))){
562     h = new TH1F("hTriggerPurity", "Trigger Purity", 10, -0.5, 9.5);
563     h->GetXaxis()->SetTitle("Trigger Cluster");
564     h->GetYaxis()->SetTitle("freq.");
565   } else h->Reset();
566   fContainer->AddAt(h, kTriggerPurity);
567
568   return fContainer;
569 }
570
571 /*
572 * Plotting Functions
573 */
574
575 //_______________________________________________________
576 TH1 *AliTRDcheckDET::PlotTrackStatus(const AliTRDtrackV1 *track)
577 {
578 //
579 // Plot the track propagation status. The following errors are defined (see AliTRDtrackV1::ETRDtrackError)
580 //   PROL - track prolongation failure
581 //   PROP - track propagation failure
582 //   AJST - crossing sectors failure
583 //   SNP  - too large bending
584 //   TINI - tracklet initialization failure
585 //   UPDT - track position/covariance update failure 
586 //
587 // Performance plot looks as below:
588 //Begin_Html
589 //<img src="TRD/trackStatus.gif">
590 //End_Html
591 //
592   if(track) fkTrack = track;
593   if(!fkTrack){
594     AliDebug(4, "No Track defined.");
595     return NULL;
596   }
597   TH1 *h = NULL;
598   if(!(h = dynamic_cast<TH1F *>(fContainer->At(kTrackStatus)))){
599     AliWarning("No Histogram defined.");
600     return NULL;
601   }
602   h->Fill(fkTrack->GetStatusTRD());
603   return h;
604 }
605
606 //_______________________________________________________
607 TH1 *AliTRDcheckDET::PlotTrackletStatus(const AliTRDtrackV1 *track)
608 {
609 //
610 // Plot the tracklet propagation status. The following errors are defined for tracklet (see AliTRDtrackV1::ETRDlayerError)
611 //   Geom   - 
612 //   Bound  - tracklet too close to chamber walls
613 //   NoCl   - no clusters in the track roads
614 //   NoAttach - fail to attach clusters
615 //   NoClTr - fail to use clusters for fit
616 //   NoFit  - tracklet fit failled
617 //   Chi2   - chi2 tracklet-track over threshold
618 //
619 // Performance plot looks as below:
620 //Begin_Html
621 //<img src="TRD/trackletStatus.gif">
622 //End_Html
623 //
624   if(track) fkTrack = track;
625   if(!fkTrack){
626     AliDebug(4, "No Track defined.");
627     return NULL;
628   }
629   TObjArray *arr =NULL;
630   if(!(arr = dynamic_cast<TObjArray*>(fContainer->At(kTrackletStatus)))){
631     AliWarning("Histograms not defined.");
632     return NULL;
633   }
634
635   TH1 *h = NULL;
636   for(Int_t ily=AliTRDgeometry::kNlayer; ily--;){
637     if(!(h = dynamic_cast<TH1F*>(arr->At(ily)))){
638       AliWarning(Form("Missing histo for layer %d.", ily));
639       continue;
640     }
641     h->Fill(fkTrack->GetStatusTRD(ily));
642   }
643   return h;
644 }
645
646 //_______________________________________________________
647 TH1 *AliTRDcheckDET::PlotNClustersTracklet(const AliTRDtrackV1 *track){
648   //
649   // Plot the mean number of clusters per tracklet
650   //
651   if(track) fkTrack = track;
652   if(!fkTrack){
653     AliDebug(4, "No Track defined.");
654     return NULL;
655   }
656   AliExternalTrackParam *par = fkTrack->GetTrackIn() ? fkTrack->GetTrackIn() : fkTrack->GetTrackOut();
657   TH1 *h = NULL;
658   TProfile2D *hlayer = NULL;
659   Double_t eta = 0., phi = 0.;
660   if(!(h = dynamic_cast<TH1F *>(fContainer->At(kNclustersTracklet)))){
661     AliWarning("No Histogram defined.");
662     return NULL;
663   }
664   AliTRDseedV1 *tracklet = NULL;
665   TObjArray *histosLayer = dynamic_cast<TObjArray *>(fContainer->At(kNclustersLayer));
666   if(!histosLayer){
667     AliWarning("No Histograms for single layer defined");
668   }
669   for(Int_t itl = 0; itl < AliTRDgeometry::kNlayer; itl++){
670     if(!(tracklet = fkTrack->GetTracklet(itl)) || !tracklet->IsOK()) continue;
671     h->Fill(tracklet->GetN2());
672     if(histosLayer && par){
673       if((hlayer = dynamic_cast<TProfile2D *>(histosLayer->At(itl)))){
674         GetEtaPhiAt(par, tracklet->GetX0(), eta, phi);
675         hlayer->Fill(eta, phi, tracklet->GetN2());
676       }
677     }
678   }
679   return h;
680 }
681
682 //_______________________________________________________
683 TH1 *AliTRDcheckDET::PlotNClustersTrack(const AliTRDtrackV1 *track){
684   //
685   // Plot the number of clusters in one track
686   //
687   if(track) fkTrack = track;
688   if(!fkTrack){
689     AliDebug(4, "No Track defined.");
690     return NULL;
691   }
692   TH1 *h = NULL;
693   if(!(h = dynamic_cast<TH1F *>(fContainer->At(kNclustersTrack)))){
694     AliWarning("No Histogram defined.");
695     return NULL;
696   }
697   
698   Int_t nclusters = 0;
699   AliTRDseedV1 *tracklet = NULL;
700   AliExternalTrackParam *par = fkTrack->GetTrackOut() ? fkTrack->GetTrackOut() : fkTrack->GetTrackIn();
701   if(!par) return NULL;
702   Double_t momentumRec = par->P();
703   for(Int_t itl = 0; itl < AliTRDgeometry::kNlayer; itl++){
704     if(!(tracklet = fkTrack->GetTracklet(itl)) || !tracklet->IsOK()) continue;
705     Int_t n(tracklet->GetN());
706     nclusters += n;
707     if(DebugLevel() > 2){
708       Int_t crossing = Int_t(tracklet->IsRowCross());
709       Int_t detector = tracklet->GetDetector();
710       Float_t theta = TMath::ATan(tracklet->GetZref(1));
711       Float_t phi = TMath::ATan(tracklet->GetYref(1));
712       Float_t momentumMC = 0.;
713       Int_t pdg = 0;
714       Int_t kinkIndex = fkESD ? fkESD->GetKinkIndex() : 0;
715       UShort_t nclsTPC = fkESD ? fkESD->GetTPCncls() : 0;
716       if(fkMC){
717         if(fkMC->GetTrackRef()) momentumMC = fkMC->GetTrackRef()->P();
718         pdg = fkMC->GetPDG();
719       }
720       (*DebugStream()) << "NClustersTrack"
721         << "Detector="  << detector
722         << "crossing="  << crossing
723         << "momentumMC="<< momentumMC
724         << "momentumRec=" << momentumRec
725         << "pdg="                               << pdg
726         << "theta="                     << theta
727         << "phi="                               << phi
728         << "kinkIndex=" << kinkIndex
729         << "TPCncls="           << nclsTPC
730         << "TRDncls="   << n
731         << "\n";
732     }
733   }
734   h->Fill(nclusters);
735   return h;
736 }
737
738
739 //_______________________________________________________
740 TH1 *AliTRDcheckDET::PlotNTrackletsTrack(const AliTRDtrackV1 *track){
741   //
742   // Plot the number of tracklets
743   //
744   if(track) fkTrack = track;
745   if(!fkTrack){
746     AliDebug(4, "No Track defined.");
747     return NULL;
748   }
749   TH1 *h = NULL, *hSta = NULL; TH2 *hBarrel = NULL;
750   if(!(h = dynamic_cast<TH1F *>(fContainer->At(kNtrackletsTrack)))){
751     AliWarning("No Histogram defined.");
752     return NULL;
753   }
754   Int_t nTracklets = fkTrack->GetNumberOfTracklets();
755   h->Fill(nTracklets);
756   if(!fkESD) return h;
757   Int_t status = fkESD->GetStatus();
758
759 /*  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);*/
760   Double_t p = 0.;
761   Int_t method = -1;    // to distinguish between stand alone and full barrel tracks in the debugging
762   if((status & AliESDtrack::kTRDin) != 0){
763     method = 1;
764     // Full Barrel Track: Save momentum dependence
765     if(!(hBarrel = dynamic_cast<TH2F *>(fContainer->At(kNtrackletsBAR)))){
766       AliWarning("Method: Barrel.  Histogram not processed!");
767       return NULL;
768     }
769     AliExternalTrackParam *par(fkTrack->GetTrackIn());
770     if(!par){
771       AliError("Input track params missing");
772       return NULL;
773     }
774     p = par->P(); // p needed later in the debug streaming
775     hBarrel->Fill(p, nTracklets);
776   } else {
777     // Stand alone Track: momentum dependence not usefull
778     method = 0;
779     if(!(hSta = dynamic_cast<TH1F *>(fContainer->At(kNtrackletsSTA)))) {
780       AliWarning("Method: StandAlone.  Histogram not processed!");
781       return NULL;
782     }
783     hSta->Fill(nTracklets);
784   }
785
786   if(DebugLevel() > 2){
787     AliTRDseedV1 *tracklet = NULL;
788     AliTRDgeometry *geo(AliTRDinfoGen::Geometry());
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 = geo->GetSector(detector);
794       stack = geo->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(AliTRDinfoGen::Reconstructor());
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   AliTRDgeometry *geo(AliTRDinfoGen::Geometry());
940   for(Int_t il = 0; il < AliTRDgeometry::kNlayer; il++){
941     y = points[il].GetY();
942     z = points[il].GetZ();
943     if((stack = geo->GetStack(z, il)) < 0) continue; // Not findable
944     pp = geo->GetPadPlane(il, stack);
945     ymin = pp->GetCol0() + epsilon;
946     ymax = pp->GetColEnd() - epsilon; 
947     zmin = pp->GetRowEnd() + epsilon; 
948     zmax = pp->GetRow0() - epsilon;
949     // ignore y-crossing (material)
950     if((z + deltaZ > zmin && z - deltaZ < zmax) && (y + deltaY > ymin && y - deltaY < ymax)) nFindable++;
951       if(DebugLevel() > 3){
952         Double_t posTracklet[2] = {tracklet ? tracklet->GetYfit(0) : 0, tracklet ? tracklet->GetZfit(0) : 0};
953         Int_t hasTracklet = tracklet ? 1 : 0;
954         (*DebugStream())   << "FindableTracklets"
955           << "layer="     << il
956           << "ytracklet=" << posTracklet[0]
957           << "ytrack="    << y
958           << "ztracklet=" << posTracklet[1]
959           << "ztrack="    << z
960           << "tracklet="  << hasTracklet
961           << "\n";
962       }
963   }
964   
965   h->Fill(nFindable > 0 ? TMath::Min(nFound/static_cast<Double_t>(nFindable), 1.) : 1);
966   AliDebug(2, Form("Findable[Found]: %d[%d|%f]", nFindable, nFound, nFound/static_cast<Float_t>(nFindable > 0 ? nFindable : 1)));
967   return h;
968 }
969
970
971 //_______________________________________________________
972 TH1 *AliTRDcheckDET::PlotChi2(const AliTRDtrackV1 *track){
973   //
974   // Plot the chi2 of the track
975   //
976   if(track) fkTrack = track;
977   if(!fkTrack){
978     AliDebug(4, "No Track defined.");
979     return NULL;
980   }
981   TH1 *h = NULL;
982   if(!(h = dynamic_cast<TH2S*>(fContainer->At(kChi2)))) {
983     AliWarning("No Histogram defined.");
984     return NULL;
985   }
986   Int_t n = fkTrack->GetNumberOfTracklets();
987   if(!n) return NULL;
988
989   h->Fill(n, fkTrack->GetChi2()/n);
990   return h;
991 }
992
993
994 //_______________________________________________________
995 TH1 *AliTRDcheckDET::PlotPHt(const AliTRDtrackV1 *track){
996   //
997   // Plot the average pulse height
998   //
999   if(track) fkTrack = track;
1000   if(!fkTrack){
1001     AliDebug(4, "No Track defined.");
1002     return NULL;
1003   }
1004   TProfile *h = NULL; TH2F *phs2D = NULL;
1005   if(!(h = dynamic_cast<TProfile *>(((TObjArray*)(fContainer->At(kPH)))->At(0)))){
1006     AliWarning("No Histogram defined.");
1007     return NULL;
1008   }
1009   if(!(phs2D = dynamic_cast<TH2F *>(((TObjArray*)(fContainer->At(kPH)))->At(2)))){
1010     AliWarning("2D Pulse Height histogram not defined. Histogramm cannot be filled");
1011   }
1012   AliTRDseedV1 *tracklet = NULL;
1013   AliTRDcluster *c = NULL;
1014   for(Int_t itl = 0; itl < AliTRDgeometry::kNlayer; itl++){
1015     if(!(tracklet = fkTrack->GetTracklet(itl)) || !tracklet->IsOK())continue;
1016     Int_t crossing = Int_t(tracklet->IsRowCross());
1017     Int_t detector = tracklet->GetDetector();
1018     tracklet->ResetClusterIter();
1019     while((c = tracklet->NextCluster())){
1020       if(!IsUsingClustersOutsideChamber() && !c->IsInChamber()) continue;
1021       Int_t localtime        = c->GetLocalTimeBin();
1022       Double_t absoluteCharge = TMath::Abs(c->GetQ());
1023       h->Fill(localtime, absoluteCharge);
1024       phs2D->Fill(localtime, absoluteCharge); 
1025       if(DebugLevel() > 3){
1026         Int_t inChamber = c->IsInChamber() ? 1 : 0;
1027         Double_t distance[2];
1028         GetDistanceToTracklet(distance, tracklet, c);
1029         Float_t theta = TMath::ATan(tracklet->GetZref(1));
1030         Float_t phi = TMath::ATan(tracklet->GetYref(1));
1031         AliExternalTrackParam *trdPar = fkTrack->GetTrackIn();
1032         Float_t momentumMC = 0, momentumRec = trdPar ? trdPar->P() : fkTrack->P(); // prefer Track Low
1033         Int_t pdg = 0;
1034         Int_t kinkIndex = fkESD ? fkESD->GetKinkIndex() : 0;
1035         UShort_t tpcNCLS = fkESD ? fkESD->GetTPCncls() : 0;
1036         if(fkMC){
1037           if(fkMC->GetTrackRef()) momentumMC = fkMC->GetTrackRef()->P();
1038           pdg = fkMC->GetPDG();
1039         }
1040         (*DebugStream()) << "PHt"
1041           << "Detector="        << detector
1042           << "crossing="        << crossing
1043           << "inChamber=" << inChamber
1044           << "Timebin="         << localtime
1045           << "Charge="          << absoluteCharge
1046           << "momentumMC="      << momentumMC
1047           << "momentumRec="     << momentumRec
1048           << "pdg="                             << pdg
1049           << "theta="                   << theta
1050           << "phi="                             << phi
1051           << "kinkIndex="       << kinkIndex
1052           << "TPCncls="         << tpcNCLS
1053           << "dy="        << distance[0]
1054           << "dz="        << distance[1]
1055           << "c.="        << c
1056           << "\n";
1057       }
1058     }
1059   }
1060   return h;
1061 }
1062
1063 //_______________________________________________________
1064 TH1 *AliTRDcheckDET::PlotPHx(const AliTRDtrackV1 *track){
1065   //
1066   // Plots the average pulse height vs the distance from the anode wire
1067   // (plus const anode wire offset)
1068   //
1069   if(track) fkTrack = track;
1070   if(!fkTrack){
1071      AliDebug(4, "No Track defined.");
1072      return NULL;
1073   }
1074   TProfile *h = NULL;
1075   if(!(h = dynamic_cast<TProfile *>(((TObjArray*)(fContainer->At(kPH)))->At(1)))){
1076     AliWarning("No Histogram defined.");
1077     return NULL;
1078   }
1079   AliTRDseedV1 *tracklet(NULL);
1080   AliTRDcluster *c(NULL);
1081   Double_t xd(0.), dqdl(0.);
1082   TVectorD vq(AliTRDseedV1::kNtb), vxd(AliTRDseedV1::kNtb), vdqdl(AliTRDseedV1::kNtb);
1083   for(Int_t itl = 0; itl < AliTRDgeometry::kNlayer; itl++){
1084     if(!(tracklet = fkTrack->GetTracklet(itl)) || !(tracklet->IsOK())) continue;
1085     Int_t det(tracklet->GetDetector());
1086     Bool_t rc(tracklet->IsRowCross());
1087     for(Int_t ic(0); ic<AliTRDseedV1::kNtb; ic++){
1088       Bool_t kFIRST(kFALSE);
1089       if(!(c = tracklet->GetClusters(ic))){
1090          if(!(c = tracklet->GetClusters(AliTRDseedV1::kNtb+ic))) continue;
1091       } else kFIRST=kTRUE;
1092       if(!IsUsingClustersOutsideChamber() && !c->IsInChamber()) continue;
1093       xd = tracklet->GetX0() - c->GetX(); vxd[ic] = xd;
1094       dqdl=tracklet->GetdQdl(ic); vdqdl[ic] = dqdl;
1095       vq[ic]=c->GetQ();
1096       if(kFIRST && (c = tracklet->GetClusters(AliTRDseedV1::kNtb+ic))) vq[ic]+=c->GetQ();
1097       h->Fill(xd, dqdl);
1098     }
1099     if(DebugLevel() > 3){
1100       (*DebugStream()) << "PHx"
1101         << "det="  << det
1102         << "rc="   << rc
1103         << "xd="   << &vxd
1104         << "q="    << &vq
1105         << "dqdl=" << &vdqdl
1106         << "\n";
1107     }
1108   }  
1109   return h;
1110 }
1111
1112 //_______________________________________________________
1113 TH1 *AliTRDcheckDET::PlotChargeCluster(const AliTRDtrackV1 *track){
1114   //
1115   // Plot the cluster charge
1116   //
1117   if(track) fkTrack = track;
1118   if(!fkTrack){
1119     AliDebug(4, "No Track defined.");
1120     return NULL;
1121   }
1122   TH1 *h = NULL;
1123   if(!(h = dynamic_cast<TH1F *>(fContainer->At(kChargeCluster)))){
1124     AliWarning("No Histogram defined.");
1125     return NULL;
1126   }
1127   AliTRDseedV1 *tracklet = NULL;
1128   AliTRDcluster *c = NULL;
1129   for(Int_t itl = 0; itl < AliTRDgeometry::kNlayer; itl++){
1130     if(!(tracklet = fkTrack->GetTracklet(itl)) || !tracklet->IsOK())continue;
1131     for(Int_t ic(0); ic < AliTRDseedV1::kNtb; ic++){
1132       Bool_t kFIRST(kFALSE);
1133       if(!(c = tracklet->GetClusters(ic))) {
1134         if(!(c = tracklet->GetClusters(AliTRDseedV1::kNtb+ic))) continue;
1135       } else kFIRST = kTRUE;
1136       Float_t q(c->GetQ());
1137       if(kFIRST && (c = tracklet->GetClusters(AliTRDseedV1::kNtb+ic))) q+=c->GetQ();
1138       h->Fill(q);
1139     }
1140   }
1141   return h;
1142 }
1143
1144 //_______________________________________________________
1145 TH1 *AliTRDcheckDET::PlotChargeTracklet(const AliTRDtrackV1 *track){
1146   //
1147   // Plot the charge deposit per chamber
1148   //
1149   if(track) fkTrack = track;
1150   if(!fkTrack){
1151     AliDebug(4, "No Track defined.");
1152     return NULL;
1153   }
1154   TH1 *h = NULL;
1155   if(!(h = dynamic_cast<TH1F *>(fContainer->At(kChargeTracklet)))){
1156     AliWarning("No Histogram defined.");
1157     return NULL;
1158   }
1159   AliTRDseedV1 *tracklet = NULL;
1160   AliTRDcluster *c = NULL;
1161   Double_t qTot = 0;
1162   Int_t nTracklets =fkTrack->GetNumberOfTracklets();
1163   for(Int_t itl(0); itl < AliTRDgeometry::kNlayer; itl++){
1164     if(!(tracklet = fkTrack->GetTracklet(itl)) || !tracklet->IsOK()) continue;
1165     qTot = 0.;
1166     for(Int_t ic = AliTRDseedV1::kNclusters; ic--;){
1167       if(!(c = tracklet->GetClusters(ic))) continue;
1168       qTot += TMath::Abs(c->GetQ());
1169     }
1170     h->Fill(qTot);
1171     if(DebugLevel() > 3){
1172       Int_t crossing = (Int_t)tracklet->IsRowCross();
1173       Int_t detector = tracklet->GetDetector();
1174       Float_t theta = TMath::ATan(tracklet->GetZfit(1));
1175       Float_t phi = TMath::ATan(tracklet->GetYfit(1));
1176       Float_t momentum = 0.;
1177       Int_t pdg = 0;
1178       Int_t kinkIndex = fkESD ? fkESD->GetKinkIndex() : 0;
1179       UShort_t nclsTPC = fkESD ? fkESD->GetTPCncls() : 0;
1180       if(fkMC){
1181               if(fkMC->GetTrackRef()) momentum = fkMC->GetTrackRef()->P();
1182         pdg = fkMC->GetPDG();
1183       }
1184       (*DebugStream()) << "ChargeTracklet"
1185         << "Detector="  << detector
1186         << "crossing="  << crossing
1187         << "momentum="  << momentum
1188         << "nTracklets="<< nTracklets
1189         << "pdg="                               << pdg
1190         << "theta="                     << theta
1191         << "phi="                               << phi
1192         << "kinkIndex=" << kinkIndex
1193         << "TPCncls="           << nclsTPC
1194         << "QT="        << qTot
1195         << "\n";
1196     }
1197   }
1198   return h;
1199 }
1200
1201 //_______________________________________________________
1202 TH1 *AliTRDcheckDET::PlotNTracksSector(const AliTRDtrackV1 *track){
1203   //
1204   // Plot the number of tracks per Sector
1205   //
1206   if(track) fkTrack = track;
1207   if(!fkTrack){
1208     AliDebug(4, "No Track defined.");
1209     return NULL;
1210   }
1211   TH1 *h = NULL;
1212   if(!(h = dynamic_cast<TH1F *>(fContainer->At(kNtracksSector)))){
1213     AliWarning("No Histogram defined.");
1214     return NULL;
1215   }
1216
1217   // TODO we should compare with
1218   // sector = Int_t(track->GetAlpha() / AliTRDgeometry::GetAlpha());
1219
1220   AliTRDseedV1 *tracklet = NULL;
1221   Int_t sector = -1;
1222   for(Int_t itl = 0; itl < AliTRDgeometry::kNlayer; itl++){
1223     if(!(tracklet = fkTrack->GetTracklet(itl)) || !tracklet->IsOK()) continue;
1224     sector = static_cast<Int_t>(tracklet->GetDetector()/AliTRDgeometry::kNdets);
1225     break;
1226   }
1227   h->Fill(sector);
1228   return h;
1229 }
1230
1231
1232 //________________________________________________________
1233 void AliTRDcheckDET::GetDistanceToTracklet(Double_t *dist, AliTRDseedV1 * const tracklet, AliTRDcluster * const c)
1234 {
1235   Float_t x = c->GetX();
1236   dist[0] = c->GetY() - tracklet->GetYat(x);
1237   dist[1] = c->GetZ() - tracklet->GetZat(x);
1238 }
1239
1240 //________________________________________________________
1241 void AliTRDcheckDET::GetEtaPhiAt(const AliExternalTrackParam *track, Double_t x, Double_t &eta, Double_t &phi){
1242   //
1243   // Get phi and eta at a given radial position
1244   // 
1245   AliExternalTrackParam workpar(*track);
1246
1247   Double_t posLocal[3];
1248   Bool_t sucPos = workpar.GetXYZAt(x, fEventInfo->GetRunInfo()->GetMagneticField(), posLocal);
1249   Double_t sagPhi = sucPos ? TMath::ATan2(posLocal[1], posLocal[0]) : 0.;
1250   phi = sagPhi;
1251   eta = workpar.Eta();
1252 }
1253
1254
1255 //_______________________________________________________
1256 TH1* AliTRDcheckDET::MakePlotChi2()
1257 {
1258 // Plot chi2/track normalized to number of degree of freedom 
1259 // (tracklets) and compare with the theoretical distribution.
1260 // 
1261 // Alex Bercuci <A.Bercuci@gsi.de>
1262
1263   return NULL;
1264
1265 /*  TH2S *h2 = (TH2S*)fContainer->At(kChi2);
1266   TF1 f("fChi2", "[0]*pow(x, [1]-1)*exp(-0.5*x)", 0., 50.);
1267   f.SetParLimits(1,1, 1e100);
1268   TLegend *leg = new TLegend(.7,.7,.95,.95);
1269   leg->SetBorderSize(1); leg->SetHeader("Tracklets per Track");
1270   TH1D *h1 = NULL;
1271   Bool_t kFIRST = kTRUE;
1272   for(Int_t il=1; il<=h2->GetNbinsX(); il++){
1273     h1 = h2->ProjectionY(Form("pyChi2%d", il), il, il);
1274     if(h1->Integral()<50) continue;
1275     h1->Scale(1./h1->Integral());
1276     h1->SetMarkerStyle(7);h1->SetMarkerColor(il);
1277     h1->SetLineColor(il);h1->SetLineStyle(2);
1278     f.SetParameter(1, .5*il);f.SetLineColor(il);
1279     h1->Fit(&f, "QW+", kFIRST ? "pc": "pcsame");
1280     leg->AddEntry(h1, Form("%d", il), "l");
1281     if(kFIRST){
1282       h1->GetXaxis()->SetRangeUser(0., 25.);
1283     }
1284     kFIRST = kFALSE;
1285   }
1286   leg->Draw();
1287   gPad->SetLogy();
1288   return h1;*/
1289 }
1290
1291
1292 //________________________________________________________
1293 TH1* AliTRDcheckDET::MakePlotNTracklets(){
1294   //
1295   // Make nice bar plot of the number of tracklets in each method
1296   //
1297   TH2F *tmp = (TH2F *)fContainer->FindObject("hNtlsBAR");
1298   TH1D *hBAR = tmp->ProjectionY();
1299   TH1F *hSTA = (TH1F *)fContainer->FindObject("hNtlsSTA");
1300   TH1F *hCON = (TH1F *)fContainer->FindObject("hNtls");
1301   TLegend *leg = new TLegend(0.13, 0.75, 0.39, 0.89);
1302   leg->SetBorderSize(1);
1303   leg->SetFillColor(0);
1304
1305   Float_t scale = hCON->Integral();
1306   if(scale) hCON->Scale(100./scale);
1307   hCON->SetFillColor(kRed);hCON->SetLineColor(kRed);
1308   hCON->SetBarWidth(0.2);
1309   hCON->SetBarOffset(0.6);
1310   hCON->SetStats(kFALSE);
1311   hCON->GetYaxis()->SetRangeUser(0.,40.);
1312   hCON->GetYaxis()->SetTitleOffset(1.2);
1313   hCON->Draw("bar1"); leg->AddEntry(hCON, "Total", "f");
1314   hCON->SetMaximum(55.);
1315
1316   if(scale) hBAR->Scale(100./scale);
1317   hBAR->SetFillColor(kGreen);hBAR->SetLineColor(kGreen);
1318   hBAR->SetBarWidth(0.2);
1319   hBAR->SetBarOffset(0.2);
1320   hBAR->SetTitle("");
1321   hBAR->SetStats(kFALSE);
1322   hBAR->GetYaxis()->SetRangeUser(0.,40.);
1323   hBAR->GetYaxis()->SetTitleOffset(1.2);
1324   hBAR->Draw("bar1same"); leg->AddEntry(hBAR, "Barrel", "f");
1325
1326   if(scale) hSTA->Scale(100./scale);
1327   hSTA->SetFillColor(kBlue);hSTA->SetLineColor(kBlue);
1328   hSTA->SetBarWidth(0.2);
1329   hSTA->SetBarOffset(0.4);
1330   hSTA->SetTitle("");
1331   hSTA->SetStats(kFALSE);
1332   hSTA->GetYaxis()->SetRangeUser(0.,40.);
1333   hSTA->GetYaxis()->SetTitleOffset(1.2);
1334   hSTA->Draw("bar1same"); leg->AddEntry(hSTA, "Stand Alone", "f");
1335   leg->Draw();
1336   gPad->Update();
1337   return hCON;
1338 }
1339
1340 //________________________________________________________
1341 void AliTRDcheckDET::MakePlotnTrackletsVsP(){
1342   //
1343   // Plot abundance of tracks with number of tracklets as function of momentum
1344   //
1345
1346
1347
1348
1349   Color_t color[AliTRDgeometry::kNlayer] = {kBlue, kOrange, kBlack, kGreen, kCyan, kRed};
1350   TH1 *h(NULL); TGraphErrors *g[AliTRDgeometry::kNlayer];
1351   for(Int_t itl(0); itl<AliTRDgeometry::kNlayer; itl++){
1352     g[itl] = new TGraphErrors();
1353     g[itl]->SetLineColor(color[itl]);
1354     g[itl]->SetMarkerColor(color[itl]);
1355     g[itl]->SetMarkerStyle(20 + itl);
1356   }
1357
1358   TH2 *hBar = (TH2F *)fContainer->FindObject("hNtlsBAR");
1359   TAxis *ax(hBar->GetXaxis());
1360   Int_t np(ax->GetNbins());
1361   for(Int_t ipBin(1); ipBin<np; ipBin++){
1362     h = hBar->ProjectionY("npBin", ipBin, ipBin);
1363     if(!Int_t(h->Integral())) continue;
1364     h->Scale(100./h->Integral());
1365     Float_t p(ax->GetBinCenter(ipBin)); 
1366     Float_t dp(ax->GetBinWidth(ipBin)); 
1367     Int_t ip(g[0]->GetN());
1368     for(Int_t itl(AliTRDgeometry::kNlayer); itl--;){
1369       g[itl]->SetPoint(ip, p, h->GetBinContent(itl+1));
1370       g[itl]->SetPointError(ip, dp/3.46, h->GetBinError(itl+1));
1371     }
1372   }
1373
1374   TLegend *leg = new TLegend(0.76, 0.6, 1., 0.9);
1375   leg->SetBorderSize(0);
1376   leg->SetHeader("Tracklet/Track");
1377   leg->SetFillStyle(0);
1378   h = hBar->ProjectionX("npxBin"); h->Reset();
1379   h->SetTitle("");
1380   h->GetYaxis()->SetRangeUser(1., 99.);
1381   h->GetYaxis()->SetMoreLogLabels();
1382   h->GetYaxis()->CenterTitle();
1383   h->GetYaxis()->SetTitleOffset(1.2);
1384   h->SetYTitle("Prob. [%]");
1385   h->GetXaxis()->SetRangeUser(0.4, 12.);
1386   h->GetXaxis()->SetMoreLogLabels();
1387   h->GetXaxis()->CenterTitle();
1388   h->Draw("p");
1389   for(Int_t itl(AliTRDgeometry::kNlayer); itl--;){
1390     g[itl]->Draw("pc");
1391     leg->AddEntry(g[itl], Form("n = %d", itl+1),"pl");
1392   }
1393
1394   leg->Draw();
1395   gPad->SetLogx();gPad->SetLogy();
1396 }
1397
1398 //________________________________________________________
1399 Bool_t AliTRDcheckDET::MakePlotPulseHeight(){
1400   //
1401   // Create Plot of the Pluse Height Spectrum
1402   //
1403   TCanvas *output = gPad->GetCanvas();
1404   output->Divide(2);
1405   output->cd(1);
1406   TH1 *h, *h1, *h2;
1407   TObjArray *arr = (TObjArray*)fContainer->FindObject("<PH>");
1408   h = (TH1F*)arr->At(0);
1409   if((Int_t)h->GetEntries()){
1410     h->SetMarkerStyle(24);
1411     h->SetMarkerColor(kBlack);
1412     h->SetLineColor(kBlack);
1413     h->GetYaxis()->SetTitleOffset(1.5);
1414     h->Draw("e1");
1415     // Trending for the pulse height: plateau value, slope and timebin of the maximum
1416     TLinearFitter fit(1,"pol1");
1417     Double_t time = 0.;
1418     for(Int_t itime = 10; itime <= 20; itime++){
1419       time = Double_t(itime);
1420       Double_t err(h->GetBinError(itime + 1));
1421       if(err>1.e-10) fit.AddPoint(&time, h->GetBinContent(itime + 1), err);
1422     }
1423     if(!fit.Eval()){
1424       Double_t plateau = fit.GetParameter(0) + 12 * fit.GetParameter(1);
1425       Double_t slope = fit.GetParameter(1);
1426       PutTrendValue("PHplateau", plateau);
1427       PutTrendValue("PHslope", slope);
1428       PutTrendValue("PHamplificationPeak", static_cast<Double_t>(h->GetMaximumBin()-1));
1429       AliDebug(1, Form("plateau %f, slope %f, MaxTime %f", plateau, slope, static_cast<Double_t>(h->GetMaximumBin()-1)));
1430     }
1431   }
1432   //   copy the second histogram in a new one with the same x-dimension as the phs with respect to time
1433   h1 = (TH1F *)arr->At(1);
1434   h2 = new TH1F("hphs1","Average PH", 31, -0.5, 30.5);
1435   for(Int_t ibin = h1->GetXaxis()->GetFirst(); ibin < h1->GetNbinsX(); ibin++) 
1436     h2->SetBinContent(ibin, h1->GetBinContent(ibin));
1437   h2->SetMarkerStyle(22);
1438   h2->SetMarkerColor(kBlue);
1439   h2->SetLineColor(kBlue);
1440   h2->Draw("e1same");
1441   gPad->Update();
1442   
1443 //   create axis according to the histogram dimensions of the original second histogram
1444   TGaxis *axis = new TGaxis(gPad->GetUxmin(),
1445                     gPad->GetUymax(),
1446                     gPad->GetUxmax(),
1447                     gPad->GetUymax(),
1448                     -0.08, 4.88, 510,"-L");
1449   axis->SetLineColor(kBlue);
1450   axis->SetLabelColor(kBlue);
1451   axis->SetTextColor(kBlue);
1452   axis->SetTitle("x_{0}-x_{c} [cm]");
1453   axis->Draw();
1454
1455   output->cd(2);
1456   TH2 *ph2d = (TH2F *)arr->At(2);
1457   ph2d->GetYaxis()->SetTitleOffset(1.8);
1458   ph2d->SetStats(kFALSE);
1459   ph2d->Draw("colz");
1460   return kTRUE;
1461 }
1462
1463 //________________________________________________________
1464 void AliTRDcheckDET::MakePlotMeanClustersLayer(){
1465   //
1466   // Create Summary plot for the mean number of clusters per layer
1467   //
1468   TCanvas *output = gPad->GetCanvas();
1469   output->Divide(3,2);
1470   TObjArray *histos = (TObjArray *)fContainer->At(kNclustersLayer);
1471   if(!histos){
1472     AliWarning("Histos for each layer not found");
1473     return;
1474   }
1475   TProfile2D *hlayer = NULL;
1476   for(Int_t ily = 0; ily < AliTRDgeometry::kNlayer; ily++){
1477     if(!(hlayer = dynamic_cast<TProfile2D *>(histos->At(ily)))) continue;
1478     output->cd(ily + 1);
1479     gPad->SetGrid(0,0);
1480     hlayer->Draw("colz");
1481   }
1482 }
1483
1484 //________________________________________________________
1485 Bool_t AliTRDcheckDET::MakeBarPlot(TH1 *histo, Int_t color){
1486   //
1487   // Draw nice bar plots
1488   //
1489   if(!histo->GetEntries()) return kFALSE;
1490   histo->Scale(100./histo->Integral());
1491   histo->SetFillColor(color);
1492   histo->SetBarOffset(.2);
1493   histo->SetBarWidth(.6);
1494   histo->Draw("bar1");
1495   return kTRUE;
1496 }