Re: How to edit the existing tree?

From: Rene Brun (Rene.Brun@cern.ch)
Date: Mon Nov 24 1997 - 15:14:59 MET


valeri@na.infn.it wrote:

> Dear Root,
>
> If I have the tree T with 1 branch "branch1" and several
> filled events "nevent" in file f.root how to:
>
> 1. Update the vertical structure if the tree
>    - add new events to the file
>    - change the existing event in the file
>    - insert a new event into the file (not in the last position)
>    - delete the existing event in the file
>

You can add new events to an existing Tree:    TFile
f("file.root","update");
    TTree *tree = (TTree*)f.Get("tree");   // for a Tree called "tree"

     // Set thebranch addresses
..
   // Continue to fill the Tree
  tree->Fill();

at the end, write a new version of the header
 tree->Write();

You cannot change an existing event in the file. There are no
technical limitations to do that. We simply want to STRONGLY
discourage people to do that !!

Events are always added to the end of the file.
You can use the TEventList mechanism to mark a list of events
that you want to process or not to process.

> 2. Update the horizontal structure
>    - add the new branch to the tree in the file
>    - fill this new branch event by event
>
> Is it possible to have branches of the same tree in the separate
> files?
>

Yes, you can add new branches to an existing Tree. The way to adda new
branch is the following.

Open the file in update mode as shown above.

Create new branches with
    TBranch *bnew1 = tree->Branch("bnew1",.........
    TBranch *bnew2 = tree->Branch("bnew2",.........

tree->SetBranchStatus("*",0);
//  activate branch status for all old branches necessary to compute the
new branches
tree->SetBranchStatus('oldb1",1);
tree->SetBranchStatus('oldbn",1);

//  loop on entries
   Int_t nentries = tree->GetEntries();
   for (Int_t i=0;i<nentries;i++) {
     tree->GetEvent(i);

     //  compute variables for new branches and fill these branches
  ...
   bnew1->Fill();
   bnew2->Fill();
}

  // save new Tree version
  tree->Write();


All branches of a Tree must currently be in  the same file.
However, you can have as many Trees as you want.
Below is an example illustrating how to loop on branches of a small Tree

and for some selected events read all branches of a big Tree.

TFile fsmall("small.root");
TFile fbig("big.root");
TTree *smalltree=(TTree*)fsmall.Get("small");
TTree *bigtree=(TTree*)fbig.Get("big");

//  set branch addresses for small and big Trees
  ......

Int_t nentries = smalltree->GetEntries();
for (Int_t i=0;i<nentries;i++) {
   smalltree->GetEvent(i);
   if (event is not accepted) continue;
  bigtree->GetEvent(i);
  ....
}

Rene Brun



This archive was generated by hypermail 2b29 : Tue Jan 04 2000 - 00:26:22 MET