26 ClassImp(TBranchClones);
39 TBranchClones::TBranchClones()
52 TBranchClones::TBranchClones(TTree *tree,
const char* name,
void* pointer, Int_t basketsize, Int_t compress, Int_t splitlevel)
60 Init(tree,0,name,pointer,basketsize,compress,splitlevel);
66 TBranchClones::TBranchClones(TBranch *parent,
const char* name,
void* pointer, Int_t basketsize, Int_t compress, Int_t splitlevel)
74 Init(0,parent,name,pointer,basketsize,compress,splitlevel);
80 void TBranchClones::Init(TTree *tree, TBranch *parent,
const char* name,
void* pointer, Int_t basketsize, Int_t compress, Int_t splitlevel)
82 if (tree==0 && parent!=0) tree = parent->GetTree();
84 fMother = parent ? parent->GetMother() :
this;
91 if ((compress == -1) && tree->GetDirectory()) {
93 if (tree->GetDirectory()) {
94 bfile = tree->GetDirectory()->GetFile();
97 compress = bfile->GetCompressionSettings();
100 char* cpointer = (
char*) pointer;
101 char** ppointer = (
char**) pointer;
102 fList = (TClonesArray*) *ppointer;
104 TClass* cl = fList->GetClass();
108 tree->BuildStreamerInfo(cl);
109 fClassName = cl->GetName();
110 fSplitLevel = splitlevel;
113 if (basketsize < 100) {
116 leaflist.Form(
"%s_/I", name);
117 branchcount.Form(
"%s_", name);
118 fBranchCount =
new TBranch(
this, branchcount, &fN, leaflist, basketsize);
119 fBranchCount->SetBit(kIsClone);
120 TLeaf* leafcount = (TLeaf*) fBranchCount->GetListOfLeaves()->UncheckedAt(0);
121 fDirectory = fTree->GetDirectory();
125 const char* itype = 0;
127 TIter next(cl->GetListOfRealData());
128 while ((rd = (TRealData *) next())) {
129 if (rd->TestBit(TRealData::kTransient))
continue;
131 if (rd->IsObject()) {
134 TDataMember* member = rd->GetDataMember();
135 if (!member->IsPersistent()) {
139 if (!member->IsBasic() || member->IsaPointer()) {
140 Warning(
"BranchClones",
"Cannot process: %s::%s", cl->GetName(), member->GetName());
144 if ((splitlevel > 1) || fList->TestBit(TClonesArray::kForgetBits) || cl->CanIgnoreTObjectStreamer()) {
145 if (!std::strcmp(member->GetName(),
"fBits")) {
148 if (!std::strcmp(member->GetName(),
"fUniqueID")) {
152 tree->BuildStreamerInfo(TClass::GetClass(member->GetFullTypeName()));
153 TDataType* membertype = member->GetDataType();
154 Int_t type = membertype->GetType();
156 Warning(
"BranchClones",
"Cannot process: %s::%s of type zero!", cl->GetName(), member->GetName());
162 }
else if (type == 2) {
164 }
else if (type == 3) {
166 }
else if (type == 5) {
168 }
else if (type == 8) {
170 }
else if (type == 9) {
172 }
else if (type == 11) {
174 }
else if (type == 12) {
176 }
else if (type == 13) {
180 leaflist.Form(
"%s[%s]/%s", member->GetName(), branchcount.Data(), itype);
181 Int_t comp = compress;
182 branchname.Form(
"%s.%s", name, rd->GetName());
183 TBranch* branch =
new TBranch(
this, branchname,
this, leaflist, basketsize, comp);
184 branch->SetBit(kIsClone);
185 TObjArray* leaves = branch->GetListOfLeaves();
186 TLeaf* leaf = (TLeaf*) leaves->UncheckedAt(0);
187 leaf->SetOffset(rd->GetThisOffset());
188 leaf->SetLeafCount(leafcount);
189 Int_t arraydim = member->GetArrayDim();
193 maxindex *= member->GetMaxIndex(--arraydim);
195 leaf->SetLen(maxindex);
197 fBranches.Add(branch);
204 TBranchClones::~TBranchClones()
216 void TBranchClones::Browse(TBrowser* b)
224 Int_t TBranchClones::FillImpl(ROOT::Internal::TBranchIMTHelper *imtHelper)
228 Int_t nbranches = fBranches.GetEntriesFast();
229 char** ppointer = (
char**) fAddress;
233 fList = (TClonesArray*) *ppointer;
234 fN = fList->GetEntriesFast();
236 if (fN > fNdataMax) {
237 fNdataMax = fList->GetSize();
239 branchcount.Form(
"%s_", GetName());
240 TLeafI* leafi = (TLeafI*) fBranchCount->GetLeaf(branchcount);
241 leafi->SetMaximum(fNdataMax);
242 for (i = 0; i < nbranches; i++) {
243 TBranch* branch = (TBranch*) fBranches.UncheckedAt(i);
244 TObjArray* leaves = branch->GetListOfLeaves();
245 TLeaf* leaf = (TLeaf*) leaves->UncheckedAt(0);
249 nbytes += fBranchCount->FillImpl(imtHelper);
250 for (i = 0; i < nbranches; i++) {
251 TBranch* branch = (TBranch*) fBranches.UncheckedAt(i);
252 TObjArray* leaves = branch->GetListOfLeaves();
253 TLeaf* leaf = (TLeaf*) leaves->UncheckedAt(0);
254 leaf->Import(fList, fN);
255 nbytes += branch->FillImpl(imtHelper);
263 Int_t TBranchClones::GetEntry(Long64_t entry, Int_t getall)
265 if (TestBit(kDoNotProcess) && !getall) {
268 Int_t nbytes = fBranchCount->GetEntry(entry, getall);
269 TLeaf* leafcount = (TLeaf*) fBranchCount->GetListOfLeaves()->UncheckedAt(0);
270 fN = Int_t(leafcount->GetValue());
278 Int_t nbranches = fBranches.GetEntriesFast();
282 fList->ExpandCreateFast(fN);
283 for (Int_t i = 0; i < nbranches; i++) {
284 branch = (TBranch*) fBranches.UncheckedAt(i);
285 if (((TLeaf*) branch->GetListOfLeaves()->UncheckedAt(0))->GetOffset() < 0) {
288 nbytes += branch->GetEntryExport(entry, getall, fList, fN);
291 for (Int_t i = 0; i < nbranches; i++) {
292 branch = (TBranch*) fBranches.UncheckedAt(i);
293 nbytes += branch->GetEntry(entry, getall);
302 void TBranchClones::Print(Option_t *option)
const
304 fBranchCount->Print(option);
305 Int_t nbranches = fBranches.GetEntriesFast();
306 for (Int_t i = 0; i < nbranches; i++) {
307 TBranch* branch = (TBranch*) fBranches.At(i);
308 branch->Print(option);
318 void TBranchClones::Reset(Option_t* option)
323 Int_t nbranches = fBranches.GetEntriesFast();
324 for (Int_t i = 0; i < nbranches; i++) {
325 TBranch* branch = (TBranch*) fBranches.At(i);
326 branch->Reset(option);
328 fBranchCount->Reset();
337 void TBranchClones::ResetAfterMerge(TFileMergeInfo *info)
342 Int_t nbranches = fBranches.GetEntriesFast();
343 for (Int_t i = 0; i < nbranches; i++) {
344 TBranch* branch = (TBranch*) fBranches.At(i);
345 branch->ResetAfterMerge(info);
347 fBranchCount->ResetAfterMerge(info);
353 void TBranchClones::SetAddress(
void* addr)
356 fAddress = (
char*) addr;
357 char** pp= (
char**) fAddress;
358 if (pp && (*pp == 0)) {
360 *pp= (
char*)
new TClonesArray(fClassName);
364 fList = (TClonesArray*) *pp;
366 fBranchCount->SetAddress(&fN);
372 void TBranchClones::SetBasketSize(Int_t buffsize)
374 TBranch::SetBasketSize(buffsize);
376 Int_t nbranches = fBranches.GetEntriesFast();
377 for (Int_t i = 0; i < nbranches; i++) {
378 TBranch* branch = (TBranch*) fBranches[i];
379 branch->SetBasketSize(fBasketSize);
386 void TBranchClones::Streamer(TBuffer& b)
390 b.ReadVersion(&R__s, &R__c);
394 b >> fEntryOffsetLen;
403 fClassName.Streamer(b);
404 fBranches.Streamer(b);
408 Int_t nbranches = fBranches.GetEntriesFast();
409 for (Int_t i = 0; i < nbranches; i++) {
410 branch = (TBranch*) fBranches[i];
411 branch->SetBit(kIsClone);
412 leaf = (TLeaf*) branch->GetListOfLeaves()->UncheckedAt(0);
416 TClass* cl = TClass::GetClass((
const char*) fClassName);
418 Warning(
"Streamer",
"Unknown class: %s. Cannot read BranchClones: %s", fClassName.Data(), GetName());
419 SetBit(kDoNotProcess);
422 if (!cl->GetListOfRealData()) {
427 TIter next(cl->GetListOfRealData());
428 while ((rd = (TRealData*) next())) {
429 if (rd->TestBit(TRealData::kTransient))
continue;
431 TDataMember* member = rd->GetDataMember();
432 if (!member || !member->IsBasic() || !member->IsPersistent()) {
435 TDataType* membertype = member->GetDataType();
436 if (!membertype->GetType()) {
439 branchname.Form(
"%s.%s", GetName(), rd->GetName());
440 branch = (TBranch*) fBranches.FindObject(branchname);
444 TObjArray* leaves = branch->GetListOfLeaves();
445 leaf = (TLeaf*) leaves->UncheckedAt(0);
446 leaf->SetOffset(rd->GetThisOffset());
448 b.CheckByteCount(R__s, R__c, TBranchClones::IsA());
450 R__c = b.WriteVersion(TBranchClones::IsA(), kTRUE);
454 b << fEntryOffsetLen;
463 fClassName.Streamer(b);
464 fBranches.Streamer(b);
465 b.SetByteCount(R__c, kTRUE);
474 void TBranchClones::UpdateFile()
476 fBranchCount->UpdateFile();
477 TBranch::UpdateFile();