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