]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCQADataMakerRec.cxx
In AddTaskPHOSPi0Flow.C set Cent. Bin past event buffers/lists,
[u/mrichter/AliRoot.git] / TPC / AliTPCQADataMakerRec.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2007, 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
17 /* $Id: $ */
18
19 /*
20   Based on AliPHOSQADataMaker
21   Produces the data needed to calculate the quality assurance. 
22   All data must be mergeable objects.
23   P. Christiansen, Lund, January 2008
24
25   Updated July 2011:
26   ==================
27
28   Major changes to accomodate updates of general DQM/QA changes to have per
29   trigger histograms (for a given event specie).
30
31   1) One instance of AliTPCdataQA only. (This also solves some old wishes by
32   offline team to use less memory because event the 2d arrays for this object
33   is not used). This now has a new flag for only keeping DQM info event by
34   event! For this reason there is no need for a special DQM reset any more
35   between runs!
36
37   2) Fill the histogram for each event. The histograms are no longer filled
38   from the AliTPCdataQA but per event.
39
40   3) Use profiles for the RAW info. By adding the profiles event by event we
41   get the correct event averages WITHOUT having to normalize in the end!
42   Results should therefore also be directly mergable when that feature will
43   come. (none of the other histograms are merged).
44
45   This means that from the DQM/QA point of view the TPC DQM is now fully
46   standard and should ease future developments.
47
48   Updated June 2010:
49   ==================
50
51   The "beautification" of the online DQM histograms have been moved to
52   an amore macro.  
53
54   The per event RAW histograms have been modified in AliTPCdataQA and
55   the copies have therefore also been modified here.
56
57   The AliTPCdataQA can now be configured a bit from here: time bin
58   range (extended default range to 1-1000, event range at start:
59   0-100000, 1000 events per bin). (At least the parameters are not
60   hardcoded:-)
61
62   Implementation:
63   ===============
64
65   For the QA of the RAW data we use the class, AliTPCdataQA, from the
66   existing TPC Calibration framework (which is more advanced than the
67   standard QA framework) and extract the histograms at the end. The
68   Analyse method of the AliTPCdataQA class is called in the method,
69   EndOfDetectorCycle, and there also: 1d histogram(s) are projected
70   and added to the QA list.
71 */
72
73 #include "AliTPCQADataMakerRec.h"
74
75 // --- ROOT system ---
76 #include <TClonesArray.h>
77 #include <TString.h>
78 #include <TSystem.h>
79 #include <TBox.h>
80 #include <TLine.h>
81 #include <TAxis.h>
82 #include <TH1.h> 
83 #include <TProfile.h> 
84 #include <TProfile2D.h> 
85
86 // --- Standard library ---
87
88 // --- AliRoot header files ---
89 #include "AliQAChecker.h"
90 #include "AliESDEvent.h"
91 #include "AliESDtrack.h"
92 #include "AliLog.h"
93 #include "AliTPCCalPad.h"
94 #include "AliTPCCalROC.h"
95 #include "AliTPCClustersRow.h"
96 #include "AliTPCclusterMI.h"
97 #include "AliSimDigits.h"
98 #include <AliDetectorRecoParam.h>
99
100
101 ClassImp(AliTPCQADataMakerRec)
102
103 //____________________________________________________________________________ 
104 AliTPCQADataMakerRec::AliTPCQADataMakerRec() : 
105 AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kTPC), 
106                   "TPC Rec Quality Assurance Data Maker"),
107 fTPCdataQA(NULL),
108 fRawFirstTimeBin(1),
109 fRawLastTimeBin(1000)
110 {
111   // ctor
112   
113   for(Int_t i = 0; i < 6; i++)
114     fMapping[i] = 0;
115 }
116
117 //____________________________________________________________________________ 
118 AliTPCQADataMakerRec::AliTPCQADataMakerRec(const AliTPCQADataMakerRec& qadm) :
119   AliQADataMakerRec(),
120   fTPCdataQA(NULL),
121   fRawFirstTimeBin(qadm.GetRawFirstTimeBin()),
122   fRawLastTimeBin(qadm.GetRawLastTimeBin())
123 {
124   //copy ctor 
125   // Does not copy the calibration object, instead InitRaws have to be
126   // called again
127   SetName((const char*)qadm.GetName()) ; 
128   SetTitle((const char*)qadm.GetTitle()); 
129
130   for(Int_t i = 0; i < 6; i++)
131     fMapping[i] = 0;
132 }
133
134 //__________________________________________________________________
135 AliTPCQADataMakerRec& AliTPCQADataMakerRec::operator = (const AliTPCQADataMakerRec& qadm )
136 {
137   // Equal operator.
138   this->~AliTPCQADataMakerRec();
139   new(this) AliTPCQADataMakerRec(qadm);
140   return *this;
141 }
142
143 //__________________________________________________________________
144 AliTPCQADataMakerRec::~AliTPCQADataMakerRec()
145 {
146   // Destructor
147   delete fTPCdataQA; 
148
149   for(Int_t i = 0; i < 6; i++) 
150     delete fMapping[i];
151 }
152  
153 //____________________________________________________________________________ 
154 void AliTPCQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
155 {
156   //Detector specific actions at end of cycle
157   ResetEventTrigClasses();
158
159   AliQAChecker::Instance()->Run(AliQAv1::kTPC, task, list) ;  
160 }
161
162
163 //____________________________________________________________________________ 
164 void AliTPCQADataMakerRec::InitESDs()
165 {
166   //create ESDs histograms in ESDs subdir  
167   const Bool_t expert   = kTRUE ; 
168   const Bool_t image    = kTRUE ; 
169   
170   TH1F * histESDclusters = 
171     new TH1F("hESDclusters", "N TPC clusters per track; N clusters; Counts",
172              160, 0, 160);
173   histESDclusters->Sumw2();
174   Add2ESDsList(histESDclusters, kClusters, !expert, image);
175
176   TH1F * histESDratio = 
177     new TH1F("hESDratio", "Ratio: TPC clusters / findable; Ratio: cluster/findable; Counts",
178              100, 0, 1);
179   histESDratio->Sumw2();
180   Add2ESDsList(histESDratio, kRatio, !expert, image);
181   
182   TH1F * histESDpt = 
183     new TH1F("hESDpt", "P_{T} distribution; p_{T} [GeV/c]; Counts",
184              50, 0, 5);
185   histESDpt->Sumw2();
186   Add2ESDsList(histESDpt, kPt, !expert, image);
187   //
188   ClonePerTrigClass(AliQAv1::kESDS); // this should be the last line
189 }
190
191 //____________________________________________________________________________ 
192 void AliTPCQADataMakerRec::InitRaws()
193 {
194   //
195   // Adding the raw 
196   //  
197
198   // Modified: 7/7 - 2008
199   // Laurent Aphecetche pointed out that the mapping was read from file
200   // for each event, so now we read in the map here and set if for 
201   // the raw data qa
202   const Bool_t expert   = kTRUE ; 
203   const Bool_t saveCorr = kTRUE ; 
204   const Bool_t image    = kTRUE ; 
205   
206   // It might happen that we will be in this method a few times
207   // (we create the dataQA at the first call to this method)
208   if(!fTPCdataQA) {
209     fTPCdataQA = new AliTPCdataQA();
210     LoadMaps(); // Load Altro maps
211     fTPCdataQA->SetAltroMapping(fMapping); // set Altro mapping
212     fTPCdataQA->SetRangeTime(fRawFirstTimeBin, fRawLastTimeBin); // set time bin interval 
213     fTPCdataQA->SetIsDQM(kTRUE);
214   }
215
216   TProfile * histRawsOccupancyVsSector = 
217     new TProfile("hRawsOccupancyVsSector", "Occupancy vs sector; Sector; Occupancy",
218              72, 0, 72);
219   histRawsOccupancyVsSector->SetMarkerStyle(20);
220   histRawsOccupancyVsSector->SetOption("P");
221   histRawsOccupancyVsSector->SetStats(kFALSE);
222   if ( GetRecoParam()->GetEventSpecie() == AliRecoParam::kCalib )
223     Add2RawsList(histRawsOccupancyVsSector, kRawsOccupancyVsSector, expert, image, !saveCorr);
224   else
225     Add2RawsList(histRawsOccupancyVsSector, kRawsOccupancyVsSector, !expert, image, !saveCorr);
226     
227   TProfile * histRawsQVsSector = 
228     new TProfile("hRawsQVsSector", "<Q> vs sector; Sector; <Q>",
229              72, 0, 72);
230   Add2RawsList(histRawsQVsSector, kRawsQVsSector, expert, !image, !saveCorr);
231
232   TProfile * histRawsQmaxVsSector = 
233     new TProfile("hRawsQmaxVsSector", "<Qmax> vs sector; Sector; <Qmax>",
234              72, 0, 72);
235   histRawsQmaxVsSector->SetMarkerStyle(20);
236   histRawsQmaxVsSector->SetOption("P");
237   histRawsQmaxVsSector->SetStats(kFALSE);
238   if ( GetRecoParam()->GetEventSpecie() == AliRecoParam::kCalib )
239     Add2RawsList(histRawsQmaxVsSector, kRawsQmaxVsSector, expert, image, !saveCorr);
240   else
241     Add2RawsList(histRawsQmaxVsSector, kRawsQmaxVsSector, !expert, image, !saveCorr);
242
243   TProfile2D * histRawsOccupancy2dVsSector = 
244     new TProfile2D("hRawsOccupancy2dVsSector", "Occupancy vs sector; Sector; Patch",
245                    72, 0, 36, 6, 0, 6);
246   histRawsOccupancy2dVsSector->SetOption("COLZ");
247   histRawsOccupancy2dVsSector->SetStats(kFALSE);
248   if ( GetRecoParam()->GetEventSpecie() == AliRecoParam::kCalib )
249     Add2RawsList(histRawsOccupancy2dVsSector, kRawsOccupancy2dVsSector, expert, image, !saveCorr);
250   else
251     Add2RawsList(histRawsOccupancy2dVsSector, kRawsOccupancy2dVsSector, !expert, image, !saveCorr);
252     
253   //
254   ClonePerTrigClass(AliQAv1::kRAWS); // this should be the last line
255 }
256
257 //____________________________________________________________________________ 
258 void AliTPCQADataMakerRec::InitDigits()
259 {
260   const Bool_t expert   = kTRUE ; 
261   const Bool_t image    = kTRUE ; 
262   TH1F * histDigitsADC = 
263     new TH1F("hDigitsADC", "Digit ADC distribution; ADC; Counts",
264              1000, 0, 1000);
265   histDigitsADC->Sumw2();
266   Add2DigitsList(histDigitsADC, kDigitsADC, !expert, image);
267   //
268   ClonePerTrigClass(AliQAv1::kDIGITS); // this should be the last line
269 }
270
271 //____________________________________________________________________________ 
272 void AliTPCQADataMakerRec::InitRecPoints()
273 {
274   const Bool_t expert   = kTRUE ; 
275   const Bool_t image    = kTRUE ; 
276   
277   TH1F * histRecPointsQmaxShort = 
278     new TH1F("hRecPointsQmaxShort", "Qmax distrbution (short pads); Qmax; Counts",
279              100, 0, 300);
280   histRecPointsQmaxShort->Sumw2();
281   Add2RecPointsList(histRecPointsQmaxShort, kQmaxShort, !expert, image);
282
283   TH1F * histRecPointsQmaxMedium = 
284     new TH1F("hRecPointsQmaxMedium", "Qmax distrbution (medium pads); Qmax; Counts",
285              100, 0, 300);
286   histRecPointsQmaxMedium->Sumw2();
287   Add2RecPointsList(histRecPointsQmaxMedium, kQmaxMedium, !expert, image);
288
289   TH1F * histRecPointsQmaxLong = 
290     new TH1F("hRecPointsQmaxLong", "Qmax distrbution (long pads); Qmax; Counts",
291              100, 0, 300);
292   histRecPointsQmaxLong->Sumw2();
293   Add2RecPointsList(histRecPointsQmaxLong, kQmaxLong, !expert, image);
294
295   TH1F * histRecPointsQShort = 
296     new TH1F("hRecPointsQShort", "Q distrbution (short pads); Q; Counts",
297              100, 0, 2000);
298   histRecPointsQShort->Sumw2();
299   Add2RecPointsList(histRecPointsQShort, kQShort, !expert, image);
300
301   TH1F * histRecPointsQMedium = 
302     new TH1F("hRecPointsQMedium", "Q distrbution (medium pads); Q; Counts",
303              100, 0, 2000);
304   histRecPointsQMedium->Sumw2();
305   Add2RecPointsList(histRecPointsQMedium, kQMedium, !expert, image);
306
307   TH1F * histRecPointsQLong = 
308     new TH1F("hRecPointsQLong", "Q distrbution (long pads); Q; Counts",
309              100, 0, 2000);
310   histRecPointsQLong->Sumw2();
311   Add2RecPointsList(histRecPointsQLong, kQLong, !expert, image);
312
313   TH1F * histRecPointsRow = 
314     new TH1F("hRecPointsRow", "Clusters per row; Row; Counts",
315              159, 0, 159);
316   histRecPointsRow->Sumw2();
317   Add2RecPointsList(histRecPointsRow, kRow, !expert, image);
318   //
319   ClonePerTrigClass(AliQAv1::kRECPOINTS); // this should be the last line
320 }
321
322 //____________________________________________________________________________
323 void AliTPCQADataMakerRec::MakeESDs(AliESDEvent * esd)
324 {
325   // make QA data from ESDs
326  
327   const Int_t nESDTracks = esd->GetNumberOfTracks();
328   Int_t nTPCtracks = 0; 
329   for(Int_t i = 0; i < nESDTracks; i++) {
330     
331     AliESDtrack * track = esd->GetTrack(i);
332     
333     if ((track->GetStatus() & AliESDtrack::kTPCrefit)==0)
334       continue;
335     
336     nTPCtracks++;
337     
338     Int_t nTPCclusters         = track->GetTPCNcls();
339     Int_t nTPCclustersFindable = track->GetTPCNclsF();
340     if ( nTPCclustersFindable<=0) continue;
341     FillESDsData(kClusters,nTPCclusters);
342     FillESDsData(kRatio,Float_t(nTPCclusters)/Float_t(nTPCclustersFindable));
343     FillESDsData(kPt,track->Pt()); 
344   }
345   //
346   IncEvCountCycleESDs();
347   IncEvCountTotalESDs();
348   //
349 }
350
351 //____________________________________________________________________________
352 void AliTPCQADataMakerRec::MakeRaws(AliRawReader* rawReader)
353 {
354   //
355   // To make QA for the RAW data we use the TPC Calibration framework 
356   // to handle the data and then in the end extract the data
357   //
358   
359   GetRawsData(0); // dummy call to init raw data
360   rawReader->Reset() ; 
361   if (! fTPCdataQA ) {
362
363     AliError("No TPC data QA (no call to InitRaws?)!!!!") ; 
364   } else {  
365
366     if(fTPCdataQA->GetIsDQM() == kFALSE)
367       AliError("Data QA has to be initialized as DQM!!!!") ; 
368
369     // Fill profile data
370     fTPCdataQA->ResetProfiles();
371     
372     if(fTPCdataQA->ProcessEvent(rawReader)) { // means that TPC data was processed  
373
374       fTPCdataQA->FillOccupancyProfile();
375       
376       // Fill histograms    
377       TObjArray *arrRW = GetMatchingRawsData(kRawsOccupancyVsSector); // all kRawsOccupancyVsSector clones matching to triggers
378       for (int ih=arrRW->GetEntriesFast();ih--;) {
379         TProfile* hRawsOccupancyVsSector = dynamic_cast<TProfile*>(arrRW->At(ih));
380         if (hRawsOccupancyVsSector) hRawsOccupancyVsSector->Add(fTPCdataQA->GetHistOccVsSector());
381       }
382       arrRW = GetMatchingRawsData(kRawsOccupancy2dVsSector);
383       for (int ih=arrRW->GetEntriesFast();ih--;) {
384         TProfile2D* hRawsOccupancy2dVsSector = dynamic_cast<TProfile2D*>(arrRW->At(ih));
385         if (hRawsOccupancy2dVsSector) hRawsOccupancy2dVsSector->Add(fTPCdataQA->GetHistOcc2dVsSector());
386       }
387       arrRW = GetMatchingRawsData(kRawsQVsSector);
388       for (int ih=arrRW->GetEntriesFast();ih--;) {
389         TProfile* hRawsQVsSector = dynamic_cast<TProfile*>(arrRW->At(ih));
390         if (hRawsQVsSector) hRawsQVsSector->Add(fTPCdataQA->GetHistQVsSector());
391       }
392       arrRW = GetMatchingRawsData(kRawsQmaxVsSector);
393       for (int ih=arrRW->GetEntriesFast();ih--;) {
394         TProfile* hRawsQmaxVsSector = dynamic_cast<TProfile*>(arrRW->At(ih));
395         if (hRawsQmaxVsSector) hRawsQmaxVsSector->Add(fTPCdataQA->GetHistQmaxVsSector());
396       }
397       //
398       IncEvCountCycleRaws();
399       IncEvCountTotalRaws();
400       //
401     }
402   }
403 }
404
405 //____________________________________________________________________________
406 void AliTPCQADataMakerRec::MakeDigits(TTree* digitTree)
407 {
408  
409   TBranch* branch = digitTree->GetBranch("Segment");
410   AliSimDigits* digArray = 0;
411   branch->SetAddress(&digArray);
412   
413   Int_t nEntries = Int_t(digitTree->GetEntries());
414   
415   for (Int_t n = 0; n < nEntries; n++) {
416     
417     digitTree->GetEvent(n);
418     
419     if (digArray->First())
420       do {
421         Float_t dig = digArray->CurrentDigit();
422         
423         FillDigitsData(kDigitsADC,dig);
424       } while (digArray->Next());    
425   }
426   //
427   IncEvCountCycleDigits();
428   IncEvCountTotalDigits();
429   //
430 }
431
432 //____________________________________________________________________________
433 void AliTPCQADataMakerRec::MakeRecPoints(TTree* recTree)
434 {
435   
436   AliTPCClustersRow* clrow = 0x0;
437   TBranch* branch = recTree->GetBranch("Segment");  
438   branch->SetAddress(&clrow);
439   TClonesArray * clarray = 0x0;
440
441   const Int_t nEntries = Int_t(recTree->GetEntries());
442   for (Int_t i = 0; i < nEntries; i++) {
443     
444     branch->GetEntry(i);
445
446     clarray = clrow->GetArray();
447
448     if (!clarray) continue;
449
450     const Int_t nClusters = clarray->GetEntriesFast();
451     for (Int_t icl=0; icl < nClusters; icl++){
452       
453       AliTPCclusterMI* cluster = 
454         (AliTPCclusterMI*)clarray->At(icl);
455       
456       Float_t Qmax = cluster->GetMax();
457       Float_t Q    = cluster->GetQ();
458       Int_t   row  = cluster->GetRow();
459
460       if(cluster->GetDetector()<36) { // IROC (short pads)
461
462         FillRecPointsData(kQmaxShort,Qmax);
463         FillRecPointsData(kQShort,Q);
464       } else { // OROC (medium and long pads)
465         row += 63;
466         if(cluster->GetRow()<64) { // medium pads
467
468           FillRecPointsData(kQmaxMedium,Qmax);
469           FillRecPointsData(kQMedium,Q);
470         } else { // long pads
471
472           FillRecPointsData(kQmaxLong,Qmax);
473           FillRecPointsData(kQLong,Q);
474         }
475       }
476       
477       FillRecPointsData(kRow,row);
478     } // end loop over clusters
479   } // end loop over tree
480   //
481   IncEvCountCycleRecPoints();
482   IncEvCountTotalRecPoints();
483   //
484 }
485
486 //____________________________________________________________________________
487 void AliTPCQADataMakerRec::LoadMaps()
488 {
489   TString path = gSystem->Getenv("ALICE_ROOT");
490   path += "/TPC/mapping/Patch";
491
492   for(Int_t i = 0; i < 6; i++) {
493
494     if(fMapping[i]!=0) // mapping already loaded
495       continue;
496     TString path2 = path;
497     path2 += i;
498     path2 += ".data";
499     fMapping[i] = new AliTPCAltroMapping(path2.Data());
500   }
501 }
502