]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/qaRec/AliTRDcheckDetector.cxx
be43a2d7bc6312fc0aeca3b286ec9de35d892b18
[u/mrichter/AliRoot.git] / TRD / qaRec / AliTRDcheckDetector.cxx
1 #include <TH1F.h>
2 #include <TMath.h>
3 #include <TObjArray.h>
4 #include <TProfile.h>
5
6 #include "AliTRDcluster.h"
7 #include "AliTRDseedV1.h"
8 #include "AliTRDtrackV1.h"
9 #include "TTreeStream.h"
10
11 #include "AliTRDtrackInfo/AliTRDtrackInfo.h"
12 #include "AliTRDcheckDetector.h"
13
14 ////////////////////////////////////////////////////////////////////////////
15 //                                                                        //
16 //  Reconstruction QA                                                     //
17 //                                                                        //
18 //  Task doing basic checks for tracking and detector performance         //
19 //                                                                        //
20 //  Authors:                                                              //
21 //    Anton Andronic <A.Andronic@gsi.de>                                  //
22 //    Markus Fasel <M.Fasel@gsi.de>                                       //
23 //                                                                        //
24 ////////////////////////////////////////////////////////////////////////////
25
26 //_______________________________________________________
27 AliTRDcheckDetector::AliTRDcheckDetector():
28   AliTRDrecoTask("DetChecker", "Basic Detector Checker")
29   ,fPHSdetector(0x0)
30   ,fPHSsector(0x0)
31   ,fQCLdetector(0x0)
32   ,fQCLsector(0x0)
33   ,fQTdetector(0x0)
34   ,fQTsector(0x0)
35 {
36   //
37   // Default constructor
38   //
39 }
40
41 //_______________________________________________________
42 AliTRDcheckDetector::~AliTRDcheckDetector(){
43   //
44   // Destructor
45   // 
46   if(fPHSdetector) delete fPHSdetector;
47   if(fPHSsector) delete fPHSsector;
48   if(fQCLdetector) delete fQCLdetector;
49   if(fQCLsector) delete fQCLsector;
50   if(fQTdetector) delete fQTdetector;
51   if(fQTsector) delete fQTsector;
52 }
53
54 //_______________________________________________________
55 void AliTRDcheckDetector::CreateOutputObjects(){
56   //
57   // Create Output Objects
58   //
59   fContainer = new TObjArray();
60   TH1F *hNtrks = new TH1F("hNtrks", "Number of Tracks per event", 100, 0, 100);
61   TH1F *hNcls = new TH1F("hNcls", "Nr. of clusters per track", 30, 0, 180);
62   TH1F *hNtls = new TH1F("hNtls", "Nr. tracklets per track", 10, 0, 10);
63   TH1F *hNclTls = new TH1F("hNclTls","Mean Number of clusters per tracklet", 30, 0, 30);
64   TH1F *hChi2 = new TH1F("hChi2", "Chi2", 200, 0, 20);
65   TH1F *hChi2n = new TH1F("hChi2n", "Norm. Chi2 (tracklets)", 50, 0, 5);
66   TH1F *hSM = new TH1F("hSM", "Track Counts in Supermodule", 18, 0., 18.);
67   TH1F *hQcl = new TH1F("hQcl1", "Cluster charge", 200, -1200, 1200);
68   TProfile *hPH = new TProfile("hPH", "Average PH", 31, -0.5, 30.5);
69   TH1F * hQT = new TH1F("hQT", "Total Charge Deposit", 6000, 0, 6000);
70   // Register Histograms
71   fContainer->Add(hNtrks);
72   fContainer->Add(hNcls);
73   fContainer->Add(hNtls);
74   fContainer->Add(hNclTls);
75   fContainer->Add(hChi2);
76   fContainer->Add(hChi2n);
77   fContainer->Add(hSM);
78   fContainer->Add(hPH);
79   fContainer->Add(hQcl);
80   fContainer->Add(hQT);
81
82   // Initialise the PHS, cluster charge and total charge deposit histograms for each detector
83   fPHSdetector = new TObjArray();
84   fPHSdetector->SetName("fPHSdetector");
85   fQCLdetector = new TObjArray();
86   fQCLdetector->SetName("fQCLdetector");
87   fQTdetector = new TObjArray();
88   fQTdetector->SetName("fQTdetector");
89   for(Int_t idet = 0; idet < kNDetectors; idet++){
90     fPHSdetector->Add(new TProfile(Form("hPHSdet%d", idet), Form("Average Pulse Height in Detector %d", idet), 31, -0.5, 30.5));
91     fQCLdetector->Add(new TH1F(Form("hQd%d",idet),Form("Qd%d",idet), 100,0,5000));
92     fQTdetector->Add(new TH1F(Form("hQTd%d", idet), Form("Qtotal%d", idet), 6000, 0, 6000));
93   }
94   fContainer->Add(fPHSdetector);
95   fContainer->Add(fQCLdetector);
96   fContainer->Add(fQTdetector);
97   
98   // Initialise the PHS, cluster charge and total charge deposit histograms for each sector
99   fPHSsector = new TObjArray();
100   fPHSsector->SetName("fPHSsector");
101   fQCLsector = new TObjArray();
102   fQCLsector->SetName("fQCLsector");
103   fQTsector = new TObjArray();
104   fQTsector->SetName("fQTsector");
105   for(Int_t isec = 0; isec < kNSectors; isec++){
106     fPHSsector->Add(new TProfile(Form("hPHSsec%d", isec), Form("Average Pulse Height in Sector %d", isec), 31, -0.5, 30.5));
107     fQCLsector->Add(new TH1F(Form("hQd%d",isec),Form("Qd%d",isec), 100,0,5000));
108     fQTsector->Add(new TH1F(Form("hQTd%d", isec), Form("Qtotal%d", isec), 6000, 0, 6000));
109   }
110   fContainer->Add(fPHSsector);
111   fContainer->Add(fQCLsector);
112   fContainer->Add(fQTsector);
113 }
114
115 //_______________________________________________________
116 void AliTRDcheckDetector::Exec(Option_t *){
117   //
118   // Execution function
119   // Filling TRD quality histos
120   //
121   Int_t nTracks = 0;            // Count the number of tracks per event
122   AliTRDtrackInfo *fTrackInfo = 0x0;
123   AliTRDtrackV1 *fTRDtrack = 0x0;
124   AliTRDseedV1 *fTracklet = 0x0;
125   AliTRDcluster *fTRDcluster = 0x0;
126   for(Int_t iti = 0; iti < fTracks->GetEntriesFast(); iti++){
127     fTrackInfo = dynamic_cast<AliTRDtrackInfo *>(fTracks->UncheckedAt(iti));
128     if(!fTrackInfo || !(fTRDtrack = fTrackInfo->GetTRDtrack())) continue;
129     Int_t nclusters = fTRDtrack->GetNumberOfClusters();
130     Int_t ntracklets = fTRDtrack->GetNumberOfTracklets();
131     Float_t chi2 = fTRDtrack->GetChi2();
132     // Fill Histograms
133     dynamic_cast<TH1F *>(fContainer->UncheckedAt(kNclustersHist))->Fill(nclusters);
134     dynamic_cast<TH1F *>(fContainer->UncheckedAt(kNtrackletsHist))->Fill(ntracklets);
135     dynamic_cast<TH1F *>(fContainer->UncheckedAt(kChi2))->Fill(chi2);
136     dynamic_cast<TH1F *>(fContainer->UncheckedAt(kChi2Normalized))->Fill(chi2/static_cast<Float_t>(ntracklets));
137     // now loop over single tracklets
138     if(ntracklets > 2){
139       Int_t sector = -1;
140       for(Int_t ilayer = 0; ilayer < kNLayers; ilayer++){
141         if(!(fTracklet = fTRDtrack->GetTracklet(ilayer)) || !fTracklet->IsOK()) continue;
142         Float_t nClustersTracklet = fTracklet->GetN();
143         dynamic_cast<TH1F *>(fContainer->UncheckedAt(kNclusterTrackletHist))->Fill(nClustersTracklet);
144         if(nClustersTracklet){
145           Float_t Qtot = 0;
146           Int_t detector = -1;
147           for(Int_t itb = 0; itb < kNTimebins; itb++){
148             if(!(fTRDcluster = fTracklet->GetClusters(itb))) continue;
149             Int_t localtime        = fTRDcluster->GetLocalTimeBin();
150             Double_t clusterCharge = fTRDcluster->GetQ();
151             detector               = fTRDcluster->GetDetector();
152             sector                 = static_cast<Int_t>(detector/kNDetectorsSector);
153             Double_t absolute_charge = TMath::Abs(clusterCharge);
154             Qtot += absolute_charge;
155             // Fill Histograms
156             // 1st. PHS
157             dynamic_cast<TProfile *>(fContainer->UncheckedAt(kPulseHeight))->Fill(localtime, absolute_charge);
158             dynamic_cast<TProfile *>(fPHSdetector->UncheckedAt(detector))->Fill(localtime, absolute_charge);
159             dynamic_cast<TProfile *>(fPHSsector->UncheckedAt(sector))->Fill(localtime, absolute_charge);
160             // 2nd. cluster charge
161             dynamic_cast<TH1F *>(fContainer->UncheckedAt(kClusterCharge))->Fill(absolute_charge);
162             dynamic_cast<TH1F *>(fQCLdetector->UncheckedAt(detector))->Fill(absolute_charge);
163             dynamic_cast<TH1F *>(fQCLsector->UncheckedAt(sector))->Fill(absolute_charge);
164           }
165           // Fill total cluster charge histograms
166           dynamic_cast<TH1F *>(fContainer->UncheckedAt(kChargeDeposit))->Fill(Qtot);
167           dynamic_cast<TH1F *>(fQTdetector->UncheckedAt(detector))->Fill(Qtot);
168           dynamic_cast<TH1F *>(fQTsector->UncheckedAt(sector))->Fill(Qtot);
169         }
170       }
171       dynamic_cast<TH1F *>(fContainer->UncheckedAt(kNTracksSectorHist))->Fill(sector);
172     }
173     nTracks++;
174   }
175   dynamic_cast<TH1F *>(fContainer->UncheckedAt(kNTracksEventHist))->Fill(nTracks);
176   PostData(0, fContainer);
177 }
178
179 //_______________________________________________________
180 void AliTRDcheckDetector::Terminate(Option_t *){
181   //
182   // Terminate function
183   //
184 }
185