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