]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCclusterer.cxx
cluster information
[u/mrichter/AliRoot.git] / TPC / AliTPCclusterer.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16 /* $Id$ */
17
18 //-------------------------------------------------------
19 //          Implementation of the TPC clusterer
20 //
21 //   Origin: Jouri Belikov, CERN, Jouri.Belikov@cern.ch 
22 //-------------------------------------------------------
23
24 #include "AliTPCclusterer.h"
25 #include "AliTPCcluster.h"
26 #include <TObjArray.h>
27 #include <TFile.h>
28 #include "AliTPCClustersRow.h"
29 #include "AliDigits.h"
30 #include "AliSimDigits.h"
31 #include "AliTPCParam.h"
32 #include <TTree.h>
33
34 ClassImp(AliTPCclusterer)
35
36 AliTPCclusterer::AliTPCclusterer(const AliTPCParam *par) { 
37 //-------------------------------------------------------
38 //  The main constructor
39 //-------------------------------------------------------
40   fClusterArray.Setup(par);
41   fClusterArray.SetClusterType("AliTPCcluster");
42 }
43
44 void AliTPCclusterer::FindPeaks(Int_t k,Int_t max,
45 AliBin *b,Int_t *idx,UInt_t *msk,Int_t& n) {
46   //find local maxima
47   if (n<31)
48   if (IsMaximum(k,max,b)) {
49     idx[n]=k; msk[n]=(2<<n);
50     n++;
51   }
52   b[k].SetMask(0);
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);
57 }
58
59 void AliTPCclusterer::MarkPeak(Int_t k, Int_t max, AliBin *bins, UInt_t m) {
60   //mark this peak
61   UShort_t q=bins[k].GetQ();
62
63   bins[k].SetMask(bins[k].GetMask()|m); 
64
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);
73 }
74
75 void AliTPCclusterer::MakeCluster(Int_t k,Int_t max,AliBin *bins,UInt_t m,
76 AliTPCcluster &c) {
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;
80
81   c.SetQ(c.GetQ()+q);
82   c.SetY(c.GetY()+i*q); 
83   c.SetZ(c.GetZ()+j*q); 
84   c.SetSigmaY2(c.GetSigmaY2()+i*i*q);
85   c.SetSigmaZ2(c.GetSigmaZ2()+j*j*q);
86
87   bins[k].SetMask(0xFFFFFFFE);
88   
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);
93 }
94
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");
101   if (!branch) {
102      Error("Digits2Cluster","Can't get the branch !");
103      return 1;
104   }
105   AliSimDigits digarr, *dummy=&digarr;
106   branch->SetAddress(&dummy);
107   
108   fClusterArray.MakeTree(cTree);
109
110   AliTPCParam *par=(AliTPCParam *)fClusterArray.GetParam();
111   const Int_t kMAXZ=par->GetMaxTBin()+2;
112
113   Int_t nclusters=0;
114
115   Int_t nentries = (Int_t)dTree->GetEntries();
116   for (Int_t n=0; n<nentries; n++) {
117    
118     Int_t sec, row;
119     dTree->GetEvent(n);
120
121     if (!par->AdjustSectorRow(digarr.GetID(),sec,row)) {
122        Error("Digits2Clusters","!nvalid segment ID ! %d",digarr.GetID());
123        continue;
124     }
125
126     AliTPCClustersRow *clrow=fClusterArray.CreateRow(sec,row);
127
128     Float_t rx=par->GetPadRowRadii(sec,row);
129
130     Int_t npads, sign;
131     {
132        const Int_t kNIS=par->GetNInnerSector(), kNOS=par->GetNOuterSector();
133        if (sec < kNIS) {
134           npads = par->GetNPadsLow(row);
135           sign = (sec < kNIS/2) ? 1 : -1;
136        } else {
137           npads = par->GetNPadsUp(row);
138           sign = ((sec-kNIS) < kNOS/2) ? 1 : -1;
139        }
140     }
141
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);
146     }
147
148     digarr.First();
149     do {
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());
156
157     Int_t ncl=0;
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);
162
163       if (npeaks>30) continue;
164
165       Int_t k,l;
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()) {
176               msk[l]=msk[k];
177               idx[l]*=-1;
178            } else {
179               msk[k]=msk[l];
180               idx[k]*=-1;
181               break;
182            } 
183         }
184       }
185
186       for (k=0; k<npeaks; k++) {
187         MarkPeak(TMath::Abs(idx[k]), kMAXZ, bins, msk[k]);
188       }
189         
190       for (k=0; k<npeaks; k++) {
191          if (idx[k] < 0) continue; //removed peak
192          AliTPCcluster c;
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());
197
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);
201          if (s2 != 0.) {
202            c.SetSigmaY2(c.GetSigmaY2()*0.108);
203            if (sec<par->GetNInnerSector()) c.SetSigmaY2(c.GetSigmaY2()*2.07);
204          }
205
206          s2 = c.GetSigmaZ2()/c.GetQ() - c.GetZ()*c.GetZ();
207          w=par->GetZWidth();
208          c.SetSigmaZ2((s2 + 1./12.)*w*w);
209          if (s2 != 0.) {
210            c.SetSigmaZ2(c.GetSigmaZ2()*0.169);
211            if (sec<par->GetNInnerSector()) c.SetSigmaZ2(c.GetSigmaZ2()*1.77);
212          }
213
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()));
218
219          if (rx<230./250.*TMath::Abs(c.GetZ())) continue;
220
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);
225
226          c.SetQ(bins[idx[k]].GetQ());
227
228          if (ki==1 || ki==npads || kj==1 || kj==kMAXZ-2) {
229            c.SetSigmaY2(c.GetSigmaY2()*25.);
230            c.SetSigmaZ2(c.GetSigmaZ2()*4.);
231          }
232          clrow->InsertCluster(&c); ncl++;
233       }
234     }
235     fClusterArray.StoreRow(sec,row);
236     fClusterArray.ClearRow(sec,row);
237
238     nclusters+=ncl;
239
240     delete[] bins;  
241   }
242
243   Info("Digits2Cluster","Number of found clusters : %d",nclusters);
244
245   return 0;
246 }
247