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