Major changes in the AliL3Transform class. The class has been made completely
[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   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   
234   if(!fDigitsTree)
235     GetDigitsTree(event);
236     
237   UShort_t dig;
238   Int_t time,pad,sector,row;
239   Int_t nrows=0;
240   Int_t ndigitcount=0;
241   Int_t entries = (Int_t)fDigitsTree->GetEntries();
242   Int_t ndigits[entries];
243   Int_t lslice,lrow;
244   
245   for(Int_t n=fLastIndex; n<fDigitsTree->GetEntries(); n++)
246     {
247       fDigitsTree->GetEvent(n);
248       fParam->AdjustSectorRow(fDigits->GetID(),sector,row);
249       AliL3Transform::Sector2Slice(lslice,lrow,sector,row);
250       //if(fSlice != lslice || lrow<fRowMin || lrow>fRowMax) continue;
251       if(lslice < fSlice) continue;
252       if(lslice != fSlice) break;
253       if(lrow < fRowMin) continue;
254       if(lrow > fRowMax) break;
255
256       Float_t xyz[3];
257       ndigits[lrow] = 0;
258       fDigits->First();
259       do {
260         time=fDigits->CurrentRow();
261         pad=fDigits->CurrentColumn();
262         dig = fDigits->GetDigit(time,pad);
263         if(dig<=fParam->GetZeroSup()) continue;
264         if(time < fParam->GetMaxTBin()-1 && time > 0)
265           if(fDigits->GetDigit(time+1,pad) <= fParam->GetZeroSup()
266              && fDigits->GetDigit(time-1,pad) <= fParam->GetZeroSup())
267             continue;
268
269         AliL3Transform::Raw2Local(xyz,sector,row,pad,time);
270         if(fParam->GetPadRowRadii(sector,row)<230./250.*fabs(xyz[2]))
271           continue;
272
273         ndigits[lrow]++; //for this row only
274         ndigitcount++;  //total number of digits to be published
275
276       } while (fDigits->Next());
277
278       nrows++;
279     }
280   Int_t size = sizeof(AliL3DigitData)*ndigitcount
281     + nrows*sizeof(AliL3DigitRowData);
282
283   LOG(AliL3Log::kDebug,"AliL3FileHandler::AliDigits2Memory","Digits")
284     <<AliL3Log::kDec<<"Found "<<ndigitcount<<" Digits"<<ENDLOG;
285   
286   data=(AliL3DigitRowData*) Allocate(size);
287   nrow = (UInt_t)nrows;
288   AliL3DigitRowData *tempPt = data;
289   for(Int_t n=fLastIndex; n<fDigitsTree->GetEntries(); n++)
290     {
291       fDigitsTree->GetEvent(n);
292       Float_t xyz[3];
293       fParam->AdjustSectorRow(fDigits->GetID(),sector,row);
294       AliL3Transform::Sector2Slice(lslice,lrow,sector,row);
295       //if(fSlice != lslice || lrow<fRowMin || lrow>fRowMax) continue;
296       if(lslice < fSlice) continue;
297       if(lslice != fSlice) break;
298       if(lrow < fRowMin) continue;
299       if(lrow > fRowMax) break;
300
301       tempPt->fRow = lrow;
302       tempPt->fNDigit = ndigits[lrow];
303
304       Int_t localcount=0;
305       fDigits->First();
306       do {
307         //dig=fDigits->CurrentDigit();
308         //if (dig<=fParam->GetZeroSup()) continue;
309         time=fDigits->CurrentRow();
310         pad=fDigits->CurrentColumn();
311         dig = fDigits->GetDigit(time,pad);
312         if (dig <= fParam->GetZeroSup()) continue;
313         if(time < fParam->GetMaxTBin()-1 && time > 0)
314           if(fDigits->GetDigit(time-1,pad) <= fParam->GetZeroSup() &&
315              fDigits->GetDigit(time+1,pad) <= fParam->GetZeroSup()) continue;
316
317         //Exclude data outside cone:
318         AliL3Transform::Raw2Local(xyz,sector,row,pad,time);
319         if(fParam->GetPadRowRadii(sector,row)<230./250.*fabs(xyz[2]))
320           continue;
321
322         if(localcount >= ndigits[lrow])
323           LOG(AliL3Log::kFatal,"AliL3FileHandler::AliDigits2Binary","Memory")
324             <<AliL3Log::kDec<<"Mismatch: localcount "<<localcount<<" ndigits "
325             <<ndigits[lrow]<<ENDLOG;
326         
327         tempPt->fDigitData[localcount].fCharge=dig;
328         tempPt->fDigitData[localcount].fPad=pad;
329         tempPt->fDigitData[localcount].fTime=time;
330 #ifdef do_mc
331         tempPt->fDigitData[localcount].fTrackID[0] = fDigits->GetTrackID(time,pad,0);
332         tempPt->fDigitData[localcount].fTrackID[1] = fDigits->GetTrackID(time,pad,1);
333         tempPt->fDigitData[localcount].fTrackID[2] = fDigits->GetTrackID(time,pad,2);
334 #endif
335         localcount++;
336       } while (fDigits->Next());
337
338       Byte_t *tmp = (Byte_t*)tempPt;
339       Int_t size = sizeof(AliL3DigitRowData)
340                                       + ndigits[lrow]*sizeof(AliL3DigitData);
341       tmp += size;
342       tempPt = (AliL3DigitRowData*)tmp;
343       //fLastIndex=n;
344     }
345   //fLastIndex++;
346   return data;
347 }
348
349 Bool_t AliL3FileHandler::GetDigitsTree(Int_t event)
350 {
351   //Connects to the TPC digit tree in the AliROOT file.
352   
353   fInAli->cd();
354   Char_t dname[100];
355   sprintf(dname,"TreeD_75x40_100x60_%d",event);
356   fDigitsTree = (TTree*)fInAli->Get("TreeD_75x40_100x60_0");
357   if(!fDigitsTree) 
358     {
359       LOG(AliL3Log::kError,"AliL3FileHandler::GetDigitsTree","Digits Tree")
360         <<AliL3Log::kHex<<"Error getting digitstree "<<(Int_t)fDigitsTree<<ENDLOG;
361       return kFALSE;
362     }
363   fDigitsTree->GetBranch("Segment")->SetAddress(&fDigits);
364   return kTRUE;
365 }
366
367 void AliL3FileHandler::AliDigits2RootFile(AliL3DigitRowData *rowPt,Char_t *new_digitsfile)
368 {
369   //Write the data stored in rowPt, into a new AliROOT file.
370   //The data is stored in the AliROOT format 
371   //This is specially a nice thing if you have modified data, and wants to run it  
372   //through the offline reconstruction chain.
373   //The arguments is a pointer to the data, and the name of the new AliROOT file.
374   //Remember to pass the original AliROOT file (the one that contains the original
375   //simulated data) to this object, in order to retrieve the MC id's of the digits.
376   
377   if(!fInAli)
378     {
379       printf("AliL3FileHandler::AliDigits2RootFile : No rootfile\n");
380       return;
381     }
382   if(!fParam)
383     {
384       printf("AliL3FileHandler::AliDigits2RootFile : No parameter object. Run on rootfile\n");
385       return;
386     }
387   
388   //Get the original digitstree:
389   fInAli->cd();
390   AliTPCDigitsArray *old_array = new AliTPCDigitsArray();
391   old_array->Setup(fParam);
392   old_array->SetClass("AliSimDigits");
393   Bool_t ok = old_array->ConnectTree("TreeD_75x40_100x60_0");
394   if(!ok)
395     {
396       printf("AliL3FileHandler::AliDigits2RootFile : No digits tree object\n");
397       return;
398     }
399   
400   Bool_t create=kFALSE;
401   TFile *digFile;
402   
403   digFile = TFile::Open(new_digitsfile,"NEW");
404   if(digFile->IsOpen())
405     {    
406       create = kTRUE;
407       fParam->Write(fParam->GetTitle());
408     }
409   else
410     {
411       LOG(AliL3Log::kDebug,"AliL3FileHandler::AliDigits2RootFile","Rootfile")
412         <<"Rootfile did already exist, so I will just open it for updates"<<ENDLOG;
413       digFile = TFile::Open(new_digitsfile,"UPDATE");
414       create=kFALSE;
415     }
416   if(!digFile->IsOpen())
417     {
418       LOG(AliL3Log::kError,"AliL3FileHandler::AliDigits2RootFile","Rootfile")
419         <<"Error opening rootfile "<<new_digitsfile<<ENDLOG;
420       return;
421     }
422   
423   digFile->cd();
424     
425   //setup a new one, or connect it to the existing one:
426   AliTPCDigitsArray *arr = new AliTPCDigitsArray; 
427   arr->SetClass("AliSimDigits");
428   arr->Setup(fParam);
429   if(create)
430     arr->MakeTree();
431   else
432     {
433       Bool_t ok = arr->ConnectTree("TreeD_75x40_100x60_0");
434       if(!ok)
435         {
436           printf("AliL3FileHandler::AliDigits2RootFile : No digits tree object in existing file\n");
437           return;
438         }
439     }
440   Int_t digcounter=0;
441
442   for(Int_t i=fRowMin; i<=fRowMax; i++)
443     {
444       
445       if((Int_t)rowPt->fRow != i) printf("AliL3FileHandler::AliDigits2RootFile : Mismatching row numbering!!!\n");
446             
447       Int_t sector,row;
448       AliL3Transform::Slice2Sector(fSlice,i,sector,row);
449       AliSimDigits * dig = (AliSimDigits*)arr->CreateRow(sector,row);
450       AliSimDigits *old_dig = (AliSimDigits*)old_array->LoadRow(sector,row);
451       if(!old_dig)
452         printf("AliL3FileHandler::AliDigits2RootFile : No padrow %d %d\n",sector,row);
453       
454       AliL3DigitData *digPt = rowPt->fDigitData;
455       digcounter=0;
456       for(UInt_t j=0; j<rowPt->fNDigit; j++)
457         {
458           Int_t charge = (Int_t)digPt[j].fCharge;
459           Int_t pad = (Int_t)digPt[j].fPad;
460           Int_t time = (Int_t)digPt[j].fTime;
461           
462           if(charge == 0) //Only write the digits that has not been removed
463             continue;
464           digcounter++;
465           dig->SetDigitFast(charge,time,pad);
466           
467           Int_t trackID[3] = {old_dig->GetTrackID(time,pad,0),old_dig->GetTrackID(time,pad,1),old_dig->GetTrackID(time,pad,2)};
468           Int_t s_pad = pad;
469           Int_t s_time = time - 1;
470           while(trackID[0] < 0)
471             {
472               if(s_time >= 0 && s_time < AliL3Transform::GetNTimeBins() && s_pad >= 0 && s_pad < AliL3Transform::GetNPads(i))
473                 {
474                   if(old_dig->GetTrackID(s_time,s_pad,0) > 0)
475                     {
476                       trackID[0]=old_dig->GetTrackID(s_time,s_pad,0); 
477                       trackID[1]=old_dig->GetTrackID(s_time,s_pad,1); 
478                       trackID[2]=old_dig->GetTrackID(s_time,s_pad,2); 
479                     }
480                 }
481               if(s_pad == pad && s_time == time - 1)
482                 s_time = time + 1;
483               else if(s_pad == pad && s_time == time + 1)
484                 {s_pad = pad - 1; s_time = time;}
485               else if(s_pad == pad - 1 && s_time == time)
486                 s_time = time - 1;
487               else if(s_pad == pad - 1 && s_time == time - 1)
488                 s_time = time + 1;
489               else if(s_pad == pad - 1 && s_time == time + 1)
490                 {s_pad = pad + 1; s_time = time;}
491               else if(s_pad == pad + 1 && s_time == time)
492                 s_time = time - 1;
493               else if(s_pad == pad + 1 && s_time == time - 1)
494                 s_time = time + 1;
495               else 
496                 break;
497             }
498           
499           dig->SetTrackIDFast(trackID[0],time,pad,0);
500           dig->SetTrackIDFast(trackID[1],time,pad,1);
501           dig->SetTrackIDFast(trackID[2],time,pad,2);
502           
503         }
504       //cout<<"Wrote "<<digcounter<<" on row "<<i<<endl;
505       UpdateRowPointer(rowPt);
506       arr->StoreRow(sector,row);
507       arr->ClearRow(sector,row);  
508       old_array->ClearRow(sector,row);
509     }
510   digFile->cd();
511   char treeName[100];
512   sprintf(treeName,"TreeD_%s_0",fParam->GetTitle());
513   printf("Writing tree to file.....");
514   arr->GetTree()->Write(treeName,TObject::kOverwrite);
515   printf("done\n");
516   digFile->Close();
517   //arr->GetTree()->Delete();
518   //delete arr;
519 }
520
521 ///////////////////////////////////////// Point IO  
522 Bool_t AliL3FileHandler::AliPoints2Binary(){
523   Bool_t out = kTRUE;
524   UInt_t npoint;
525   AliL3SpacePointData *data = AliPoints2Memory(npoint);
526   out = Memory2Binary(npoint,data);
527   Free();
528   return out;
529 }
530
531 AliL3SpacePointData * AliL3FileHandler::AliPoints2Memory(UInt_t & npoint){
532   AliL3SpacePointData *data = 0;
533   npoint=0;
534   if(!fInAli){
535     LOG(AliL3Log::kWarning,"AliL3FileHandler::AliPoints2Memory","File")
536     <<"No Input avalible: no object TFile"<<ENDLOG;
537     return 0;
538   }
539   if(!fInAli->IsOpen()){
540     LOG(AliL3Log::kWarning,"AliL3FileHandler::AliPoints2Memory","File")
541     <<"No Input avalible: TFile not opend"<<ENDLOG;
542     return 0;
543   }
544
545   TDirectory *savedir = gDirectory;
546   fInAli->cd();
547   
548   Char_t cname[100];
549   Int_t eventn = 0;
550   sprintf(cname,"TreeC_TPC_%d",eventn);
551   AliTPCClustersArray carray;
552   carray.Setup(fParam);
553   carray.SetClusterType("AliTPCcluster");
554   Bool_t clusterok = carray.ConnectTree(cname);
555   if(!clusterok) return 0;
556
557   AliTPCClustersRow ** clusterrow = 
558                new AliTPCClustersRow*[ (int)carray.GetTree()->GetEntries()];
559   Int_t *rows = new int[ (int)carray.GetTree()->GetEntries()];
560   Int_t *sects = new int[  (int)carray.GetTree()->GetEntries()];
561   Int_t sum=0;
562
563   Int_t lslice,lrow;
564   for(Int_t i=0; i<carray.GetTree()->GetEntries(); i++){
565     AliSegmentID *s = carray.LoadEntry(i);
566     Int_t sector,row;
567     fParam->AdjustSectorRow(s->GetID(),sector,row);
568     rows[i] = row;
569     sects[i] = sector;
570     clusterrow[i] = 0;
571     AliL3Transform::Sector2Slice(lslice,lrow,sector,row);
572     if(fSlice != lslice || lrow<fRowMin || lrow>fRowMax) continue;
573     clusterrow[i] = carray.GetRow(sector,row);
574     if(clusterrow[i])
575       sum+=clusterrow[i]->GetArray()->GetEntriesFast();
576   }
577   UInt_t size = sum*sizeof(AliL3SpacePointData);
578
579   LOG(AliL3Log::kDebug,"AliL3FileHandler::AliPoints2Memory","File")
580   <<AliL3Log::kDec<<"Found "<<sum<<" SpacePoints"<<ENDLOG;
581
582   data = (AliL3SpacePointData *) Allocate(size);
583   npoint = sum;
584   UInt_t n=0; 
585   for(Int_t i=0; i<carray.GetTree()->GetEntries(); i++){
586     if(!clusterrow[i]) continue;
587     Int_t row = rows[i];
588     Int_t sector = sects[i];
589     AliL3Transform::Sector2Slice(lslice,lrow,sector,row);
590     Int_t entries_in_row = clusterrow[i]->GetArray()->GetEntriesFast();
591     for(Int_t j = 0;j<entries_in_row;j++){
592       AliTPCcluster *c = (AliTPCcluster*)(*clusterrow[i])[j];
593       data[n].fZ = c->GetZ();
594       data[n].fY = c->GetY();
595       data[n].fX = fParam->GetPadRowRadii(sector,row);
596       data[n].fID = n+((fSlice&0x7f)<<25)+((fPatch&0x7)<<22);//uli
597       data[n].fPadRow = lrow;
598       data[n].fXYErr = c->GetSigmaY2();
599       data[n].fZErr = c->GetSigmaZ2();
600       if(fMC) fprintf(fMC,"%d %d\n",data[n].fID,c->GetLabel(0));
601       n++;
602     }
603   }
604   for(Int_t i=0;i<carray.GetTree()->GetEntries();i++){
605     Int_t row = rows[i];
606     Int_t sector = sects[i];
607     if(carray.GetRow(sector,row))
608       carray.ClearRow(sector,row);
609   }
610
611   delete [] clusterrow;
612   delete [] rows;
613   delete [] sects;
614   savedir->cd();   
615
616   return data;
617 }
618