]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliQADataMakerSim.cxx
bring back some functionality which was deleted by accident
[u/mrichter/AliRoot.git] / STEER / AliQADataMakerSim.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
17 /* $Id$ */
18
19 //
20 //  Base Class
21 //  Produces the data needed to calculate the quality assurance. 
22 //  All data must be mergeable objects.
23 //  Y. Schutz CERN July 2007
24 //
25
26 // --- ROOT system ---
27 #include <TCanvas.h>
28 #include <TFile.h>
29 #include <TTree.h>
30 #include <TClonesArray.h>
31
32 // --- Standard library ---
33
34 // --- AliRoot header files ---
35 #include "AliLog.h"
36 #include "AliQADataMakerSim.h"
37
38 ClassImp(AliQADataMakerSim)
39              
40 //____________________________________________________________________________ 
41 AliQADataMakerSim::AliQADataMakerSim(const char * name, const char * title) : 
42   AliQADataMaker(name, title), 
43   fDigitsQAList(NULL), 
44   fHitsQAList(NULL),
45   fSDigitsQAList(NULL),  
46   fHitsArray(NULL),
47   fSDigitsArray(NULL)
48 {
49         // ctor
50         fDetectorDirName = GetName() ; 
51 }
52
53 //____________________________________________________________________________ 
54 AliQADataMakerSim::AliQADataMakerSim(const AliQADataMakerSim& qadm) :
55   AliQADataMaker(qadm.GetName(), qadm.GetTitle()), 
56   fDigitsQAList(qadm.fDigitsQAList),
57   fHitsQAList(qadm.fHitsQAList),
58   fSDigitsQAList(qadm.fSDigitsQAList),  
59   fHitsArray(NULL),
60   fSDigitsArray(NULL)
61 {
62   //copy ctor
63   fDetectorDirName = GetName() ; 
64 }
65
66 //____________________________________________________________________________ 
67 AliQADataMakerSim::~AliQADataMakerSim()
68 {
69         //dtor: delete the TObjArray and thei content
70         if ( fDigitsQAList ) { 
71     for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
72       if ( fDigitsQAList[specie]->IsOwner() )
73                         fDigitsQAList[specie]->Delete() ;
74     }
75                 delete[] fDigitsQAList ;
76   }
77         if ( fHitsQAList ) {
78     for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
79       if ( fHitsQAList[specie]->IsOwner() ) 
80                         fHitsQAList[specie]->Delete() ;
81     }
82         delete[] fHitsQAList ;
83   }
84         if ( fSDigitsQAList ) { 
85     for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
86       if ( fSDigitsQAList[specie]->IsOwner() ) 
87                         fSDigitsQAList[specie]->Delete() ; 
88     }
89                 delete[] fSDigitsQAList ;
90   }
91   if (fHitsArray) {
92     fHitsArray->Clear() ; 
93     delete fHitsArray ;
94   }
95   if (fSDigitsArray) {
96     fSDigitsArray->Clear() ; 
97     delete fSDigitsArray ;
98   }  
99 }
100
101 //__________________________________________________________________
102 AliQADataMakerSim& AliQADataMakerSim::operator = (const AliQADataMakerSim& qadm )
103 {
104   // Assignment operator.
105   this->~AliQADataMakerSim();
106   new(this) AliQADataMakerSim(qadm);
107   return *this;
108 }
109
110 //____________________________________________________________________________
111 void AliQADataMakerSim::EndOfCycle() 
112
113   // Finishes a cycle of QA for all tasks
114   EndOfCycle(AliQAv1::kHITS) ; 
115   EndOfCycle(AliQAv1::kSDIGITS) ; 
116   EndOfCycle(AliQAv1::kDIGITS) ;
117   ResetCycle() ; 
118 }
119
120 //____________________________________________________________________________
121 void AliQADataMakerSim::EndOfCycle(AliQAv1::TASKINDEX_t task) 
122
123   // Finishes a cycle of QA data acquistion
124         TObjArray ** list = NULL ; 
125         
126         if ( task == AliQAv1::kHITS ) 
127                 list = fHitsQAList ; 
128         else if ( task == AliQAv1::kSDIGITS )
129                 list = fSDigitsQAList ; 
130         else if ( task == AliQAv1::kDIGITS ) 
131                 list = fDigitsQAList ; 
132   
133   if ( ! list ) 
134     return ; 
135         EndOfDetectorCycle(task, list) ; 
136   fDetectorDir = fOutput->GetDirectory(GetDetectorDirName()) ;
137         if (!fDetectorDir) 
138     fDetectorDir = fOutput->mkdir(GetDetectorDirName()) ; 
139   TDirectory * subDir = fDetectorDir->GetDirectory(AliQAv1::GetTaskName(task)) ; 
140   if (!subDir)
141     subDir = fDetectorDir->mkdir(AliQAv1::GetTaskName(task)) ;  
142   subDir->cd() ; 
143   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
144     if (! AliQAv1::Instance(AliQAv1::GetDetIndex(GetName()))->IsEventSpecieSet(AliRecoParam::ConvertIndex(specie)) ) 
145       continue ;
146     TDirectory * eventSpecieDir = subDir->GetDirectory(AliRecoParam::GetEventSpecieName(specie)) ;
147     if (!eventSpecieDir) 
148       eventSpecieDir = subDir->mkdir(AliRecoParam::GetEventSpecieName(specie)) ; 
149     eventSpecieDir->cd() ; 
150     TIter next(list[specie]) ; 
151     TObject * obj ; 
152     while ( (obj = next()) )  {
153       if (!obj->TestBit(AliQAv1::GetExpertBit()))
154         obj->Write() ;
155     }
156     if (WriteExpert()) {
157       TDirectory * expertDir = eventSpecieDir->GetDirectory(AliQAv1::GetExpert()) ; 
158       if (!expertDir) 
159         expertDir = eventSpecieDir->mkdir(AliQAv1::GetExpert()) ; 
160       expertDir->cd() ;
161       next.Reset() ; 
162       while ( (obj = next()) ) {
163         if (!obj->TestBit(AliQAv1::GetExpertBit()))
164           continue ; 
165         obj->Write() ;
166       }      
167     }
168     fOutput->Save() ; 
169   }
170   if (fPrintImage) 
171     MakeImage(task) ; 
172 }
173
174 //____________________________________________________________________________
175 void AliQADataMakerSim::Exec(AliQAv1::TASKINDEX_t task, TObject * data) 
176
177   // creates the quality assurance data for the various tasks (Hits, SDigits, Digits, ESDs)
178   
179         if ( task == AliQAv1::kHITS ) {  
180                 AliDebug(AliQAv1::GetQADebugLevel(), "Processing Hits QA") ; 
181                 if (strcmp(data->ClassName(), "TClonesArray") == 0) { 
182       fHitsArray = static_cast<TClonesArray *>(data) ; 
183                         MakeHits() ;
184                 } else if (strcmp(data->ClassName(), "TTree") == 0) {
185                         TTree * tree = static_cast<TTree *>(data) ; 
186       MakeHits(tree) ; 
187     } else {
188       AliWarning("data are neither a TClonesArray nor a TTree") ; 
189     }
190         } else if ( task == AliQAv1::kSDIGITS ) {
191                 AliDebug(AliQAv1::GetQADebugLevel(), "Processing SDigits QA") ; 
192                 if (strcmp(data->ClassName(), "TClonesArray") == 0) { 
193       fSDigitsArray = static_cast<TClonesArray *>(data) ; 
194                         MakeSDigits() ;
195                 } else if (strcmp(data->ClassName(), "TTree") == 0) {
196                         TTree * tree = static_cast<TTree *>(data) ; 
197       MakeSDigits(tree) ; 
198     } else {
199       AliWarning("data are neither a TClonesArray nor a TTree") ; 
200     }
201         } else if ( task == AliQAv1::kDIGITS ) {
202                 AliDebug(AliQAv1::GetQADebugLevel(), "Processing Digits QA") ; 
203                 if (strcmp(data->ClassName(), "TClonesArray") == 0) { 
204       fDigitsArray = static_cast<TClonesArray *>(data) ; 
205                         MakeDigits() ;
206                 } else if (strcmp(data->ClassName(), "TTree") == 0)  {
207                         TTree * tree = static_cast<TTree *>(data) ; 
208       MakeDigits(tree) ; 
209     } else {
210       AliWarning("data are neither a TClonesArray nor a TTree") ; 
211     }
212   }
213 }
214
215 //____________________________________________________________________________ 
216 void AliQADataMakerSim::MakeImage(AliQAv1::TASKINDEX_t task)
217 {
218   // create a drawing of detetor defined histograms
219   TObjArray ** list = NULL ;  
220   switch (task) {
221     case AliQAv1::kRAWS:
222       break;
223     case AliQAv1::kHITS:
224       list = fHitsQAList ;
225       break;
226     case AliQAv1::kSDIGITS:
227       list = fSDigitsQAList ;
228       break;  
229     case AliQAv1::kDIGITS:
230       list = fDigitsQAList ;
231       break;  
232     case AliQAv1::kDIGITSR:
233       break;
234     case AliQAv1::kRECPOINTS:
235       break;
236     case AliQAv1::kTRACKSEGMENTS:
237       break;
238     case AliQAv1::kRECPARTICLES:
239       break;
240     case AliQAv1::kESDS:
241       break;
242     case AliQAv1::kNTASKINDEX:
243       break;
244     default:
245     break;
246   }
247   if ( !list) {
248     AliFatal("data not initialized, call AliQADataMaker::Init"); 
249   return ; 
250   }
251   MakeTheImage(list, task, "Sim") ; 
252 }
253
254 //____________________________________________________________________________ 
255 TObjArray **  AliQADataMakerSim::Init(AliQAv1::TASKINDEX_t task, Int_t cycles)
256 {
257   // general intialisation
258         
259         if (cycles > 0)
260                 SetCycle(cycles) ;  
261         TObjArray ** rv = NULL ; 
262         if ( task == AliQAv1::kHITS ) {
263                 if ( ! fHitsQAList ) {
264       fHitsQAList = new TObjArray *[AliRecoParam::kNSpecies] ; 
265       for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
266         fHitsQAList[specie] = new TObjArray(AliQAv1::GetMaxQAObj()) ;    
267         fHitsQAList[specie]->SetName(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(task).Data(), AliRecoParam::GetEventSpecieName(specie))) ;
268       }
269                 }
270                 rv = fHitsQAList ;
271         } else if ( task == AliQAv1::kSDIGITS ) {
272                 if ( ! fSDigitsQAList ) {
273       fSDigitsQAList = new TObjArray *[AliRecoParam::kNSpecies] ; 
274       for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
275         fSDigitsQAList[specie] = new TObjArray(AliQAv1::GetMaxQAObj()) ; 
276         fSDigitsQAList[specie]->SetName(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(task).Data(), AliRecoParam::GetEventSpecieName(specie))) ; 
277       }
278                 }
279                 rv = fSDigitsQAList ;
280    } else if ( task == AliQAv1::kDIGITS ) {
281            if ( ! fDigitsQAList ) {
282        fDigitsQAList = new TObjArray *[AliRecoParam::kNSpecies] ; 
283        for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {    
284          fDigitsQAList[specie] = new TObjArray(AliQAv1::GetMaxQAObj()) ;
285          fDigitsQAList[specie]->SetName(Form("%s_%s_%s", GetName(), AliQAv1::GetTaskName(task).Data(), AliRecoParam::GetEventSpecieName(specie))) ;
286        }
287            }
288            rv =  fDigitsQAList ;
289    }
290   
291         return rv ; 
292
293
294 //____________________________________________________________________________ 
295 void AliQADataMakerSim::Init(AliQAv1::TASKINDEX_t task, TObjArray ** list, Int_t run, Int_t cycles)
296 {
297   // Intialisation by passing the list of QA data booked elsewhere
298   
299         fRun = run ;
300         if (cycles > 0)
301                 SetCycle(cycles) ;  
302         
303         if ( task == AliQAv1::kHITS ) {
304                 fHitsQAList = list ;     
305         } else if ( task == AliQAv1::kSDIGITS) {
306                 fSDigitsQAList = list ; 
307         } else if ( task == AliQAv1::kDIGITS ) {
308                 fDigitsQAList = list ; 
309         } 
310 }
311
312 //____________________________________________________________________________
313 void AliQADataMakerSim::StartOfCycle(Int_t run) 
314
315   // Finishes a cycle of QA for all tasks
316   Bool_t samecycle = kFALSE ; 
317   StartOfCycle(AliQAv1::kHITS,    run, samecycle) ;
318   samecycle = kTRUE ; 
319   StartOfCycle(AliQAv1::kSDIGITS, run, samecycle) ;
320   StartOfCycle(AliQAv1::kDIGITS,  run, samecycle) ;
321 }
322
323 //____________________________________________________________________________
324 void AliQADataMakerSim::StartOfCycle(AliQAv1::TASKINDEX_t task, Int_t run, const Bool_t sameCycle) 
325
326   // Finishes a cycle of QA data acquistion
327   if ( run > 0 ) 
328     fRun = run ; 
329         if ( !sameCycle || fCurrentCycle == -1) {
330                 ResetCycle() ;
331         if (fOutput) 
332                 fOutput->Close() ; 
333         fOutput = AliQAv1::GetQADataFile(GetName(), fRun) ;     
334         }       
335
336         AliDebug(AliQAv1::GetQADebugLevel(), Form(" Run %d Cycle %d task %s file %s", 
337                                  fRun, fCurrentCycle, AliQAv1::GetTaskName(task).Data(), fOutput->GetName() )) ;
338
339         //fDetectorDir = fOutput->GetDirectory(GetDetectorDirName()) ; 
340 //      if (!fDetectorDir)
341 //              fDetectorDir = fOutput->mkdir(GetDetectorDirName()) ; 
342 //
343 //      TDirectory * subDir = fDetectorDir->GetDirectory(AliQAv1::GetTaskName(task)) ; 
344 //      if (!subDir)
345 //              subDir = fDetectorDir->mkdir(AliQAv1::GetTaskName(task)) ;  
346 //  
347 //  for ( Int_t index = AliRecoParam::kDefault ; index < AliRecoParam::kNSpecies ; index++ ) {
348 //    TDirectory * eventSpecieDir = subDir->GetDirectory(AliRecoParam::GetEventSpecieName(index)) ; 
349 //    if (!eventSpecieDir) 
350 //      eventSpecieDir = subDir->mkdir(AliRecoParam::GetEventSpecieName(index)) ; 
351 //    TDirectory * expertDir = eventSpecieDir->GetDirectory(AliQAv1::GetExpert()) ; 
352 //    if (!expertDir) 
353 //      expertDir = eventSpecieDir->mkdir(AliQAv1::GetExpert()) ; 
354 //   }   
355         StartOfDetectorCycle() ; 
356 }