]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/PHOSbase/AliPHOSCpvPedProducer.h
74c6ee28e1e46481788ef40fb8cf0480aba85b3c
[u/mrichter/AliRoot.git] / PHOS / PHOSbase / AliPHOSCpvPedProducer.h
1 #ifndef AliPHOSCpvPedProducer_h
2 #define AliPHOSCpvPedProducer_h
3
4 // Authors: Mikhail.Stolpovskiy@cern.ch, Sergey.Evdokimov@cern.ch
5 // The AliPHOSCpvPedProducer class calculates pedestals using AliPHOSCpvRawStream
6 // Also writes pedestals to files
7 // And creates a ROOT file with some histograms.
8 // this class supposed to be used in Cpv DA programm
9
10 #include <TH1.h>
11 #include <TH2.h>
12 #include <TMath.h>
13 #include <TFile.h>
14 #include <TString.h>
15 #include <THnSparse.h>
16 #include "AliPHOSCpvParam.h"
17 #include "AliPHOSCpvRawStream.h"
18 #include <TClonesArray.h>
19
20 class TFile;
21 class AliPHOSCpvPedProducer: public TObject {
22
23
24 public:
25   AliPHOSCpvPedProducer(Int_t sigcut = 3);
26   virtual ~AliPHOSCpvPedProducer();
27   void   SetSigCut(Int_t sigcut = 3) {fSigCut = sigcut;} //set n. of pedestal distribution sigmas used to create zero suppresion table
28   Bool_t LoadNewEvent(AliRawReader *& rawReader); // returns true, if ok
29   void   SetTurbo(Bool_t turbo);                  // if turbo==true then do read without error checking
30   Bool_t GetTurbo() const {return fTurbo;}
31
32   Bool_t FillPedestal(Int_t pad,Float_t q);   // pad - absolute pad number; q - charge of the pad
33   Bool_t FillPedestal();  // analyse event and fill pedestals
34   Bool_t CalcPedestal(Int_t iDDL);        //  analyse pedestals when all events processed for indicated DDL
35
36   TH2F* GetPedMeanMap(Int_t iDDL) const {return fPedMeanMap[iDDL];}  //Get the pedestal mean map for a DDL to send to AMORE
37   TH2F* GetPedSigMap(Int_t iDDL)  const {return fPedSigMap[iDDL];}   //Get the pedestal sigma map for a DDL to send to AMORE
38   TH1F* GetPedMean(Int_t iChFee)  const {return f1DPedMean[iChFee];} //Get the pedestal mean map for a FEE channel to send to AMORE
39   TH1F* GetPedSigma(Int_t iChFee) const {return f1DPedSigma[iChFee];}//Get the pedestal Sigma map for a FEE channel to send to AMORE
40   void  WritePedFiles(Int_t iDDL) const;                             // write pedestals to load(?) to RCB for indicated DDL
41   void WriteAllHistsToFile(const char* name) const;                  // create and write a new root file with hists
42
43   void SetErrorsHist(TH1I * pHist) {fhErrors = new TH1I(*pHist);}    //Set the histogram of errors, taken from AliPHOSCpvRawDigiProdicer
44   static Int_t GetMaxThr() {return fMaxThr;}                         // Get maximal threshold
45
46 protected:
47   void CreateErrHist();    // initialize histogram of errors
48   void CreateDDLHistos(Int_t iDDL);  // initialize histograms for pedestal calculation and representation
49
50   TH1F       *fPadAdc[AliPHOSCpvParam::kNDDL][AliPHOSCpvParam::kPadPcX][AliPHOSCpvParam::kPadPcY];        //Charge distribution for pads
51   Int_t       fSigCut;                         //n. of pedestal distribution sigmas used to create zero suppresion table
52   static const Int_t fMaxThr = 511;            //maximal threshold (9 bits all with 1)
53   Bool_t fTurbo;           // if true, then read without error checking
54
55   TH2F       *fPedMeanMap[AliPHOSCpvParam::kNDDL]; //2D mean pedestal map to export to AMORE
56   TH2F       *fPedSigMap [AliPHOSCpvParam::kNDDL]; //2D pedestal sigma map to export to AMORE
57   TH1F       *f1DPedMean [AliPHOSCpvParam::kNDDL]; //1D mean pedestal map to export to AMORE
58   TH1F       *f1DPedSigma[AliPHOSCpvParam::kNDDL]; //1D pedestal sigma map to export to AMORE
59   TH1I       *fhErrors;                        //histogram of errors from AliPHOSCpvRawDigiProducer
60   AliPHOSCpvRawStream         * fRawStream;       //! Raw data stream
61 private:
62   ClassDef(AliPHOSCpvPedProducer,1);                                                  //Cpv calibration and pedestal class
63 };
64 #endif