]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/src/AliL3FileHandler.cxx
c54ce76237842febf78208442d9358547b761654
[u/mrichter/AliRoot.git] / HLT / src / AliL3FileHandler.cxx
1 //$Id$
2
3 // Author: Uli Frankenfeld <mailto:franken@fi.uib.no>, Anders Vestbo <mailto:vestbo$fi.uib.no>
4 //*-- Copyright &copy Uli 
5
6 #include <math.h>
7 #include <iostream.h>
8
9 #include "AliL3Transform.h"
10 #include "AliL3Logging.h"
11 #include "AliL3MemHandler.h"
12 #include "AliL3FileHandler.h"
13
14 #include "AliTPCDigitsArray.h"
15 #include "AliTPCClustersArray.h"
16 #include "AliTPCcluster.h"
17 #include "AliTPCClustersRow.h"
18
19 #include "AliL3DigitData.h"
20 #include "AliL3TrackSegmentData.h"
21 #include "AliL3SpacePointData.h"
22 #include "AliL3TrackArray.h"
23
24 //_____________________________________________________________
25 // AliL3FileHandler
26 //
27 // The HLT ROOT <-> binary files handling class
28 //
29 // This class provides the interface between AliROOT files,
30 // and HLT binary files. It should be used for converting 
31 // TPC data stored in AliROOT format (outputfile from a simulation),
32 // into the data format currently used by in the HLT framework. 
33 // This enables the possibility to always use the same data format, 
34 // whether you are using a binary file as an input, or a AliROOT file.
35 //
36 // For example on how to create binary files from a AliROOT simulation,
37 // see example macro exa/Binary.C.
38 //
39 // For reading a AliROOT file into HLT format in memory, do the following:
40 //
41 // AliL3FileHandler file;
42 // file.SetAliInput("galice.root");
43 // AliL3DigitRowData *dataPt = (AliL3DigitRowData*)file.AliDigits2Memory(nrows,eventnr);
44 // 
45 // All the data are then stored in memory and accessible via the pointer dataPt.
46 // Accesing the data is then identical to the example 1) showed in AliL3MemHandler class.
47 //
48 // For converting the data back, and writing it to a new AliROOT file do:
49 //
50 // AliL3FileHandler file;
51 // file.SetAliInput("galice.root");
52 // file.Init(slice,patch,NumberOfRowsInPatch);
53 // file.AliDigits2RootFile(dataPt,"new_galice.root");
54 // file.CloseAliInput();
55
56 ClassImp(AliL3FileHandler)
57
58 AliL3FileHandler::AliL3FileHandler()
59 {
60   //Default constructor
61   fInAli = 0;
62   fParam = 0;
63   fMC =0;
64   fLastIndex=0;
65   fDigits=0;
66   fDigitsTree=0;
67 }
68
69 AliL3FileHandler::~AliL3FileHandler()
70 {
71   //Destructor
72   if(fMC) CloseMCOutput();
73   FreeDigitsTree();
74   if(fInAli) CloseAliInput();
75   
76 }
77
78 void AliL3FileHandler::FreeDigitsTree()
79 {
80   if(!fDigitsTree)
81     {
82       LOG(AliL3Log::kWarning,"AliL3FileHandler::FreeDigitsTree()","Pointer")
83         <<"Cannot free digitstree, it is not present"<<ENDLOG;
84       return;
85     }
86   fDigits=0;
87   fDigitsTree->Delete();
88   fDigitsTree=0;
89 }
90
91 Bool_t AliL3FileHandler::SetMCOutput(char *name)
92 {
93   fMC = fopen(name,"w");
94   if(!fMC){
95     LOG(AliL3Log::kWarning,"AliL3FileHandler::SetMCOutput","File Open")
96       <<"Pointer to File = 0x0 "<<ENDLOG;
97     return kFALSE;
98   }
99   return kTRUE;
100 }
101
102 Bool_t AliL3FileHandler::SetMCOutput(FILE *file)
103 {
104   fMC = file;
105   if(!fMC){
106     LOG(AliL3Log::kWarning,"AliL3FileHandler::SetMCOutput","File Open")
107       <<"Pointer to File = 0x0 "<<ENDLOG;
108     return kFALSE;
109   }
110   return kTRUE;
111 }
112
113 void AliL3FileHandler::CloseMCOutput()
114 {
115   if(!fMC){
116     LOG(AliL3Log::kWarning,"AliL3FileHandler::CloseMCOutPut","File Close")
117       <<"Nothing to Close"<<ENDLOG;
118     return;
119   }
120   fclose(fMC);
121   fMC =0;
122 }
123
124 Bool_t AliL3FileHandler::SetAliInput()
125 {
126   if(!fInAli->IsOpen()){
127     LOG(AliL3Log::kError,"AliL3FileHandler::SetAliInput","File Open")
128       <<"Ali File "<<fInAli->GetName()<<" does not exist"<<ENDLOG;
129     return kFALSE;
130   }
131   fParam = (AliTPCParam*)fInAli->Get("75x40_100x60");
132   if(!fParam){ 
133     LOG(AliL3Log::kError,"AliL3FileHandler::SetAliInput","File Open")
134       <<"No AliTPCParam 75x40_100x60 in File "<<fInAli->GetName()<<ENDLOG;
135     return kFALSE;
136   }
137   return kTRUE;
138 }
139
140 Bool_t AliL3FileHandler::SetAliInput(char *name)
141 {
142   //Open the AliROOT file with name.
143   
144   fInAli= new TFile(name,"READ");
145   if(!fInAli){
146     LOG(AliL3Log::kWarning,"AliL3FileHandler::SetAliInput","File Open")
147     <<"Pointer to TFile = 0x0 "<<ENDLOG;
148     return kFALSE;
149   }
150   return SetAliInput();
151 }
152
153 Bool_t AliL3FileHandler::SetAliInput(TFile *file)
154 {
155   //Specify already opened AliROOT file to use as an input.
156   
157   fInAli=file;
158   if(!fInAli){
159     LOG(AliL3Log::kWarning,"AliL3FileHandler::SetAliInput","File Open")
160     <<"Pointer to TFile = 0x0 "<<ENDLOG;
161     return kFALSE;
162   }
163   return SetAliInput();
164 }
165
166 void AliL3FileHandler::CloseAliInput()
167 {
168   if(!fInAli){
169     LOG(AliL3Log::kWarning,"AliL3FileHandler::CloseAliInput","File Close")
170       <<"Nothing to Close"<<ENDLOG;
171     return;
172   }
173   if(fInAli->IsOpen()) fInAli->Close();
174   delete fInAli;
175   fInAli = 0;
176   
177 }
178
179 Bool_t AliL3FileHandler::IsDigit()
180 {
181   //Check if there is a TPC digit tree in the current file.
182   //Return kTRUE if tree was found, and kFALSE if not found.
183   
184   if(!fInAli){
185     LOG(AliL3Log::kWarning,"AliL3FileHandler::IsDigit","File")
186     <<"Pointer to TFile = 0x0 "<<ENDLOG;
187     return kTRUE;  //may you are use binary input which is Digits!!
188   }
189   TTree *t=(TTree*)fInAli->Get("TreeD_75x40_100x60_0");
190   if(t){
191     LOG(AliL3Log::kInformational,"AliL3FileHandler::IsDigit","File Type")
192     <<"Found Digit Tree -> Use Fast Cluster Finder"<<ENDLOG;
193     return kTRUE;
194   }
195   else{
196     LOG(AliL3Log::kInformational,"AliL3FileHandler::IsDigit","File Type")
197     <<"No Digit Tree -> Use Cluster Tree"<<ENDLOG;
198     return kFALSE;
199   }
200 }
201
202 ///////////////////////////////////////// Digit IO  
203 Bool_t AliL3FileHandler::AliDigits2Binary(Int_t event)
204 {
205   
206   Bool_t out = kTRUE;
207   UInt_t nrow;
208   AliL3DigitRowData* data = AliDigits2Memory(nrow,event);
209   out = Memory2Binary(nrow,data);
210   Free();
211   return out;
212 }
213
214 Bool_t AliL3FileHandler::AliDigits2CompBinary(Int_t event)
215 {
216   //Convert AliROOT TPC data, into HLT data format.
217   //event specifies the event you want in the aliroot file.
218   
219   Bool_t out = kTRUE;
220   UInt_t ndigits=0;
221   AliL3DigitRowData* digits = AliDigits2Memory(ndigits,event);
222   out = Memory2CompBinary(ndigits,digits);
223   Free();
224   return out;
225 }
226
227 AliL3DigitRowData * AliL3FileHandler::AliDigits2Memory(UInt_t & nrow,Int_t event)
228 {
229   //Read data from AliROOT file into memory, and store it in the HLT data format.
230   //Returns a pointer to the data.
231
232   AliL3DigitRowData *data = 0;
233   nrow=0;
234   
235   if(!fInAli){
236     LOG(AliL3Log::kWarning,"AliL3FileHandler::AliDigits2Memory","File")
237     <<"No Input avalible: no object TFile"<<ENDLOG;
238     return 0; 
239   }
240   if(!fInAli->IsOpen()){
241     LOG(AliL3Log::kWarning,"AliL3FileHandler::AliDigits2Memory","File")
242     <<"No Input avalible: TFile not opend"<<ENDLOG;
243     return 0;
244   }
245
246   if(!fDigitsTree)
247     GetDigitsTree(event);
248   
249   UShort_t dig;
250   Int_t time,pad,sector,row;
251   Int_t nrows=0;
252   Int_t ndigitcount=0;
253   Int_t entries = (Int_t)fDigitsTree->GetEntries();
254   Int_t ndigits[entries];
255   Int_t lslice,lrow;
256   Float_t xyz[3];
257   
258   for(Int_t n=fLastIndex; n<fDigitsTree->GetEntries(); n++)
259     {
260       fDigitsTree->GetEvent(n);
261       fParam->AdjustSectorRow(fDigits->GetID(),sector,row);
262       AliL3Transform::Sector2Slice(lslice,lrow,sector,row);
263       //if(fSlice != lslice || lrow<fRowMin || lrow>fRowMax) continue;
264       if(lslice < fSlice) continue;
265       if(lslice != fSlice) break;
266       if(lrow < fRowMin) continue;
267       if(lrow > fRowMax) break;
268
269       ndigits[lrow] = 0;
270       fDigits->First();
271       do {
272         time=fDigits->CurrentRow();
273         pad=fDigits->CurrentColumn();
274         dig = fDigits->GetDigit(time,pad);
275         if(dig<=fParam->GetZeroSup()) continue;
276         
277         /*
278           if(time < fParam->GetMaxTBin()-1 && time > 0)
279           if(fDigits->GetDigit(time+1,pad) <= fParam->GetZeroSup()
280           && fDigits->GetDigit(time-1,pad) <= fParam->GetZeroSup())
281           continue;
282         */
283         
284         AliL3Transform::Raw2Local(xyz,sector,row,pad,time);
285         if(fParam->GetPadRowRadii(sector,row)<230./250.*fabs(xyz[2]))
286           continue;
287
288         ndigits[lrow]++; //for this row only
289         ndigitcount++;   //total number of digits to be published
290
291       } while (fDigits->Next());
292       //cout << lrow << " " << ndigits[lrow] << " - " << ndigitcount << endl;
293       nrows++;
294     }
295
296   Int_t size = sizeof(AliL3DigitData)*ndigitcount
297     + nrows*sizeof(AliL3DigitRowData);
298
299   LOG(AliL3Log::kDebug,"AliL3FileHandler::AliDigits2Memory","Digits")
300     <<AliL3Log::kDec<<"Found "<<ndigitcount<<" Digits"<<ENDLOG;
301   
302   data=(AliL3DigitRowData*) Allocate(size);
303   nrow = (UInt_t)nrows;
304   AliL3DigitRowData *tempPt = data;
305   for(Int_t n=fLastIndex; n<fDigitsTree->GetEntries(); n++)
306     {
307       fDigitsTree->GetEvent(n);
308       fParam->AdjustSectorRow(fDigits->GetID(),sector,row);
309       AliL3Transform::Sector2Slice(lslice,lrow,sector,row);
310       //if(fSlice != lslice || lrow<fRowMin || lrow>fRowMax) continue;
311       if(lslice < fSlice) continue;
312       if(lslice != fSlice) break;
313       if(lrow < fRowMin) continue;
314       if(lrow > fRowMax) break;
315
316       tempPt->fRow = lrow;
317       tempPt->fNDigit = ndigits[lrow];
318
319       Int_t localcount=0;
320       fDigits->First();
321       do {
322         //dig=fDigits->CurrentDigit();
323         //if (dig<=fParam->GetZeroSup()) continue;
324         time=fDigits->CurrentRow();
325         pad=fDigits->CurrentColumn();
326         dig = fDigits->GetDigit(time,pad);
327         if (dig <= fParam->GetZeroSup()) continue;
328         
329         /*
330           if(time < fParam->GetMaxTBin()-1 && time > 0)
331           if(fDigits->GetDigit(time-1,pad) <= fParam->GetZeroSup() &&
332           fDigits->GetDigit(time+1,pad) <= fParam->GetZeroSup()) continue;
333         */
334         
335         //Exclude data outside cone:
336         AliL3Transform::Raw2Local(xyz,sector,row,pad,time);
337         if(fParam->GetPadRowRadii(sector,row)<230./250.*fabs(xyz[2]))
338           continue;
339
340         if(localcount >= ndigits[lrow])
341           LOG(AliL3Log::kFatal,"AliL3FileHandler::AliDigits2Binary","Memory")
342             <<AliL3Log::kDec<<"Mismatch: localcount "<<localcount<<" ndigits "
343             <<ndigits[lrow]<<ENDLOG;
344         
345         tempPt->fDigitData[localcount].fCharge=dig;
346         tempPt->fDigitData[localcount].fPad=pad;
347         tempPt->fDigitData[localcount].fTime=time;
348 #ifdef do_mc
349         tempPt->fDigitData[localcount].fTrackID[0] = fDigits->GetTrackID(time,pad,0);
350         tempPt->fDigitData[localcount].fTrackID[1] = fDigits->GetTrackID(time,pad,1);
351         tempPt->fDigitData[localcount].fTrackID[2] = fDigits->GetTrackID(time,pad,2);
352 #endif
353         localcount++;
354       } while (fDigits->Next());
355
356       Byte_t *tmp = (Byte_t*)tempPt;
357       Int_t size = sizeof(AliL3DigitRowData)
358                                       + ndigits[lrow]*sizeof(AliL3DigitData);
359       tmp += size;
360       tempPt = (AliL3DigitRowData*)tmp;
361       //fLastIndex=n;
362     }
363   //fLastIndex++;
364   return data;
365 }
366
367 Bool_t AliL3FileHandler::GetDigitsTree(Int_t event)
368 {
369   //Connects to the TPC digit tree in the AliROOT file.
370   
371   fInAli->cd();
372   Char_t dname[100];
373   sprintf(dname,"TreeD_75x40_100x60_%d",event);
374   fDigitsTree = (TTree*)fInAli->Get(dname);
375   if(!fDigitsTree) 
376     {
377       LOG(AliL3Log::kError,"AliL3FileHandler::GetDigitsTree","Digits Tree")
378         <<AliL3Log::kHex<<"Error getting digitstree "<<(Int_t)fDigitsTree<<ENDLOG;
379       return kFALSE;
380     }
381   fDigitsTree->GetBranch("Segment")->SetAddress(&fDigits);
382   return kTRUE;
383 }
384
385 void AliL3FileHandler::AliDigits2RootFile(AliL3DigitRowData *rowPt,Char_t *new_digitsfile)
386 {
387   //Write the data stored in rowPt, into a new AliROOT file.
388   //The data is stored in the AliROOT format 
389   //This is specially a nice thing if you have modified data, and wants to run it  
390   //through the offline reconstruction chain.
391   //The arguments is a pointer to the data, and the name of the new AliROOT file.
392   //Remember to pass the original AliROOT file (the one that contains the original
393   //simulated data) to this object, in order to retrieve the MC id's of the digits.
394   
395   if(!fInAli)
396     {
397       printf("AliL3FileHandler::AliDigits2RootFile : No rootfile\n");
398       return;
399     }
400   if(!fParam)
401     {
402       printf("AliL3FileHandler::AliDigits2RootFile : No parameter object. Run on rootfile\n");
403       return;
404     }
405   
406   //Get the original digitstree:
407   fInAli->cd();
408   AliTPCDigitsArray *old_array = new AliTPCDigitsArray();
409   old_array->Setup(fParam);
410   old_array->SetClass("AliSimDigits");
411   Bool_t ok = old_array->ConnectTree("TreeD_75x40_100x60_0");
412   if(!ok)
413     {
414       printf("AliL3FileHandler::AliDigits2RootFile : No digits tree object\n");
415       return;
416     }
417   
418   Bool_t create=kFALSE;
419   TFile *digFile;
420   
421   digFile = TFile::Open(new_digitsfile,"NEW");
422   if(digFile->IsOpen())
423     {    
424       create = kTRUE;
425       fParam->Write(fParam->GetTitle());
426     }
427   else
428     {
429       LOG(AliL3Log::kDebug,"AliL3FileHandler::AliDigits2RootFile","Rootfile")
430         <<"Rootfile did already exist, so I will just open it for updates"<<ENDLOG;
431       digFile = TFile::Open(new_digitsfile,"UPDATE");
432       create=kFALSE;
433     }
434   if(!digFile->IsOpen())
435     {
436       LOG(AliL3Log::kError,"AliL3FileHandler::AliDigits2RootFile","Rootfile")
437         <<"Error opening rootfile "<<new_digitsfile<<ENDLOG;
438       return;
439     }
440   
441   digFile->cd();
442     
443   //setup a new one, or connect it to the existing one:
444   AliTPCDigitsArray *arr = new AliTPCDigitsArray; 
445   arr->SetClass("AliSimDigits");
446   arr->Setup(fParam);
447   if(create)
448     arr->MakeTree();
449   else
450     {
451       Bool_t ok = arr->ConnectTree("TreeD_75x40_100x60_0");
452       if(!ok)
453         {
454           printf("AliL3FileHandler::AliDigits2RootFile : No digits tree object in existing file\n");
455           return;
456         }
457     }
458   Int_t digcounter=0;
459
460   for(Int_t i=fRowMin; i<=fRowMax; i++)
461     {
462       
463       if((Int_t)rowPt->fRow != i) printf("AliL3FileHandler::AliDigits2RootFile : Mismatching row numbering!!!\n");
464             
465       Int_t sector,row;
466       AliL3Transform::Slice2Sector(fSlice,i,sector,row);
467       AliSimDigits * dig = (AliSimDigits*)arr->CreateRow(sector,row);
468       AliSimDigits *old_dig = (AliSimDigits*)old_array->LoadRow(sector,row);
469       if(!old_dig)
470         printf("AliL3FileHandler::AliDigits2RootFile : No padrow %d %d\n",sector,row);
471       
472       AliL3DigitData *digPt = rowPt->fDigitData;
473       digcounter=0;
474       for(UInt_t j=0; j<rowPt->fNDigit; j++)
475         {
476           Int_t charge = (Int_t)digPt[j].fCharge;
477           Int_t pad = (Int_t)digPt[j].fPad;
478           Int_t time = (Int_t)digPt[j].fTime;
479           
480           if(charge == 0) //Only write the digits that has not been removed
481             continue;
482           digcounter++;
483           dig->SetDigitFast(charge,time,pad);
484           
485           Int_t trackID[3] = {old_dig->GetTrackID(time,pad,0),old_dig->GetTrackID(time,pad,1),old_dig->GetTrackID(time,pad,2)};
486           Int_t s_pad = pad;
487           Int_t s_time = time - 1;
488           while(trackID[0] < 0)
489             {
490               if(s_time >= 0 && s_time < AliL3Transform::GetNTimeBins() && s_pad >= 0 && s_pad < AliL3Transform::GetNPads(i))
491                 {
492                   if(old_dig->GetTrackID(s_time,s_pad,0) > 0)
493                     {
494                       trackID[0]=old_dig->GetTrackID(s_time,s_pad,0); 
495                       trackID[1]=old_dig->GetTrackID(s_time,s_pad,1); 
496                       trackID[2]=old_dig->GetTrackID(s_time,s_pad,2); 
497                     }
498                 }
499               if(s_pad == pad && s_time == time - 1)
500                 s_time = time + 1;
501               else if(s_pad == pad && s_time == time + 1)
502                 {s_pad = pad - 1; s_time = time;}
503               else if(s_pad == pad - 1 && s_time == time)
504                 s_time = time - 1;
505               else if(s_pad == pad - 1 && s_time == time - 1)
506                 s_time = time + 1;
507               else if(s_pad == pad - 1 && s_time == time + 1)
508                 {s_pad = pad + 1; s_time = time;}
509               else if(s_pad == pad + 1 && s_time == time)
510                 s_time = time - 1;
511               else if(s_pad == pad + 1 && s_time == time - 1)
512                 s_time = time + 1;
513               else 
514                 break;
515             }
516           
517           dig->SetTrackIDFast(trackID[0],time,pad,0);
518           dig->SetTrackIDFast(trackID[1],time,pad,1);
519           dig->SetTrackIDFast(trackID[2],time,pad,2);
520           
521         }
522       //cout<<"Wrote "<<digcounter<<" on row "<<i<<endl;
523       UpdateRowPointer(rowPt);
524       arr->StoreRow(sector,row);
525       arr->ClearRow(sector,row);  
526       old_array->ClearRow(sector,row);
527     }
528   digFile->cd();
529   char treeName[100];
530   sprintf(treeName,"TreeD_%s_0",fParam->GetTitle());
531   printf("Writing tree to file.....");
532   arr->GetTree()->Write(treeName,TObject::kOverwrite);
533   printf("done\n");
534   digFile->Close();
535   //arr->GetTree()->Delete();
536   //delete arr;
537 }
538
539 ///////////////////////////////////////// Point IO  
540 Bool_t AliL3FileHandler::AliPoints2Binary(){
541   Bool_t out = kTRUE;
542   UInt_t npoint;
543   AliL3SpacePointData *data = AliPoints2Memory(npoint);
544   out = Memory2Binary(npoint,data);
545   Free();
546   return out;
547 }
548
549 AliL3SpacePointData * AliL3FileHandler::AliPoints2Memory(UInt_t & npoint){
550   AliL3SpacePointData *data = 0;
551   npoint=0;
552   if(!fInAli){
553     LOG(AliL3Log::kWarning,"AliL3FileHandler::AliPoints2Memory","File")
554     <<"No Input avalible: no object TFile"<<ENDLOG;
555     return 0;
556   }
557   if(!fInAli->IsOpen()){
558     LOG(AliL3Log::kWarning,"AliL3FileHandler::AliPoints2Memory","File")
559     <<"No Input avalible: TFile not opend"<<ENDLOG;
560     return 0;
561   }
562
563   TDirectory *savedir = gDirectory;
564   fInAli->cd();
565   
566   Char_t cname[100];
567   Int_t eventn = 0;
568   sprintf(cname,"TreeC_TPC_%d",eventn);
569   AliTPCClustersArray carray;
570   carray.Setup(fParam);
571   carray.SetClusterType("AliTPCcluster");
572   Bool_t clusterok = carray.ConnectTree(cname);
573   if(!clusterok) return 0;
574
575   AliTPCClustersRow ** clusterrow = 
576                new AliTPCClustersRow*[ (int)carray.GetTree()->GetEntries()];
577   Int_t *rows = new int[ (int)carray.GetTree()->GetEntries()];
578   Int_t *sects = new int[  (int)carray.GetTree()->GetEntries()];
579   Int_t sum=0;
580
581   Int_t lslice,lrow;
582   for(Int_t i=0; i<carray.GetTree()->GetEntries(); i++){
583     AliSegmentID *s = carray.LoadEntry(i);
584     Int_t sector,row;
585     fParam->AdjustSectorRow(s->GetID(),sector,row);
586     rows[i] = row;
587     sects[i] = sector;
588     clusterrow[i] = 0;
589     AliL3Transform::Sector2Slice(lslice,lrow,sector,row);
590     if(fSlice != lslice || lrow<fRowMin || lrow>fRowMax) continue;
591     clusterrow[i] = carray.GetRow(sector,row);
592     if(clusterrow[i])
593       sum+=clusterrow[i]->GetArray()->GetEntriesFast();
594   }
595   UInt_t size = sum*sizeof(AliL3SpacePointData);
596
597   LOG(AliL3Log::kDebug,"AliL3FileHandler::AliPoints2Memory","File")
598   <<AliL3Log::kDec<<"Found "<<sum<<" SpacePoints"<<ENDLOG;
599
600   data = (AliL3SpacePointData *) Allocate(size);
601   npoint = sum;
602   UInt_t n=0; 
603   for(Int_t i=0; i<carray.GetTree()->GetEntries(); i++){
604     if(!clusterrow[i]) continue;
605     Int_t row = rows[i];
606     Int_t sector = sects[i];
607     AliL3Transform::Sector2Slice(lslice,lrow,sector,row);
608     Int_t entries_in_row = clusterrow[i]->GetArray()->GetEntriesFast();
609     for(Int_t j = 0;j<entries_in_row;j++){
610       AliTPCcluster *c = (AliTPCcluster*)(*clusterrow[i])[j];
611       data[n].fZ = c->GetZ();
612       data[n].fY = c->GetY();
613       data[n].fX = fParam->GetPadRowRadii(sector,row);
614       data[n].fID = n+((fSlice&0x7f)<<25)+((fPatch&0x7)<<22);//uli
615       data[n].fPadRow = lrow;
616       data[n].fXYErr = c->GetSigmaY2();
617       data[n].fZErr = c->GetSigmaZ2();
618       if(fMC) fprintf(fMC,"%d %d\n",data[n].fID,c->GetLabel(0));
619       n++;
620     }
621   }
622   for(Int_t i=0;i<carray.GetTree()->GetEntries();i++){
623     Int_t row = rows[i];
624     Int_t sector = sects[i];
625     if(carray.GetRow(sector,row))
626       carray.ClearRow(sector,row);
627   }
628
629   delete [] clusterrow;
630   delete [] rows;
631   delete [] sects;
632   savedir->cd();   
633
634   return data;
635 }
636