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