]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliQADataMaker.h
New revision of the CTP configuration and simulation. For more details look in the...
[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   const Int_t         Add2DigitsList(TH1 * hist, const Int_t index)    { return Add2List(hist, index, fDigitsQAList) ; }
44   const Int_t         Add2ESDsList(TH1 * hist, const Int_t index)      { return Add2List(hist, index, fESDsQAList) ; }
45   const Int_t         Add2HitsList(TH1 * hist, const Int_t index)      { return Add2List(hist, index, fHitsQAList) ; }
46   const Int_t         Add2RecPointsList(TH1 * hist, const Int_t index) { return Add2List(hist, index, fRecPointsQAList) ; }
47   const Int_t         Add2RawsList(TH1 * hist, const Int_t index)      { return Add2List(hist, index, fRawsQAList) ; }
48   const Int_t         Add2SDigitsList(TH1 * hist, const Int_t index)   { return Add2List(hist, index, fSDigitsQAList) ; }
49   virtual void        Exec(AliQA::TASKINDEX, TObject * data) ;
50   void                EndOfCycle(AliQA::TASKINDEX) ;
51   void                Finish() const ; 
52   TH1 *               GetDigitsData(const Int_t index)    { return dynamic_cast<TH1 *>(GetData(fDigitsQAList, index)) ; }
53   TH1 *               GetESDsData(const Int_t index)      { return dynamic_cast<TH1 *>(GetData(fESDsQAList, index)) ; }
54   TH1 *               GetHitsData(const Int_t index)      { return dynamic_cast<TH1 *>(GetData(fHitsQAList, index)) ; }
55   TH1 *               GetRecPointsData(const Int_t index) { return dynamic_cast<TH1 *>(GetData(fRecPointsQAList, index)) ; }
56   TH1 *               GetRawsData(const Int_t index)      { return dynamic_cast<TH1 *>(GetData(fRawsQAList, index)) ; }
57   TH1 *               GetSDigitsData(const Int_t index)   { return dynamic_cast<TH1 *>(GetData(fSDigitsQAList, index)) ; }
58   const char *        GetDetectorDirName() { return fDetectorDirName.Data() ; }
59   const Int_t         Increment() { return ++fCycleCounter ; } 
60   TList *             Init(AliQA::TASKINDEX, Int_t run, Int_t cycles = -1) ;
61   void                Init(AliQA::TASKINDEX, TList * list, Int_t run, Int_t cycles = -1) ;
62   const Bool_t        IsCycleDone() const { return fCycleCounter > fCycle ? kTRUE : kFALSE ; }
63   void                Reset() ;         
64   void                SetCycle(Int_t nevts) { fCycle = nevts ; } 
65   void                StartOfCycle(AliQA::TASKINDEX, const Bool_t sameCycle = kFALSE) ;
66
67 protected: 
68
69   Int_t          Add2List(TH1 * hist, const Int_t index, TList * list) { list->AddAt(hist, index) ; return list->LastIndex() ; }
70   virtual void   EndOfDetectorCycle(AliQA::TASKINDEX, TList * ) {AliInfo("To be implemented by detectors");} 
71   TObject *      GetData(TList * list, const Int_t index)  { return list->At(index) ; } 
72   virtual void   InitDigits()        {AliInfo("To be implemented by detectors");}
73   virtual void   InitESDs()          {AliInfo("To be implemented by detectors");}
74   virtual void   InitHits()          {AliInfo("To be implemented by detectors");}
75   //virtual void   InitRecParticles()  {AliInfo("To be implemented by detectors");}
76   virtual void   InitRecPoints()     {AliInfo("To be implemented by detectors");}
77   virtual void   InitRaws()          {AliInfo("To be implemented by detectors");}
78   virtual void   InitSDigits()       {AliInfo("To be implemented by detectors");}
79   //virtual void   InitTrackSegments() {AliInfo("To ne implemented by detectors");}
80   virtual void   MakeESDs(AliESDEvent * )          {AliInfo("To be implemented by detectors");} 
81   virtual void   MakeHits(TClonesArray * )         {AliInfo("To be implemented by detectors");} 
82   virtual void   MakeHits(TTree * )                {AliInfo("To be implemented by detectors");} 
83   virtual void   MakeDigits(TClonesArray * )       {AliInfo("To be implemented by detectors");} 
84   virtual void   MakeDigits(TTree * )              {AliInfo("To be implemented by detectors");} 
85   //  virtual void   MakeRecParticles(TClonesArray * ) {AliInfo("To be implemented by detectors");} 
86   virtual void   MakeRaws(AliRawReader *)          {AliInfo("To be implemented by detectors");} 
87   virtual void   MakeRecPoints(TTree * )           {AliInfo("To be implemented by detectors");} 
88   virtual void   MakeSDigits(TClonesArray * )      {AliInfo("To be implemented by detectors");} 
89   virtual void   MakeSDigits(TTree * )             {AliInfo("To be implemented by detectors");} 
90   //virtual void   MakeTrackSegments(TTree * )       {AliInfo("To be implemented by detectors");} 
91   void           ResetCycle() { fCurrentCycle++ ; fCycleCounter = 0 ; } 
92   virtual void   StartOfDetectorCycle() {AliInfo("To be implemented by detectors");} 
93
94   TFile *        fOutput ;          //! output root file
95   TDirectory *   fDetectorDir ;     //! directory for the given detector in the file
96   TString        fDetectorDirName ; //! detector directory name in the quality assurance data file
97   TList *        fDigitsQAList ;    //! list of the digits QA data objects
98   TList *        fESDsQAList ;      //! list of the ESDs QA data objects
99   TList *        fHitsQAList ;      //! list of the hits QA data objects
100   TList *        fRawsQAList ;      //! list of the raws QA data objects
101   TList *        fRecPointsQAList ; //! list of the recpoints QA data objects
102   TList *        fSDigitsQAList ;   //! list of the sdigits QA data objects
103   Int_t          fCurrentCycle ;    //! current cycle number
104   Int_t          fCycle ;           //! length (# events) of the QA data acquisition cycle  
105   Int_t          fCycleCounter ;    //! cycle counter
106   Int_t          fRun ;             //! run number
107   
108  ClassDef(AliQADataMaker,1)  // description 
109
110 };
111
112 #endif // AliQADataMaker_H