3980b9f22fb3b5f16a0a62a90267808ede03abe2
[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 /*
17 $Log$
18 Revision 1.3  2000/10/05 16:14:01  kowal2
19 Forward declarations.
20
21 Revision 1.2  2000/06/30 12:07:50  kowal2
22 Updated from the TPC-PreRelease branch
23
24 Revision 1.1.2.1  2000/06/25 08:53:55  kowal2
25 Splitted from AliTPCtracking
26
27 */
28
29 //-------------------------------------------------------
30 //          Implementation of the TPC clusterer
31 //
32 //   Origin: Jouri Belikov, CERN, Jouri.Belikov@cern.ch 
33 //-------------------------------------------------------
34
35 #include "AliTPCclusterer.h"
36 #include "AliTPCcluster.h"
37 #include <TObjArray.h>
38 #include <TFile.h>
39 #include "AliTPCClustersArray.h"
40 #include "AliTPCClustersRow.h"
41 #include "AliDigits.h"
42 #include "AliSimDigits.h"
43 #include "AliTPCParam.h"
44 #include <iostream.h>
45 #include <TTree.h>
46
47 void AliTPCclusterer::FindPeaks(Int_t k,Int_t max,
48 AliBin *b,Int_t *idx,UInt_t *msk,Int_t& n) {
49   //find local maxima
50   if (n<31)
51   if (IsMaximum(k,max,b)) {
52     idx[n]=k; msk[n]=(2<<n);
53     n++;
54   }
55   b[k].SetMask(0);
56   if (b[k-max].GetMask()&1) FindPeaks(k-max,max,b,idx,msk,n);
57   if (b[k-1  ].GetMask()&1) FindPeaks(k-1  ,max,b,idx,msk,n);
58   if (b[k+max].GetMask()&1) FindPeaks(k+max,max,b,idx,msk,n);
59   if (b[k+1  ].GetMask()&1) FindPeaks(k+1  ,max,b,idx,msk,n);
60 }
61
62 void AliTPCclusterer::MarkPeak(Int_t k, Int_t max, AliBin *bins, UInt_t m) {
63   //mark this peak
64   UShort_t q=bins[k].GetQ();
65
66   bins[k].SetMask(bins[k].GetMask()|m); 
67
68   if (bins[k-max].GetQ() <= q)
69      if ((bins[k-max].GetMask()&m) == 0) MarkPeak(k-max,max,bins,m);
70   if (bins[k-1  ].GetQ() <= q)
71      if ((bins[k-1  ].GetMask()&m) == 0) MarkPeak(k-1  ,max,bins,m);
72   if (bins[k+max].GetQ() <= q)
73      if ((bins[k+max].GetMask()&m) == 0) MarkPeak(k+max,max,bins,m);
74   if (bins[k+1  ].GetQ() <= q)
75      if ((bins[k+1  ].GetMask()&m) == 0) MarkPeak(k+1  ,max,bins,m);
76 }
77
78 void AliTPCclusterer::MakeCluster(Int_t k,Int_t max,AliBin *bins,UInt_t m,
79 AliTPCcluster &c) {
80   //make cluster using digits of this peak
81   Float_t q=(Float_t)bins[k].GetQ();
82   Int_t i=k/max, j=k-i*max;
83
84   c.SetQ(c.GetQ()+q);
85   c.SetY(c.GetY()+i*q); 
86   c.SetZ(c.GetZ()+j*q); 
87   c.SetSigmaY2(c.GetSigmaY2()+i*i*q);
88   c.SetSigmaZ2(c.GetSigmaZ2()+j*j*q);
89
90   bins[k].SetMask(0xFFFFFFFE);
91   
92   if (bins[k-max].GetMask() == m) MakeCluster(k-max,max,bins,m,c);
93   if (bins[k-1  ].GetMask() == m) MakeCluster(k-1  ,max,bins,m,c);
94   if (bins[k+max].GetMask() == m) MakeCluster(k+max,max,bins,m,c);
95   if (bins[k+1  ].GetMask() == m) MakeCluster(k+1  ,max,bins,m,c);
96 }
97
98 //_____________________________________________________________________________
99 void AliTPCclusterer::Digits2Clusters(const AliTPCParam *par, TFile *of)
100 {
101   //-----------------------------------------------------------------
102   // This is a simple cluster finder.
103   //-----------------------------------------------------------------
104   TDirectory *savedir=gDirectory; 
105
106   if (!of->IsOpen()) {
107      cerr<<"AliTPC::Digits2Clusters(): output file not open !\n";
108      return;
109   }
110
111   const Int_t kMAXZ=par->GetMaxTBin()+2;
112
113   TTree *t = (TTree *)gDirectory->Get("TreeD_75x40_100x60");
114   AliSimDigits digarr, *dummy=&digarr;
115   t->GetBranch("Segment")->SetAddress(&dummy);
116   Stat_t nentries = t->GetEntries();
117
118   of->cd();
119
120   ((AliTPCParam*)par)->Write(par->GetTitle());
121   AliTPCClustersArray carray;
122   carray.Setup(par);
123   carray.SetClusterType("AliTPCcluster");
124   carray.MakeTree();
125
126   Int_t nclusters=0;
127
128   for (Int_t n=0; n<nentries; n++) {
129     t->GetEvent(n);
130     Int_t sec, row;
131     if (!par->AdjustSectorRow(digarr.GetID(),sec,row)) {
132        cerr<<"AliTPC warning: invalid segment ID ! "<<digarr.GetID()<<endl;
133        continue;
134     }
135
136     AliTPCClustersRow *clrow=carray.CreateRow(sec,row);
137
138     Float_t rx=par->GetPadRowRadii(sec,row);
139
140     Int_t npads, sign;
141     {
142        const Int_t kNIS=par->GetNInnerSector(), kNOS=par->GetNOuterSector();
143        if (sec < kNIS) {
144           npads = par->GetNPadsLow(row);
145           sign = (sec < kNIS/2) ? 1 : -1;
146        } else {
147           npads = par->GetNPadsUp(row);
148           sign = ((sec-kNIS) < kNOS/2) ? 1 : -1;
149        }
150     }
151
152     const Int_t kMAXBIN=kMAXZ*(npads+2);
153     AliBin *bins=new AliBin[kMAXBIN];
154     for (Int_t ii=0;ii<kMAXBIN;ii++) {
155        bins[ii].SetQ(0); bins[ii].SetMask(0xFFFFFFFE);
156     }
157
158     digarr.First();
159     do {
160        Short_t dig=digarr.CurrentDigit();
161        if (dig<=par->GetZeroSup()) continue;
162        Int_t j=digarr.CurrentRow()+1, i=digarr.CurrentColumn()+1;
163        bins[i*kMAXZ+j].SetQ(dig);
164        bins[i*kMAXZ+j].SetMask(1);
165     } while (digarr.Next());
166
167     Int_t ncl=0;
168     for (Int_t i=0; i<kMAXBIN; i++) {
169       if ((bins[i].GetMask()&1) == 0) continue;
170       Int_t idx[32]; UInt_t msk[32]; Int_t npeaks=0;
171       FindPeaks(i, kMAXZ, bins, idx, msk, npeaks);
172
173       if (npeaks>30) continue;
174
175       Int_t k,l;
176       for (k=0; k<npeaks-1; k++){//mark adjacent peaks
177         if (idx[k] < 0) continue; //this peak is already removed
178         for (l=k+1; l<npeaks; l++) {
179            if (idx[l] < 0) continue; //this peak is already removed
180            Int_t ki=idx[k]/kMAXZ, kj=idx[k] - ki*kMAXZ;
181            Int_t li=idx[l]/kMAXZ, lj=idx[l] - li*kMAXZ;
182            Int_t di=TMath::Abs(ki - li);
183            Int_t dj=TMath::Abs(kj - lj);
184            if (di>1 || dj>1) continue;
185            if (bins[idx[k]].GetQ() > bins[idx[l]].GetQ()) {
186               msk[l]=msk[k];
187               idx[l]*=-1;
188            } else {
189               msk[k]=msk[l];
190               idx[k]*=-1;
191               break;
192            } 
193         }
194       }
195
196       for (k=0; k<npeaks; k++) {
197         MarkPeak(TMath::Abs(idx[k]), kMAXZ, bins, msk[k]);
198       }
199         
200       for (k=0; k<npeaks; k++) {
201          if (idx[k] < 0) continue; //removed peak
202          AliTPCcluster c;
203          MakeCluster(idx[k], kMAXZ, bins, msk[k], c);
204          if (c.GetQ() < 5) continue; //noise cluster
205          c.SetY(c.GetY()/c.GetQ());
206          c.SetZ(c.GetZ()/c.GetQ());
207
208          Float_t s2 = c.GetSigmaY2()/c.GetQ() - c.GetY()*c.GetY();
209          Float_t w=par->GetPadPitchWidth(sec);
210          c.SetSigmaY2((s2 + 1./12.)*w*w);
211          if (s2 != 0.) {
212            c.SetSigmaY2(c.GetSigmaY2()*0.108);
213            if (sec<par->GetNInnerSector()) c.SetSigmaY2(c.GetSigmaY2()*2.07);
214          }
215
216          s2 = c.GetSigmaZ2()/c.GetQ() - c.GetZ()*c.GetZ();
217          w=par->GetZWidth();
218          c.SetSigmaZ2((s2 + 1./12.)*w*w);
219          if (s2 != 0.) {
220            c.SetSigmaZ2(c.GetSigmaZ2()*0.169);
221            if (sec<par->GetNInnerSector()) c.SetSigmaZ2(c.GetSigmaZ2()*1.77);
222          }
223
224          c.SetY((c.GetY() - 0.5 - 0.5*npads)*par->GetPadPitchWidth(sec));
225          c.SetZ(par->GetZWidth()*(c.GetZ()-1)); 
226          c.SetZ(c.GetZ() - 3.*par->GetZSigma()); // PASA delay 
227          c.SetZ(sign*(par->GetZLength() - c.GetZ()));
228
229          if (rx<230./250.*TMath::Abs(c.GetZ())) continue;
230
231          Int_t ki=idx[k]/kMAXZ, kj=idx[k] - ki*kMAXZ;
232          c.SetLabel(digarr.GetTrackID(kj-1,ki-1,0),0);
233          c.SetLabel(digarr.GetTrackID(kj-1,ki-1,1),1);
234          c.SetLabel(digarr.GetTrackID(kj-1,ki-1,2),2);
235
236          c.SetQ(bins[idx[k]].GetQ());
237
238          if (ki==1 || ki==npads || kj==1 || kj==kMAXZ-2) {
239            c.SetSigmaY2(c.GetSigmaY2()*25.);
240            c.SetSigmaZ2(c.GetSigmaZ2()*4.);
241          }
242          clrow->InsertCluster(&c); ncl++;
243       }
244     }
245     carray.StoreRow(sec,row);
246     carray.ClearRow(sec,row);
247
248     //cerr<<"sector, row, compressed digits, clusters: "
249     //<<sec<<' '<<row<<' '<<digarr.GetSize()<<' '<<ncl<<"                  \r";
250
251     nclusters+=ncl;
252
253     delete[] bins;  
254   }
255
256   cerr<<"Number of found clusters : "<<nclusters<<"                        \n";
257
258   carray.GetTree()->Write();
259   savedir->cd();
260 }
261