]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliQADataMakerSim.cxx
Getter for 4-momentum added.
[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::operator = (const AliQADataMakerSim& qadm )
63 {
64   // Assignment operator.
65   this->~AliQADataMakerSim();
66   new(this) AliQADataMakerSim(qadm);
67   return *this;
68 }
69
70 //____________________________________________________________________________
71 void AliQADataMakerSim::EndOfCycle(AliQA::TASKINDEX_t task) 
72
73   // Finishes a cycle of QA data acquistion
74         TObjArray * list = 0x0 ; 
75         
76         if ( task == AliQA::kHITS ) 
77                 list = fHitsQAList ; 
78         else if ( task == AliQA::kSDIGITS )
79                 list = fSDigitsQAList ; 
80         else if ( task == AliQA::kDIGITS ) 
81                 list = fDigitsQAList ; 
82         
83         EndOfDetectorCycle(task, list) ; 
84         TDirectory * subDir = fDetectorDir->GetDirectory(AliQA::GetTaskName(task)) ; 
85         if (subDir) { 
86                 subDir->cd() ; 
87                 list->Write() ; 
88         }
89 }
90  
91 //____________________________________________________________________________
92 void AliQADataMakerSim::Exec(AliQA::TASKINDEX_t task, TObject * data) 
93
94   // creates the quality assurance data for the various tasks (Hits, SDigits, Digits, ESDs)
95     
96         if ( task == AliQA::kHITS ) {  
97                 AliDebug(1, "Processing Hits QA") ; 
98                 TClonesArray * arr = dynamic_cast<TClonesArray *>(data) ; 
99                 if (arr) { 
100                         MakeHits(arr) ;
101                 } else {
102                         TTree * tree = dynamic_cast<TTree *>(data) ; 
103                         if (tree) {
104                                 MakeHits(tree) ; 
105                         } else {
106                                 AliWarning("data are neither a TClonesArray nor a TTree") ; 
107                         }
108                 }
109         } else if ( task == AliQA::kSDIGITS ) {
110                 AliDebug(1, "Processing SDigits QA") ; 
111                 TClonesArray * arr = dynamic_cast<TClonesArray *>(data) ; 
112                 if (arr) { 
113                         MakeSDigits(arr) ;
114                 } else {
115                         TTree * tree = dynamic_cast<TTree *>(data) ; 
116                         if (tree) {
117                                 MakeSDigits(tree) ; 
118                         } else {
119                                 AliWarning("data are neither a TClonesArray nor a TTree") ; 
120                         }
121                 }
122         } else if ( task == AliQA::kDIGITS ) {
123                 AliDebug(1, "Processing Digits QA") ; 
124                 TClonesArray * arr = dynamic_cast<TClonesArray *>(data) ; 
125                 if (arr) { 
126                         MakeDigits(arr) ;
127                 } else {
128                         TTree * tree = dynamic_cast<TTree *>(data) ; 
129                         if (tree) {
130                                 MakeDigits(tree) ; 
131                         } else {
132                                 AliWarning("data are neither a TClonesArray nor a TTree") ; 
133                         }
134                 }
135         }
136 }
137
138 //____________________________________________________________________________ 
139 TObjArray *  AliQADataMakerSim::Init(AliQA::TASKINDEX_t task, Int_t run, Int_t cycles)
140 {
141   // general intialisation
142         
143         fRun = run ;
144         if (cycles > 0)
145                 SetCycle(cycles) ;  
146         TObjArray * rv = NULL ; 
147         if ( task == AliQA::kHITS ) {
148                 fHitsQAList = new TObjArray(100) ;       
149                 InitHits() ;
150                 rv = fHitsQAList ;
151         } else if ( task == AliQA::kSDIGITS ) {
152                 fSDigitsQAList = new TObjArray(100) ; 
153                 InitSDigits() ;
154                 rv = fSDigitsQAList ;
155    } else if ( task == AliQA::kDIGITS ) {
156            fDigitsQAList = new TObjArray(100); 
157            InitDigits() ;
158            rv =  fDigitsQAList ;
159    }
160   
161         return rv ; 
162 }
163
164 //____________________________________________________________________________ 
165 void AliQADataMakerSim::Init(AliQA::TASKINDEX_t task, TObjArray * list, Int_t run, Int_t cycles)
166 {
167   // Intialisation by passing the list of QA data booked elsewhere
168   
169         fRun = run ;
170         if (cycles > 0)
171                 SetCycle(cycles) ;  
172         
173         if ( task == AliQA::kHITS ) {
174                 fHitsQAList = list ;     
175         } else if ( task == AliQA::kSDIGITS) {
176                 fSDigitsQAList = list ; 
177         } else if ( task == AliQA::kDIGITS ) {
178                 fDigitsQAList = list ; 
179         } 
180 }
181
182 //____________________________________________________________________________
183 void AliQADataMakerSim::StartOfCycle(AliQA::TASKINDEX_t task, const Bool_t sameCycle) 
184
185   // Finishes a cycle of QA data acquistion
186         if ( !sameCycle || fCurrentCycle == -1) {
187                 ResetCycle() ;
188         if (fOutput) 
189                 fOutput->Close() ; 
190         fOutput = AliQA::GetQADataFile(GetName(), fRun, fCurrentCycle) ;        
191         }       
192
193         AliInfo(Form(" Run %d Cycle %d task %s file %s", 
194                                  fRun, fCurrentCycle, AliQA::GetTaskName(task).Data(), fOutput->GetName() )) ;
195
196         fDetectorDir = fOutput->GetDirectory(GetDetectorDirName()) ; 
197         if (!fDetectorDir)
198                 fDetectorDir = fOutput->mkdir(GetDetectorDirName()) ; 
199
200         TDirectory * subDir = fDetectorDir->GetDirectory(AliQA::GetTaskName(task)) ; 
201         if (!subDir)
202                 subDir = fDetectorDir->mkdir(AliQA::GetTaskName(task)) ;  
203         subDir->cd() ; 
204         
205         TObjArray * list = 0x0 ; 
206   
207         if ( task == AliQA::kHITS )  
208                 list = fHitsQAList ;
209         else if ( task == AliQA::kSDIGITS )  
210                 list = fSDigitsQAList ;
211         else  if ( task == AliQA::kDIGITS ) 
212                 list = fDigitsQAList ;
213         
214 // Should be the choice of detectors
215 //      TIter next(list) ;
216 //      TH1 * h ; 
217 //      while ( (h = dynamic_cast<TH1 *>(next())) )
218 //              h->Reset() ;  
219 //
220         StartOfDetectorCycle() ; 
221 }