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