ROOT logo
#ifndef ALITPCCLUSTERER_H
#define ALITPCCLUSTERER_H

/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
 * See cxx source for full Copyright notice                               */

/* $Id$ */



//-------------------------------------------------------
//                       TPC clusterer
//
//   Origin: Marian Ivanov  
//-------------------------------------------------------
#include <Rtypes.h>
#include <TObject.h>
#include <AliTPCRecoParam.h>
#define kMAXCLUSTER 2500

class TFile;
class TClonesArray;
class AliTPCParam;
class AliTPCRecoParam;
class AliTPCclusterMI;
class AliTPCClustersRow;
class AliRawReader;
class AliSimDigits;
class TTree;
class TTreeSRedirector;
class  AliRawEventHeaderBase;
class AliTPCCalROC;

class AliTPCclusterer : public TObject{
public:
  AliTPCclusterer(const AliTPCParam* par, const AliTPCRecoParam * recoParam = 0);
  virtual ~AliTPCclusterer();
  virtual void Digits2Clusters();
  virtual void Digits2Clusters(AliRawReader* rawReader);
  virtual void SetInput(TTree * tree);  // set input tree with digits
  virtual void SetOutput(TTree * tree); // set output tree with 
  virtual void FillRow();               // fill the output container - Tree or TObjArray
  TObjArray * GetOutputArray(){return fOutputArray;}
  TClonesArray * GetOutputClonesArray(){return fOutputClonesArray;}

  void StoreInClonesArray(Bool_t bOutput = kTRUE) {fBClonesArray = bOutput;} // store output clusters in TClonesArray

  void SetUseHLTClusters(Int_t useHLTClusters) {fUseHLTClusters = useHLTClusters;} // set usage from HLT clusters from rec.C options

private:
  AliTPCclusterer(const AliTPCclusterer &param); // copy constructor
  AliTPCclusterer &operator = (const AliTPCclusterer & param); //assignment

  Bool_t IsMaximum(Float_t k, Int_t max, const Float_t *bins) const; 
  void MakeCluster2(Int_t k,Int_t max,Float_t *bins,UInt_t m,
   AliTPCclusterMI &c);  
  void MakeCluster(Int_t k,Int_t max,Float_t *bins,UInt_t m,
   AliTPCclusterMI &c); 
  Float_t  GetSigmaY2(Int_t iz);
  Float_t  GetSigmaZ2(Int_t iz);
  Float_t  FitMax(Float_t vmatrix[5][5], Float_t y, Float_t z, Float_t sigmay, Float_t sigmaz);
  void AddCluster(AliTPCclusterMI &c, Float_t *matrix, Int_t pos);  // add the cluster to the array
  void UnfoldCluster(Float_t * matrix[7], Float_t recmatrix[5][5], 
		     Float_t & meani, Float_t & meanj, Float_t & sum, Float_t &overlap );
  void FindClusters(AliTPCCalROC * noiseROC);
  Bool_t AcceptCluster(AliTPCclusterMI*c);
  Double_t  ProcesSignal(Float_t * signal, Int_t nchannels, Int_t id[3], Double_t &rms, Double_t &pedestalCalib);
  void ProcessSectorData();
  Int_t ReadHLTClusters();
  
  Float_t * fBins;       //!digits array
  Int_t   * fSigBins; //!digits array containg only timebins above threshold
  Int_t     fNSigBins;//!size of fSigBins
  Int_t fLoop;         //loop - cf in 2 loops
  Int_t fMaxBin;       //current ( for current sector)  maximal bin
  Int_t fMaxTime;      //current ( for current sector)  maximal time
  Int_t fMaxPad;       //current ( for current sector)  maximal pad
  Int_t fSector;      //!current sector
  Int_t fRow;         //!current row
  Float_t fSign;      //!current sign 
  Float_t fRx;        // current radius
  Float_t fPadWidth;  // the width of the pad
  Float_t fPadLength;  // the width of the pad
  Float_t fZWidth;     //the z bin width
  Bool_t  fPedSubtraction; // perform pedestal subtraction or not
  AliRawEventHeaderBase *fEventHeader; //! event header information
  UInt_t  fTimeStamp;   // Time Stamp
  UInt_t  fEventType;   // Event Type
  TTree * fInput;   //!input  tree with digits - object not owner
  TTree * fOutput;   //!output tree with digits - object not owner
  TObjArray *fOutputArray;     //! output TObjArray with pointers arrays of cluster
  TClonesArray *fOutputClonesArray; //! output TClonesArray with clusters
  AliTPCClustersRow * fRowCl;  //! current cluster row
  AliSimDigits * fRowDig;      //! current digits row
  const AliTPCParam * fParam;        //! tpc parameters
  Int_t fNcluster;             // number of clusters - for given row
  Int_t fNclusters;            // tot number of clusters
  TTreeSRedirector *fDebugStreamer;     //!debug streamer
  const AliTPCRecoParam  * fRecoParam;        //! reconstruction parameters
  Bool_t  fBDumpSignal; // dump signal flag
  Bool_t  fBClonesArray; // output clusters stored in TClonesArray 
  Int_t  fUseHLTClusters; // use HLT clusters instead of offline clusters

  // Non-persistent arrays

  Float_t** fAllBins; //! All sector bins
  Int_t** fAllSigBins;//! All signal bins in a sector
  Int_t*  fAllNSigBins;//! Number of signal bins in a sector
  TObject* fHLTClusterAccess;// interface to HLT clusters

  ClassDef(AliTPCclusterer,0)  // TPC cluster finder
};

inline Bool_t AliTPCclusterer::IsMaximum(Float_t q,Int_t max,const Float_t *bins) const {
  //is this a local maximum ?
  if (bins[-max] >= q) return kFALSE;
  if (bins[-1  ] >= q) return kFALSE; 
  if (bins[+max] > q) return kFALSE; 
  if (bins[+1  ] > q) return kFALSE; 
  if (bins[-max-1] >= q) return kFALSE;
  if (bins[+max-1] >= q) return kFALSE; 
  if (bins[+max+1] > q) return kFALSE; 
  if (bins[-max+1] >= q) return kFALSE;
  //
  //
  if (fRecoParam->GetClusterMaxRange(1)>0){  //local maxima in z on more than 1 time bin 
    if (bins[-2] > q) return kFALSE;
    if (bins[ 2] > q) return kFALSE;
  }
  if (fRecoParam->GetClusterMaxRange(0)>0){  //local maxima in y on more than pad 
    if (bins[-2*max] > q) return kFALSE;
    if (bins[ 2*max] > q) return kFALSE;
  }

  return kTRUE; 
}



//-----------------------------------------------------------------

#endif


 AliTPCclusterer.h:1
 AliTPCclusterer.h:2
 AliTPCclusterer.h:3
 AliTPCclusterer.h:4
 AliTPCclusterer.h:5
 AliTPCclusterer.h:6
 AliTPCclusterer.h:7
 AliTPCclusterer.h:8
 AliTPCclusterer.h:9
 AliTPCclusterer.h:10
 AliTPCclusterer.h:11
 AliTPCclusterer.h:12
 AliTPCclusterer.h:13
 AliTPCclusterer.h:14
 AliTPCclusterer.h:15
 AliTPCclusterer.h:16
 AliTPCclusterer.h:17
 AliTPCclusterer.h:18
 AliTPCclusterer.h:19
 AliTPCclusterer.h:20
 AliTPCclusterer.h:21
 AliTPCclusterer.h:22
 AliTPCclusterer.h:23
 AliTPCclusterer.h:24
 AliTPCclusterer.h:25
 AliTPCclusterer.h:26
 AliTPCclusterer.h:27
 AliTPCclusterer.h:28
 AliTPCclusterer.h:29
 AliTPCclusterer.h:30
 AliTPCclusterer.h:31
 AliTPCclusterer.h:32
 AliTPCclusterer.h:33
 AliTPCclusterer.h:34
 AliTPCclusterer.h:35
 AliTPCclusterer.h:36
 AliTPCclusterer.h:37
 AliTPCclusterer.h:38
 AliTPCclusterer.h:39
 AliTPCclusterer.h:40
 AliTPCclusterer.h:41
 AliTPCclusterer.h:42
 AliTPCclusterer.h:43
 AliTPCclusterer.h:44
 AliTPCclusterer.h:45
 AliTPCclusterer.h:46
 AliTPCclusterer.h:47
 AliTPCclusterer.h:48
 AliTPCclusterer.h:49
 AliTPCclusterer.h:50
 AliTPCclusterer.h:51
 AliTPCclusterer.h:52
 AliTPCclusterer.h:53
 AliTPCclusterer.h:54
 AliTPCclusterer.h:55
 AliTPCclusterer.h:56
 AliTPCclusterer.h:57
 AliTPCclusterer.h:58
 AliTPCclusterer.h:59
 AliTPCclusterer.h:60
 AliTPCclusterer.h:61
 AliTPCclusterer.h:62
 AliTPCclusterer.h:63
 AliTPCclusterer.h:64
 AliTPCclusterer.h:65
 AliTPCclusterer.h:66
 AliTPCclusterer.h:67
 AliTPCclusterer.h:68
 AliTPCclusterer.h:69
 AliTPCclusterer.h:70
 AliTPCclusterer.h:71
 AliTPCclusterer.h:72
 AliTPCclusterer.h:73
 AliTPCclusterer.h:74
 AliTPCclusterer.h:75
 AliTPCclusterer.h:76
 AliTPCclusterer.h:77
 AliTPCclusterer.h:78
 AliTPCclusterer.h:79
 AliTPCclusterer.h:80
 AliTPCclusterer.h:81
 AliTPCclusterer.h:82
 AliTPCclusterer.h:83
 AliTPCclusterer.h:84
 AliTPCclusterer.h:85
 AliTPCclusterer.h:86
 AliTPCclusterer.h:87
 AliTPCclusterer.h:88
 AliTPCclusterer.h:89
 AliTPCclusterer.h:90
 AliTPCclusterer.h:91
 AliTPCclusterer.h:92
 AliTPCclusterer.h:93
 AliTPCclusterer.h:94
 AliTPCclusterer.h:95
 AliTPCclusterer.h:96
 AliTPCclusterer.h:97
 AliTPCclusterer.h:98
 AliTPCclusterer.h:99
 AliTPCclusterer.h:100
 AliTPCclusterer.h:101
 AliTPCclusterer.h:102
 AliTPCclusterer.h:103
 AliTPCclusterer.h:104
 AliTPCclusterer.h:105
 AliTPCclusterer.h:106
 AliTPCclusterer.h:107
 AliTPCclusterer.h:108
 AliTPCclusterer.h:109
 AliTPCclusterer.h:110
 AliTPCclusterer.h:111
 AliTPCclusterer.h:112
 AliTPCclusterer.h:113
 AliTPCclusterer.h:114
 AliTPCclusterer.h:115
 AliTPCclusterer.h:116
 AliTPCclusterer.h:117
 AliTPCclusterer.h:118
 AliTPCclusterer.h:119
 AliTPCclusterer.h:120
 AliTPCclusterer.h:121
 AliTPCclusterer.h:122
 AliTPCclusterer.h:123
 AliTPCclusterer.h:124
 AliTPCclusterer.h:125
 AliTPCclusterer.h:126
 AliTPCclusterer.h:127
 AliTPCclusterer.h:128
 AliTPCclusterer.h:129
 AliTPCclusterer.h:130
 AliTPCclusterer.h:131
 AliTPCclusterer.h:132
 AliTPCclusterer.h:133
 AliTPCclusterer.h:134
 AliTPCclusterer.h:135
 AliTPCclusterer.h:136
 AliTPCclusterer.h:137
 AliTPCclusterer.h:138
 AliTPCclusterer.h:139
 AliTPCclusterer.h:140
 AliTPCclusterer.h:141
 AliTPCclusterer.h:142
 AliTPCclusterer.h:143
 AliTPCclusterer.h:144