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