]> git.uio.no Git - u/mrichter/AliRoot.git/blob - CONTAINERS/AliArrayBranch.cxx
Removing extra namespace (alpha)
[u/mrichter/AliRoot.git] / CONTAINERS / AliArrayBranch.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 //-----------------------------------------------------
19 // A Branch for the case of an array of clone objects. 
20 //-----------------------------------------------------
21 // Author M.Ivanov
22 //*KEND.
23
24 #include "TROOT.h"
25 #include "AliArrayBranch.h"
26 #include "TFile.h"
27 #include "TTree.h" 
28 #include "TBasket.h"
29 #include "TClass.h"
30 #include "TRealData.h"
31 #include "TDataType.h"
32 #include "TDataMember.h"
33
34 #include "TBranch.h"
35 #include "TBranchClones.h"
36 #include "TLeaf.h"
37 #include "TLeafB.h"
38 #include "TLeafC.h"
39 #include "TLeafF.h"
40 #include "TLeafD.h"
41 #include "TLeafI.h"
42 #include "TLeafS.h"
43 #include "TLeafObject.h"
44
45 #include "AliObjectArray.h"
46 #include "AliDataType.h"
47
48
49 R__EXTERN TTree *gTree;
50
51 ClassImp(AliArraySubBranch)
52 ClassImp(AliArrayBranch)
53 ClassImp(AliObjectBranch)
54 ClassImp(AliTree)
55
56
57
58   Int_t AliArraySubBranch::GetEntryExport(Int_t entry, Int_t /*getall*/, AliObjectArray *list, Int_t nentries)
59 {
60 //*-*-*-*-*-*Read all leaves of entry and return total number of bytes*-*-*
61 //*-* export buffers to real objects in the AliObjectArray list.
62 //*-*
63
64    if (TestBit(kDoNotProcess)) return 0;
65    if (fReadEntry == entry) return 1;
66    if (entry < 0 || entry >= fEntryNumber) return 0;
67    Int_t nbytes;
68    Int_t first  = fBasketEntry[fReadBasket];
69    Int_t last;
70    if (fReadBasket == fWriteBasket) last = fEntryNumber - 1;
71    else                             last = fBasketEntry[fReadBasket+1] - 1;
72 //
73 //      Are we still in the same ReadBasket?
74    if (entry < first || entry > last) {
75       fReadBasket = TMath::BinarySearch(fWriteBasket+1, fBasketEntry, entry);
76       first       = fBasketEntry[fReadBasket];
77    }
78
79 //     We have found the basket containing this entry.
80 //     make sure basket buffers are in memory.
81    TBasket *basket = GetBasket(fReadBasket);
82    if (!basket) return 0;
83    TBuffer *buf    = basket->GetBufferRef();
84 //     Set entry offset in buffer and read data from all leaves
85    if (!buf->IsReading()) {
86       basket->SetReadMode();
87    }
88 //   Int_t bufbegin = basket->GetEntryPointer(entry-first);
89    Int_t bufbegin;
90    Int_t *entryOffset = basket->GetEntryOffset();
91    if (entryOffset) bufbegin = entryOffset[entry-first];
92    else             bufbegin = basket->GetKeylen() + (entry-first)*basket->GetNevBufSize();
93    buf->SetBufferOffset(bufbegin);
94
95    TLeaf *leaf = (TLeaf*)fLeaves.UncheckedAt(0);
96    //   leaf->ReadBasketExport(*buf,list,nentries);  //!!! MI
97    ReadBasketExport(*buf,leaf, list,nentries);
98    nbytes = buf->Length() - bufbegin;
99    fReadEntry = entry;
100
101    return nbytes;
102 }
103
104 void AliArraySubBranch::ReadBasketExport(TBuffer &b, TLeaf *leaf, AliObjectArray *list, Int_t n)
105 {
106   //
107   // 
108   Int_t len    = leaf->GetLenStatic();
109   Int_t offset = leaf->GetOffset();
110   void *value  = leaf->GetValuePointer();
111
112   //8bit integer
113   if (leaf->IsA()==TLeafB::Class()){   
114     Int_t j = 0;
115     for (Int_t i=0;i<n;i++) {
116       memcpy((char*)list->UncheckedAt(i) + offset,&((Char_t*)value)[j], len);
117       j += len;
118     } 
119   } 
120   //variable length string.
121   if (leaf->IsA()==TLeafC::Class()){  
122     UChar_t len;
123     b >> len;
124     if (len) {
125       if (len >= len) len = len-1;
126       b.ReadFastArray((Char_t*)value,len);
127       ((Char_t*)value)[len] = 0;
128     } else {
129       value = 0;
130     }    
131     Int_t j = 0;
132     for (Int_t i=0;i<n;i++) {
133       memcpy((char*)list->UncheckedAt(i) + offset,&((Char_t*)value)[j], 1);
134       j += len;
135     }
136   }
137   //double
138   if (leaf->IsA()==TLeafD::Class()){   
139     b.ReadFastArray(((Double_t*)value),n*len);    
140     Int_t j = 0;
141     for (Int_t i=0;i<n;i++) {
142       memcpy((char*)list->UncheckedAt(i) + offset,&((Double_t*)value)[j], 8*len);
143       j += len;
144    }
145   }
146   //float
147   if (leaf->IsA()==TLeafF::Class()){   
148     if (n*len == 1) {
149       b >> ((Float_t*)value)[0];
150     } else {
151       b.ReadFastArray(((Float_t*)value),n*len);
152     }
153     
154     Float_t *val = (Float_t*)value;
155     for (Int_t i=0;i<n;i++) {
156       char *first = (char*)list->UncheckedAt(i);
157       Float_t *ff = (Float_t*)&first[offset];
158       for (Int_t j=0;j<len;j++) {
159         ff[j] = val[j];
160       }
161       val += len;
162     }
163     return;
164   }
165   //int2
166   if (leaf->IsA()==TLeafS::Class()){       
167     if (n*len == 1) {
168       b >> ((Short_t*)value)[0];
169     } else {
170       b.ReadFastArray(((Short_t*)value),n*len);
171     }
172     Short_t *val = (Short_t*)value;
173     for (Int_t i=0;i<n;i++) {
174       char *first = (char*)list->UncheckedAt(i);
175       Short_t *ii = (Short_t*)&first[offset];
176       for (Int_t j=0;j<len;j++) {
177         ii[j] = val[j];
178       }
179       val += len;
180     }    
181     return;
182   }     
183   //int4
184   if (leaf->IsA()==TLeafI::Class()){       
185     if (n*len == 1) {
186       b >> ((Int_t*)value)[0];
187     } else {
188       b.ReadFastArray(((Int_t*)value),n*len);
189     }
190     Int_t *val = (Int_t*)value;
191     for (Int_t i=0;i<n;i++) {
192       char *first = (char*)list->UncheckedAt(i);
193       Int_t *ii = (Int_t*)&first[offset];
194       for (Int_t j=0;j<len;j++) {
195         ii[j] = val[j];
196       }
197       val += len;
198     }    
199     return;
200   }      
201 }
202
203
204 //______________________________________________________________________________
205 AliArrayBranch::AliArrayBranch(): TBranch()
206 {
207 //*-*-*-*-*-*Default constructor for BranchClones*-*-*-*-*-*-*-*-*-*
208 //*-*        ====================================
209
210    fList        = 0;
211    fRead        = 0;
212    fN           = 0;
213    fNdataMax    = 0;
214    fBranchCount = 0;
215 }
216
217
218 //______________________________________________________________________________
219 AliArrayBranch::AliArrayBranch(const Text_t *name, void *pointer, TTree * tree,  Int_t basketsize, Int_t compress)
220     :TBranch()
221 {
222 //*-*-*-*-*-*-*-*-*-*-*-*-*Create a BranchClones*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
223 //*-*                      =====================
224 //
225    char leaflist[80];
226    char branchname[80];
227    char branchcount[64];
228    fTree       = tree;
229    gTree       = tree; // MI because some bug in ROOT I didn't obtain proper gTree
230    // it is necesary to set gTree because oder subranchces  defined below need it
231    SetName(name);
232    if (compress == -1) {
233       TFile *bfile = fTree->GetDirectory()->GetFile();
234       if (bfile) compress = bfile->GetCompressionLevel();
235    }
236    char *cpointer  = (char*)pointer;
237    char **ppointer = (char**)(cpointer);
238    fList     = (AliObjectArray*)(*ppointer);
239    fAddress  = cpointer;
240    fRead     = 0;
241    fN        = 0;
242    fNdataMax = 0;
243
244    AliClassInfo *clinfo = fList->GetClassInfo();
245    if (!clinfo) return;
246    fClassName = clinfo->GetName();   
247    
248 //*-*- Create a branch to store the array count
249    if (basketsize < 100) basketsize = 100;
250    sprintf(leaflist,"%s_/I",name);
251    sprintf(branchcount,"%s_",name);
252    fBranchCount = new TBranch(branchcount,&fN,leaflist,basketsize);
253    fBranchCount->SetBit(kIsClone);
254    TLeaf *leafcount = (TLeaf*)fBranchCount->GetListOfLeaves()->UncheckedAt(0);
255
256 //*-*-  Create the first basket
257
258    fDirectory  = fTree->GetDirectory();
259    fFileName   = "";
260
261    TBasket *basket = new TBasket(branchcount,fTree->GetName(),this);
262    fBaskets.Add(basket);
263
264 //*-*- Loop on all public data members of the class and its base classes
265
266    TClass *cl = fList->GetClass();
267    if (cl){
268      if (!cl->GetListOfRealData())  cl->BuildRealData();
269      
270      const char *itype = 0;
271      TRealData *rd;
272      TIter      next(cl->GetListOfRealData());
273      while ((rd = (TRealData *) next())) {
274        TDataMember *member = rd->GetDataMember();
275        if (!member->IsBasic()) {
276          Warning("BranchClones","Cannot process member:%s",member->GetName());
277          continue;
278        }
279        if (!member->IsPersistent()) continue; //do not process members with a ! as the first
280        // character in the comment field
281        TDataType *membertype = member->GetDataType();
282        Int_t type = membertype->GetType();
283        if (type == 0) {
284          Warning("BranchClones","Cannot process member:%s",member->GetName());
285          continue;
286        }
287        if (type == 1)  itype = "B";
288        if (type == 11) itype = "b";
289        if (type == 3)  itype = "I";
290        if (type == 5)  itype = "F";
291        if (type == 8)  itype = "D";
292        if (type == 13) itype = "i";
293        if (type == 2)  itype = "S";
294        if (type == 12) itype = "s";
295        
296        
297        Int_t arraydim = member->GetArrayDim();
298        if (arraydim!=1){
299          //   OLD Version 
300          sprintf(leaflist,"%s[%s]/%s",member->GetName(),branchcount,itype);
301          Int_t comp = compress;
302          if (type == 5) comp--;
303          sprintf(branchname,"%s.%s",name,rd->GetName());
304          TBranch *branch  = new AliArraySubBranch(branchname,this,leaflist,basketsize,comp);
305          branch->SetBit(kIsClone);
306          TObjArray *leaves = branch->GetListOfLeaves();
307          TLeaf *leaf = (TLeaf*)leaves->UncheckedAt(0);
308          leaf->SetOffset(rd->GetThisOffset());
309          leaf->SetLeafCount(leafcount);
310          Int_t arraydim = member->GetArrayDim();
311          if (arraydim) {
312            Int_t maxindex = member->GetMaxIndex(arraydim-1);
313            leaf->SetLen(maxindex);
314          }
315          fBranches.Add(branch);
316        }
317        else
318          for (Int_t i=0;i< member->GetMaxIndex(0);i++){
319            const char * dmname = member->GetName() ;
320            char  bname[200];
321            Int_t j=0;
322            while ( (dmname[j]!='[') && (dmname[j]!=0) ){
323              bname[j]=dmname[j];
324              j++;
325            }
326            bname[j]=0;
327            sprintf(leaflist,"%s(%d)[%s]/%s",bname,i,branchcount,itype);
328            Int_t comp = compress;
329            if (type == 5) comp--;
330            sprintf(branchname,"%s.%s(%d)",name,bname,i);
331            TBranch *branch  = new AliArraySubBranch(branchname,this,leaflist,basketsize,comp);
332            branch->SetBit(kIsClone);
333            TObjArray *leaves = branch->GetListOfLeaves();
334            TLeaf *leaf = (TLeaf*)leaves->UncheckedAt(0);
335            leaf->SetOffset(rd->GetThisOffset()+membertype->Size()*i);
336            leaf->SetLeafCount(leafcount);
337            fBranches.Add(branch);
338          }                
339        
340      }     
341    }
342    else if (clinfo->IsA()->InheritsFrom("AliDataType")){ //branch for basic type
343      Int_t type = (((AliDataType*)clinfo)->GetDataType())->GetType();
344      const char *itype = 0;
345      if (type <=0) 
346        Warning("BranchClones","Cannot process member:%s",clinfo->GetName());       
347      else{
348      
349        if (type == 1)  itype = "B";
350        if (type == 11) itype = "b";
351        if (type == 3)  itype = "I";
352        if (type == 5)  itype = "F";
353        if (type == 8)  itype = "D";
354        if (type == 13) itype = "i";
355        if (type == 2)  itype = "S";
356        if (type == 12) itype = "s";
357        sprintf(leaflist,"%s[%s]/%s",name,branchcount,itype);
358        Int_t comp = compress;
359        if (type == 5) comp--;
360        sprintf(branchname,"%s",clinfo->GetName());       
361        TBranch *branch  = new AliArraySubBranch(branchname,this,leaflist,basketsize,comp);
362        branch->SetBit(kIsClone);
363        TObjArray *leaves = branch->GetListOfLeaves();
364        TLeaf *leaf = (TLeaf*)leaves->UncheckedAt(0);
365        leaf->SetOffset(0);
366        leaf->SetLeafCount(leafcount);
367        fBranches.Add(branch);  
368      }
369    }
370
371 }
372
373
374 //______________________________________________________________________________
375 AliArrayBranch::~AliArrayBranch()
376 {
377 //*-*-*-*-*-*Default destructor for a BranchClones*-*-*-*-*-*-*-*-*-*-*-*
378 //*-*        =====================================
379
380    delete fBranchCount;
381    fBranchCount = 0;
382    fBranches.Delete();
383    fList = 0;
384 }
385
386
387 //______________________________________________________________________________
388 void AliArrayBranch::Browse(TBrowser *b)
389 {
390    fBranches.Browse( b );
391 }
392
393 //______________________________________________________________________________
394 Int_t AliArrayBranch::Fill()
395 {
396 //*-*-*-*-*Loop on all Branches of this BranchClones to fill Basket buffer*-*
397 //*-*      ===============================================================
398
399    Int_t i;
400    Int_t nbytes = 0;
401    Int_t nbranches = fBranches.GetEntriesFast();
402    char **ppointer = (char**)(fAddress);
403    if (ppointer == 0) return 0;
404    fList = (AliObjectArray*)(*ppointer);
405    //   fN    = fList->GetEntriesFast();
406    fN    = fList->GetSize();
407    fEntries++;
408
409    if (fN > fNdataMax) {
410       fNdataMax = fList->GetSize();
411       char branchcount[64];
412       sprintf(branchcount,"%s_",GetName());
413       TLeafI *leafi = (TLeafI*)fBranchCount->GetLeaf(branchcount);
414       leafi->SetMaximum(fNdataMax);
415       for (i=0;i<nbranches;i++)  {
416          TBranch *branch = (TBranch*)fBranches.UncheckedAt(i);
417          TObjArray *leaves = branch->GetListOfLeaves();
418          TLeaf *leaf = (TLeaf*)leaves->UncheckedAt(0);
419          leaf->SetAddress();
420       }
421    }
422    nbytes += fBranchCount->Fill();
423    for (i=0;i<nbranches;i++)  {
424       TBranch *branch = (TBranch*)fBranches.UncheckedAt(i);
425       TObjArray *leaves = branch->GetListOfLeaves();
426       TLeaf *leaf = (TLeaf*)leaves->UncheckedAt(0);
427       // leaf->Import(fList, fN);   // MI
428       Import(leaf,fN);              // MI change
429       nbytes += branch->Fill();
430    }
431    return nbytes;
432 }
433
434 //______________________________________________________________________________
435 Int_t AliArrayBranch::GetEntry(Long64_t entry, Int_t getall)
436 {
437 //*-*-*-*-*Read all branches of a BranchClones and return total number of bytes
438 //*-*      ====================================================================
439
440    if (TestBit(kDoNotProcess) && !getall) return 0;
441    Int_t nbytes = fBranchCount->GetEntry(entry);
442    TLeaf *leafcount = (TLeaf*)fBranchCount->GetListOfLeaves()->UncheckedAt(0);
443    fN = Int_t(leafcount->GetValue());
444    if (fN <= 0) return 0;
445
446    TBranch *branch;
447    Int_t nbranches = fBranches.GetEntriesFast();
448
449      // if fList exists, create clonesarray objects
450    if (fList) {
451      //fList->ExpandCreateFast(fN);   //MI  
452      fList->Resize(fN);    //MI change 
453      for (Int_t i=0;i<nbranches;i++)  {
454        branch = (TBranch*)fBranches.UncheckedAt(i);  
455        nbytes += ((AliArraySubBranch*)branch)->GetEntryExport(entry, getall, fList, fN);  // !!!MI
456       }
457    } else {
458       for (Int_t i=0;i<nbranches;i++)  {
459          branch = (TBranch*)fBranches.UncheckedAt(i);
460          nbytes += branch->GetEntry(entry, getall);
461       }
462    }
463   return nbytes;
464 }
465
466 //______________________________________________________________________________
467 void AliArrayBranch::Print(Option_t *option) const
468 {
469 //*-*-*-*-*-*-*-*-*-*-*-*Print TBranch parameters*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
470 //*-*                    ========================
471
472    fBranchCount->Print(option);
473    Int_t i;
474    Int_t nbranches = fBranches.GetEntriesFast();
475    for (i=0;i<nbranches;i++)  {
476       TBranch *branch = (TBranch*)fBranches[i];
477       branch->Print(option);
478    }
479 }
480
481 //______________________________________________________________________________
482 void AliArrayBranch::Reset(Option_t *option)
483 {
484 //*-*-*-*-*-*-*-*Reset a Branch*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
485 //*-*            ====================
486 //
487 //    Existing buffers are deleted
488 //    Entries, max and min are reset
489 //
490
491    fEntries        = 0;
492    fTotBytes       = 0;
493    fZipBytes       = 0;
494    Int_t i;
495    Int_t nbranches = fBranches.GetEntriesFast();
496    for (i=0;i<nbranches;i++)  {
497       TBranch *branch = (TBranch*)fBranches[i];
498       branch->Reset(option);
499    }
500    fBranchCount->Reset();
501
502
503 //______________________________________________________________________________
504 void  AliArrayBranch::SetAddress(void *add)  
505 {
506 //*-*-*-*-*-*-*-*Set address of this branch*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
507 //*-*            ====================
508 //*-*
509
510    fReadEntry = -1;
511    fAddress = (char*)add;
512    char **ppointer = (char**)(fAddress);
513    if ( (*ppointer)==0 ) {  //MI change 
514      *ppointer = (char*) new AliObjectArray(fClassName);
515      fAddress = (char*)ppointer;
516    }
517    fList = (AliObjectArray*)(*ppointer);
518    fBranchCount->SetAddress(&fN);
519
520 }
521
522 //______________________________________________________________________________
523 void AliArrayBranch::SetBasketSize(Int_t buffsize)
524 {
525 //*-*-*-*-*-*-*-*Reset basket size for all subbranches of this branchclones
526 //*-*            ==========================================================
527 //
528
529    fBasketSize = buffsize;
530    Int_t i;
531    Int_t nbranches = fBranches.GetEntriesFast();
532    for (i=0;i<nbranches;i++)  {
533       TBranch *branch = (TBranch*)fBranches[i];
534       branch->SetBasketSize(buffsize);
535    }
536 }
537
538 //_______________________________________________________________________
539 void AliArrayBranch::Streamer(TBuffer &b)
540 {
541 //*-*-*-*-*-*-*-*-*Stream a class object*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
542 //*-*              =========================================
543    if (b.IsReading()) {
544       b.ReadVersion();  //Version_t v = b.ReadVersion();
545       TNamed::Streamer(b);
546       b >> fCompress;
547       b >> fBasketSize;
548       b >> fEntryOffsetLen;
549       b >> fMaxBaskets;
550       b >> fWriteBasket;
551       b >> fEntryNumber;
552       b >> fEntries;
553       b >> fTotBytes;
554       b >> fZipBytes;
555       b >> fOffset;
556       b >> fBranchCount;
557       fClassName.Streamer(b);
558       fBranches.Streamer(b);
559       fTree = gTree;
560       TBranch *branch;
561       TLeaf *leaf;
562       Int_t nbranches = fBranches.GetEntriesFast();
563       for (Int_t i=0;i<nbranches;i++)  {
564          branch = (TBranch*)fBranches[i];
565          branch->SetBit(kIsClone);
566          leaf = (TLeaf*)branch->GetListOfLeaves()->UncheckedAt(0);
567          leaf->SetOffset(0);
568       }
569       fRead = 1;
570
571       AliClassInfo *clinfo = AliClassInfo::FindClassInfo(fClassName);
572       if (!clinfo) {
573         AliObjectArray tmp(fClassName); 
574         //MI change - object array to construct class description
575         clinfo = AliClassInfo::FindClassInfo(fClassName);
576       }
577       if (!clinfo) return;
578       
579       TClass *cl = clinfo->GetClass();
580       //      if (!cl) {
581       //   Warning("Streamer","Unknow class: %s. Cannot read BranchClones: %s",
582       //      fClassName.Data(),GetName());
583       //   return;
584       //}
585       if (cl){
586
587         if (!cl->GetListOfRealData())  cl->BuildRealData();
588         char branchname[80];
589         TRealData *rd;
590         TIter      next(cl->GetListOfRealData());
591         while ((rd = (TRealData *) next())) {
592           TDataMember *member = rd->GetDataMember();
593           if (!member->IsBasic())      continue;
594           if (!member->IsPersistent()) continue;
595           TDataType *membertype = member->GetDataType();
596           if (membertype->GetType() == 0) continue; 
597          //MI change - for array spliting
598           Int_t arraydim = member->GetArrayDim();
599           if (arraydim==1){
600             for (Int_t i=0;i< member->GetMaxIndex(0);i++){
601               const char * dmname = member->GetName() ;
602               char  bname[200];
603               Int_t j=0;
604               while ( (dmname[j]!='[') && (dmname[j]!=0) ){
605                 bname[j]=dmname[j];
606                 j++;
607               }
608               bname[j]=0;
609               sprintf(branchname,"%s.%s(%d)",GetName(),bname,i);
610               branch  = (TBranch*)fBranches.FindObject(branchname);
611               if (!branch) continue;
612               TObjArray *leaves = branch->GetListOfLeaves();
613               leaf = (TLeaf*)leaves->UncheckedAt(0);
614               leaf->SetOffset(rd->GetThisOffset()+membertype->Size()*i);
615             }
616           }
617           sprintf(branchname,"%s.%s",GetName(),rd->GetName());
618           branch  = (TBranch*)fBranches.FindObject(branchname);
619           if (!branch) continue;
620           TObjArray *leaves = branch->GetListOfLeaves();
621           leaf = (TLeaf*)leaves->UncheckedAt(0);
622           leaf->SetOffset(rd->GetThisOffset());         
623         }
624
625       }
626       else if (clinfo->IsA()->InheritsFrom("AliDataType")){ //branch for basic type
627         char branchname[100];
628         sprintf(branchname,"%s",clinfo->GetName());  
629         branch  = (TBranch*)fBranches.FindObject(branchname);
630         if (branch){
631           TObjArray *leaves = branch->GetListOfLeaves();
632           leaf = (TLeaf*)leaves->UncheckedAt(0);
633           leaf->SetOffset(0);
634         }   
635       }
636    }
637    else{
638         
639       b.WriteVersion(AliArrayBranch::IsA());
640       TNamed::Streamer(b);
641       b << fCompress;
642       b << fBasketSize;
643       b << fEntryOffsetLen;
644       b << fMaxBaskets;
645       b << fWriteBasket;
646       b << fEntryNumber;
647       b << fEntries;
648       b << fTotBytes;
649       b << fZipBytes;
650       b << fOffset;
651       b << fBranchCount;
652       fClassName.Streamer(b);
653       fBranches.Streamer(b);
654    }
655 }
656
657
658
659 void AliArrayBranch::Import(TLeaf * leaf, Int_t n)
660 {
661   //
662   // Import leaf n from the branch
663   //
664
665   const Int_t kIntUndefined = -9999;
666   Int_t j = 0;
667   char *clone;
668   Int_t len    = leaf->GetLenStatic();
669   Int_t offset = leaf->GetOffset();
670   void *value  = leaf->GetValuePointer();
671   //
672   for (Int_t i=0;i<n;i++) {
673     clone = (char*)fList->UncheckedAt(i);
674     //8bit int
675     if (leaf->IsA()==TLeafB::Class()){  
676       memcpy(&((Char_t*)value)[j],clone + offset, len);
677     }
678     //var size
679     if (leaf->IsA()==TLeafC::Class()){  
680       memcpy(&((Char_t*)value)[j],clone + offset, 1);
681     }
682     //double
683     if (leaf->IsA()==TLeafD::Class()){  
684       if (clone) memcpy(&((Double_t*)value)[j],clone + offset, 8*len);
685       else       memcpy(&((Double_t*)value)[j],&kIntUndefined,  8*len);      
686     }
687     //float
688     if (leaf->IsA()==TLeafF::Class()){   
689       if (clone) memcpy(&((Float_t*)value)[j],clone + offset, 4*len);
690       else       memcpy(&((Float_t*)value)[j],&kIntUndefined,  4*len);     
691     }    
692     //int
693     if (leaf->IsA()==TLeafI::Class()){  
694       if (clone) memcpy(&((Int_t*)value)[j],clone + offset, 4*len);
695       else       memcpy(&((Int_t*)value)[j],&kIntUndefined,  4*len);
696     }
697    //short
698     if (leaf->IsA()==TLeafS::Class()){  
699       if (clone) memcpy(&((Short_t*)value)[j],clone + offset, 2*len);
700       else       memcpy(&((Short_t*)value)[j],&kIntUndefined,  2*len);
701     }     
702     j += len;  
703   }
704   //
705 }
706
707 AliObjectBranch::AliObjectBranch(const Text_t *name, const Text_t *classname, void *addobj,
708                                  TTree * tree, 
709                                  Int_t basketsize, Int_t splitlevel, Int_t compress): TBranchObject()
710 {
711 //*-*-*-*-*-*-*-*-*-*-*-*-*Create a BranchObject*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
712 //*-*                      =====================
713 //
714    TClass *cl      = gROOT->GetClass(classname);
715    fTree =tree;  //MI change
716    if (!cl) {
717       Error("TBranchObject","Cannot find class:%s",classname);
718       return;
719    }
720    if (!cl->GetListOfRealData())  cl->BuildRealData();
721    Int_t bufsize=basketsize; //MI ?
722    SetName(name);
723    SetTitle(name);
724    fCompress = compress;
725    if (compress == -1) {
726      TFile *bfile = fTree->GetDirectory()->GetFile(); //MI chnge fTrre - gTree
727       if (bfile) fCompress = bfile->GetCompressionLevel();
728    }
729    if (basketsize < 100) basketsize = 100;
730    fBasketSize     = basketsize;
731    fAddress        = (char*)addobj;
732    fClassName      = classname;
733 #if ROOT_VERSION_CODE >= 262401
734    fBasketEntry    = new Long64_t[fMaxBaskets];
735 #else
736    fBasketEntry    = new Int_t[fMaxBaskets];
737 #endif
738    fBasketBytes    = new Int_t[fMaxBaskets];
739 #if ROOT_VERSION_CODE >= 262146
740    fBasketSeek     = new Long64_t[fMaxBaskets];
741 #else
742    fBasketSeek     = new Seek_t[fMaxBaskets];
743 #endif
744    fOldObject      = 0;
745
746    fBasketEntry[0] = fEntryNumber;
747    fBasketBytes[0] = 0;
748
749    TLeaf *leaf     = new TLeafObject(name,classname);
750    leaf->SetBranch(this);
751    leaf->SetAddress(addobj);
752    fNleaves = 1;
753    fLeaves.Add(leaf);
754    fTree->GetListOfLeaves()->Add(leaf);  //MI change fTree-gTree
755
756 // Set the bit kAutoDelete to specify that when reading
757 // in TLeafObject::ReadBasket, the object should be deleted
758 // before calling Streamer.
759 // It is foreseen to not set this bit in a future version.
760    SetAutoDelete(kTRUE);
761
762 //*-*-  Create the first basket
763    //   fTree       = gTree;  //MI change - no need anymore 
764    fDirectory  = fTree->GetDirectory();
765    fFileName   = "";
766
767    if (!splitlevel){
768      TBasket *basket = new TBasket(name,fTree->GetName(),this);
769      fBaskets.Add(basket);
770      return;
771    }
772
773    //
774    // 
775    TBranch * branch =this;
776    TObjArray *blist = branch->GetListOfBranches();
777    const char *rdname;
778    const char *dname;
779    char branchname[64];
780    if (!cl->GetListOfRealData()) cl->BuildRealData();
781    char **apointer = (char**)(addobj);
782    TObject *obj = (TObject*)(*apointer);
783    Bool_t delobj = kFALSE;
784    if (!obj) {
785       obj = (TObject*)cl->New();
786       delobj = kTRUE;
787    }
788 //*-*- Loop on all public data members of the class and its base classes
789    Int_t lenName = strlen(name);
790    Int_t isDot = 0;
791    if (name[lenName-1] == '.') isDot = 1;
792    TBranch *branch1 = 0;
793    TRealData *rd;
794    TIter      next(cl->GetListOfRealData());
795    while ((rd = (TRealData *) next())) {
796       TDataMember *dm = rd->GetDataMember();
797       if (!dm->IsPersistent()) continue; //do not process members with a ! as the first
798                                          // character in the comment field
799       rdname = rd->GetName();
800       dname  = dm->GetName();
801
802   //  Next line now commented, functionality to process arrays is now implemented
803   //  the statement is left to show how to use Property() and kIsArray
804   //     if (dm->Property() & kIsArray) continue;
805
806       TDataType *dtype = dm->GetDataType();
807       Int_t code = 0;
808       if (dtype) code = dm->GetDataType()->GetType();
809
810 //*-*- Encode branch name. Use real data member name
811       sprintf(branchname,"%s",rdname);
812       if (isDot) {
813          if (dm->IsaPointer()) sprintf(branchname,"%s%s",name,&rdname[1]);
814          else                  sprintf(branchname,"%s%s",name,&rdname[0]);
815       }
816       char leaflist[64];
817       Int_t offset    = rd->GetThisOffset();
818       char *pointer   = (char*)obj + offset;
819       if (dm->IsaPointer()) {
820          TClass *clobj = 0;
821          if (!dm->IsBasic()) clobj = gROOT->GetClass(dm->GetTypeName());
822          if (clobj && !strcmp("TClonesArray",clobj->GetName())) {
823             char *cpointer  =(char*)pointer;
824             char **ppointer =(char**)cpointer;
825             TClonesArray *list = (TClonesArray*)(*ppointer);
826             if (splitlevel != 100) {
827                if (isDot) branch1 = new TBranchClones(&branchname[0],pointer,bufsize);
828                else       branch1 = new TBranchClones(&branchname[1],pointer,bufsize);
829                blist->Add(branch1);
830             } else {
831                if (isDot) branch1 = new TBranchObject(&branchname[0],list->ClassName(),pointer,bufsize);
832                else       branch1 = new TBranchObject(&branchname[1],list->ClassName(),pointer,bufsize);
833                blist->Add(branch1);
834             }
835          }
836          else
837            if (clobj && !strcmp("AliObjectArray",clobj->GetName())) {
838              char *cpointer  =(char*)pointer;
839              char **ppointer =(char**)cpointer;
840              TClonesArray *list = (TClonesArray*)(*ppointer);
841              if (splitlevel != 100) {
842                if (isDot) branch1 = new AliArrayBranch(&branchname[0],pointer,fTree,bufsize,compress);
843                else       branch1 = new AliArrayBranch(&branchname[1],pointer,fTree,bufsize,compress);
844                blist->Add(branch1);
845              } else {
846                if (isDot) branch1 = new AliObjectBranch(&branchname[0],list->ClassName(),pointer,fTree,bufsize);
847                else       branch1 = new AliObjectBranch(&branchname[1],list->ClassName(),pointer,fTree,bufsize);
848                blist->Add(branch1);
849              }
850            }
851            else {
852              if (!clobj) {
853                if (code != 1) continue;
854                sprintf(leaflist,"%s/%s",dname,"C");
855                branch1 = new TBranch(branchname,pointer,leaflist,bufsize);
856                branch1->SetTitle(dname);
857                blist->Add(branch1);
858              } else {
859                if (!clobj->InheritsFrom(TObject::Class())) continue;
860                //branch1 = new TBranchObject(dname,clobj->GetName(),pointer,bufsize,0); //MI change
861                branch1 = new AliObjectBranch(dname,clobj->GetName(),pointer,fTree,bufsize,splitlevel);
862                if (isDot) branch1->SetName(&branchname[0]);
863                else       branch1->SetName(&branchname[1]);  //do not use the first character (*)
864                blist->Add(branch1);
865              }
866            }
867       }else {
868 //*-*-------------Data Member is a basic data type----------
869         if (dm->IsBasic()) {
870           if      (code ==  1) sprintf(leaflist,"%s/%s",rdname,"B");
871           else if (code == 11) sprintf(leaflist,"%s/%s",rdname,"b");
872           else if (code ==  2) sprintf(leaflist,"%s/%s",rdname,"S");
873           else if (code == 12) sprintf(leaflist,"%s/%s",rdname,"s");
874           else if (code ==  3) sprintf(leaflist,"%s/%s",rdname,"I");
875           else if (code == 13) sprintf(leaflist,"%s/%s",rdname,"i");
876           else if (code ==  5) sprintf(leaflist,"%s/%s",rdname,"F");
877           else if (code ==  8) sprintf(leaflist,"%s/%s",rdname,"D");
878           else {
879             printf("Cannot create branch for rdname=%s, code=%d\n",branchname, code);
880             leaflist[0] = 0;
881           }
882           branch1 = new TBranch(branchname,pointer,leaflist,bufsize);
883           branch1->SetTitle(rdname);
884           blist->Add(branch1);
885         }
886       }
887       if (branch1) branch1->SetOffset(offset);
888       else Warning("Branch","Cannot process member:%s",rdname);      
889    }
890    if (delobj) delete obj;
891 }
892
893
894
895 //______________________________________________________________________________
896 void AliObjectBranch::SetAddress(void *add)
897 {
898 //*-*-*-*-*-*-*-*Set address of this branch*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
899 //*-*            ====================
900 //
901
902    //special case when called from code generated by TTree::MakeClass
903    if (Long_t(add) == -1) {
904       SetBit(kWarn);
905       return;
906    }
907    fReadEntry = -1;
908    Int_t nbranches = fBranches.GetEntriesFast();
909    TLeaf *leaf = (TLeaf*)fLeaves.UncheckedAt(0);
910    if (leaf) leaf->SetAddress(add);
911    TBranch *branch;
912    fAddress = (char*)add;
913    char *pointer   = fAddress;
914    void **ppointer = (void**)add;
915    TObject *obj = (TObject*)(*ppointer);
916    TClass *cl = gROOT->GetClass(fClassName.Data());
917    if (!obj && cl) {
918       obj = (TObject*)cl->New();
919       *ppointer = (void*)obj;
920    }
921    Int_t i, offset;
922    if (!cl) {
923       for (i=0;i<nbranches;i++)  {
924          branch  = (TBranch*)fBranches[i];
925          pointer = (char*)obj;
926          branch->SetAddress(pointer);
927       }
928       return;
929    }
930    if (!cl->GetListOfRealData())  cl->BuildRealData();
931    char *fullname = new char[200];
932    const char *bname = GetName();
933    Int_t lenName = strlen(bname);
934    Int_t isDot = 0;
935    if (bname[lenName-1] == '.') isDot = 1;
936    const char *rdname;
937    TRealData *rd;
938    TIter      next(cl->GetListOfRealData());
939    while ((rd = (TRealData *) next())) {
940       TDataMember *dm = rd->GetDataMember();
941       if (!dm->IsPersistent()) continue;
942       rdname = rd->GetName();
943       TDataType *dtype = dm->GetDataType();
944       Int_t code = 0;
945       if (dtype) code = dm->GetDataType()->GetType();
946       offset  = rd->GetThisOffset();
947       pointer = (char*)obj + offset;
948       branch  = 0;
949       if (dm->IsaPointer()) {
950          TClass *clobj = 0;
951          if (!dm->IsBasic()) clobj = gROOT->GetClass(dm->GetTypeName());
952          if (clobj && !strcmp("TClonesArray",clobj->GetName())) {
953             if (isDot) sprintf(fullname,"%s%s",bname,&rdname[1]);
954             else       sprintf(fullname,"%s",&rdname[1]);
955             branch = (TBranch*)fBranches.FindObject(fullname);
956          }
957          else
958            if (clobj && !strcmp("AliObjectArray",clobj->GetName())) {
959              if (isDot) sprintf(fullname,"%s%s",bname,&rdname[1]);
960              else       sprintf(fullname,"%s",&rdname[1]);
961              branch = (TBranch*)fBranches.FindObject(fullname);
962            }
963          else {
964             if (!clobj) {
965                if (code != 1) continue;
966                if (isDot) sprintf(fullname,"%s%s",bname,&rdname[0]);
967                else       sprintf(fullname,"%s",&rdname[0]);
968                branch = (TBranch*)fBranches.FindObject(fullname);
969             } else {
970                if (!clobj->InheritsFrom(TObject::Class())) continue;
971                if (isDot) sprintf(fullname,"%s%s",bname,&rdname[1]);
972                else       sprintf(fullname,"%s",&rdname[1]);
973                branch = (TBranch*)fBranches.FindObject(fullname);
974             }
975          }
976       } else {
977          if (dm->IsBasic()) {
978             if (isDot) sprintf(fullname,"%s%s",bname,&rdname[0]);
979             else       sprintf(fullname,"%s",&rdname[0]);
980             branch = (TBranch*)fBranches.FindObject(fullname);
981          }
982       }
983       if(branch) branch->SetAddress(pointer);
984    }
985    delete [] fullname;
986 }
987
988
989
990
991
992 AliTree::AliTree(const char *name,const char *title, Int_t maxvirtualsize):
993   TTree(name,title,maxvirtualsize)
994 {
995   //
996   //default constructor for AliTree
997   gTree =this;
998 }
999
1000 TBranch * AliTree::AliBranch(const char *name, void *clonesaddress, Int_t bufsize, Int_t splitlevel,
1001                              Int_t compres)
1002 {
1003   //
1004   // Create an AliBranch and returns a reference to the TBranch
1005   //
1006   if (clonesaddress == 0) return 0;
1007   char *cpointer =(char*)clonesaddress;
1008   char **ppointer =(char**)cpointer;
1009   AliObjectArray *list = (AliObjectArray*)(*ppointer);
1010   if (list == 0) return 0;
1011   gTree = this;
1012   if (splitlevel) {
1013     TBranch *branch = new AliArrayBranch(name,clonesaddress,this,bufsize, compres);
1014     fBranches.Add(branch);
1015     return branch;
1016   } else {
1017     TBranchObject *branch = new TBranchObject(name,list->ClassName(),clonesaddress,bufsize,0);
1018     fBranches.Add(branch);
1019     return branch;
1020   }
1021 }
1022
1023 TBranch* AliTree::AliBranch(const char *name, const char *classname, void *addobj, 
1024                      Int_t bufsize, Int_t splitlevel)
1025 {
1026   //
1027   // Create an AliBranch and returns a reference to the TBranch
1028   //
1029   gTree = this;
1030   TClass *cl = gROOT->GetClass(classname);
1031   if (!cl) {
1032     Error("BranchObject","Cannot find class:%s",classname);
1033     return 0;
1034   }
1035   TBranch * branch = new AliObjectBranch(name,classname,addobj, this, bufsize,splitlevel);
1036   fBranches.Add(branch);
1037   return branch;
1038 }
1039
1040