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