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