]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliCorrQADataMakerRec.cxx
added the notion of expert QA data
[u/mrichter/AliRoot.git] / STEER / AliCorrQADataMakerRec.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: AliCorrQADataMakerRec.cxx 27570 2008-07-24 21:49:27Z cvetan $ */
18
19 /*
20   Produces the data needed to calculate the quality assurance. 
21   All data must be mergeable objects.
22   Y. Schutz CERN July 2007
23 */
24
25 // --- ROOT system ---
26 #include <TClonesArray.h>
27 #include <TFile.h> 
28 #include <TH1F.h> 
29 #include <TH1I.h> 
30 #include <TH2F.h> 
31 #include <TNtupleD.h>
32 #include <TParameter.h>
33
34 // --- Standard library ---
35
36 // --- AliRoot header files ---
37 #include "AliLog.h"
38 #include "AliCorrQADataMakerRec.h"
39 #include "AliQAChecker.h"
40
41 ClassImp(AliCorrQADataMakerRec)
42            
43 //____________________________________________________________________________ 
44 AliCorrQADataMakerRec::AliCorrQADataMakerRec(AliQADataMaker ** qadm ) : 
45 AliQADataMakerRec(AliQA::GetDetName(AliQA::kCORR), "Corr Quality Assurance Data Maker"),
46   fMaxRawVar(0),  
47   fqadm(qadm)
48 {
49   // ctor
50
51 }
52
53 //____________________________________________________________________________ 
54 AliCorrQADataMakerRec::AliCorrQADataMakerRec(const AliCorrQADataMakerRec& qadm) :
55   AliQADataMakerRec(),
56   fMaxRawVar(qadm.fMaxRawVar), 
57   fqadm(qadm.fqadm)
58 {
59   //copy ctor 
60   SetName((const char*)qadm.GetName()) ; 
61   SetTitle((const char*)qadm.GetTitle()); 
62 }
63
64 //__________________________________________________________________
65 AliCorrQADataMakerRec& AliCorrQADataMakerRec::operator = (const AliCorrQADataMakerRec& qadm )
66 {
67   // Equal operator.
68   this->~AliCorrQADataMakerRec();
69   new(this) AliCorrQADataMakerRec(qadm);
70   return *this;
71 }
72  
73 //____________________________________________________________________________ 
74 void AliCorrQADataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray * /*list*/)
75 {
76   //Detector specific actions at end of cycle
77   // do the QA checking
78   if (task == AliQA::kRAWS) {
79      AliQAChecker::Instance()->Run(AliQA::kCORR, task, fObject) ; 
80   }
81 }
82
83 //____________________________________________________________________________ 
84 void AliCorrQADataMakerRec::InitESDs()
85 {
86   //Create histograms to controll ESD
87
88   AliInfo("TO BE IMPLEMENTED") ; 
89 }
90
91 //____________________________________________________________________________ 
92 void AliCorrQADataMakerRec::InitRecPoints()
93 {
94   // create Reconstructed Points histograms in RecPoints subdir
95
96   AliInfo("TO BE IMPLEMENTED") ; 
97 }
98
99 //____________________________________________________________________________ 
100 void AliCorrQADataMakerRec::InitRaws()
101 {
102   // createa ntuple taking all the parameters declared by detectors
103   if (fObject) 
104     return ; 
105   delete fRawsQAList ; // not needed for the time being 
106   fRawsQAList = NULL ; 
107   TString varlist("") ;
108   for ( Int_t detIndex = 0 ; detIndex < AliQA::kNDET ; detIndex++ ) {
109     AliQADataMaker * qadm = fqadm[detIndex] ; 
110     if ( ! qadm ) 
111       continue ;
112     TList * list = qadm->GetParameterList() ; 
113     if (list) {
114       TIter next(list) ; 
115       TParameter<double> * p ; 
116       while ( (p = dynamic_cast<TParameter<double>*>(next()) ) ) {
117         varlist.Append(p->GetName()) ; 
118         varlist.Append(":") ; 
119         fMaxRawVar++ ; 
120       }
121     }
122   }
123   varlist = varlist.Strip(TString::kTrailing, ':') ; 
124   if (fMaxRawVar == 0) { 
125     AliWarning("NTUPLE not created") ; 
126   } else {
127     fObject = new TNtupleD(AliQA::GetQACorrName(), "Raws data correlation among detectors", varlist.Data()) ;  
128   }  
129 }
130
131 //____________________________________________________________________________
132 void AliCorrQADataMakerRec::MakeESDs(AliESDEvent * /*esd*/)
133 {
134   // make QA data from ESDs
135
136   AliInfo("TO BE IMPLEMENTED") ; 
137
138 }
139
140 //____________________________________________________________________________
141 void AliCorrQADataMakerRec::MakeRaws(AliRawReader *)
142 {
143   //Fill prepared histograms with Raw digit properties
144   if ( fMaxRawVar > 0 ) {
145     const Int_t kSize = fMaxRawVar ; 
146     Double_t  *varvalue = new Double_t[kSize] ;
147     Int_t index = 0 ;
148     for ( Int_t detIndex = 0 ; detIndex < AliQA::kNDET ; detIndex++ ) {
149       AliQADataMaker * qadm = fqadm[detIndex] ; 
150       if ( ! qadm ) 
151         continue ;
152       TList * list = qadm->GetParameterList() ; 
153       TIter next(list) ; 
154       TParameter<double> * p ; 
155       while ( (p = dynamic_cast<TParameter<double>*>(next()) ) ) {
156         varvalue[index++] = p->GetVal() ; 
157       }
158     }
159     (dynamic_cast<TNtupleD*>(fObject))->Fill(varvalue);
160     delete [] varvalue;
161   }
162 }
163
164 //____________________________________________________________________________
165 void AliCorrQADataMakerRec::MakeRecPoints(TTree * /*clustersTree*/)
166 {
167   AliInfo("TO BE IMPLEMENTED") ; 
168 }
169
170 //____________________________________________________________________________ 
171 void AliCorrQADataMakerRec::StartOfDetectorCycle()
172 {
173   //Detector specific actions at start of cycle  
174
175 }