]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCclustererMI.h
CsI QE maps file in OCDB
[u/mrichter/AliRoot.git] / TPC / AliTPCclustererMI.h
1 #ifndef ALITPCCLUSTERERMI_H
2 #define ALITPCCLUSTERERMI_H
3
4 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5  * See cxx source for full Copyright notice                               */
6
7 /* $Id$ */
8
9
10
11 //-------------------------------------------------------
12 //                       TPC clusterer
13 //
14 //   Origin: Marian Ivanov  
15 //-------------------------------------------------------
16 #include <Rtypes.h>
17 #include <TObject.h>
18 #define kMAXCLUSTER 2500
19
20 class TFile;
21 class TClonesArray;
22 class AliTPCParam;
23 class AliTPCRecoParam;
24 class AliTPCclusterMI;
25 class AliTPCClustersRow;
26 class AliRawReader;
27 class AliSimDigits;
28 class TTree;
29 class TTreeSRedirector;
30 class  AliRawEventHeaderBase;
31 class AliTPCCalROC;
32
33 class AliTPCclustererMI : public TObject{
34 public:
35   AliTPCclustererMI(const AliTPCParam* par, const AliTPCRecoParam * recoParam = 0);
36   AliTPCclustererMI(const AliTPCclustererMI &param); // copy constructor
37   AliTPCclustererMI &operator = (const AliTPCclustererMI & param); //assignment
38   virtual ~AliTPCclustererMI();
39   virtual void Digits2Clusters();
40   virtual void Digits2ClustersOld(AliRawReader* rawReader);
41   virtual void Digits2Clusters(AliRawReader* rawReader);
42   virtual void SetInput(TTree * tree);  // set input tree with digits
43   virtual void SetOutput(TTree * tree); // set output tree with 
44   virtual void FillRow();               // fill the output container - Tree or TObjArray
45   TObjArray * GetOutputArray(){return fOutputArray;}
46   TClonesArray * GetOutputClonesArray(){return fOutputClonesArray;}
47
48   void StoreInClonesArray(Bool_t bOutput = kTRUE) {fBClonesArray = bOutput;} // store output clusters in TClonesArray
49
50 private:
51   Bool_t IsMaximum(Float_t k, Int_t max, const Float_t *bins) const; 
52   void MakeCluster2(Int_t k,Int_t max,Float_t *bins,UInt_t m,
53    AliTPCclusterMI &c);  
54   void MakeCluster(Int_t k,Int_t max,Float_t *bins,UInt_t m,
55    AliTPCclusterMI &c); 
56   Float_t  GetSigmaY2(Int_t iz);
57   Float_t  GetSigmaZ2(Int_t iz);
58   Float_t  FitMax(Float_t vmatrix[5][5], Float_t y, Float_t z, Float_t sigmay, Float_t sigmaz);
59   void AddCluster(AliTPCclusterMI &c, Float_t *matrix, Int_t pos);  // add the cluster to the array
60   void UnfoldCluster(Float_t * matrix[7], Float_t recmatrix[5][5], 
61                      Float_t & meani, Float_t & meanj, Float_t & sum, Float_t &overlap );
62   void FindClusters(AliTPCCalROC * noiseROC);
63   Bool_t AcceptCluster(AliTPCclusterMI*c);
64   Double_t  ProcesSignal(Float_t * signal, Int_t nchannels, Int_t id[3], Double_t &rms, Double_t &pedestalCalib);
65   void ProcessSectorData(Float_t** allBins, Int_t** allSigBins, Int_t*  allNSigBins);
66   
67   Float_t * fBins;       //!digits array
68   Int_t   * fSigBins; //!digits array containg only timebins above threshold
69   Int_t     fNSigBins;//!size of fSigBins
70   Int_t fLoop;         //loop - cf in 2 loops
71   Int_t fMaxBin;       //current ( for current sector)  maximal bin
72   Int_t fMaxTime;      //current ( for current sector)  maximal time
73   Int_t fMaxPad;       //current ( for current sector)  maximal pad
74   Int_t fSector;      //!current sector
75   Int_t fRow;         //!current row
76   Float_t fSign;      //!current sign 
77   Float_t fRx;        // current radius
78   Float_t fPadWidth;  // the width of the pad
79   Float_t fPadLength;  // the width of the pad
80   Float_t fZWidth;     //the z bin width
81   Bool_t  fPedSubtraction; // perform pedestal subtraction or not
82   AliRawEventHeaderBase *fEventHeader; //! event header information
83   UInt_t  fTimeStamp;   // Time Stamp
84   UInt_t  fEventType;   // Event Type
85   TTree * fInput;   //!input  tree with digits - object not owner
86   TTree * fOutput;   //!output tree with digits - object not owner
87   TObjArray *fOutputArray;     //! output TObjArray with pointers arrays of cluster
88   TClonesArray *fOutputClonesArray; //! output TClonesArray with clusters
89   AliTPCClustersRow * fRowCl;  //! current cluster row
90   AliSimDigits * fRowDig;      //! current digits row
91   const AliTPCParam * fParam;        //! tpc parameters
92   Int_t fNcluster;             // number of clusters - for given row
93   Int_t fNclusters;            // tot number of clusters
94   TTreeSRedirector *fDebugStreamer;     //!debug streamer
95   const AliTPCRecoParam  * fRecoParam;        //! reconstruction parameters
96   Bool_t  fBDumpSignal; // dump signal flag
97   Bool_t  fBClonesArray; // output clusters stored in TClonesArray 
98
99   ClassDef(AliTPCclustererMI,2)  // Time Projection Chamber digits
100 };
101
102 inline Bool_t AliTPCclustererMI::IsMaximum(Float_t q,Int_t max,const Float_t *bins) const {
103   //is this a local maximum ?
104   if (bins[-max] >= q) return kFALSE;
105   if (bins[-1  ] >= q) return kFALSE; 
106   if (bins[+max] > q) return kFALSE; 
107   if (bins[+1  ] > q) return kFALSE; 
108   if (bins[-max-1] >= q) return kFALSE;
109   if (bins[+max-1] >= q) return kFALSE; 
110   if (bins[+max+1] > q) return kFALSE; 
111   if (bins[-max+1] >= q) return kFALSE;
112   return kTRUE; 
113 }
114
115
116
117 //-----------------------------------------------------------------
118
119 #endif
120
121