*SPAM* Re: subbranches

From: Axel Naumann <axel_at_fnal.gov>
Date: Wed, 23 Mar 2005 11:56:14 -0600


Hi Chiara,

your macro has several problems (no class def, no pragma to generate your classes' dictionary needed for the TClonesArray, no default constructor needed to be able to read the objects from file, illegal c++, wrong allocation, ...). I've attached a (hopefully) corrected version.

To retrieve your histograms, just create a TClonesArray as you did in your macro, call TTree::SetBranchAddress on it, and after a TTree::GetEntry call all your objects will be filled just as they were when you wrote them, i.e. the lists will be filled with the proper objects, which you can access using (yet to be written) member functions of your classes.

Axel.

chiara zampolli wrote:
> Dear Axel,
>
> Thank you for your answer.
> Francesco and I have tried to write a short macro (see attachment)
> following your suggestions. After running it, what we find in the
> prova.root file (exploring it via the TBrowser) is a leaf (fStrips) with
> two entries, which we suppose to belong to the two plates we have
> defined in the macro itself. Now, is this correct? Are the histograms
> stored as we would like them to be, even if we do not see them? If this
> is the case, could you suggest us how to read/retrieve them?
>
> Thank you. Cheers,
> chiara & nofero.
>
> Axel Naumann wrote:
>

>> Hi Chiara,
>> one way of doing that is by creating an class which represents a
>> certain part of your hierarchy. I assume you have several entries of
>> (sector1, sector2) pairs, i.e. you'd fill your tree more than once
>> (otherwise you probably shouldn't use a TTree), and that you have
>> several histos per strip object (otherwise get rid of the TString
>> class, and store the list of histos directly in TPlate). In that case
>> this should work:
>>
>> class TStrip: public: TObject {
>> ...
>> private:
>>   std::list<TH1*> fHists;
>> };
>>
>> class TPlate: public TObject {
>> ...
>> private:
>> std::list<TStrip> fStrips;
>> };
>>
>> TClonesArray* ca=new TClonesArray("TPlate");
>>
>> TTree* t=new TTree("ch_nofero", "tree of sectors, plates, and strips");
>> t->Branch("sectors", &ca);
>>
>> You could also create one branch for each strip, named e.g.
>> "Sector1Plate4Strip2". This is pretty efficient, but a bit ugly.
>>
>> Axel.
>>
>> chiara zampolli wrote:
>>
>>> Dear rooters,
>>>
>>>     we are trying to create a tree with four levels of subbranches,
>>> with the last one made up of histograms (TH1F):
>>>
>>>   LEVEL1:               sector1                                     
>>> sector2
>>>                                        
>>> |                                               |
>>>                                    /       
>>> \                                     /         \
>>>   LEVEL 2:           plate1   plate2                          
>>> plate1   plate2
>>>                                   |
>>>                               /
>>>   LEVEL3:        strip1  ....         .....                        
>>> .......        ......
>>>                             |
>>>   LEVEL4:      histo1
>>>
>>> and so on. Could you suggest me how to create such a tree? Since we
>>> are not able to split more than once the tree we have defined.
>>> Till now, we have tried using the TClonesArray class and the
>>> TTree->Branch(....) method.
>>>
>>>      Thank you.
>>>      Cheers,
>>>          ch & nofero
>>>

>
>
>
> ------------------------------------------------------------------------
>
> #include "TFile.h"
> #include "TObject.h"
> #include "TH1F.h"
> #include "TTree.h"
> #include <list>
> //_____________________________________________________________________________
>
> class TStrip: public TObject {
>
> public:
> TStrip(const char *name){};
> void pushback(TH1F *histog){fHists.push_back(histog);};
> private:
> std:list<TH1F*> fHists;
> };
>
> class TPlate: public TObject {
>
> public:
> TPlate(const char *name){};
> void pushback(TStrip *strip){fStrips.push_back(strip);};
>
> private:
> std:list<TStrip*> fStrips;
>
> };
>
> void lista(){
>
> // Char_t hnameMEAN[100];
> TH1F *hist[5];
> TTree *fTreeEq;
> fTreeEq = new TTree("TreeEq","Tree of Plates and Sectors");
>
> TStrip *Strips[3];
> TPlate *Plate[2];
> TH1F *hist[5];
> TH1F *hist[0] = new TH1F("hist0", "hist0", 100,-1.,1.);
> TH1F *hist[1] = new TH1F("hist1", "hist1", 100,-1.,1.);
> TH1F *hist[2] = new TH1F("hist2", "hist2", 100,-1.,1.);
> TH1F *hist[3] = new TH1F("hist3", "hist3", 100,-1.,1.);
> TH1F *hist[4] = new TH1F("hist4", "hist4", 100,-1.,1.);
>
> char stripname[10];
> char platename[10];
>
> TClonesArray *ca=new TClonesArray("TPlate");
> TClonesArray &aca=*ca;
>
> for (Int_t k=0; k<2; k++){
> sprintf(platename,"plate%i",k);
> new(aca[k]) TPlate (platename);
> Plate[k] = new TPlate(platename);
> for (Int_t i=0; i<1; i++){
> sprintf(stripname,"strip%i",i);
> Strips[i] = new TStrip(stripname);
> for (Int_t j=0; j<5; j++){
> Strips[i]->pushback(hist[j]);
> }
> Plate[k]->pushback(Strips[i]);
> }
> }
>
> TFile *fileout = new TFile("prova.root","recreate");
> fTreeEq->Branch("sectors", &ca);
> fTreeEq->Fill();
> fTreeEq->Write();
> fileout->Close();
>
> }
>

#include "TFile.h"
#include "TObject.h"
#include "TH1F.h"
#include "TTree.h"
#include <list>

//_____________________________________________________________________________

class TStrip: public TObject {   

public:
  TStrip(const char *name=0){};
  ~TStrip() { // delete all hists here - well, not in your example case where you have 5 global ones   }
  void pushback(TH1F *histog){fHists.push_back(histog);}; private:
  std::list<TH1F*> fHists;
  ClassDef(TStrip,1); // a strip
};

class TPlate: public TObject {   

public:
  TPlate(const char *name=0){};
  ~TPlate() { // delete all entries in fStrips here   }
  void pushback(TStrip *strip){fStrips.push_back(strip);};   

private:
  std::list<TStrip*> fStrips;
  ClassDef(TPlate,1); // a plate
};

#ifdef __MAKECINT__
#pragma link C++ class TStrip+;
#pragma link C++ class TPlate+;
#endif


void lista(){      

  TH1F *hist[5];

  hist[0] = new TH1F("hist0", "hist0", 100,-1.,1.);
  hist[1] = new TH1F("hist1", "hist1", 100,-1.,1.);
  hist[2] = new TH1F("hist2", "hist2", 100,-1.,1.);
  hist[3] = new TH1F("hist3", "hist3", 100,-1.,1.);
  hist[4] = new TH1F("hist4", "hist4", 100,-1.,1.);
  

  char stripname[10];
  char platename[10];

  TFile *fileout = new TFile("prova.root","recreate");

  TClonesArray *ca=new TClonesArray("TPlate");   TClonesArray &aca=*ca;

  TTree *fTreeEq = new TTree("TreeEq","Tree of Plates and Sectors");   fTreeEq->Branch("sectors", &ca);

  // for (Int_t entry=0; entry<10; entry++) {   for (Int_t k=0; k<2; k++){
    sprintf(platename,"plate%i",k);
    TPlate *Plate = new(aca[k]) TPlate (platename);     for (Int_t i=0; i<1; i++){

      sprintf(stripname,"strip%i",i);
      TStrip *Strip = new TStrip(stripname);
      for (Int_t j=0; j<5; j++){
	Strip->pushback(hist[j]);
      }
      Plate->pushback(Strip);

    }
  }
  fTreeEq->Fill();
  // for (Int_t k=0; k<2; k++){
  //   TPlate *Plate = (TPlate*)aca[k];
  //   Plate->Clear() // deletes all entries in the list
  // }
  // ca->Clear();
  

  fTreeEq->Write();
  delete fileout;

} Received on Wed Mar 23 2005 - 18:56:35 MET

This archive was generated by hypermail 2.2.0 : Tue Jan 02 2007 - 14:45:06 MET