TRef between Branch

From: Stefano Dusini <stefano.dusini_at_pd.infn.it>
Date: Mon, 12 Oct 2009 15:47:35 +0200


Hi

I would like to write a tree with two branch one containg the vertecies (Vertex class) and the other the tracks of the event (Track class). Clarly the number of entries in the brach containing the vertex is smaler then the number of entries in the track branch.

To link the tracks to the verticies I use TRefArray.

here below an example where I create some dummy vertexes and tracks.  

   TFile *hFile = TFile::Open("myFile.root","recreate");     TTree *myTree = new TTree("myTree","myTree");     myTree->BranchRef();
    Vertex *vertex = new Vertex();
    Track *track = new Track();
    TBranch *vertexBranch =
myTree->Branch("VertexBranch","Vertex",&vertex,16000,2);

    TBranch *trackBranch = myTree->Branch("TrackBranch","Track",&track,16000,2);     

    for(int iv=0;iv<3;iv++) {

      vertex->SetID(iv); // I set the id of the vertex
      int itt=0;
      for(int it=0;it<iv+1;it++) {
	track->SetID(itt); //I set the id of the Track
	vertex->AddTrack(track); //This metod add the pointer track to a TRefArray of
the class VERTEX
	trackBranch->Fill();
	track->Clear();
      }
      vertexBranch->Fill();
      vertex->Clear();

    }
    hFile->Write();
    myTree->Print();
    hFile->Close();

The code to read the file is    

    TFile *hFile = TFile::Open("emuHist.root");     TTree *myTree = (TTree*) hFile->Get("myTree");     myTree->BranchRef();
    myTree->Print();     

    Vertex *vertex = new Vertex();
    Track *track = new Track();
    TBranch *vertexBranch = myTree->GetBranch("VertexBranch");     TBranch *trackBranch = myTree->GetBranch("TrackBranch");     

    vertexBranch->SetAddress(&vertex);
    trackBranch->SetAddress(&track);     

    int nv = vertexBranch->GetEntries();     int nt = trackBranch->GetEntries();     

    for(int iv=0;iv<nv;iv++) {

      vertexBranch->GetEntry(iv);
      cout << "Vertex [" << iv << "] = " << vertex->GetID() << endl;
      TRefArray *aTrk = vertex->GetTrackArray();
      cout << "Number of tracks at the vertex " << aTrk->GetEntries() << endl;
      for(int it=0;it<aTrk->GetEntries();it++) {
	Track *track = vertex->GetTrack(it); //return (Track*)(eTrackArray->At(index))
	cout << "Track [" << it << "] ID " << track->GetID() << endl;
      }
      vertex->Clear();

    }
  }

I'm able to read the vertexes but the pointers at the track attached to the vertex, i.e. the elements of the TRefArray aTrk are NULL pointers.

In the TTree::Print() of the tree that I have create is I see the correct number of entries for the Vertex and Track branch but zero entry for the TRefTable. Do you have any explanation?

Thanks a lot
Stefano



*Tree :myTree : myTree *
*Entries : 0 : Total = 15243 bytes File Size = 2890 *
--
Stefano Dusin
INFN - sezione di Padova
via Marzolo 8
Office: +39-049-8277323
Received on Mon Oct 12 2009 - 15:47:59 CEST

This archive was generated by hypermail 2.2.0 : Mon Oct 12 2009 - 17:50:04 CEST