#include "TBasket.h"
#include "TBranch.h"
#include "TBranchClones.h"
#include "TBranchElement.h"
#include "TStreamerInfo.h"
#include "TBranchRef.h"
#include "TError.h"
#include "TProcessID.h"
#include "TMath.h"
#include "TTree.h"
#include "TTreeCloner.h"
#include "TFile.h"
#include "TLeafB.h"
#include "TLeafI.h"
#include "TLeafL.h"
TTreeCloner::TTreeCloner(TTree *from, TTree *to, Option_t *method) :
fIsValid(kTRUE),
fFromTree(from),
fToTree(to),
fMethod(method),
fFromBranches( from ? from->GetListOfLeaves()->GetEntries()+1 : 0),
fToBranches( to ? from->GetListOfLeaves()->GetEntries()+1 : 0),
fMaxBaskets(CollectBranches()),
fBasketBranchNum(new UInt_t[fMaxBaskets]),
fBasketNum(new UInt_t[fMaxBaskets]),
fBasketSeek(new Long64_t[fMaxBaskets]),
fBasketEntry(new Long64_t[fMaxBaskets]),
fBasketIndex(new UInt_t[fMaxBaskets]),
fCloneMethod(TTreeCloner::kDefault),
fToStartEntries(0)
{
TString opt(method);
opt.ToLower();
if (opt.Contains("sortbasketsbybranch")) {
fCloneMethod = TTreeCloner::kSortBasketsByBranch;
} else if (opt.Contains("sortbasketsbyentry")) {
fCloneMethod = TTreeCloner::kSortBasketsByEntry;
} else {
fCloneMethod = TTreeCloner::kSortBasketsByOffset;
}
if (fToTree) fToStartEntries = fToTree->GetEntries();
}
Bool_t TTreeCloner::Exec()
{
CopyStreamerInfos();
CopyProcessIds();
CloseOutWriteBaskets();
CollectBaskets();
SortBaskets();
WriteBaskets();
CopyMemoryBaskets();
return kTRUE;
}
TTreeCloner::~TTreeCloner()
{
delete [] fBasketBranchNum;
delete [] fBasketNum;
delete [] fBasketSeek;
delete [] fBasketEntry;
delete [] fBasketIndex;
}
void TTreeCloner::CloseOutWriteBaskets()
{
for(Int_t i=0; i<fToBranches.GetEntries(); ++i) {
TBranch *to = (TBranch*)fToBranches.UncheckedAt(i);
to->FlushOneBasket(to->GetWriteBasket());
}
}
UInt_t TTreeCloner::CollectBranches(TBranch *from, TBranch *to)
{
UInt_t numBaskets = 0;
if (from->InheritsFrom(TBranchClones::Class())) {
TBranchClones *fromclones = (TBranchClones*)from;
TBranchClones *toclones = (TBranchClones*)to;
numBaskets += CollectBranches(fromclones->fBranchCount,toclones->fBranchCount);
} else if (from->InheritsFrom(TBranchElement::Class())) {
Int_t nb = from->GetListOfLeaves()->GetEntries();
Int_t fnb= to->GetListOfLeaves()->GetEntries();
if (nb!=fnb && (nb==0 || fnb==0)) {
Error("TTreeCloner::CollectBranches",
"The export branch and the import branch do not have the same split level. (The branch name is %s.)",
from->GetName());
fIsValid = kFALSE;
return 0;
}
TBranchElement *fromelem = (TBranchElement*)from;
TBranchElement *toelem = (TBranchElement*)to;
if (fromelem->fMaximum > toelem->fMaximum) toelem->fMaximum = fromelem->fMaximum;
} else {
Int_t nb = from->GetListOfLeaves()->GetEntries();
Int_t fnb= to->GetListOfLeaves()->GetEntries();
if (nb!=fnb) {
Error("TTreeCloner::CollectBranches",
"The export branch and the import branch (%s) do not have the same number of leaves (%d vs %d)",
from->GetName(), fnb,nb);
fIsValid = kFALSE;
return 0;
}
for (Int_t i=0;i<nb;i++) {
TLeaf *fromleaf_gen = (TLeaf*)from->GetListOfLeaves()->At(i);
if (fromleaf_gen->IsA()==TLeafI::Class()) {
TLeafI *fromleaf = (TLeafI*)from->GetListOfLeaves()->At(i);
TLeafI *toleaf = (TLeafI*)to->GetListOfLeaves()->At(i);
if (fromleaf->GetMaximum() > toleaf->GetMaximum())
toleaf->SetMaximum( fromleaf->GetMaximum() );
if (fromleaf->GetMinimum() < toleaf->GetMinimum())
toleaf->SetMinimum( fromleaf->GetMinimum() );
} else if (fromleaf_gen->IsA()==TLeafL::Class()) {
TLeafL *fromleaf = (TLeafL*)from->GetListOfLeaves()->At(i);
TLeafL *toleaf = (TLeafL*)to->GetListOfLeaves()->At(i);
if (fromleaf->GetMaximum() > toleaf->GetMaximum())
toleaf->SetMaximum( fromleaf->GetMaximum() );
if (fromleaf->GetMinimum() < toleaf->GetMinimum())
toleaf->SetMinimum( fromleaf->GetMinimum() );
} else if (fromleaf_gen->IsA()==TLeafB::Class()) {
TLeafB *fromleaf = (TLeafB*)from->GetListOfLeaves()->At(i);
TLeafB *toleaf = (TLeafB*)to->GetListOfLeaves()->At(i);
if (fromleaf->GetMaximum() > toleaf->GetMaximum())
toleaf->SetMaximum( fromleaf->GetMaximum() );
if (fromleaf->GetMinimum() < toleaf->GetMinimum())
toleaf->SetMinimum( fromleaf->GetMinimum() );
}
}
}
fFromBranches.AddLast(from);
fToBranches.AddLast(to);
numBaskets += from->GetWriteBasket();
numBaskets += CollectBranches(from->GetListOfBranches(),to->GetListOfBranches());
return numBaskets;
}
UInt_t TTreeCloner::CollectBranches(TObjArray *from, TObjArray *to)
{
Int_t fnb = from->GetEntries();
Int_t tnb = to->GetEntries();
if (!fnb || !tnb) {
return 0;
}
UInt_t numBasket = 0;
Int_t fi = 0;
Int_t ti = 0;
while (ti < tnb) {
TBranch* fb = (TBranch*) from->UncheckedAt(fi);
TBranch* tb = (TBranch*) to->UncheckedAt(ti);
Int_t firstfi = fi;
while (strcmp(fb->GetName(), tb->GetName())) {
++fi;
if (fi==firstfi) {
fb = 0;
break;
}
if (fi >= fnb) {
fi = 0;
}
fb = (TBranch*) from->UncheckedAt(fi);
}
if (fb) {
numBasket += CollectBranches(fb, tb);
++fi;
if (fi >= fnb) {
fi = 0;
}
}
++ti;
}
return numBasket;
}
UInt_t TTreeCloner::CollectBranches()
{
UInt_t numBasket = CollectBranches(fFromTree->GetListOfBranches(),
fToTree->GetListOfBranches());
if (fFromTree->GetBranchRef()) {
fToTree->BranchRef();
numBasket += CollectBranches(fFromTree->GetBranchRef(),fToTree->GetBranchRef());
}
return numBasket;
}
void TTreeCloner::CollectBaskets()
{
UInt_t len = fFromBranches.GetEntries();
for(UInt_t i=0,bi=0; i<len; ++i) {
TBranch *from = (TBranch*)fFromBranches.UncheckedAt(i);
for(Int_t b=0; b<from->GetWriteBasket(); ++b,++bi) {
fBasketBranchNum[bi] = i;
fBasketNum[bi] = b;
fBasketSeek[bi] = from->GetBasketSeek(b);
fBasketEntry[bi] = from->GetBasketEntry()[b];
fBasketIndex[bi] = bi;
}
}
}
void TTreeCloner::CopyStreamerInfos()
{
TFile *fromFile = fFromTree->GetDirectory()->GetFile();
TFile *toFile = fToTree->GetDirectory()->GetFile();
TList *l = fromFile->GetStreamerInfoList();
TIter next(l);
TStreamerInfo *oldInfo;
while ( (oldInfo = (TStreamerInfo*)next()) ) {
TStreamerInfo *curInfo = 0;
TClass *cl = TClass::GetClass(oldInfo->GetName());
if ((cl->IsLoaded() && (cl->GetNew()!=0 || cl->HasDefaultConstructor()))
|| !cl->IsLoaded()) {
curInfo = (TStreamerInfo*)cl->GetStreamerInfo(oldInfo->GetClassVersion());
if (oldInfo->GetClassVersion()==1) {
TStreamerInfo *matchInfo = (TStreamerInfo*)cl->FindStreamerInfo(oldInfo->GetCheckSum());
if (matchInfo) {
curInfo = matchInfo;
}
}
curInfo->TagFile(toFile);
} else {
oldInfo->TagFile(toFile);
}
}
delete l;
}
void TTreeCloner::CopyMemoryBaskets()
{
TBasket *basket = 0;
for(Int_t i=0; i<fToBranches.GetEntries(); ++i) {
TBranch *from = (TBranch*)fFromBranches.UncheckedAt( i );
TBranch *to = (TBranch*)fToBranches.UncheckedAt( i );
basket = from->GetListOfBaskets()->GetEntries() ? from->GetBasket(from->GetWriteBasket()) : 0;
if (basket) {
basket = (TBasket*)basket->Clone();
basket->SetBranch(to);
to->AddBasket(*basket, kFALSE, fToStartEntries+from->GetBasketEntry()[from->GetWriteBasket()]);
} else {
to->AddLastBasket( fToStartEntries+from->GetBasketEntry()[from->GetWriteBasket()] );
}
if (from->GetEntries()!=0 && from->GetWriteBasket()==0 && (basket==0 || basket->GetNevBuf()==0)) {
to->SetEntries(to->GetEntries()+from->GetEntries());
}
}
}
void TTreeCloner::CopyProcessIds()
{
TFile *fromfile = fFromTree->GetDirectory()->GetFile();
TFile *tofile = fToTree->GetDirectory()->GetFile();
fPidOffset = tofile->GetNProcessIDs();
TIter next(fromfile->GetListOfKeys());
TKey *key;
TDirectory::TContext cur(gDirectory,fromfile);
while ((key = (TKey*)next())) {
if (!strcmp(key->GetClassName(),"TProcessID")) {
TProcessID *pid = (TProcessID*)key->ReadObjectAny(0);
UShort_t out = 0;
TObjArray *pids = tofile->GetListOfProcessIDs();
Int_t npids = tofile->GetNProcessIDs();
Bool_t wasIn = kFALSE;
for (Int_t i=0;i<npids;i++) {
if (pids->At(i) == pid) {out = (UShort_t)i; wasIn = kTRUE; break;}
}
if (!wasIn) {
TDirectory *dirsav = gDirectory;
tofile->cd();
tofile->SetBit(TFile::kHasReferences);
pids->AddAtAndExpand(pid,npids);
pid->IncrementCount();
char name[32];
sprintf(name,"ProcessID%d",npids);
pid->Write(name);
tofile->IncrementProcessIDs();
if (gDebug > 0) {
printf("WriteProcessID, name=%s, file=%s\n",name,tofile->GetName());
}
if (dirsav) dirsav->cd();
out = (UShort_t)npids;
}
if (out<fPidOffset) {
Error("CopyProcessIDs","Copied %s from %s might already exist!\n",
pid->GetName(),fromfile->GetName());
}
}
}
}
void TTreeCloner::SortBaskets()
{
switch (fCloneMethod) {
case kSortBasketsByBranch:
break;
case kSortBasketsByEntry:
TMath::Sort( fMaxBaskets, fBasketEntry, fBasketIndex, kFALSE );
break;
case kSortBasketsByOffset:
default:
TMath::Sort( fMaxBaskets, fBasketSeek, fBasketIndex, kFALSE );
break;
}
}
void TTreeCloner::WriteBaskets()
{
TBasket *basket = new TBasket();
for(UInt_t j=0; j<fMaxBaskets; ++j) {
TBranch *from = (TBranch*)fFromBranches.UncheckedAt( fBasketBranchNum[ fBasketIndex[j] ] );
TBranch *to = (TBranch*)fToBranches.UncheckedAt( fBasketBranchNum[ fBasketIndex[j] ] );
TFile *tofile = to->GetFile(0);
TFile *fromfile = from->GetFile(0);
Int_t index = fBasketNum[ fBasketIndex[j] ];
Long64_t pos = from->GetBasketSeek(index);
if (pos!=0) {
if (from->GetBasketBytes()[index] == 0) {
from->GetBasketBytes()[index] = basket->ReadBasketBytes(pos, fromfile);
}
Int_t len = from->GetBasketBytes()[index];
basket->LoadBasketBuffers(pos,len,fromfile);
basket->IncrementPidOffset(fPidOffset);
basket->CopyTo(tofile);
to->AddBasket(*basket,kTRUE,fToStartEntries + from->GetBasketEntry()[index]);
} else {
TBasket *frombasket = from->GetBasket( index );
if (frombasket->GetNevBuf()>0) {
TBasket *tobasket = (TBasket*)frombasket->Clone();
tobasket->SetBranch(to);
to->AddBasket(*tobasket, kFALSE, fToStartEntries+from->GetBasketEntry()[index]);
to->FlushOneBasket(to->GetWriteBasket());
}
}
}
delete basket;
}
Last change: Fri Oct 31 14:49:44 2008
Last generated: 2008-10-31 14:49
This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.