]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDQADataMakerSim.cxx
adding Mihaela's macros to be used for the analysis train...
[u/mrichter/AliRoot.git] / TRD / AliTRDQADataMakerSim.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 ////////////////////////////////////////////////////////////////////////////
19 //                                                                        //
20 //  Produces the data needed to calculate the quality assurance.          //
21 //  All data must be mergeable objects.                                   //
22 //                                                                        //
23 //  Author:                                                               //
24 //    Sylwester Radomski (radomski@physi.uni-heidelberg.de)               //
25 //                                                                        //
26 ////////////////////////////////////////////////////////////////////////////
27
28 // --- ROOT system ---
29 #include <TClonesArray.h>
30 #include <TFile.h> 
31 #include <TH1D.h> 
32 #include <TH2D.h>
33 #include <TH3D.h>
34 #include <TProfile.h>
35 #include <TF1.h>
36 #include <TCanvas.h>
37 #include <TTree.h>
38
39 // --- AliRoot header files ---
40 //#include "AliESDEvent.h"
41 #include "AliLog.h"
42 #include "AliTRDdigit.h"
43 #include "AliTRDhit.h"
44 //#include "AliTRDcluster.h"
45 #include "AliTRDQADataMakerSim.h"
46 #include "AliTRDdigitsManager.h"
47 #include "AliTRDgeometry.h"
48 #include "AliTRDarrayADC.h"
49 #include "AliTRDarraySignal.h"
50 //#include "AliTRDrawStream.h"
51
52 #include "AliQAChecker.h"
53
54 ClassImp(AliTRDQADataMakerSim)
55
56 //____________________________________________________________________________ 
57   AliTRDQADataMakerSim::AliTRDQADataMakerSim() : 
58   AliQADataMakerSim(AliQA::GetDetName(AliQA::kTRD), "TRD Quality Assurance Data Maker")
59 {
60   //
61   // Default constructor
62 }
63
64 //____________________________________________________________________________ 
65 AliTRDQADataMakerSim::AliTRDQADataMakerSim(const AliTRDQADataMakerSim& qadm) :
66   AliQADataMakerSim()
67 {
68   //
69   // Copy constructor 
70   //
71
72   SetName((const char*)qadm.GetName()) ; 
73   SetTitle((const char*)qadm.GetTitle()); 
74
75 }
76
77 //__________________________________________________________________
78 AliTRDQADataMakerSim& AliTRDQADataMakerSim::operator=(const AliTRDQADataMakerSim& qadm)
79 {
80   //
81   // Equal operator.
82   //
83
84   this->~AliTRDQADataMakerSim();
85   new(this) AliTRDQADataMakerSim(qadm);
86   return *this;
87
88 }
89
90 //____________________________________________________________________________ 
91 void AliTRDQADataMakerSim::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray ** list)
92 {
93   //
94   // Detector specific actions at end of cycle
95   //
96
97   //AliInfo(Form("EndOfCycle", "Fitting RecPoints %d", task));
98
99   // call the checker
100   AliQAChecker::Instance()->Run(AliQA::kTRD, task, list) ;    
101
102
103 }
104
105 //____________________________________________________________________________ 
106 void AliTRDQADataMakerSim::InitHits()
107 {
108   //
109   // Create Hits histograms in Hits subdir
110   //
111
112   const Int_t kNhist = 4;
113   TH1D *hist[kNhist];
114
115   hist[0] = new TH1D("qaTRD_hits_det", ";Detector Id of the hit", 540, -0.5, 539.5) ; 
116
117   hist[1] = new TH1D("qaTRD_hist_Qdrift", ";Charge from tracks", 100, 0, 100);
118   hist[2] = new TH1D("qaTRD_hist_Qamp", ";Charge from TRD photon", 100, 0, 100);
119   hist[3] = new TH1D("qaTRD_hist_Qphoton", ";Charge from TRD photon", 100, 0, 100);
120
121   for(Int_t i=0; i<kNhist; i++) {
122     //hist[i]->Sumw2();
123     Add2HitsList(hist[i], i);
124   }
125
126 }
127
128 //____________________________________________________________________________ 
129 void AliTRDQADataMakerSim::InitDigits()
130 {
131   //
132   // Create Digits histograms in Digits subdir
133   //
134
135   const Int_t kNhist = 3;
136   TH1D *hist[kNhist];
137
138   hist[0] = new TH1D("qaTRD_digits_det", ";Detector Id of the digit", 540, -0.5, 539.5);
139   hist[1] = new TH1D("qaTRD_digits_time", ";Time bin", 40, -0.5, 39.5);
140   hist[2] = new TH1D("qaTRD_digits_amp", ";Amplitude", 100, -5.5, 94.5);
141
142   for(Int_t i=0; i<kNhist; i++) {
143     hist[i]->Sumw2();
144     Add2DigitsList(hist[i], i);
145   }
146
147 }
148
149 //____________________________________________________________________________ 
150 void AliTRDQADataMakerSim::InitSDigits()
151 {
152   //
153   // Create SDigits histograms in SDigits subdir
154   //
155
156   const Int_t kNhist = 3;
157   TH1D *hist[kNhist];
158
159   hist[0] = new TH1D("qaTRD_sdigits_det", ";Detector Id of the digit", 540, -0.5, 539.5);
160   hist[1] = new TH1D("qaTRD_sdigits_time", ";Time bin", 40, -0.5, 39.5);
161   hist[2] = new TH1D("qaTRD_sdigits_amp", ";Amplitude", 100, 0, 1e7);
162
163   for(Int_t i=0; i<kNhist; i++) {
164     hist[i]->Sumw2();
165     Add2SDigitsList(hist[i], i);
166   }
167
168 }
169
170 //____________________________________________________________________________
171 void AliTRDQADataMakerSim::MakeHits(TClonesArray * hits)
172 {
173   //
174   // Make QA data from Hits
175   //
176
177   TIter next(hits); 
178   AliTRDhit * hit; 
179
180   while ( (hit = dynamic_cast<AliTRDhit *>(next())) ) {
181     GetHitsData(0)->Fill(hit->GetDetector());
182     Double_t q = TMath::Abs(hit->GetCharge());
183
184     if (hit->FromDrift()) GetHitsData(1)->Fill(q);
185     if (hit->FromAmplification()) GetHitsData(2)->Fill(q);
186     if (hit->FromTRphoton()) GetHitsData(3)->Fill(q);
187   }
188
189 }
190
191 //____________________________________________________________________________
192 void AliTRDQADataMakerSim::MakeHits(TTree * hitTree)
193 {
194   //
195   // Make QA data from Hits
196   //
197
198   if (!CheckPointer(hitTree, "TRD hits tree")) return;
199
200   TBranch *branch = hitTree->GetBranch("TRD");
201   if (!CheckPointer(branch, "TRD hits branch")) return;
202
203   Int_t nhits = (Int_t)(hitTree->GetTotBytes()/sizeof(AliTRDhit));
204   TClonesArray *hits = new TClonesArray("AliTRDhit", nhits+1000);
205   TClonesArray *tmp = new TClonesArray("AliTRDhit", 1000);
206   branch->SetAddress(&tmp);
207
208   Int_t index = 0;
209   Int_t nEntries = (Int_t)branch->GetEntries();
210   for(Int_t i = 0; i < nEntries; i++) {
211     branch->GetEntry(i);
212     Int_t nHits = (Int_t)tmp->GetEntries();
213     for(Int_t j=0; j<nHits; j++) {
214       AliTRDhit *hit = (AliTRDhit*)tmp->At(j);
215       new((*hits)[index++]) AliTRDhit(*hit);
216     }
217   }
218
219   tmp->Delete();
220   delete tmp;
221   MakeHits(hits);
222   hits->Delete();
223   delete hits;
224
225 }
226
227 //____________________________________________________________________________
228 void AliTRDQADataMakerSim::MakeDigits(TClonesArray * digits)
229 {
230   //
231   // Makes data from Digits
232   //
233
234   TIter next(digits) ; 
235   AliTRDdigit * digit ; 
236   
237   // Info("Make digits", "From the arrya");
238   
239   while ( (digit = dynamic_cast<AliTRDdigit *>(next())) ) {
240     if (digit->GetAmp() < 1) continue;
241     GetDigitsData(0)->Fill(digit->GetDetector());
242     GetDigitsData(1)->Fill(digit->GetTime());
243     GetDigitsData(2)->Fill(digit->GetAmp());
244   }
245
246 }
247
248 //____________________________________________________________________________
249 void AliTRDQADataMakerSim::MakeDigits(TTree * digits)
250 {
251   //
252   // Makes data from digits tree
253   //
254
255   // Info("Make digits", "From a tree");
256
257   AliTRDdigitsManager *digitsManager = new AliTRDdigitsManager();
258   digitsManager->CreateArrays();
259   digitsManager->ReadDigits(digits);
260
261   TH1D *histDet = (TH1D*)GetDigitsData(0);
262   TH1D *histTime = (TH1D*)GetDigitsData(1);
263   TH1D *histSignal = (TH1D*)GetDigitsData(2);
264
265   for (Int_t i = 0; i < AliTRDgeometry::kNdet; i++) 
266     {
267       AliTRDarrayADC *digitsIn = (AliTRDarrayADC *) digitsManager->GetDigits(i);      
268
269       // This is to take care of switched off super modules
270       if (digitsIn->GetNtime() == 0) continue;
271       
272       digitsIn->Expand();
273       
274       //AliTRDSignalIndex* indexes = digitsManager->GetIndexes(i);
275       //if (indexes->IsAllocated() == kFALSE) digitsManager->BuildIndexes(i);
276       
277       Int_t nRows = digitsIn->GetNrow();
278       Int_t nCols = digitsIn->GetNcol();
279       Int_t nTbins = digitsIn->GetNtime();
280       
281       for(Int_t row = 0; row < nRows; row++) 
282         for(Int_t col = 0; col < nCols; col++) 
283           for(Int_t time = 0; time < nTbins; time++) 
284             {   
285               Float_t signal = digitsIn->GetData(row,col,time);
286               if (signal < 1) continue;
287               histDet->Fill(i);
288               histTime->Fill(time);
289               histSignal->Fill(signal);
290             }
291
292     //delete digitsIn;
293   }
294   delete digitsManager;
295 }
296
297 //____________________________________________________________________________
298 void AliTRDQADataMakerSim::MakeSDigits(TClonesArray * sdigits)
299 {
300   //
301   // Makes data from Digits
302   //
303
304   TIter next(sdigits) ; 
305   AliTRDdigit * digit ; 
306   while ( (digit = dynamic_cast<AliTRDdigit *>(next())) ) {
307     GetDigitsData(0)->Fill(digit->GetDetector());
308     GetDigitsData(1)->Fill(digit->GetTime());
309     GetDigitsData(2)->Fill(digit->GetAmp());
310   }
311
312 }
313
314 //____________________________________________________________________________
315 void AliTRDQADataMakerSim::MakeSDigits(TTree * digits)
316 {
317   //
318   // Makes data from SDigits
319   //
320
321   AliTRDdigitsManager *digitsManager = new AliTRDdigitsManager();
322   digitsManager->SetSDigits(kTRUE);
323   digitsManager->CreateArrays();
324   digitsManager->ReadDigits(digits);
325
326   TH1D *histDet = (TH1D*)GetSDigitsData(0);
327   TH1D *histTime = (TH1D*)GetSDigitsData(1);
328   TH1D *histSignal = (TH1D*)GetSDigitsData(2);
329
330   for (Int_t i = 0; i < AliTRDgeometry::kNdet; i++) 
331     {
332       AliTRDarraySignal *digitsIn = (AliTRDarraySignal *) digitsManager->GetSDigits(i);
333   
334     // This is to take care of switched off super modules
335       if (digitsIn->GetNtime() == 0) continue;
336       
337     digitsIn->Expand(); 
338       
339       //AliTRDSignalIndex* indexes = digitsManager->GetIndexes(i);
340       //if (indexes->IsAllocated() == kFALSE) digitsManager->BuildIndexes(i);
341       Int_t nRows = digitsIn->GetNrow();
342       Int_t nCols = digitsIn->GetNcol();
343       Int_t nTbins = digitsIn->GetNtime();
344
345       for(Int_t row = 0; row < nRows; row++) 
346         {
347           for(Int_t col = 0; col < nCols; col++) 
348             {
349               for(Int_t time = 0; time < nTbins; time++) 
350                 {           
351                   Float_t signal = digitsIn->GetData(row,col,time);
352                   if (signal < 1) continue;
353                   histDet->Fill(i);
354                   histTime->Fill(time);
355                   histSignal->Fill(signal);
356                 }
357             }
358         }
359       // delete digitsIn;
360     }
361   delete digitsManager;
362 }
363
364 //____________________________________________________________________________ 
365 void AliTRDQADataMakerSim::StartOfDetectorCycle()
366 {
367   //
368   // Detector specific actions at start of cycle
369   //
370
371 }
372
373 //__________________________________________________________________________
374 Int_t AliTRDQADataMakerSim::CheckPointer(TObject *obj, const char *name) 
375 {
376   //
377   // Checks initialization of pointers
378   //
379
380   if (!obj) AliWarning(Form("null pointer: %s", name));
381   return !!obj;
382
383 }