Re: Branch: I want easy access to data

From: Rene Brun (Rene.Brun@cern.ch)
Date: Wed Mar 25 1998 - 12:01:09 MET


Apologies again if you experience long delays to answer your mails.
I am now in a totally divergent process with my mail folders.

Rene Brun


Witold Przygoda wrote:
> 
> Hallo,
> 
> excuse me if I repeat the well known question but this is
> a very important problem to me.
> 
> I will show it on an easy example.
> Imagine that I have hit class like this:
> 
> class MyHits : public TObject {
> 
> public:
> 
> int SimpleType;
> MyDataClass finalData[3];
> 
> // of course Streamer function defined in .cc file
> 
> MyHits();
> virtual ~MyHits();
> ClassDef(MyHits,1)
> 
> };
> 
> MyDataClass : public TObject {
> 
> int SimpleType;
> int someData[10];
> float *myArrayPtr;
> 
> // again I provide Streamer which takes care of pointer also
> 
> MyDataClass();
> virtual ~MyDataClass();
> ClassDef(MyDataClass,1)
> };
> 
> What I know:
> I can write it to one branch (just dumping data into one solid
> block without spliting). This works fine and I can read so
> created file, fill back data structure, make loops over events
> and fill some histograms.
> 
> What I would like and donīt know:
> It is demanded by some people to provide them much easier (and faster!)
> way to access all data. I do not want to define the detailed (leaf by
> leaf) structure of my branch(es) but to stay with Streamer-mechanism.
> I know I can set up split mode and each variable will create another
> branch.
> 
> My question: I failed with having access to array myData (of type
> MyDataClass). The branch seems to be not created. Is it possible
> (and how) finally to have the access like
> T.Draw("simpleType");
> 
> also to T.Draw("myData[0].SimpleType");
>         T.Draw("myData[1].someData");
>         T.Draw("myData[1].someData[3]");
> (I know it is wrong what I write above but I hope it explains
> my wish). In other words how to use Streamer method of storing
> of my Class(es) into root-file (branch) and the pointers (keys) to all
> fields of my Class (and also fields of component which is of another
> class type).
> 
> Thank you a lot in advance,
> 
> Witek Przygoda

The rules for the split mode are described at URL:
  http://root.cern.ch/root/HowtoWritTree.html
In particular, the split mode cannot support yet the case
of an array of objects. In your case, you also have a pointer
to an array of floats.
I have modified your classes, such that , at least, you can run
rootcint and use the non-split mode.
I also show an example of macro to generate a file with a Tree (wit.C)
and another macro to analyze this tree (wit2.C)

//file witold.h
#ifndef ROOT_event
#define ROOT_event

#include <TObject.h>


class MyDataClass : public TObject {

public:
int event;
int SimpleType;
int someData[10];
float *myArrayPtr;

// again I provide Streamer which takes care of pointer also

MyDataClass();
MyDataClass(const MyDataClass &dc);
virtual ~MyDataClass();
void SetType(int typ);
ClassDef(MyDataClass,1)
};

class MyHits : public TObject {

public:

int SimpleType;
MyDataClass finalData[3];

// of course Streamer function defined in .cc file

MyHits();
MyHits(const MyHits &mh);
virtual ~MyHits();
int Type(int n) {return finalData[n].SimpleType;}
int Data(int n, int i=0) {return finalData[n].someData[i];}
float ArrayPtr(int n, int i = 0) {return finalData[n].myArrayPtr[i];}
ClassDef(MyHits,1)

};

#endif


//--------file witold.C
#include "witold.h"

Int_t event = 0;
ClassImp(MyDataClass)

MyDataClass::MyDataClass()
{
   SimpleType = 1;
   for (int i=0;i<10;i++) someData[i]=i;
   myArrayPtr = new float[8];
   for (int j=0;j<8;j++) myArrayPtr[j] = 1.1+10.*j;
}
MyDataClass::MyDataClass(const MyDataClass &dc)
{
}

MyDataClass::~MyDataClass()
{
   delete [] myArrayPtr;
   myArrayPtr = 0;
}

void MyDataClass::SetType(int typ)
{
   SimpleType = typ;
}


ClassImp(MyHits)

MyHits::MyHits()
{
   SimpleType = 2;
}

MyHits::MyHits(const MyHits &mh)
{
}

MyHits::~MyHits()
{
}


//______________________________________________________________________________
void MyDataClass::Streamer(TBuffer &R__b)
{
   // Stream an object of class MyDataClass.

   if (R__b.IsReading()) {
      Version_t R__v = R__b.ReadVersion(); if (R__v) { }
      TObject::Streamer(R__b);
      R__b >> SimpleType;
      R__b.ReadFastArray(someData,10);
      if (!myArrayPtr) myArrayPtr = new float[8];
      R__b.ReadFastArray(myArrayPtr,8);
   } else {
      R__b.WriteVersion(MyDataClass::IsA());
      TObject::Streamer(R__b);
      R__b << SimpleType;
      R__b.WriteFastArray(someData, 10);
      R__b.WriteFastArray(myArrayPtr, 8);  // I put a constant here
   }
}

//______________________________________________________________________________
void MyHits::Streamer(TBuffer &R__b)
{
   // Stream an object of class MyHits.

   if (R__b.IsReading()) {
      Version_t R__v = R__b.ReadVersion(); if (R__v) { }
      TObject::Streamer(R__b);
      R__b >> SimpleType;
      int R__i;
      for (R__i=0; R__i < 3; R__i++)
         finalData[R__i].Streamer(R__b);
   } else {
      R__b.WriteVersion(MyHits::IsA());
      TObject::Streamer(R__b);
      R__b << SimpleType;
      int R__i;
      for (R__i=0; R__i < 3; R__i++)
         finalData[R__i].Streamer(R__b);
   }
}


//--------file wLinkDef.h
#ifdef __CINT__

#pragma link off all globals;
#pragma link off all classes;
#pragma link off all functions;


#pragma link C++ class MyDataClass-;
#pragma link C++ class MyHits-;

#endif


//--example of commands to run rootcint/compile on hp-ux
rootcint -f wcint.C -c witold.h wLinkDef.h
CC +a1 +z -I$ROOTSYS/include -c wcint.C witold.C
CC -b  -g +a1 -z wcint.o witold.o -o witold.so



//------example of macro to generate a Tree  (wit.C)
{
   gROOT->Reset();

   gSystem->Load("witold.so");
   TFile f("witold.root","recreate");
   MyHits *mh = new MyHits();
   TTree *T = new TTree("T","witold tree");
   T->Branch("hits","MyHits",&mh,32000,0);
   
   for (Int_t i=0;i<100;i++) {
     mh->finalData[0].SetType(i);
     mh->finalData[1].SetType(1000+i);
     mh->finalData[2].SetType(2000+i);
     T->Fill();
   }
   T->Print();
   T->Write();
}



//------example of macro toanalyze the Tree  (wit2.C)
{
   gROOT->Reset();
   gSystem->Load("witold.so");
   TFile f("witold.root");
   MyHits *mh = new MyHits();
   TTree *T = (TTree*)f.Get("T");
   T->SetBranchAddress("hits",&mh);

   TH1F *h1 = new TH1F("h1","Types",100,0,3000);
   Int_t nentries = T->GetEntries();
   for (Int_t i=0;i<nentries;i++) {
      T->GetEvent(i);
      h1->Fill(mh->Type(0));
      h1->Fill(mh->Type(1));
      h1->Fill(mh->Type(2));
   }
   h1->Draw();
}



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