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