[Brahms-dev-l] Bdst Looping - a HOWTO

From: Bjorn H Samset <bjornhs@rcf2.rhic.bnl.gov>
Date: Wed Mar 24 2004 - 05:12:22 EST
Dear dev'ils.

Last night I got an email from Radek asking how to read and analyze our
official dsts. This brought home the fact that it was high time we added
some example code to bdst to help people get started on this - it's one of
those things that's easy once you know the tricks... So I've done the
following:

In bdst (get it by saying first 'klog' and then 'cvs checkout bdst')
you'll find a new directory
bdst/looper
Here I've put an example looper class
bdst/looper/BdstExampleLooper.C
bdst/looper/BdstExampleLooper.h
and a file
bdst/looper/BdstLooperHowto.txt
that details the following:
1) Drawing and projecting - the simplest way
2) Using a looper class - the way to do proper analysis
3) How to make a general looper class

Part 1 is not meant to replace the chapter on trees in the ROOT users
guide (http://root.cern.ch/root/doc/RootDoc.html) but it includes some
local information. Part 2 describes how you use my example class, which
can be copied and modified to suit your own analysis. Part 3 shows how I
used ROOTs builting functions to make the class - use it as a tutorial of
you want to make your own or want to loop over another type of ROOT tree.

I'm sure the Howto file could be better and more detailed - Radek, could
you check it out and let me know where you run into problems? I'll expand
the file as necesseary.

I hope this can encoureage more people to analyze the official dsts :-)

A note to bdst developers: To keep this example class operational, any new
data members in bdst should also be added here. It's backwards compatible,
but not forwards...

For completeness (and to have it stored in one more place...) I add here
the full text of BdstLooperHowto.txt:

*************************************************************************
*************************************************************************
How to read a DST!
A short guide by Bjorn H. Samset <bjornhs@fys.uio.no>
*************************************************************************
*************************************************************************
Contents:
1) Drawing and projecting - the simplest way
2) Using a looper class - the way to do proper analysis
3) How to make a general looper class
*************************************************************************

NB: This file assumes that you have a file with a dst generated by the
bdst code. If you're looking at an official dst, this should always be
the case.

*************************************************************************
1 - DRAWING AND PROJECTING
*************************************************************************

Say you have a file with a dst, called dst009999v02.root, i.e. a dst
from run 9999, version 02. To draw variables from it, open a root
session and do the following:

Root> TFile *f = new TFile("dst009999v02.root");
Root> TTree *t = f->Get("T1");

'T1' is our default tree name! Now you can look at it using a
TBrowser, or plot from the command line saying e.g.

Root> t->Draw("G.fBbVtxZ","G.fCent<20");
Root> TH2F *h1 = new TH2F("h1","h1",200,-10,10,200,0.8,1.5);
Root>
t->Project("h1","1/MRS.fTofw.fBeta:MRS.fD5.fP","MRS.fD5.fStatus==1");

For details on this, refer to the chapter on trees in the ROOT Users
Guide:
ftp://root.cern.ch/root/doc/chapter12.pdf

*************************************************************************
2 - USING A LOOPER CLASS
*************************************************************************

The most general way of reading data from a tree is to loop over its
entries using a pre-generated looper class. There is an example that
you can modify in CVS, in the same dir. as this file:
bdst/looper/BdstExampleLooper.C
bdst/looper/BdstExampleLooper.h

To use it, do the following:
* Start a root session and load the class:
  brat [0] .L BdstExampleLooper.C

* Make an instance of the class BdstExampleLooper, giving the filename
  as argument:
  brat [1] BdstExampleLooper *b = new
BdstExampleLooper("dst009999v02.root");
  You may get a number of messages, some of the type
    Warning in <TClass::TClass>: no dictionary for class BdstRun is
      available
  and some like this:
    Error in <TTree::SetBranchStatus>: unknown branch ->
FS.fD1EntranceLine.fOrigin.fX
  These can safely be ignored.

* Loop over it by calling its Loop() function:
  brat [2] b->Loop();

You can also add a number of dsts to the same looper and analyze them
together:
   brat [1] BdstExampleLooper *b = new
BdstExampleLooper("/my/dst/dir/dst*v02.root");
This would make a looper for all dsts of version 2 in the directory
/my/dst/dir. Note that the dsts must all be in the same dir. for this
to work - the usual trick to do this is to make a dir with softlinks
to all the dsts you want to analyze - e.g.
> mkdir /my/dst/dir
> cd /my/dst/dir
> ln -s /brahms/data07/data/run03/dau/200/r*/dst/dst*v09.root .
This gives you links to all the official d+Au dsts of version 9 in
/my/dst/dir.

The function Loop() is defined in the file BdstExampleLooper.C - this
is where you will code your analysis. It's just a skeleton method that
reads information from the tree per event, so it should have this
general structure:
1) Initialize variables, make histograms etc.
2) The loop, where you fill histograms, calculate rapidity, make cuts
   etc.
3) Finalizing part, where you save the histograms to a file or just
   plot them.

NB: All the members in a bdst tree are arrays! In a looper class, you
access global and run info like this:
   cout << G_fBbVtxZ[0] << endl;
i.e. as the first entry in an array of length one. To loop over
information from each track in the FS for this event, do the
following:
   for (Int_t i = 0; i<FS_; i++)
     cout << FS_fVtxZ[i] - G_fBbVtxZ[0] << endl;
which gives you the track projection distance to the BB vertex for
each track.



*************************************************************************
3 - HOW TO MAKE A GENERAL LOOPER CLASS
*************************************************************************

This section details how I went about in creating the class
BdstExampleLooper. It's very easy since ROOT does all the hard stuff -
just follow my lead if you want a general looper for another type of tree.

1) Manually load a tree of the type you want a looper for:

brat [0] TFile *f = new TFile("dst006411.root");
Warning in <TClass::TClass>: no dictionary for class BdstRun is available
Warning in <TClass::TClass>: no dictionary for class BdstGlobal is
available
Warning in <TClass::TClass>: no dictionary for class BdstFsTrack is
available
Warning in <TClass::TClass>: no dictionary for class BdstTrack is
available
Warning in <TClass::TClass>: no dictionary for class BdstMagnet is
available
Warning in <TClass::TClass>: no dictionary for class BdstTof is available
Warning in <TClass::TClass>: no dictionary for class BdstC1 is available
Warning in <TClass::TClass>: no dictionary for class BdstRich is available
Warning in <TClass::TClass>: no dictionary for class BdstMrsTrack is
available
Warning in <TClass::TClass>: no dictionary for class BdstC4 is available
Warning in <TClass::TClass>: no dictionary for class BdstTpm1Track is
available
brat [1] TTree *t = f->Get("T1");


2) Call the function MakeClass() with the name of the looper as argument:
brat [2] t->MakeClass("BdstExampleLooper");
Info in <TTreePlayer::MakeClass>: Files: BdstExampleLooper.h and
BdstExampleLooper.C generated from Tree: T1

3) You now have an almost complete looper, except that by default it's
   tied to the tree you loaded. Make the following minor changes in
   the header file (BdstExampleLooper.h in my case) to fix this:

a) Increase the array limits, in case other files have more entries. Cange

   const Int_t kMaxR = 1;
   const Int_t kMaxG = 1;
   const Int_t kMaxFS = 4;
   const Int_t kMaxMRS = 4;
   const Int_t kMaxTPM1 = 1;

to e.g.

   const Int_t kMaxR = 1;
   const Int_t kMaxG = 1;
   const Int_t kMaxFS = 20;
   const Int_t kMaxMRS = 20;
   const Int_t kMaxTPM1 = 50;

(You have to think a bit here - how many tracks to you forsee per
event?)

b) Query-replace 'TTree' with 'TChain' (6 places) all through the header
file.

c) Change the class constructor to take the filename as an argument.
Change

   BdstExampleLooper(TChain *tree=0);

to

   BdstExampleLooper(Char_t *input, TChain *tree=0);

in the declarations and change

   BdstExampleLooper::BdstExampleLooper(TChain *tree)

to

   BdstExampleLooper::BdstExampleLooper(Char_t *input, TChain *tree)

in the function itself (still in the header file!)

d) In the creator, make the following change:

  if (tree == 0) {
      TFile *f =
(TFile*)gROOT->GetListOfFiles()->FindObject("dst006411.root");
      if (!f) {
         f = new TFile("dst006411.root");
      }
      tree = (TChain*)gDirectory->Get("T1");

   }

should be simplified to this:

  if (tree == 0) {
    tree = new TChain("T1");
    tree->Add(input);
  }


That's it :-) You now have a general looper class that you can get to
loop over one or more dsts by giving an appropriate input
argument. Add any analysis code you want in BdstExampleLooper.C and
use it as detailed in part 2 above.

Happy analyzing :-)

--Bjrn


--
Bjorn H. Samset                           Phone: 22856465/92051998
PhD student, heavy ion physics            Adr:   Schouterrassen 6
Inst. of Physics, University of Oslo             0573 Oslo
                              \|/
----------------------------> -*- <-----------------------------
                              /|\

_______________________________________________
Brahms-dev-l mailing list
Brahms-dev-l@lists.bnl.gov
http://lists.bnl.gov/mailman/listinfo/brahms-dev-l
Received on Wed Mar 24 06:20:44 2004

This archive was generated by hypermail 2.1.8 : Wed Mar 24 2004 - 06:21:05 EST