Re: [ROOT] Re: Root Questions

From: Valeri Fine (Faine) (fine@bnl.gov)
Date: Fri Feb 23 2001 - 20:25:58 MET


Hello, Rene

We found it is very common task to iterate TFile object.
At the moment I am debugging my TFileIter class. It works.
If one finds it is a general interest I can spend some extra time
to make this class free of the particular experiment specific to contribute.

I am taking a liberty to post the very preliminary version of such iterator interface
expecting the useful feedback from the rest ROOT people.

At present the this  iterates  the TFile keys by  the key names. It can return things
by the class name those keys point to or by the object names as well. Something like:

    TObject *NextClassObjectGet(const char *className);

    TObject *NextObjectByNameGet(const char *objectName);



#ifndef ROOT_TFileIter
#define ROOT_TFileIter

///////////////////////////////////////////////////////////////////////////
//                                                                       //
// Class to iterate (read / write ) the events written to /from TFile.   //
// The event is supposed to assing an unique ID in form of               //
//                                                                       //
// TKey name ::= event Id ::= eventName "." run_number "." event_number  //
//                                                                       //
// and stored as the TKey name of the object written                     //
//                                                                       //
//  author Valeri Fine                                                   //
//                                                                       //
///////////////////////////////////////////////////////////////////////////

#include "TString.h"
#include "TIterator.h"

class TObject;
class TDirectory;
class TFile;
class TIter;

class TFileIter : public TIterator {

  private:

    TFile      *fFileBackUp;
    TDirectory *fDirectoryBackUp;

  protected:
    TFile   *fRootFile;
    TObject *fObj;
    TString  fEventName;
    Int_t    fRunNumber;
    Int_t    fEventNumber;
    TIter   *fKeysIterator;

  protected:
    void     SaveFileScope();
    void     RestoreFileScope();
    TKey    *NextEventKey(Int_t eventNumber=-1, Int_t runNumber=-1, const char *name="*");

  public:

    TFileIter(TFile *file=0);
    virtual ~TFileIter();
    virtual TFile *GetTFile(){ return fRootFile; }
    virtual TObject *NextEventGet(Int_t eventNumber=-1, Int_t runNumber=-1, const char *name="*");
    // virtual Int_t    NextEventPut(TObject *obj, Int_t eventNum,  Int_t runNumber, const char *name);
    TObject *operator++();
    TObject *operator*() const;

  public:  // abstract TIterator methods implementations:

    virtual const TCollection *GetCollection() const;
    virtual TObject *Next();
    virtual void Reset();
    TObject *operator()();

};

const TCollection *TFileIter::GetCollection() const
{ return fKeysIterator ? fKeysIterator->GetCollection() : 0; }

TObject *TFileIter::Next()
{ return NextEventGet(); }

void TFileIter::Reset()
{ if (fKeysIterator) fKeysIterator->Reset(); }

TObject *TFileIter::operator()(){ return Next(); }

inline void  TFileIter::SaveFileScope()
{ fFileBackUp = gFile; fDirectoryBackUp = gDirectory; }

inline void TFileIter::RestoreFileScope()
{  gFile = fFileBackUp; gDirectory = fDirectoryBackUp; }

inline TObject *TFileIter::operator*() const { return fObj; }

inline TObject *TFileIter::operator++()
{ return Next(); }

#endif


Any idea  ?

 Send me directly fine@bnl.gov if appropriated.


                  Valeri

----- Original Message -----
From: Rene Brun <Rene.Brun@cern.ch>
To: corey a stambaugh <stambaug@cis.ohio-state.edu>; <roottalk@pcroot.cern.ch>; <rootdoc@pcroot.cern.ch>
Sent: 23 февраля 2001 г. 4:07
Subject: [ROOT] Re: Root Questions


> Hi Casey,
> This is a message of public interest, so I post it to roottalk.
> In the attachment you will find a file work.C that fills the missing pieces
> in your code. I have added this example to the tutorials.
>
> Suzanne, could you find a place in the Users Guide to add this example?
>
> Rene Brun
>
> corey a stambaugh wrote:
> >
> > Hello,
> >         My name is Casey Stambaugh and I am using Root at Ohio State
> > Univ. for the CDF group.
> >
> >         I apologize for using this email address, I could not figure out
> > how to use the roottalk digest.
> >
> >         I am trying to do a simple root program to run through a .hbk file
> > and print out all the graphs.
> >
> >         This is my current code
> >
> > void work(int x, int y, int page){
> >
> >    int page_cnt, hist_per_page, i;
> >    TFile *f1 = new TFile("xftRoot_out.hbk");
> >    TCanvas *c1 = new TCanvas("c1");
> >    TPostScript *ps = new TPostScript("file.ps",112);
> >    c1->Divide(x,y);
> >    hist_per_page = x*y;
> >
> >         while (page_cnt < page)
> >         {
> >                 ps->NewPage();
> >                 i=1;
> >                         while (hist_per_page >= i)
> >                         {
> >                                 c1->cd(i);
> >                                 ??????->Draw();
> >                                 i++;
> >                         }
> >                 c1->Update();
> >         page_cnt++;
> >    }
> >    ps->Close();
> >    fgSystem->Exec("lpr -Psmith2079 file.ps");
> >    gSystem->Exec("gv file.ps");
> >
> >    }
> >
> > I was wondering if it was possible to somehow extract individual
> > histograms  to the ????? with a simple call or somehow run through
> > histograms in this manner.
> >
> > If this is not possible, that will suffice but I have spent many days
> > reading the root page finding nothing to support or contradict my plans.
> >
> > Thank you for your time and I apologize again for using this email
> >
> > Casey Stambaugh
> > stambaug@cis.ohio-state.edu


--------------------------------------------------------------------------------


> void work(int x, int y, int
e){ 
> file://example of script to plot all histograms in a Root file
> file://on a Postcript file (x times y per page)
> file://The following logic can be used to loop on all the keys of a file:
> // TFile f("myfile.root");
> // TIter next(f1->GetListOfKeys());
> // TKey *key;
> // while ((key = (TKey*)next())) {
> //    file://key->GetClassName() returns the name of the object class
> //    TObject *obj = key->ReadObj(); file://read object from file
> //    file://obj->ClassName() should be equal to key->GetClassName()
> //    file://obj->InheritsFrom("someclass") test if obj inherits from someclass
> // }    
>  int page_cnt, hist_per_page, i;  
>  TFile *f1 = new TFile("hsimple.root");
>  TCanvas *c1 = new TCanvas("c1");
>  TPostScript *ps = new TPostScript("file.ps",112);
>  c1->Divide(x,y);
>  hist_per_page = x*y; 
>  TIter next(f1->GetListOfKeys());       file://make an iterator on list of keys
>  TKey *key;
>  while (page_cnt < page) {
>     ps->NewPage();
>     i=1;
>     while (hist_per_page >= i) {
>        c1->cd(i);
>        key = (TKey*)next();             file://get next key on the file
>        if (!key) break;                 file://if no more keys, key=0
>        TObject *obj = key->ReadObj();   file://read object associated to key
>        if (obj->InheritsFrom("TH1")) {  file://interested by histograms only
>           obj->Draw();                  file://draw histogram with default option
>           i++;
>        }
>     }
>     c1->Update();
>     if (!key) break;
>     page_cnt++;     
>  }
>  ps->Close();
>  file://gSystem->Exec("lpr -Psmith2079 file.ps");
>  file://gSystem->Exec("gv file.ps");
> }
> 



This archive was generated by hypermail 2b29 : Tue Jan 01 2002 - 17:50:37 MET