]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/src/AliL3FileHandler.cxx
Some more comments.
[u/mrichter/AliRoot.git] / HLT / src / AliL3FileHandler.cxx
1 //$Id$
2
3 // Author: Uli Frankenfeld <mailto:franken@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   if(fDigitsTree) delete fDigitsTree;
74   if(fInAli) CloseAliInput();
75   
76 }
77
78 Bool_t AliL3FileHandler::SetMCOutput(char *name)
79 {
80   fMC = fopen(name,"w");
81   if(!fMC){
82     LOG(AliL3Log::kWarning,"AliL3FileHandler::SetMCOutput","File Open")
83       <<"Pointer to File = 0x0 "<<ENDLOG;
84     return kFALSE;
85   }
86   return kTRUE;
87 }
88
89 Bool_t AliL3FileHandler::SetMCOutput(FILE *file)
90 {
91   fMC = file;
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 void AliL3FileHandler::CloseMCOutput()
101 {
102   if(!fMC){
103     LOG(AliL3Log::kWarning,"AliL3FileHandler::CloseMCOutPut","File Close")
104       <<"Nothing to Close"<<ENDLOG;
105     return;
106   }
107   fclose(fMC);
108   fMC =0;
109 }
110
111 Bool_t AliL3FileHandler::SetAliInput()
112 {
113   if(!fInAli->IsOpen()){
114     LOG(AliL3Log::kError,"AliL3FileHandler::SetAliInput","File Open")
115       <<"Ali File "<<fInAli->GetName()<<" does not exist"<<ENDLOG;
116     return kFALSE;
117   }
118   fParam = (AliTPCParam*)fInAli->Get("75x40_100x60");
119   if(!fParam){ 
120     LOG(AliL3Log::kError,"AliL3FileHandler::SetAliInput","File Open")
121       <<"No AliTPCParam 75x40_100x60 in File "<<fInAli->GetName()<<ENDLOG;
122     return kFALSE;
123   }
124   return kTRUE;
125 }
126
127 Bool_t AliL3FileHandler::SetAliInput(char *name)
128 {
129   //Open the AliROOT file with name.
130   
131   fInAli= new TFile(name,"READ");
132   if(!fInAli){
133     LOG(AliL3Log::kWarning,"AliL3FileHandler::SetAliInput","File Open")
134     <<"Pointer to TFile = 0x0 "<<ENDLOG;
135     return kFALSE;
136   }
137   return SetAliInput();
138 }
139
140 Bool_t AliL3FileHandler::SetAliInput(TFile *file)
141 {
142   //Specify already opened AliROOT file to use as an input.
143   
144   fInAli=file;
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 void AliL3FileHandler::CloseAliInput()
154 {
155   if(!fInAli){
156     LOG(AliL3Log::kWarning,"AliL3FileHandler::CloseAliInput","File Close")
157       <<"Nothing to Close"<<ENDLOG;
158     return;
159   }
160   if(fInAli->IsOpen()) fInAli->Close();
161   delete fInAli;
162   fInAli = 0;
163   
164 }
165
166 Bool_t AliL3FileHandler::IsDigit()
167 {
168   //Check if there is a TPC digit tree in the current file.
169   //Return kTRUE if tree was found, and kFALSE if not found.
170   
171   if(!fInAli){
172     LOG(AliL3Log::kWarning,"AliL3FileHandler::IsDigit","File")
173     <<"Pointer to TFile = 0x0 "<<ENDLOG;
174     return kTRUE;  //may you are use binary input which is Digits!!
175   }
176   TTree *t=(TTree*)fInAli->Get("TreeD_75x40_100x60_0");
177   if(t){
178     LOG(AliL3Log::kInformational,"AliL3FileHandler::IsDigit","File Type")
179     <<"Found Digit Tree -> Use Fast Cluster Finder"<<ENDLOG;
180     return kTRUE;
181   }
182   else{
183     LOG(AliL3Log::kInformational,"AliL3FileHandler::IsDigit","File Type")
184     <<"No Digit Tree -> Use Cluster Tree"<<ENDLOG;
185     return kFALSE;
186   }
187 }
188
189 ///////////////////////////////////////// Digit IO  
190 Bool_t AliL3FileHandler::AliDigits2Binary(Int_t event)
191 {
192   
193   Bool_t out = kTRUE;
194   UInt_t nrow;
195   AliL3DigitRowData* data = AliDigits2Memory(nrow,event);
196   out = Memory2Binary(nrow,data);
197   Free();
198   return out;
199 }
200
201
202 Bool_t AliL3FileHandler::AliDigits2CompBinary(Int_t event)
203 {
204   //Convert AliROOT TPC data, into HLT data format.
205   //event specifies the event you want in the aliroot file.
206   
207   Bool_t out = kTRUE;
208   UInt_t ndigits=0;
209   AliL3DigitRowData* digits = AliDigits2Memory(ndigits,event);
210   out = Memory2CompBinary(ndigits,digits);
211   Free();
212   return out;
213 }
214
215
216 AliL3DigitRowData * AliL3FileHandler::AliDigits2Memory(UInt_t & nrow,Int_t event)
217 {
218   //Read data from AliROOT file into memory, and store it in the HLT data format.
219   //Returns a pointer to the data.
220
221   AliL3DigitRowData *data = 0;
222   nrow=0;
223   if(!fInAli){
224     LOG(AliL3Log::kWarning,"AliL3FileHandler::AliDigits2Memory","File")
225     <<"No Input avalible: no object TFile"<<ENDLOG;
226     return 0; 
227   }
228   if(!fInAli->IsOpen()){
229     LOG(AliL3Log::kWarning,"AliL3FileHandler::AliDigits2Memory","File")
230     <<"No Input avalible: TFile not opend"<<ENDLOG;
231     return 0;
232   }
233   if(!fTransformer)
234     {
235       LOG(AliL3Log::kWarning,"AliL3FileHandler::AliDigits2Memory","Transformer")
236         <<"No transformer object"<<ENDLOG;
237       return 0;
238     }
239   
240   if(!fDigitsTree)
241     GetDigitsTree(event);
242     
243   UShort_t dig;
244   Int_t time,pad,sector,row;
245   Int_t nrows=0;
246   Int_t ndigitcount=0;
247   Int_t entries = (Int_t)fDigitsTree->GetEntries();
248   Int_t ndigits[entries];
249   Int_t lslice,lrow;
250   
251   for(Int_t n=fLastIndex; n<fDigitsTree->GetEntries(); n++)
252     {
253       fDigitsTree->GetEvent(n);
254       fParam->AdjustSectorRow(fDigits->GetID(),sector,row);
255       fTransformer->Sector2Slice(lslice,lrow,sector,row);
256       //if(fSlice != lslice || lrow<fRowMin || lrow>fRowMax) continue;
257       if(lslice < fSlice) continue;
258       if(lslice != fSlice) break;
259       if(lrow < fRowMin) continue;
260       if(lrow > fRowMax) break;
261
262       Float_t xyz[3];
263       ndigits[lrow] = 0;
264       fDigits->First();
265       do {
266         time=fDigits->CurrentRow();
267         pad=fDigits->CurrentColumn();
268         dig = fDigits->GetDigit(time,pad);
269         if(dig<=fParam->GetZeroSup()) continue;
270         if(time < fParam->GetMaxTBin()-1 && time > 0)
271           if(fDigits->GetDigit(time+1,pad) <= fParam->GetZeroSup()
272              && fDigits->GetDigit(time-1,pad) <= fParam->GetZeroSup())
273             continue;
274
275         fTransformer->Raw2Local(xyz,sector,row,pad,time);
276         if(fParam->GetPadRowRadii(sector,row)<230./250.*fabs(xyz[2]))
277           continue;
278
279         ndigits[lrow]++; //for this row only
280         ndigitcount++;  //total number of digits to be published
281
282       } while (fDigits->Next());
283
284       nrows++;
285     }
286   Int_t size = sizeof(AliL3DigitData)*ndigitcount
287     + nrows*sizeof(AliL3DigitRowData);
288
289   LOG(AliL3Log::kDebug,"AliL3FileHandler::AliDigits2Memory","Digits")
290     <<AliL3Log::kDec<<"Found "<<ndigitcount<<" Digits"<<ENDLOG;
291   
292   data=(AliL3DigitRowData*) Allocate(size);
293   nrow = (UInt_t)nrows;
294   AliL3DigitRowData *tempPt = data;
295   for(Int_t n=fLastIndex; n<fDigitsTree->GetEntries(); n++)
296     {
297       fDigitsTree->GetEvent(n);
298       Float_t xyz[3];
299       fParam->AdjustSectorRow(fDigits->GetID(),sector,row);
300       fTransformer->Sector2Slice(lslice,lrow,sector,row);
301       //if(fSlice != lslice || lrow<fRowMin || lrow>fRowMax) continue;
302       if(lslice < fSlice) continue;
303       if(lslice != fSlice) break;
304       if(lrow < fRowMin) continue;
305       if(lrow > fRowMax) break;
306
307       tempPt->fRow = lrow;
308       tempPt->fNDigit = ndigits[lrow];
309
310       Int_t localcount=0;
311       fDigits->First();
312       do {
313         //dig=fDigits->CurrentDigit();
314         //if (dig<=fParam->GetZeroSup()) continue;
315         time=fDigits->CurrentRow();
316         pad=fDigits->CurrentColumn();
317         dig = fDigits->GetDigit(time,pad);
318         if (dig <= fParam->GetZeroSup()) continue;
319         if(time < fParam->GetMaxTBin()-1 && time > 0)
320           if(fDigits->GetDigit(time-1,pad) <= fParam->GetZeroSup() &&
321              fDigits->GetDigit(time+1,pad) <= fParam->GetZeroSup()) continue;
322
323         //Exclude data outside cone:
324         fTransformer->Raw2Local(xyz,sector,row,pad,time);
325         if(fParam->GetPadRowRadii(sector,row)<230./250.*fabs(xyz[2]))
326           continue;
327
328         if(localcount >= ndigits[lrow])
329           LOG(AliL3Log::kFatal,"AliL3FileHandler::AliDigits2Binary","Memory")
330             <<AliL3Log::kDec<<"Mismatch: localcount "<<localcount<<" ndigits "
331             <<ndigits[lrow]<<ENDLOG;
332         
333         tempPt->fDigitData[localcount].fCharge=dig;
334         tempPt->fDigitData[localcount].fPad=pad;
335         tempPt->fDigitData[localcount].fTime=time;
336 #ifdef do_mc
337         tempPt->fDigitData[localcount].fTrackID[0] = fDigits->GetTrackID(time,pad,0);
338         tempPt->fDigitData[localcount].fTrackID[1] = fDigits->GetTrackID(time,pad,1);
339         tempPt->fDigitData[localcount].fTrackID[2] = fDigits->GetTrackID(time,pad,2);
340 #endif
341         localcount++;
342       } while (fDigits->Next());
343
344       Byte_t *tmp = (Byte_t*)tempPt;
345       Int_t size = sizeof(AliL3DigitRowData)
346                                       + ndigits[lrow]*sizeof(AliL3DigitData);
347       tmp += size;
348       tempPt = (AliL3DigitRowData*)tmp;
349       //fLastIndex=n;
350     }
351   //fLastIndex++;
352   return data;
353 }
354
355 Bool_t AliL3FileHandler::GetDigitsTree(Int_t event)
356 {
357   //Connects to the TPC digit tree in the AliROOT file.
358   
359   fInAli->cd();
360   Char_t dname[100];
361   sprintf(dname,"TreeD_75x40_100x60_%d",event);
362   fDigitsTree = (TTree*)fInAli->Get("TreeD_75x40_100x60_0");
363   if(!fDigitsTree) 
364     {
365       LOG(AliL3Log::kError,"AliL3FileHandler::GetDigitsTree","Digits Tree")
366         <<AliL3Log::kHex<<"Error getting digitstree "<<(Int_t)fDigitsTree<<ENDLOG;
367       return kFALSE;
368     }
369   fDigitsTree->GetBranch("Segment")->SetAddress(&fDigits);
370   return kTRUE;
371 }
372
373 void AliL3FileHandler::AliDigits2RootFile(AliL3DigitRowData *rowPt,Char_t *new_digitsfile)
374 {
375   //Write the data stored in rowPt, into a new AliROOT file.
376   //The data is stored in the AliROOT format 
377   //This is specially a nice thing if you have modified data, and wants to run it  
378   //through the offline reconstruction chain.
379   //The arguments is a pointer to the data, and the name of the new AliROOT file.
380   //Remember to pass the original AliROOT file (the one that contains the original
381   //simulated data) to this object, in order to retrieve the MC id's of the digits.
382   
383   if(!fInAli)
384     {
385       printf("AliL3FileHandler::AliDigits2RootFile : No rootfile\n");
386       return;
387     }
388   if(!fParam)
389     {
390       printf("AliL3FileHandler::AliDigits2RootFile : No parameter object. Run on rootfile\n");
391       return;
392     }
393   if(!fTransformer)
394     {
395       printf("AliL3FileHandler::AliDigits2RootFile : No transform object\n");
396       return;
397     }
398   
399   //Get the original digitstree:
400   fInAli->cd();
401   AliTPCDigitsArray *old_array = new AliTPCDigitsArray();
402   old_array->Setup(fParam);
403   old_array->SetClass("AliSimDigits");
404   Bool_t ok = old_array->ConnectTree("TreeD_75x40_100x60_0");
405   if(!ok)
406     {
407       printf("AliL3FileHandler::AliDigits2RootFile : No digits tree object\n");
408       return;
409     }
410   
411   Bool_t create=kFALSE;
412   TFile *digFile;
413   
414   digFile = TFile::Open(new_digitsfile,"NEW");
415   if(digFile->IsOpen())
416     {    
417       create = kTRUE;
418       fParam->Write(fParam->GetTitle());
419     }
420   else
421     {
422       LOG(AliL3Log::kDebug,"AliL3FileHandler::AliDigits2RootFile","Rootfile")
423         <<"Rootfile did already exist, so I will just open it for updates"<<ENDLOG;
424       digFile = TFile::Open(new_digitsfile,"UPDATE");
425       create=kFALSE;
426     }
427   if(!digFile->IsOpen())
428     {
429       LOG(AliL3Log::kError,"AliL3FileHandler::AliDigits2RootFile","Rootfile")
430         <<"Error opening rootfile "<<new_digitsfile<<ENDLOG;
431       return;
432     }
433   
434   digFile->cd();
435     
436   //setup a new one, or connect it to the existing one:
437   AliTPCDigitsArray *arr = new AliTPCDigitsArray; 
438   arr->SetClass("AliSimDigits");
439   arr->Setup(fParam);
440   if(create)
441     arr->MakeTree();
442   else
443     {
444       Bool_t ok = arr->ConnectTree("TreeD_75x40_100x60_0");
445       if(!ok)
446         {
447           printf("AliL3FileHandler::AliDigits2RootFile : No digits tree object in existing file\n");
448           return;
449         }
450     }
451   Int_t digcounter=0;
452
453   for(Int_t i=fRowMin; i<=fRowMax; i++)
454     {
455       
456       if((Int_t)rowPt->fRow != i) printf("AliL3FileHandler::AliDigits2RootFile : Mismatching row numbering!!!\n");
457             
458       Int_t sector,row;
459       fTransformer->Slice2Sector(fSlice,i,sector,row);
460       AliSimDigits * dig = (AliSimDigits*)arr->CreateRow(sector,row);
461       AliSimDigits *old_dig = (AliSimDigits*)old_array->LoadRow(sector,row);
462       if(!old_dig)
463         printf("AliL3FileHandler::AliDigits2RootFile : No padrow %d %d\n",sector,row);
464       
465       AliL3DigitData *digPt = rowPt->fDigitData;
466       digcounter=0;
467       for(UInt_t j=0; j<rowPt->fNDigit; j++)
468         {
469           Int_t charge = (Int_t)digPt[j].fCharge;
470           Int_t pad = (Int_t)digPt[j].fPad;
471           Int_t time = (Int_t)digPt[j].fTime;
472           
473           if(charge == 0) //Only write the digits that has not been removed
474             continue;
475           digcounter++;
476           dig->SetDigitFast(charge,time,pad);
477           
478           Int_t trackID[3] = {old_dig->GetTrackID(time,pad,0),old_dig->GetTrackID(time,pad,1),old_dig->GetTrackID(time,pad,2)};
479           Int_t s_pad = pad;
480           Int_t s_time = time - 1;
481           while(trackID[0] < 0)
482             {
483               if(s_time >= 0 && s_time < fTransformer->GetNTimeBins() && s_pad >= 0 && s_pad < fTransformer->GetNPads(i))
484                 {
485                   if(old_dig->GetTrackID(s_time,s_pad,0) > 0)
486                     {
487                       trackID[0]=old_dig->GetTrackID(s_time,s_pad,0); 
488                       trackID[1]=old_dig->GetTrackID(s_time,s_pad,1); 
489                       trackID[2]=old_dig->GetTrackID(s_time,s_pad,2); 
490                     }
491                 }
492               if(s_pad == pad && s_time == time - 1)
493                 s_time = time + 1;
494               else if(s_pad == pad && s_time == time + 1)
495                 {s_pad = pad - 1; s_time = time;}
496               else if(s_pad == pad - 1 && s_time == time)
497                 s_time = time - 1;
498               else if(s_pad == pad - 1 && s_time == time - 1)
499                 s_time = time + 1;
500               else if(s_pad == pad - 1 && s_time == time + 1)
501                 {s_pad = pad + 1; s_time = time;}
502               else if(s_pad == pad + 1 && s_time == time)
503                 s_time = time - 1;
504               else if(s_pad == pad + 1 && s_time == time - 1)
505                 s_time = time + 1;
506               else 
507                 break;
508             }
509           
510           dig->SetTrackIDFast(trackID[0],time,pad,0);
511           dig->SetTrackIDFast(trackID[1],time,pad,1);
512           dig->SetTrackIDFast(trackID[2],time,pad,2);
513           
514         }
515       //cout<<"Wrote "<<digcounter<<" on row "<<i<<endl;
516       UpdateRowPointer(rowPt);
517       arr->StoreRow(sector,row);
518       arr->ClearRow(sector,row);  
519       old_array->ClearRow(sector,row);
520     }
521   digFile->cd();
522   char treeName[100];
523   sprintf(treeName,"TreeD_%s_0",fParam->GetTitle());
524   printf("Writing tree to file.....");
525   arr->GetTree()->Write(treeName,TObject::kOverwrite);
526   printf("done\n");
527   digFile->Close();
528   //arr->GetTree()->Delete();
529   //delete arr;
530 }
531
532 ///////////////////////////////////////// Point IO  
533 Bool_t AliL3FileHandler::AliPoints2Binary(){
534   Bool_t out = kTRUE;
535   UInt_t npoint;
536   AliL3SpacePointData *data = AliPoints2Memory(npoint);
537   out = Memory2Binary(npoint,data);
538   Free();
539   return out;
540 }
541
542 AliL3SpacePointData * AliL3FileHandler::AliPoints2Memory(UInt_t & npoint){
543   AliL3SpacePointData *data = 0;
544   npoint=0;
545   if(!fInAli){
546     LOG(AliL3Log::kWarning,"AliL3FileHandler::AliPoints2Memory","File")
547     <<"No Input avalible: no object TFile"<<ENDLOG;
548     return 0;
549   }
550   if(!fInAli->IsOpen()){
551     LOG(AliL3Log::kWarning,"AliL3FileHandler::AliPoints2Memory","File")
552     <<"No Input avalible: TFile not opend"<<ENDLOG;
553     return 0;
554   }
555   if(!fTransformer)
556     {
557       LOG(AliL3Log::kWarning,"AliL3FileHandler::AliPoints2Memory","Transformer")
558         <<"No transformer object"<<ENDLOG;
559       return 0;
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     fTransformer->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     fTransformer->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