]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGGA/PHOSTasks/CaloCellQA/AliCaloCellsQA.h
Time histo range adjusted
[u/mrichter/AliRoot.git] / PWGGA / PHOSTasks / CaloCellQA / AliCaloCellsQA.h
1 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
2  * See cxx source for full Copyright notice                               */
3
4 //_________________________________________________________________________
5 // Class for bad channels & bad runs identification
6 // Author: Olga Driga (SUBATECH)
7
8 #ifndef ALICALOCELLSQA_H
9 #define ALICALOCELLSQA_H
10
11 // --- ROOT system ---
12 #include <TObjArray.h>
13 #include <TH1D.h>
14 #include <TH1F.h>
15 #include <TH2F.h>
16 class AliLog;
17
18 // --- AliRoot header files ---
19 #include <AliVCaloCells.h>
20 #include <AliVCluster.h>
21
22 class AliCaloCellsQA : public TObject {
23
24 public:
25
26   // detectors
27   enum {
28     kEMCAL = 0,
29     kPHOS  = 1
30 // ,kDCAL  = 2      // not yet implemented
31   };
32
33   AliCaloCellsQA();
34   AliCaloCellsQA(Int_t nmods, Int_t det = kEMCAL, Int_t startRunNumber = 100000, Int_t endRunNumber = 200000);
35
36   virtual ~AliCaloCellsQA();
37
38   virtual void ActivateFullAnalysis();
39   virtual void InitSummaryHistograms(Int_t nbins = 400, Double_t emax = 4.,
40                                      Int_t nbinsh = 100, Double_t emaxh = 300.,
41                                      Int_t nbinst = 250, Double_t tmin =-0.1e-6, Double_t tmax = 0.15e-6);
42   virtual void InitTransientFindCurrentRun(Int_t runNumber);
43   virtual void Fill(Int_t runNumber, TObjArray *clusArray, AliVCaloCells *cells, Double_t vertexXYZ[3]);  // main method
44
45   // getters
46   virtual Int_t      GetDetector()     { return fDetector; }
47   virtual Int_t      GetNMods()        { return fNMods; }
48   virtual Double_t   GetClusElowMin()  { return fClusElowMin; }
49   virtual Double_t   GetClusEhighMin() { return fClusEhighMin; }
50   virtual Double_t   GetPi0EClusMin()  { return fPi0EClusMin; }
51   virtual Bool_t     GetFullAnalysis() { return fkFullAnalysis; }
52
53   virtual TObjArray* GetListOfHistos() { return fListOfHistos; }
54
55   // setters
56   virtual void SetClusterEnergyCuts(Double_t pi0EClusMin = 0.5, Double_t ElowMin = 0.3, Double_t EhighMin = 1.0);
57   virtual void SetBinningParameters(Int_t nbins1 = 100, Int_t nbins2 = 250, Int_t nbins3x = 200, Int_t nbins3y = 30,
58                                     Double_t xmax1 = 5., Double_t xmax2 = 0.5, Double_t xmax3 = 20.);
59
60 protected:
61
62   virtual void    Init(Int_t nmods, Int_t det, Int_t startRunNumber, Int_t endRunNumber);
63   virtual void    InitTransientMembers(Int_t run);
64   virtual void    InitHistosForRun(Int_t run);
65
66   virtual void    FillCellsInCluster(TObjArray *clusArray, AliVCaloCells *cells);
67   virtual void    FillJustCells(AliVCaloCells *cells);
68   virtual void    FillPi0Mass(TObjArray *clusArray, Double_t vertexXYZ[3]);
69
70   virtual Int_t   CheckClusterGetSM(AliVCluster* clus);
71   virtual Int_t   GetSM(Int_t absId);
72   virtual Bool_t  IsCellLocalMaximum(Int_t c, AliVCluster* clus, AliVCaloCells* cells);
73   virtual Bool_t  IsCellLocalMaximum(Int_t absId, AliVCaloCells* cells);
74   virtual void    AbsIdToSMEtaPhi(Int_t absId, Int_t &sm, Int_t &eta, Int_t &phi);
75
76 private:
77
78   AliCaloCellsQA(const AliCaloCellsQA &);
79   AliCaloCellsQA & operator = (const AliCaloCellsQA &);
80
81 private:
82
83   // changeable parameters
84   Int_t     fDetector;                  // kEMCAL or kPHOS
85   Int_t     fNMods;                     // maximum supermodule number + 1 (4 or 10)
86   Double_t  fClusElowMin;               // minimum cluster energy cut for low energy per run histograms
87   Double_t  fClusEhighMin;              // minimum cluster energy cut for high energy per run histograms
88   Double_t  fPi0EClusMin;               // minimum cluster energy cut for pi0 mass histograms
89   Bool_t    fkFullAnalysis;             // flag to activate all available histograms
90
91   // changeable binning parameters
92   Int_t     fNBinsECells;               // number of bins in hECells
93   Int_t     fNBinsPi0Mass;              // number of bins in hPi0Mass
94   Int_t     fNBinsXNCellsInCluster;     // number of bins in hNCellsInCluster, X axis
95   Int_t     fNBinsYNCellsInCluster;     // number of bins in hNCellsInCluster, Y axis
96   Double_t  fXMaxECells;                // X axis maximum in hECells
97   Double_t  fXMaxPi0Mass;               // X axis maximum in hPi0Mass
98   Double_t  fXMaxNCellsInCluster;       // X axis maximum in hNCellsInCluster
99
100   // internal parameters
101   Int_t     fRunNumbers[1000];          // already encountered runs
102   Int_t     fNRuns;                     // number of encountered runs
103   Int_t     fRI;                        //! current (cached) run index
104
105   // internal parameters, used for coding convenience
106   Int_t     fAbsIdMin;                  // minimum absId number (0/EMCAL, 1/PHOS)
107   Int_t     fAbsIdMax;                  // maximum absId number + 1
108
109   TObjArray  *fListOfHistos;            // array with all the histograms
110
111
112   /* All the histograms below are present in fListOfHistos.
113    *
114    * NOTE: strictly speaking, the transient members below are not necessary, but the
115    * caching mechanism used in this class boosts the analysis speed considerably
116    * (no searching for a histogram by its name is performed in each event).
117    */
118
119   TH1D *fhNEventsProcessedPerRun;             //! number of processed events per run
120
121   // current run histograms; X axis -- cell absId number;
122   TH1F *fhCellLocMaxNTimesInClusterElow;      //! number of times cell was local maximum in a low energy cluster
123   TH1F *fhCellLocMaxNTimesInClusterEhigh;     //! number of times cell was local maximum in a high energy cluster
124   TH1F *fhCellLocMaxETotalClusterElow;        //! total cluster energy for local maximum cell, low energy
125   TH1F *fhCellLocMaxETotalClusterEhigh;       //! total cluster energy for local maximum cell, high energy
126   TH1F *fhCellNonLocMaxNTimesInClusterElow;   //! number of times cell wasn't local maximum in a low energy cluster
127   TH1F *fhCellNonLocMaxNTimesInClusterEhigh;  //! number of times cell wasn't local maximum in a high energy cluster
128   TH1F *fhCellNonLocMaxETotalClusterElow;     //! total cluster energy for not local maximum cell, low energy
129   TH1F *fhCellNonLocMaxETotalClusterEhigh;    //! total cluster energy for not local maximum cell, high energy
130
131   // current run, per supermodule histograms; the maximum number of supermodules is 10
132   TH1F *fhECells[10];                         //! cell amplitude distribution
133   TH1F *fhPi0Mass[10][10];                    //! pi0 mass spectrum
134   TH2F *fhNCellsInCluster[10];                //! distribution of number of cells in cluster vs cluster energy
135
136   // summary histograms: cells spectra at different conditions; X axis -- cell absId number
137   TH2F *fhCellAmplitude;                      //! amplitude distribution per cell
138   TH2F *fhCellAmplitudeEhigh;                 //! amplitude distribution per cell, high energies
139   TH2F *fhCellAmplitudeNonLocMax;             //! amplitude distribution per not a local maximum cell
140   TH2F *fhCellAmplitudeEhighNonLocMax;        //! amplitude distribution per not a local maximum cell, high energies
141   TH2F *fhCellTime;                           //! time distribution per cell
142
143   ClassDef(AliCaloCellsQA,2)
144 };
145
146 #endif