swh:1:snp:af87cd67498ef4fe47c76ed3e7caffe5b61facaf
Raw File
Tip revision: db2b2a8ff8a518aa377b920053f2e7c59e8ff605 authored by Pere Mato on 18 May 2017, 14:33:47 UTC
Update ROOT version files to v6.09/05.
Tip revision: db2b2a8
copytree2.C
/// \file
/// \ingroup tutorial_tree
/// \notebook -nodraw
/// Copy a subset of a Tree to a new Tree, one branch in a separate file.
///
/// One branch of the new Tree is written to a separate file
/// The input file has been generated by the program in `$ROOTSYS/test/Event`
/// with the command `Event 1000 1 1 1`
///
/// \macro_code
///
/// \author Rene Brun

R__LOAD_LIBRARY($ROOTSYS/test/libEvent.so)

void copytree2() {

   //Get old file, old tree and set top branch address
   TFile *oldfile;
   TString dir = "$ROOTSYS/test/Event.root";
   gSystem->ExpandPathName(dir);
   if (!gSystem->AccessPathName(dir))
       {oldfile = new TFile("$ROOTSYS/test/Event.root");}
   else {oldfile = new TFile("./Event.root");}   TTree *oldtree = (TTree*)oldfile->Get("T");
   Event *event   = new Event();
   oldtree->SetBranchAddress("event",&event);
   oldtree->SetBranchStatus("*",0);
   oldtree->SetBranchStatus("event",1);
   oldtree->SetBranchStatus("fNtrack",1);
   oldtree->SetBranchStatus("fNseg",1);
   oldtree->SetBranchStatus("fH",1);


   //Create a new file + a clone of old tree header. Do not copy events
   TFile *newfile = new TFile("small.root","recreate");
   TTree *newtree = oldtree->CloneTree(0);

   //Divert branch fH to a separate file and copy all events
   newtree->GetBranch("fH")->SetFile("small_fH.root");
   newtree->CopyEntries(oldtree);

   newtree->Print();
   newfile->Write();
   delete oldfile;
   delete newfile;
}
back to top