]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/src/AliL3FileHandler.cxx
Load digitstree only once
[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 #include <TObject.h>
9 #include <TFile.h>
10 #include <TTree.h>
11 #include <stdio.h>
12
13 #include "AliL3Transform.h"
14 #include "AliL3Logging.h"
15 #include "AliL3MemHandler.h"
16 #include "AliL3FileHandler.h"
17
18 #include "AliTPCDigitsArray.h"
19 #include "AliTPCClustersArray.h"
20 #include "AliTPCcluster.h"
21 #include "AliTPCClustersRow.h"
22 #include "AliTPCParam.h"
23 #include "AliSimDigits.h"
24
25 #include "AliL3DigitData.h"
26 #include "AliL3TrackSegmentData.h"
27 #include "AliL3SpacePointData.h"
28 #include "AliL3TrackArray.h"
29
30 //_____________________________________________________________
31 // AliL3FileHandler
32 //
33 // The L3 Binary File handler 
34 //
35
36 ClassImp(AliL3FileHandler)
37
38 AliL3FileHandler::AliL3FileHandler(){
39   //Default constructor
40   fInAli = 0;
41   fParam = 0;
42   fTransformer = 0;
43   fMC =0;
44   fLastIndex=0;
45   fDigits=0;
46   fDigitsTree=0;
47 }
48
49
50 AliL3FileHandler::~AliL3FileHandler(){
51   //Destructor
52   if(fTransformer) delete fTransformer;
53   if(fMC) CloseMCOutput();
54   if(fDigitsTree) delete fDigitsTree;
55   if(fInAli) CloseAliInput();
56   
57 }
58
59 Bool_t AliL3FileHandler::SetMCOutput(char *name){
60   fMC = fopen(name,"w");
61   if(!fMC){
62     LOG(AliL3Log::kWarning,"AliL3FileHandler::SetMCOutput","File Open")
63     <<"Pointer to File = 0x0 "<<ENDLOG;
64     return kFALSE;
65   }
66   return kTRUE;
67 }
68
69 Bool_t AliL3FileHandler::SetMCOutput(FILE *file){
70   fMC = file;
71   if(!fMC){
72     LOG(AliL3Log::kWarning,"AliL3FileHandler::SetMCOutput","File Open")
73     <<"Pointer to File = 0x0 "<<ENDLOG;
74     return kFALSE;
75   }
76   return kTRUE;
77 }
78
79 void AliL3FileHandler::CloseMCOutput(){
80   if(!fMC){
81     LOG(AliL3Log::kWarning,"AliL3FileHandler::CloseMCOutPut","File Close")
82     <<"Nothing to Close"<<ENDLOG;
83     return;
84   }
85   fclose(fMC);
86   fMC =0;
87 }
88
89 Bool_t AliL3FileHandler::SetAliInput(){
90   if(!fInAli->IsOpen()){
91     LOG(AliL3Log::kError,"AliL3FileHandler::SetAliInput","File Open")
92     <<"Ali File "<<fInAli->GetName()<<" does not exist"<<ENDLOG;
93     return kFALSE;
94   }
95   fParam = (AliTPCParam*)fInAli->Get("75x40_100x60");
96   if(!fParam){ 
97     LOG(AliL3Log::kError,"AliL3FileHandler::SetAliInput","File Open")
98     <<"No AliTPCParam 75x40_100x60 in File "<<fInAli->GetName()<<ENDLOG;
99      return kFALSE;
100   }
101   fTransformer = new AliL3Transform();
102   return kTRUE;
103 }
104
105 Bool_t AliL3FileHandler::SetAliInput(char *name){
106   fInAli= new TFile(name,"READ");
107   if(!fInAli){
108     LOG(AliL3Log::kWarning,"AliL3FileHandler::SetAliInput","File Open")
109     <<"Pointer to TFile = 0x0 "<<ENDLOG;
110     return kFALSE;
111   }
112   return SetAliInput();
113 }
114
115 Bool_t AliL3FileHandler::SetAliInput(TFile *file){
116   fInAli=file;
117   if(!fInAli){
118     LOG(AliL3Log::kWarning,"AliL3FileHandler::SetAliInput","File Open")
119     <<"Pointer to TFile = 0x0 "<<ENDLOG;
120     return kFALSE;
121   }
122   return SetAliInput();
123 }
124
125 void AliL3FileHandler::CloseAliInput(){
126   if(!fInAli){
127     LOG(AliL3Log::kWarning,"AliL3FileHandler::CloseAliInput","File Close")
128     <<"Nothing to Close"<<ENDLOG;
129      return;
130   }
131   if(fInAli->IsOpen()) fInAli->Close();
132   delete fInAli;
133   fInAli = 0;
134     
135 }
136
137 Bool_t AliL3FileHandler::IsDigit(){
138   if(!fInAli){
139     LOG(AliL3Log::kWarning,"AliL3FileHandler::IsDigit","File")
140     <<"Pointer to TFile = 0x0 "<<ENDLOG;
141     return kTRUE;  //may you are use binary input which is Digits!!
142   }
143   TTree *t=(TTree*)fInAli->Get("TreeD_75x40_100x60_0");
144   if(t){
145     LOG(AliL3Log::kInformational,"AliL3FileHandler::IsDigit","File Type")
146     <<"Found Digit Tree -> Use Fast Cluster Finder"<<ENDLOG;
147     return kTRUE;
148   }
149   else{
150     LOG(AliL3Log::kInformational,"AliL3FileHandler::IsDigit","File Type")
151     <<"No Digit Tree -> Use Cluster Tree"<<ENDLOG;
152     return kFALSE;
153   }
154 }
155
156 ///////////////////////////////////////// Digit IO  
157 Bool_t AliL3FileHandler::AliDigits2Binary(Int_t event){
158   Bool_t out = kTRUE;
159   UInt_t nrow;
160   AliL3DigitRowData* data = AliDigits2Memory(nrow,event);
161   out = Memory2Binary(nrow,data);
162   Free();
163   return out;
164 }
165
166
167 Bool_t AliL3FileHandler::AliDigits2CompBinary(Int_t event){
168   Bool_t out = kTRUE;
169   UInt_t ndigits=0;
170   AliL3DigitRowData* digits = AliDigits2Memory(ndigits,event);
171   out = Memory2CompBinary(ndigits,digits);
172   Free();
173   return out;
174 }
175
176
177 AliL3DigitRowData * AliL3FileHandler::AliDigits2Memory(UInt_t & nrow,Int_t event){
178   AliL3DigitRowData *data = 0;
179   nrow=0;
180   if(!fInAli){
181     LOG(AliL3Log::kWarning,"AliL3FileHandler::AliDigits2Memory","File")
182     <<"No Input avalible: no object TFile"<<ENDLOG;
183     return 0; 
184   }
185   if(!fInAli->IsOpen()){
186     LOG(AliL3Log::kWarning,"AliL3FileHandler::AliDigits2Memory","File")
187     <<"No Input avalible: TFile not opend"<<ENDLOG;
188     return 0;
189   }
190
191   if(!fDigitsTree)
192     GetDigitsTree(event);
193     
194   UShort_t dig;
195   Int_t time,pad,sector,row;
196   Int_t nrows=0;
197   Int_t ndigitcount=0;
198   Int_t entries = (Int_t)fDigitsTree->GetEntries();
199   Int_t ndigits[entries];
200   Int_t lslice,lrow;
201   
202   for(Int_t n=fLastIndex; n<fDigitsTree->GetEntries(); n++)
203     {
204       fDigitsTree->GetEvent(n);
205       fParam->AdjustSectorRow(fDigits->GetID(),sector,row);
206       fTransformer->Sector2Slice(lslice,lrow,sector,row);
207       //if(fSlice != lslice || lrow<fRowMin || lrow>fRowMax) continue;
208       if(lslice < fSlice) continue;
209       if(lslice != fSlice) break;
210       if(lrow < fRowMin) continue;
211       if(lrow > fRowMax) break;
212
213       Float_t xyz[3];
214       ndigits[lrow] = 0;
215       fDigits->First();
216       do {
217         time=fDigits->CurrentRow();
218         pad=fDigits->CurrentColumn();
219         dig = fDigits->GetDigit(time,pad);
220         if(dig<=fParam->GetZeroSup()) continue;
221         if(time < fParam->GetMaxTBin()-1 && time > 0)
222           if(fDigits->GetDigit(time+1,pad) <= fParam->GetZeroSup()
223              && fDigits->GetDigit(time-1,pad) <= fParam->GetZeroSup())
224             continue;
225
226         fTransformer->Raw2Local(xyz,sector,row,pad,time);
227         if(fParam->GetPadRowRadii(sector,row)<230./250.*fabs(xyz[2]))
228           continue;
229
230         ndigits[lrow]++; //for this row only
231         ndigitcount++;  //total number of digits to be published
232
233       } while (fDigits->Next());
234
235       nrows++;
236     }
237   Int_t size = sizeof(AliL3DigitData)*ndigitcount
238     + nrows*sizeof(AliL3DigitRowData);
239
240   LOG(AliL3Log::kDebug,"AliL3FileHandler::AliDigits2Memory","Digits")
241     <<AliL3Log::kDec<<"Found "<<ndigitcount<<" Digits"<<ENDLOG;
242   
243   data=(AliL3DigitRowData*) Allocate(size);
244   nrow = (UInt_t)nrows;
245   AliL3DigitRowData *tempPt = data;
246   for(Int_t n=fLastIndex; n<fDigitsTree->GetEntries(); n++)
247     {
248       fDigitsTree->GetEvent(n);
249       Float_t xyz[3];
250       fParam->AdjustSectorRow(fDigits->GetID(),sector,row);
251       fTransformer->Sector2Slice(lslice,lrow,sector,row);
252       //if(fSlice != lslice || lrow<fRowMin || lrow>fRowMax) continue;
253       if(lslice < fSlice) continue;
254       if(lslice != fSlice) break;
255       if(lrow < fRowMin) continue;
256       if(lrow > fRowMax) break;
257
258       tempPt->fRow = lrow;
259       tempPt->fNDigit = ndigits[lrow];
260
261       Int_t localcount=0;
262       fDigits->First();
263       do {
264         //dig=fDigits->CurrentDigit();
265         //if (dig<=fParam->GetZeroSup()) continue;
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()) continue;
273
274         //Exclude data outside cone:
275         fTransformer->Raw2Local(xyz,sector,row,pad,time);
276         if(fParam->GetPadRowRadii(sector,row)<230./250.*fabs(xyz[2]))
277           continue;
278
279         if(localcount >= ndigits[lrow])
280           LOG(AliL3Log::kFatal,"AliL3FileHandler::AliDigits2Binary","Memory")
281             <<AliL3Log::kDec<<"Mismatch: localcount "<<localcount<<" ndigits "
282             <<ndigits[lrow]<<ENDLOG;
283         
284         tempPt->fDigitData[localcount].fCharge=dig;
285         tempPt->fDigitData[localcount].fPad=pad;
286         tempPt->fDigitData[localcount].fTime=time;
287         localcount++;
288       } while (fDigits->Next());
289
290       Byte_t *tmp = (Byte_t*)tempPt;
291       Int_t size = sizeof(AliL3DigitRowData)
292                                       + ndigits[lrow]*sizeof(AliL3DigitData);
293       tmp += size;
294       tempPt = (AliL3DigitRowData*)tmp;
295       //fLastIndex=n;
296     }
297   //fLastIndex++;
298   return data;
299 }
300
301 Bool_t AliL3FileHandler::GetDigitsTree(Int_t event)
302 {
303   
304   fInAli->cd();
305   Char_t dname[100];
306   sprintf(dname,"TreeD_75x40_100x60_%d",event);
307   fDigitsTree = (TTree*)fInAli->Get("TreeD_75x40_100x60_0");
308   if(!fDigitsTree) 
309     {
310       LOG(AliL3Log::kError,"AliL3FileHandler::GetDigitsTree","Digits Tree")
311         <<AliL3Log::kHex<<"Error getting digitstree "<<(Int_t)fDigitsTree<<ENDLOG;
312       return kFALSE;
313     }
314   fDigitsTree->GetBranch("Segment")->SetAddress(&fDigits);
315   return kTRUE;
316 }
317
318 void AliL3FileHandler::AliDigits2RootFile(AliL3DigitRowData *rowPt,Char_t *new_digitsfile)
319 {
320   //Write digits to a new alirootfile.
321
322   if(!fInAli)
323     {
324       printf("AliL3FileHandler::AliDigits2RootFile : No rootfile\n");
325       return;
326     }
327   if(!fParam)
328     {
329       printf("AliL3FileHandler::AliDigits2RootFile : No parameter object. Run on rootfile\n");
330       return;
331     }
332   if(!fTransformer)
333     {
334       printf("AliL3FileHandler::AliDigits2RootFile : No transform object\n");
335       return;
336     }
337   
338   //Get the original digitstree:
339   fInAli->cd();
340   AliTPCDigitsArray *old_array = new AliTPCDigitsArray();
341   old_array->Setup(fParam);
342   old_array->SetClass("AliSimDigits");
343   Bool_t ok = old_array->ConnectTree("TreeD_75x40_100x60_0");
344   if(!ok)
345     {
346       printf("AliL3FileHandler::AliDigits2RootFile : No digits tree object\n");
347       return;
348     }
349   
350   Bool_t create=kFALSE;
351   TFile *digFile;
352   
353   digFile = TFile::Open(new_digitsfile,"NEW");
354   if(digFile->IsOpen())
355     {    
356       create = kTRUE;
357       fParam->Write(fParam->GetTitle());
358     }
359   else
360     {
361       LOG(AliL3Log::kDebug,"AliL3FileHandler::AliDigits2RootFile","Rootfile")
362         <<"Rootfile did already exist, so I will just open it for updates"<<ENDLOG;
363       digFile = TFile::Open(new_digitsfile,"UPDATE");
364       create=kFALSE;
365     }
366   if(!digFile->IsOpen())
367     {
368       LOG(AliL3Log::kError,"AliL3FileHandler::AliDigits2RootFile","Rootfile")
369         <<"Error opening rootfile "<<new_digitsfile<<ENDLOG;
370       return;
371     }
372   
373   digFile->cd();
374     
375   //setup a new one, or connect it to the existing one:
376   AliTPCDigitsArray *arr = new AliTPCDigitsArray; 
377   arr->SetClass("AliSimDigits");
378   arr->Setup(fParam);
379   if(create)
380     arr->MakeTree();
381   else
382     {
383       Bool_t ok = arr->ConnectTree("TreeD_75x40_100x60_0");
384       if(!ok)
385         {
386           printf("AliL3FileHandler::AliDigits2RootFile : No digits tree object in existing file\n");
387           return;
388         }
389     }
390   for(Int_t i=fRowMin; i<=fRowMax; i++)
391     {
392       
393       if(rowPt->fRow != i) printf("AliL3FileHandler::AliDigits2RootFile : Mismatching row numbering!!!\n");
394             
395       Int_t sector,row;
396       fTransformer->Slice2Sector(fSlice,i,sector,row);
397       AliDigits * dig = arr->CreateRow(sector,row);
398       AliDigits *old_dig = old_array->LoadRow(sector,row);
399       if(!old_dig)
400         printf("AliL3FileHandler::AliDigits2RootFile : No padrow %d %d\n",sector,row);
401       
402       AliL3DigitData *digPt = rowPt->fDigitData;
403       for(UInt_t j=0; j<rowPt->fNDigit; j++)
404         {
405           UShort_t charge = digPt[j].fCharge;
406           UChar_t pad = digPt[j].fPad;
407           UShort_t time = digPt[j].fTime;
408           
409           if(charge == 0) //Only write the digits that has not been removed
410             continue;
411           dig->SetDigitFast(old_dig->GetDigit(time,pad),time,pad);
412           ((AliSimDigits*)dig)->SetTrackIDFast(((AliSimDigits*)old_dig)->GetTrackID((Int_t)time,(Int_t)pad,0),time,pad,0);
413           ((AliSimDigits*)dig)->SetTrackIDFast(((AliSimDigits*)old_dig)->GetTrackID((Int_t)time,(Int_t)pad,1),time,pad,1);
414           ((AliSimDigits*)dig)->SetTrackIDFast(((AliSimDigits*)old_dig)->GetTrackID((Int_t)time,(Int_t)pad,2),time,pad,2);
415           
416         }
417       UpdateRowPointer(rowPt);
418       arr->StoreRow(sector,row);
419       arr->ClearRow(sector,row);  
420       old_array->ClearRow(sector,row);
421     }
422   digFile->cd();
423   char treeName[100];
424   sprintf(treeName,"TreeD_%s_0",fParam->GetTitle());
425   printf("Writing tree to file.....");
426   arr->GetTree()->Write(treeName,TObject::kOverwrite);
427   printf("done\n");
428   digFile->Close();
429   //arr->GetTree()->Delete();
430   //delete arr;
431 }
432
433 ///////////////////////////////////////// Point IO  
434 Bool_t AliL3FileHandler::AliPoints2Binary(){
435   Bool_t out = kTRUE;
436   UInt_t npoint;
437   AliL3SpacePointData *data = AliPoints2Memory(npoint);
438   out = Memory2Binary(npoint,data);
439   Free();
440   return out;
441 }
442
443 AliL3SpacePointData * AliL3FileHandler::AliPoints2Memory(UInt_t & npoint){
444   AliL3SpacePointData *data = 0;
445   npoint=0;
446   if(!fInAli){
447     LOG(AliL3Log::kWarning,"AliL3FileHandler::AliPoints2Memory","File")
448     <<"No Input avalible: no object TFile"<<ENDLOG;
449     return 0;
450   }
451   if(!fInAli->IsOpen()){
452     LOG(AliL3Log::kWarning,"AliL3FileHandler::AliPoints2Memory","File")
453     <<"No Input avalible: TFile not opend"<<ENDLOG;
454     return 0;
455   }
456   TDirectory *savedir = gDirectory;
457   fInAli->cd();
458   
459   Char_t cname[100];
460   Int_t eventn = 0;
461   sprintf(cname,"TreeC_TPC_%d",eventn);
462   AliTPCClustersArray carray;
463   carray.Setup(fParam);
464   carray.SetClusterType("AliTPCcluster");
465   Bool_t clusterok = carray.ConnectTree(cname);
466   if(!clusterok) return 0;
467
468   AliTPCClustersRow ** clusterrow = 
469                new AliTPCClustersRow*[ (int)carray.GetTree()->GetEntries()];
470   Int_t *rows = new int[ (int)carray.GetTree()->GetEntries()];
471   Int_t *sects = new int[  (int)carray.GetTree()->GetEntries()];
472   Int_t sum=0;
473
474   Int_t lslice,lrow;
475   for(Int_t i=0; i<carray.GetTree()->GetEntries(); i++){
476     AliSegmentID *s = carray.LoadEntry(i);
477     Int_t sector,row;
478     fParam->AdjustSectorRow(s->GetID(),sector,row);
479     rows[i] = row;
480     sects[i] = sector;
481     clusterrow[i] = 0;
482     fTransformer->Sector2Slice(lslice,lrow,sector,row);
483     if(fSlice != lslice || lrow<fRowMin || lrow>fRowMax) continue;
484     clusterrow[i] = carray.GetRow(sector,row);
485     if(clusterrow[i])
486       sum+=clusterrow[i]->GetArray()->GetEntriesFast();
487   }
488   UInt_t size = sum*sizeof(AliL3SpacePointData);
489
490   LOG(AliL3Log::kDebug,"AliL3FileHandler::AliPoints2Memory","File")
491   <<AliL3Log::kDec<<"Found "<<sum<<" SpacePoints"<<ENDLOG;
492
493   data = (AliL3SpacePointData *) Allocate(size);
494   npoint = sum;
495   UInt_t n=0; 
496   for(Int_t i=0; i<carray.GetTree()->GetEntries(); i++){
497     if(!clusterrow[i]) continue;
498     Int_t row = rows[i];
499     Int_t sector = sects[i];
500     fTransformer->Sector2Slice(lslice,lrow,sector,row);
501     Int_t entries_in_row = clusterrow[i]->GetArray()->GetEntriesFast();
502     for(Int_t j = 0;j<entries_in_row;j++){
503       AliTPCcluster *c = (AliTPCcluster*)(*clusterrow[i])[j];
504       data[n].fZ = c->GetZ();
505       data[n].fY = c->GetY();
506       data[n].fX = fParam->GetPadRowRadii(sector,row);
507       data[n].fID = n+((fSlice&0x7f)<<25)+((fPatch&0x7)<<22);//uli
508       data[n].fPadRow = lrow;
509       data[n].fXYErr = c->GetSigmaY2();
510       data[n].fZErr = c->GetSigmaZ2();
511       if(fMC) fprintf(fMC,"%d %d\n",data[n].fID,c->GetLabel(0));
512       n++;
513     }
514   }
515   for(Int_t i=0;i<carray.GetTree()->GetEntries();i++){
516     Int_t row = rows[i];
517     Int_t sector = sects[i];
518     if(carray.GetRow(sector,row))
519       carray.ClearRow(sector,row);
520   }
521
522   delete [] clusterrow;
523   delete [] rows;
524   delete [] sects;
525   savedir->cd();   
526
527   return data;
528 }
529