THbookFile


class description - source file - inheritance tree

class THbookFile : public TNamed


    public:
THbookFile THbookFile() THbookFile THbookFile(const char* fname, Int_t lrecl = 1024) THbookFile THbookFile(const THbookFile&) virtual void ~THbookFile() virtual void Browse(TBrowser* b) virtual Bool_t cd(const char* dirname) static TClass* Class() virtual void Close(Option_t* option) virtual TObject* Convert1D(Int_t id) virtual TObject* Convert2D(Int_t id) virtual TFile* Convert2root(const char* rootname, Int_t lrecl = 0, Option_t* option) virtual TObject* ConvertCWN(Int_t id) virtual TObject* ConvertProfile(Int_t id) virtual TObject* ConvertRWN(Int_t id) void DeleteID(Int_t id) virtual TObject* FindObject(const char* name) const virtual TObject* FindObject(const TObject* obj) const TObject* Get(Int_t id) const char* GetCurDir() const Int_t GetEntry(Int_t entry, Int_t id, Int_t atype, Float_t* x) Int_t GetEntryBranch(Int_t entry, Int_t id) TList* GetList() const TList* GetListOfKeys() const Seek_t GetSize() const void InitLeaves(Int_t id, Int_t var, TTreeFormula* formula) virtual TClass* IsA() const virtual Bool_t IsFolder() const virtual Bool_t IsOpen() const virtual void ls(const char* path) const virtual void SetBranchAddress(Int_t id, const char* bname, void* add) virtual void ShowMembers(TMemberInspector& insp, char* parent) virtual void Streamer(TBuffer& b) void StreamerNVirtual(TBuffer& b)

Data Members


    protected:
Int_t fLun Fortran logical unit for this file Int_t fLrecl Record length in Hbook machine words TList* fList list of objects in memory TList* fKeys list of Hbook keys (Ids) on disk TString fCurDir name of current directory static Bool_t fgPawInit static Int_t* fgLuns

Class Description

  This class is an interface to the Hbook objects in Hbook files
  Any Hbook object (1-D, 2-D, Profile, RWN or CWN can be read
  NB: a THbookFile can only be used in READ mode
      Use the utility in $ROOTSYS/bin/h2root to convert Hbook to Root

 Example of use:
  gSystem->Load("libHbook");
  THbookFile f("myfile.hbook");
  f.ls();
  TH1F *h1 = (TH1F*)f.Get(1);  //import histogram ID=1 in h1
  h1->Fit("gaus");
  THbookTree *T = (THbookTree*)f.Get(111); //import ntuple header
  T->Print();  //show the Hbook ntuple variables
  T->Draw("x","y<0"); // as in normal TTree::Draw

  THbookFile can be browsed via TBrowser.



THbookFile() : TNamed()

THbookFile(const char *fname, Int_t lrecl) :TNamed(fname,"")
  Constructor for an HBook file object

~THbookFile()

void Browse(TBrowser *b)
 to be implemented

Bool_t cd(const char *dirname)
 change directory to dirname

void Close(Option_t *)
 Close the Hbook file

void DeleteID(Int_t id)

TObject* FindObject(const char *name) const
 return object with name in fList in memory

TObject* FindObject(const TObject *obj) const
 return object with pointer obj in fList in memory

TObject* Get(Int_t idd)
 import Hbook object with identifier idd in memory

Int_t GetEntry(Int_t entry, Int_t id, Int_t atype, Float_t *x)
 Read in memory all columns of entry number of ntuple id from the Hbook file

Int_t GetEntryBranch(Int_t entry, Int_t id)
 Read in memory only the branch bname

void InitLeaves(Int_t id, Int_t var, TTreeFormula *formula)
 This function is called from the first entry in TTreePlayer::InitLoop
 It analyzes the list of variables involved in the current query
 and pre-process the internal Hbook tables to speed-up the search
 at the next entries.

Bool_t IsOpen() const
 Returns kTRUE in case file is open and kFALSE if file is not open.

void SetBranchAddress(Int_t id, const char *bname, void *add)

TFile* Convert2root(const char *rootname, Int_t /*lrecl*/, Option_t *option)
 Convert this Hbook file to a Root file with name rootname.
 if rootname="', rootname = hbook file name with .root instead of .hbook
 By default, the Root file is connected and returned
 option:
       - "NO" do not connect the Root file
       - "C"  do not compress file (default is to compress)
       - "L"  do not convert names to lower case (default is to convert)

TObject* ConvertCWN(Int_t id)
 Convert the Column-Wise-Ntuple id to a Root Tree

TObject* ConvertRWN(Int_t id)
 Convert the Row-Wise-Ntuple id to a Root Tree

TObject* ConvertProfile(Int_t id)
 Convert an Hbook profile histogram into a Root TProfile

 the following structure is used in Hbook
    lcid points to the profile in array iq
    lcont = lq(lcid-1)
    lw    = lq(lcont)
    ln    = lq(lw)
      if option S jbyt(iq(lw),1,2) = 1
      if option I jbyt(iq(lw),1,2) = 2

TObject* Convert1D(Int_t id)
 Convert an Hbook 1-d histogram into a Root TH1F

TObject* Convert2D(Int_t id)
 Convert an Hbook 2-d histogram into a Root TH2F

void ls(const char *path) const
 List contents of Hbook directory



Inline Functions


        const char* GetCurDir() const
             Seek_t GetSize() const
             TList* GetList() const
             TList* GetListOfKeys() const
             Bool_t IsFolder() const
            TClass* Class()
            TClass* IsA() const
               void ShowMembers(TMemberInspector& insp, char* parent)
               void Streamer(TBuffer& b)
               void StreamerNVirtual(TBuffer& b)
         THbookFile THbookFile(const THbookFile&)


Author: Rene Brun 18/02/2002
Last update: root/hbook:$Name: $:$Id: THbookFile.cxx,v 1.17 2003/01/20 10:25:57 brun Exp $
Copyright (C) 1995-2002, Rene Brun and Fons Rademakers. *


ROOT page - Class index - Top of the page

This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.