]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliQADataMakerSim.cxx
Change according to Federico's suggestions (Massimo)
[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 <TFile.h>
28 #include <TTree.h>
29 #include <TClonesArray.h>
30
31 // --- Standard library ---
32
33 // --- AliRoot header files ---
34 #include "AliLog.h"
35 #include "AliQADataMakerSim.h"
36
37 ClassImp(AliQADataMakerSim)
38              
39 //____________________________________________________________________________ 
40 AliQADataMakerSim::AliQADataMakerSim(const char * name, const char * title) : 
41   AliQADataMaker(name, title), 
42   fDigitsQAList(0x0), 
43   fHitsQAList(0x0),
44   fSDigitsQAList(0x0)
45 {
46         // ctor
47         fDetectorDirName = GetName() ; 
48 }
49
50 //____________________________________________________________________________ 
51 AliQADataMakerSim::AliQADataMakerSim(const AliQADataMakerSim& qadm) :
52   AliQADataMaker(qadm.GetName(), qadm.GetTitle()), 
53   fDigitsQAList(qadm.fDigitsQAList),
54   fHitsQAList(qadm.fHitsQAList),
55   fSDigitsQAList(qadm.fSDigitsQAList) 
56 {
57   //copy ctor
58   fDetectorDirName = GetName() ; 
59 }
60
61 //____________________________________________________________________________ 
62 AliQADataMakerSim::~AliQADataMakerSim()
63 {
64         //dtor: delete the TObjArray and thei content
65         if ( fDigitsQAList ) { 
66                 if ( fDigitsQAList->IsOwner() )
67                         fDigitsQAList->Delete() ;     
68                 delete fDigitsQAList ;     
69         }
70         if ( fHitsQAList ) {
71                 if ( fHitsQAList->IsOwner() ) 
72                         fHitsQAList->Delete() ;
73                 delete fHitsQAList ;
74         }
75         if ( fSDigitsQAList ) { 
76                 if ( fSDigitsQAList->IsOwner() ) 
77                         fSDigitsQAList->Delete() ; 
78                 delete fSDigitsQAList ; 
79         }
80 }
81
82 //__________________________________________________________________
83 AliQADataMakerSim& AliQADataMakerSim::operator = (const AliQADataMakerSim& qadm )
84 {
85   // Assignment operator.
86   this->~AliQADataMakerSim();
87   new(this) AliQADataMakerSim(qadm);
88   return *this;
89 }
90
91 //____________________________________________________________________________
92 void AliQADataMakerSim::EndOfCycle(AliQA::TASKINDEX_t task) 
93
94   // Finishes a cycle of QA data acquistion
95         TObjArray * list = 0x0 ; 
96         
97         if ( task == AliQA::kHITS ) 
98                 list = fHitsQAList ; 
99         else if ( task == AliQA::kSDIGITS )
100                 list = fSDigitsQAList ; 
101         else if ( task == AliQA::kDIGITS ) 
102                 list = fDigitsQAList ; 
103         
104         EndOfDetectorCycle(task, list) ; 
105         TDirectory * subDir = fDetectorDir->GetDirectory(AliQA::GetTaskName(task)) ; 
106         if (subDir) { 
107                 subDir->cd() ; 
108                 list->Write() ; 
109         }
110 }
111  
112 //____________________________________________________________________________
113 void AliQADataMakerSim::Exec(AliQA::TASKINDEX_t task, TObject * data) 
114
115   // creates the quality assurance data for the various tasks (Hits, SDigits, Digits, ESDs)
116     
117         if ( task == AliQA::kHITS ) {  
118                 AliDebug(1, "Processing Hits QA") ; 
119                 TClonesArray * arr = dynamic_cast<TClonesArray *>(data) ; 
120                 if (arr) { 
121                         MakeHits(arr) ;
122                 } else {
123                         TTree * tree = dynamic_cast<TTree *>(data) ; 
124                         if (tree) {
125                                 MakeHits(tree) ; 
126                         } else {
127                                 AliWarning("data are neither a TClonesArray nor a TTree") ; 
128                         }
129                 }
130         } else if ( task == AliQA::kSDIGITS ) {
131                 AliDebug(1, "Processing SDigits QA") ; 
132                 TClonesArray * arr = dynamic_cast<TClonesArray *>(data) ; 
133                 if (arr) { 
134                         MakeSDigits(arr) ;
135                 } else {
136                         TTree * tree = dynamic_cast<TTree *>(data) ; 
137                         if (tree) {
138                                 MakeSDigits(tree) ; 
139                         } else {
140                                 AliWarning("data are neither a TClonesArray nor a TTree") ; 
141                         }
142                 }
143         } else if ( task == AliQA::kDIGITS ) {
144                 AliDebug(1, "Processing Digits QA") ; 
145                 TClonesArray * arr = dynamic_cast<TClonesArray *>(data) ; 
146                 if (arr) { 
147                         MakeDigits(arr) ;
148                 } else {
149                         TTree * tree = dynamic_cast<TTree *>(data) ; 
150                         if (tree) {
151                                 MakeDigits(tree) ; 
152                         } else {
153                                 AliWarning("data are neither a TClonesArray nor a TTree") ; 
154                         }
155                 }
156         }
157 }
158
159 //____________________________________________________________________________ 
160 TObjArray *  AliQADataMakerSim::Init(AliQA::TASKINDEX_t task, Int_t run, Int_t cycles)
161 {
162   // general intialisation
163         
164         fRun = run ;
165         if (cycles > 0)
166                 SetCycle(cycles) ;  
167         TObjArray * rv = NULL ; 
168         if ( task == AliQA::kHITS ) {
169                 if ( ! fHitsQAList ) {
170                         fHitsQAList = new TObjArray(100) ;       
171                         InitHits() ;
172                 }
173                 rv = fHitsQAList ;
174         } else if ( task == AliQA::kSDIGITS ) {
175                 if ( ! fSDigitsQAList ) {
176                         fSDigitsQAList = new TObjArray(100) ; 
177                         InitSDigits() ;
178                 }
179                 rv = fSDigitsQAList ;
180    } else if ( task == AliQA::kDIGITS ) {
181            if ( ! fDigitsQAList ) {
182                    fDigitsQAList = new TObjArray(100) ;
183                    InitDigits() ;
184            }
185            rv =  fDigitsQAList ;
186    }
187   
188         return rv ; 
189 }
190
191 //____________________________________________________________________________ 
192 void AliQADataMakerSim::Init(AliQA::TASKINDEX_t task, TObjArray * list, Int_t run, Int_t cycles)
193 {
194   // Intialisation by passing the list of QA data booked elsewhere
195   
196         fRun = run ;
197         if (cycles > 0)
198                 SetCycle(cycles) ;  
199         
200         if ( task == AliQA::kHITS ) {
201                 fHitsQAList = list ;     
202         } else if ( task == AliQA::kSDIGITS) {
203                 fSDigitsQAList = list ; 
204         } else if ( task == AliQA::kDIGITS ) {
205                 fDigitsQAList = list ; 
206         } 
207 }
208
209 //____________________________________________________________________________
210 void AliQADataMakerSim::StartOfCycle(AliQA::TASKINDEX_t task, const Bool_t sameCycle) 
211
212   // Finishes a cycle of QA data acquistion
213         if ( !sameCycle || fCurrentCycle == -1) {
214                 ResetCycle() ;
215         if (fOutput) 
216                 fOutput->Close() ; 
217         fOutput = AliQA::GetQADataFile(GetName(), fRun, fCurrentCycle) ;        
218         }       
219
220         AliInfo(Form(" Run %d Cycle %d task %s file %s", 
221                                  fRun, fCurrentCycle, AliQA::GetTaskName(task).Data(), fOutput->GetName() )) ;
222
223         fDetectorDir = fOutput->GetDirectory(GetDetectorDirName()) ; 
224         if (!fDetectorDir)
225                 fDetectorDir = fOutput->mkdir(GetDetectorDirName()) ; 
226
227         TDirectory * subDir = fDetectorDir->GetDirectory(AliQA::GetTaskName(task)) ; 
228         if (!subDir)
229                 subDir = fDetectorDir->mkdir(AliQA::GetTaskName(task)) ;  
230         subDir->cd() ; 
231         
232         TObjArray * list = 0x0 ; 
233   
234         if ( task == AliQA::kHITS )  
235                 list = fHitsQAList ;
236         else if ( task == AliQA::kSDIGITS )  
237                 list = fSDigitsQAList ;
238         else  if ( task == AliQA::kDIGITS ) 
239                 list = fDigitsQAList ;
240         
241 // Should be the choice of detectors
242 //      TIter next(list) ;
243 //      TH1 * h ; 
244 //      while ( (h = dynamic_cast<TH1 *>(next())) )
245 //              h->Reset() ;  
246 //
247         StartOfDetectorCycle() ; 
248 }