1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
18 //-------------------------------------------------------
19 // Implementation of the TPC clusterer
21 // Origin: Jouri Belikov, CERN, Jouri.Belikov@cern.ch
22 //-------------------------------------------------------
24 #include "AliTPCclusterer.h"
25 #include "AliTPCcluster.h"
26 #include <TObjArray.h>
28 #include "AliTPCClustersRow.h"
29 #include "AliDigits.h"
30 #include "AliSimDigits.h"
31 #include "AliTPCParam.h"
34 ClassImp(AliTPCclusterer)
36 AliTPCclusterer::AliTPCclusterer(const AliTPCParam *par) {
37 //-------------------------------------------------------
38 // The main constructor
39 //-------------------------------------------------------
40 fClusterArray.Setup(par);
41 fClusterArray.SetClusterType("AliTPCcluster");
44 void AliTPCclusterer::FindPeaks(Int_t k,Int_t max,
45 AliBin *b,Int_t *idx,UInt_t *msk,Int_t& n) {
48 if (IsMaximum(k,max,b)) {
49 idx[n]=k; msk[n]=(2<<n);
53 if (b[k-max].GetMask()&1) FindPeaks(k-max,max,b,idx,msk,n);
54 if (b[k-1 ].GetMask()&1) FindPeaks(k-1 ,max,b,idx,msk,n);
55 if (b[k+max].GetMask()&1) FindPeaks(k+max,max,b,idx,msk,n);
56 if (b[k+1 ].GetMask()&1) FindPeaks(k+1 ,max,b,idx,msk,n);
59 void AliTPCclusterer::MarkPeak(Int_t k, Int_t max, AliBin *bins, UInt_t m) {
61 UShort_t q=bins[k].GetQ();
63 bins[k].SetMask(bins[k].GetMask()|m);
65 if (bins[k-max].GetQ() <= q)
66 if ((bins[k-max].GetMask()&m) == 0) MarkPeak(k-max,max,bins,m);
67 if (bins[k-1 ].GetQ() <= q)
68 if ((bins[k-1 ].GetMask()&m) == 0) MarkPeak(k-1 ,max,bins,m);
69 if (bins[k+max].GetQ() <= q)
70 if ((bins[k+max].GetMask()&m) == 0) MarkPeak(k+max,max,bins,m);
71 if (bins[k+1 ].GetQ() <= q)
72 if ((bins[k+1 ].GetMask()&m) == 0) MarkPeak(k+1 ,max,bins,m);
75 void AliTPCclusterer::MakeCluster(Int_t k,Int_t max,AliBin *bins,UInt_t m,
77 //make cluster using digits of this peak
78 Float_t q=(Float_t)bins[k].GetQ();
79 Int_t i=k/max, j=k-i*max;
84 c.SetSigmaY2(c.GetSigmaY2()+i*i*q);
85 c.SetSigmaZ2(c.GetSigmaZ2()+j*j*q);
87 bins[k].SetMask(0xFFFFFFFE);
89 if (bins[k-max].GetMask() == m) MakeCluster(k-max,max,bins,m,c);
90 if (bins[k-1 ].GetMask() == m) MakeCluster(k-1 ,max,bins,m,c);
91 if (bins[k+max].GetMask() == m) MakeCluster(k+max,max,bins,m,c);
92 if (bins[k+1 ].GetMask() == m) MakeCluster(k+1 ,max,bins,m,c);
95 //_____________________________________________________________________________
96 Int_t AliTPCclusterer::Digits2Clusters(TTree *dTree, TTree *cTree) {
97 //-----------------------------------------------------------------
98 // This is a simple cluster finder.
99 //-----------------------------------------------------------------
100 TBranch *branch=dTree->GetBranch("Segment");
102 Error("Digits2Cluster","Can't get the branch !");
105 AliSimDigits digarr, *dummy=&digarr;
106 branch->SetAddress(&dummy);
108 fClusterArray.MakeTree(cTree);
110 AliTPCParam *par=(AliTPCParam *)fClusterArray.GetParam();
111 const Int_t kMAXZ=par->GetMaxTBin()+2;
115 Int_t nentries = (Int_t)dTree->GetEntries();
116 for (Int_t n=0; n<nentries; n++) {
121 if (!par->AdjustSectorRow(digarr.GetID(),sec,row)) {
122 Error("Digits2Clusters","!nvalid segment ID ! %d",digarr.GetID());
126 AliTPCClustersRow *clrow=fClusterArray.CreateRow(sec,row);
128 Float_t rx=par->GetPadRowRadii(sec,row);
132 const Int_t kNIS=par->GetNInnerSector(), kNOS=par->GetNOuterSector();
134 npads = par->GetNPadsLow(row);
135 sign = (sec < kNIS/2) ? 1 : -1;
137 npads = par->GetNPadsUp(row);
138 sign = ((sec-kNIS) < kNOS/2) ? 1 : -1;
142 const Int_t kMAXBIN=kMAXZ*(npads+2);
143 AliBin *bins=new AliBin[kMAXBIN];
144 for (Int_t ii=0;ii<kMAXBIN;ii++) {
145 bins[ii].SetQ(0); bins[ii].SetMask(0xFFFFFFFE);
150 Short_t dig=digarr.CurrentDigit();
151 if (dig<=par->GetZeroSup()) continue;
152 Int_t j=digarr.CurrentRow()+1, i=digarr.CurrentColumn()+1;
153 bins[i*kMAXZ+j].SetQ(dig);
154 bins[i*kMAXZ+j].SetMask(1);
155 } while (digarr.Next());
158 for (Int_t i=0; i<kMAXBIN; i++) {
159 if ((bins[i].GetMask()&1) == 0) continue;
160 Int_t idx[32]; UInt_t msk[32]; Int_t npeaks=0;
161 FindPeaks(i, kMAXZ, bins, idx, msk, npeaks);
163 if (npeaks>30) continue;
166 for (k=0; k<npeaks-1; k++){//mark adjacent peaks
167 if (idx[k] < 0) continue; //this peak is already removed
168 for (l=k+1; l<npeaks; l++) {
169 if (idx[l] < 0) continue; //this peak is already removed
170 Int_t ki=idx[k]/kMAXZ, kj=idx[k] - ki*kMAXZ;
171 Int_t li=idx[l]/kMAXZ, lj=idx[l] - li*kMAXZ;
172 Int_t di=TMath::Abs(ki - li);
173 Int_t dj=TMath::Abs(kj - lj);
174 if (di>1 || dj>1) continue;
175 if (bins[idx[k]].GetQ() > bins[idx[l]].GetQ()) {
186 for (k=0; k<npeaks; k++) {
187 MarkPeak(TMath::Abs(idx[k]), kMAXZ, bins, msk[k]);
190 for (k=0; k<npeaks; k++) {
191 if (idx[k] < 0) continue; //removed peak
193 MakeCluster(idx[k], kMAXZ, bins, msk[k], c);
194 if (c.GetQ() < 5) continue; //noise cluster
195 c.SetY(c.GetY()/c.GetQ());
196 c.SetZ(c.GetZ()/c.GetQ());
198 Float_t s2 = c.GetSigmaY2()/c.GetQ() - c.GetY()*c.GetY();
199 Float_t w=par->GetPadPitchWidth(sec);
200 c.SetSigmaY2((s2 + 1./12.)*w*w);
202 c.SetSigmaY2(c.GetSigmaY2()*0.108);
203 if (sec<par->GetNInnerSector()) c.SetSigmaY2(c.GetSigmaY2()*2.07);
206 s2 = c.GetSigmaZ2()/c.GetQ() - c.GetZ()*c.GetZ();
208 c.SetSigmaZ2((s2 + 1./12.)*w*w);
210 c.SetSigmaZ2(c.GetSigmaZ2()*0.169);
211 if (sec<par->GetNInnerSector()) c.SetSigmaZ2(c.GetSigmaZ2()*1.77);
214 c.SetY((c.GetY() - 0.5 - 0.5*npads)*par->GetPadPitchWidth(sec));
215 c.SetZ(par->GetZWidth()*(c.GetZ()-1));
216 c.SetZ(c.GetZ() - 3.*par->GetZSigma()); // PASA delay
217 c.SetZ(sign*(par->GetZLength() - c.GetZ()));
219 if (rx<230./250.*TMath::Abs(c.GetZ())) continue;
221 Int_t ki=idx[k]/kMAXZ, kj=idx[k] - ki*kMAXZ;
222 c.SetLabel(digarr.GetTrackID(kj-1,ki-1,0),0);
223 c.SetLabel(digarr.GetTrackID(kj-1,ki-1,1),1);
224 c.SetLabel(digarr.GetTrackID(kj-1,ki-1,2),2);
226 c.SetQ(bins[idx[k]].GetQ());
228 if (ki==1 || ki==npads || kj==1 || kj==kMAXZ-2) {
229 c.SetSigmaY2(c.GetSigmaY2()*25.);
230 c.SetSigmaZ2(c.GetSigmaZ2()*4.);
232 clrow->InsertCluster(&c); ncl++;
235 fClusterArray.StoreRow(sec,row);
236 fClusterArray.ClearRow(sec,row);
243 Info("Digits2Cluster","Number of found clusters : %d",nclusters);