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