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