]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/AliZDCReconstructor.h
Now, the estimator can work with the tracklets too (R.Shahoyan)
[u/mrichter/AliRoot.git] / ZDC / AliZDCReconstructor.h
1 #ifndef ALIZDCRECONSTRUCTOR_H
2 #define ALIZDCRECONSTRUCTOR_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  * See cxx source for full Copyright notice                               */
5
6 /* $Id$ */
7
8 /////////////////////////////////////////////
9 //                                         //
10 //       class for ZDC reconstruction      //
11 //                                         //
12 /////////////////////////////////////////////
13
14 #include "AliReconstructor.h"
15 #include "AliCDBManager.h"
16 #include "AliCDBStorage.h"
17 #include "AliZDCPedestals.h"
18 #include "AliZDCEnCalib.h"
19 #include "AliZDCTowerCalib.h"
20 #include "AliZDCMBCalib.h"
21 #include "AliZDCRecoParam.h"
22 #include "AliZDCRecoParampp.h"
23 #include "AliZDCRecoParamPbPb.h"
24 #include "AliLog.h"
25
26 class AliLoader;
27
28 class AliZDCReconstructor: public AliReconstructor {
29 public:
30   AliZDCReconstructor();
31   virtual ~AliZDCReconstructor();
32
33   virtual void   Init();
34   virtual Bool_t HasDigitConversion() const {return kFALSE;};
35   
36   virtual void Reconstruct(TTree* digitsTree, TTree* clustersTree) const; 
37   virtual void Reconstruct(AliRawReader* rawReader, TTree* clustersTree) const;
38
39   virtual void FillESD(TTree* /*digitsTree*/, TTree* clustersTree, AliESDEvent* esd) const 
40                 {FillZDCintoESD(clustersTree, esd);}
41   virtual void FillESD(AliRawReader* /*rawReader*/, TTree* clustersTree, AliESDEvent* esd) const 
42                 {FillZDCintoESD(clustersTree, esd);}
43   
44   // parameter settings for reconstruction
45   void SetRecoMode(Int_t recoMode, Float_t beamEnergy) 
46                   {fRecoMode=recoMode; fBeamEnergy=beamEnergy;}
47   static void SetRecoParam(AliZDCRecoParam * param){fRecoParam = param;}
48   
49   Int_t   GetRecoMode() {return fRecoMode;}
50   Float_t GetBeamEnergy() {return fBeamEnergy;}
51   
52   static const AliZDCRecoParam* GetRecoParam() {return dynamic_cast<const AliZDCRecoParam*>(AliReconstructor::GetRecoParam(9));}
53     
54   void  SetPedSubMode(Int_t pedsubMode) {fPedSubMode=pedsubMode;}
55   Int_t GetPedSubMode() {return fPedSubMode;}
56   
57   void    SetSignalThreshold(Float_t val) {fSignalThreshold=val;}
58   Float_t GetSignalThreshold() {return fSignalThreshold;}
59   
60   // OCDB objects for reconstruction
61   AliCDBStorage       *SetStorage(const char* uri);
62   AliZDCPedestals     *GetPedestalData() const; 
63   AliZDCEnCalib       *GetEnergyCalibData() const; 
64   AliZDCTowerCalib    *GetTowerCalibData() const; 
65   AliZDCMBCalib       *GetMBCalibData() const; 
66   AliZDCRecoParampp   *GetppRecoParamFromOCDB() const;  
67   AliZDCRecoParamPbPb *GetPbPbRecoParamFromOCDB() const;  
68   
69 private:
70   AliZDCReconstructor(const AliZDCReconstructor&); //Not implemented
71   AliZDCReconstructor& operator =(const AliZDCReconstructor&); //Not implemented
72
73   void   ReconstructEventpp(TTree *clustersTree, 
74             Float_t* ZN1ADCCorr, Float_t* ZP1ADCCorr, Float_t* ZN2ADCCorr, Float_t* ZP2ADCCorr,
75             Float_t* ZEM1ADCCorr, Float_t* ZEM2ADCCorr, Float_t* PMRef1, Float_t* PMRef2,
76             Bool_t isScalerOn, UInt_t* scaler, 
77             Int_t* evQualityBlock, Int_t* triggerBlock, Int_t* chBlock, UInt_t puBits) const;
78   void   ReconstructEventPbPb(TTree *clustersTree, 
79             Float_t* ZN1ADCCorr, Float_t* ZP1ADCCorr, Float_t* ZN2ADCCorr, Float_t* ZP2ADCCorr,
80             Float_t* ZEM1ADCCorr, Float_t* ZEM2ADCCorr, Float_t* PMRef1, Float_t* PMRef2,
81             Bool_t isScalerOn, UInt_t* scaler, 
82             Int_t* evQualityBlock, Int_t* triggerBlock, Int_t* chBlock, UInt_t puBits) const;
83   
84   void   FillZDCintoESD(TTree *clustersTree, AliESDEvent*esd) const;
85
86   static AliZDCRecoParam *fRecoParam; // reconstruction parameters
87
88   AliZDCPedestals  *fPedData;       //! pedestal calibration data
89   AliZDCEnCalib    *fEnCalibData;   //! energy calibration data
90   AliZDCTowerCalib *fTowCalibData;  //! equalization calibration data
91   AliZDCMBCalib    *fMBCalibData;   //! mb calibration data
92   
93   Int_t   fRecoMode;        // =1->p-p, =2->A-A
94   Float_t fBeamEnergy;      // beam energy
95   Int_t   fNRun;            // Run Number (from raw data)
96   Bool_t  fIsCalibrationMB; // true if run type = "CALIBRATION_MB"
97   Int_t   fPedSubMode;      // =0->mean values, =1->from correlations
98   Float_t fSignalThreshold; // Threshold value for "triggering" in p-p
99
100   ClassDef(AliZDCReconstructor, 10)   // class for the ZDC reconstruction
101 };
102
103 #endif