]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliQADataMaker.h
bug fix in T00 alias names
[u/mrichter/AliRoot.git] / STEER / AliQADataMaker.h
1 #ifndef ALIQADATAMAKER_H
2 #define ALIQADATAMAKER_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  * See cxx source for full Copyright notice                               */
5
6
7 /* $Id$ */
8
9 /*
10   Base Class:
11   Produces the data needed to calculate the quality assurance. 
12   All data must be mergeable objects.
13   Y. Schutz CERN July 2007
14 */
15
16
17 // --- ROOT system ---
18 #include <TH1.h>
19 #include <TList.h>
20 #include <TNamed.h>  
21 class TFile;  
22 class TDirectory;
23 class TObject; 
24 class TTree; 
25 class AliESDEvent;
26 class AliRawReader;
27 class TClonesArray;
28
29 // --- Standard library ---
30
31 // --- AliRoot header files ---
32 #include "AliQA.h"
33
34 class AliQADataMaker: public TNamed {
35   
36 public:
37   
38   AliQADataMaker(const char * name="", const char * title="") ;          // ctor
39   AliQADataMaker(const AliQADataMaker& qadm) ;   
40   AliQADataMaker& operator = (const AliQADataMaker& qadm) ;
41   virtual ~AliQADataMaker() ; // dtor
42   
43   virtual void        Exec(AliQA::TASKINDEX, TObject * data) ;
44   void                EndOfCycle(AliQA::TASKINDEX) ;
45   void                Finish(AliQA::TASKINDEX task) const ; 
46   static const char * GetDetectorDirName() { return fDetectorDirName.Data() ; }
47   const Int_t         Increment() { return ++fCycleCounter ; } 
48   TList *             Init(AliQA::TASKINDEX, Int_t run, Int_t cycles = -1) ;
49   const Bool_t        IsCycleDone() const { return fCycleCounter > fCycle ? kTRUE : kFALSE ; }
50   const Int_t         Add2DigitsList(TH1 * hist, const Int_t index)    { return Add2List(hist, index, fDigitsQAList) ; }
51   const Int_t         Add2ESDsList(TH1 * hist, const Int_t index)      { return Add2List(hist, index, fESDsQAList) ; }
52   const Int_t         Add2HitsList(TH1 * hist, const Int_t index)      { return Add2List(hist, index, fHitsQAList) ; }
53   const Int_t         Add2RecPointsList(TH1 * hist, const Int_t index) { return Add2List(hist, index, fRecPointsQAList) ; }
54   const Int_t         Add2RawsList(TH1 * hist, const Int_t index)      { return Add2List(hist, index, fRawsQAList) ; }
55   const Int_t         Add2SDigitsList(TH1 * hist, const Int_t index)   { return Add2List(hist, index, fSDigitsQAList) ; }
56   TH1 *               GetDigitsData(const Int_t index)    { return dynamic_cast<TH1 *>(GetData(fDigitsQAList, index)) ; }
57   TH1 *               GetESDsData(const Int_t index)      { return dynamic_cast<TH1 *>(GetData(fESDsQAList, index)) ; }
58   TH1 *               GetHitsData(const Int_t index)      { return dynamic_cast<TH1 *>(GetData(fHitsQAList, index)) ; }
59   TH1 *               GetRecPointsData(const Int_t index) { return dynamic_cast<TH1 *>(GetData(fRecPointsQAList, index)) ; }
60   TH1 *               GetRawsData(const Int_t index)      { return dynamic_cast<TH1 *>(GetData(fRawsQAList, index)) ; }
61   TH1 *               GetSDigitsData(const Int_t index)   { return dynamic_cast<TH1 *>(GetData(fSDigitsQAList, index)) ; }
62   void                SetCycle(Int_t nevts) { fCycle = nevts ; } 
63   void                StartOfCycle(AliQA::TASKINDEX, Option_t * sameCycle = "new") ;
64
65 protected: 
66
67   Int_t          Add2List(TH1 * hist, const Int_t index, TList * list) { list->AddAt(hist, index) ; return list->LastIndex() ; }
68   virtual void   EndOfDetectorCycle(AliQA::TASKINDEX, TList * ) {AliInfo("To be implemented by detectors");} 
69   TObject *      GetData(TList * list, const Int_t index)  { return list->At(index) ; } 
70   virtual void   InitDigits()        {AliInfo("To be implemented by detectors");}
71   virtual void   InitESDs()          {AliInfo("To be implemented by detectors");}
72   virtual void   InitHits()          {AliInfo("To be implemented by detectors");}
73   //virtual void   InitRecParticles()  {AliInfo("To be implemented by detectors");}
74   virtual void   InitRecPoints()     {AliInfo("To be implemented by detectors");}
75   virtual void   InitRaws()          {AliInfo("To be implemented by detectors");}
76   virtual void   InitSDigits()       {AliInfo("To be implemented by detectors");}
77   //virtual void   InitTrackSegments() {AliInfo("To ne implemented by detectors");}
78   virtual void   MakeESDs(AliESDEvent * )          {AliInfo("To be implemented by detectors");} 
79   virtual void   MakeHits(TClonesArray * )         {AliInfo("To be implemented by detectors");} 
80   virtual void   MakeDigits(TClonesArray * )       {AliInfo("To be implemented by detectors");} 
81   //  virtual void   MakeRecParticles(TClonesArray * ) {AliInfo("To be implemented by detectors");} 
82   virtual void   MakeRaws(AliRawReader *)          {AliInfo("To be implemented by detectors");} 
83   virtual void   MakeRecPoints(TTree * )           {AliInfo("To be implemented by detectors");} 
84   virtual void   MakeSDigits(TClonesArray * )      {AliInfo("To be implemented by detectors");} 
85   //virtual void   MakeTrackSegments(TTree * )       {AliInfo("To be implemented by detectors");} 
86   void           ResetCycle() { fCurrentCycle++ ; fCycleCounter = 0 ; } 
87   virtual void   StartOfDetectorCycle() {AliInfo("To be implemented by detectors");} 
88
89   TFile *        fOutput ;          //! output root file
90   TDirectory *   fDetectorDir ;     //! directory for the given detector in the file
91   static TString fDetectorDirName ; //! detector directory name in the quality assurance data file
92   TList *        fDigitsQAList ;    //! list of the digits QA data objects
93   TList *        fESDsQAList ;      //! list of the ESDs QA data objects
94   TList *        fHitsQAList ;      //! list of the hits QA data objects
95   TList *        fRawsQAList ;      //! list of the raws QA data objects
96   TList *        fRecPointsQAList ; //! list of the recpoints QA data objects
97   TList *        fSDigitsQAList ;   //! list of the sdigits QA data objects
98   Int_t          fCurrentCycle ;    //! current cycle number
99   Int_t          fCycle ;           //! length (# events) of the QA data acquisition cycle  
100   Int_t          fCycleCounter ;    //! cycle counter
101   Int_t          fRun ;             //! run number
102   
103  ClassDef(AliQADataMaker,1)  // description 
104
105 };
106
107 #endif // AliQADataMaker_H