STAR   Computing   Documentation  
StStrangeMuDstMaker: the Maker
Offline Software using ROOT in STAR Maintained by G. Van Buren
Last modified Tue Apr 11 15:30:19 2000
StStrangeMuDstMaker


The Maker

The Code
Specifying Content
Creating
Reading
Filtering
Tools
Buffer Sizes
Adding Content


The Code

The self-documented code can be viewed through doxygen here.


Specifying Content for Micro-DSTs

After instantiating a StStrangeMuDstMaker, the user should tell the maker which components they want to work with. No components are turned on by default. Component content is controlled via the Do*() member functions: For example, calling the DoV0() member function will cause V0 components to be created/read for the micro-DST associated with this maker. In read mode, the user need not turn on components they will not be using even if they are on the micro-DST on file. Here is a brief example:
  StStrangeMuDstMaker strangeDst();
  strangeDst.DoV0();

Note that DoMc() turns on processing of the Monte Carlo data and association branches for each of the content types turned on by any of the first three Do*() member functions.


Creating Micro-DSTs

Create mode is the default mode for the StStrangeMuDstMaker maker. The maker and supporting code libraries must be loaded, and the maker should be used in a standard STAR chain which provides an StEvent structure from which the micro-DST can be made. Currently, persistent StEvent is not supported, so a standard DST should be read in, and StEventMaker should be used to build the StEvent before the StStrangeMuDstMaker. The user need not do anything during the event loop of the chain as the filling of the micro-DST from StEvent is done automatically. A good example is provided in the macro:
$STAR/StRoot/macros/analysis/makeStrangeMuDst.C,
where something like the following is used:
  StStrangeMuDstMaker strangeDst();
  strangeDst.DoV0();     // Selects V0 vertices for micro-DST
  strangeDst.DoXi();     // Selects Xi vertices for micro-DST
  strangeDst.DoMc();     // Selects MC info for both V0 and Xi for micro-DST
  strangeDst.SetWrite(); // Sets "write" mode (using default filenames)

Creating on File

The above macro demonstrates creating a micro-DST from StEvent and saving it to a file. Saving to file is specified with the SetWrite() member function of the maker. This member function should be called after the maker is instantiated, but before the Init() member function of the chain is called. SetWrite() takes up to four paramters, but only the first is important at present. These four parameters are the filenames associated with the four component branches of the micro-DST (described in the section on Structure of the Micro-DST). Each filename has a default, and not entering parameters leads to the use of the default values. The root of the TTree itself is kept in the Event component file, which defaults to "evMuDst.root".

Creating in Memory

This is the default creation mode for the maker. If no Set*() member functions are called, a micro-DST is created, and stored in memory. It is not written to file. It can be accessed via the maker's GetTree() member function. See the Tools from StStrangeMuDstMaker section below. A word of caution: the amount of memory that the micro-DST uses as it grows from event to event is not limited. It is possible to use all the available memory on the machine, leading to a crash and loss of data. This option is thus good for relatively small micro-DSTs; larger ones should either be written to disk (see Creating on File above), or erased after each event (see Persistence in Memory below).

Persistence in Memory

When creating a micro-DST, it may be that the user does not want to retain information from event to event. For example, the HBT maker might wish to use the strangeness micro-DST and its interface, but might need it only on an event-by-event basis. There is no need in this case to let the micro-DST grow, consuming more and more memory as each event progresses. The SetNoKeep() member function implements this functionality, causing the micro-DST to be erased after each event. Usage is as follows:
  StStrangeMuDstMaker strangeDst();
  strangeDst.DoV0();
  strangeDst.SetNoKeep();

Reading Micro-DSTs

Read mode is set by using the SetRead() member function. As with the other Set*() member functions, this should be called between the instantiation of the maker and the Init() member function of the chain. Also, as with the SetWrite() member function, up to four parameters can be entered specifying the files on which the micro-DST is located. The user should, in general, need only specify one file: that of the event component branch, where the root of the TTree used in the micro-DST resides. This file defaults to "evMuDst.root" if no parameters are given to SetRead().

If the user is interested in reading in several micro-DSTs, the use of the StFile class is employed. StFile can take as an argument in its constructor an array of file names, any of which may also include wild cards (as in "ev*Dst.root"). Passing SetRead() a pointer to the StFile class then teaches StStrangeMuDstMaker what files to read. Here is an example:

  char* files[3];
  files[0] = "/bigfiles/run471.evMuDst.root";
  files[1] = "/bigfiles/run475.evMuDst.root";
  files[2] = "/oldfiles/run*.evMuDst.root";
  StFiles* flist = new StFile(files);
  StStrangeMuDstMaker strangeDst();
  strangeDst.SetRead(flist);

Access to the micro-DST is then available inside the event loop. Member functions outlined in the section below on Tools from StStrangeMuDstMaker allow the user to find out how many V0's, Xi's, or Kink's, as well as their respective Monte Carlo counterparts, are in the given event and access individual instances of the micro-DST component classes. The user can then loop over all vertices from an event, and can call member functions from the micro-DST component classes. A good example is provided in:
$STAR/StRoot/macros/analysis/readStrangeMuDst.C,
with an event loop looking something like this:

  Int_t Nevents = 10;
  // Loop over events
  for( Int_t i=0; i < Nevents; i++ ) {
    if( chain.Make(i) ) break;
    // Loop over V0 vertices
    for( Int_t j=0; j < strangeDst.GetNV0(); j++ ) {
      StV0MuDst *v0m = strangeDst.GetV0(j);
      hist->Fill(v0m->massLambda());
    }
    if( i != Nevents) chain.Clear();
  }
  chain.Finish();

Additionally, the user need not examine every event in the micro-DST. If the user specifies an integer parameter in the Make() member function call for the chain, then the event associated with that integer is read from the micro-DST. For example, if the event loop runs from i = 5...9, and the chain is run with chain.Make(i), then the 6th through 10th events on the file will be read in (the first event on file is event zero). This allows the user to skip events or read in only a portion of the events on the micro-DST. If no parameter is given to Make(), the micro-DST file is read from the beginning.

Monte Carlo Data

Monte Carlo data is accessed similarly using functions with "Mc" in their names. Associations between Monte Carlo and reconstructed data can be accessed via functions with "Assoc" in their names. Examining associated vertices may look something like this:
  Int_t Nevents = 10;
  // Loop over events
  for( Int_t i=0; i < Nevents; i++ ) {
    if( chain.Make(i) ) break;
    // Loop over V0 vertices
    for( Int_t j=0; j < strangeDst.GetNV0Assoc(); j++ ) {
      StStrangeAssoc* v0a = strangeDst.GetV0Assoc(j);
      StV0MuDst* v0m = strangeDst.GetV0(v0a->indexRecoArray());
      StV0Mc* v0mc = strangeDst.GetV0Mc(v0a->indexMcArray());
      hist->Fill(v0m->ptV0() - v0mc->ptV0());
    }
    if( i != Nevents) chain.Clear();
  }
  chain.Finish();

Code can be written that works identically for reconstructed and Monte Carlo data using the Get*I() functions and interface classes:

  Bool_t MC = kTRUE; // kFALSE for reconstructed data
  Int_t Nevents = 10;
  // Loop over events
  for( Int_t i=0; i < Nevents; i++ ) {
    if( chain.Make(i) ) break;
    // Loop over V0 vertices
    for( Int_t j=0; j < strangeDst.GetNV0I(MC); j++ ) {
      StV0I *v0i = strangeDst.GetV0I(j,MC);
      hist->Fill(v0i->pseudoRapV0());
    }
    if( i != Nevents) chain.Clear();
  }
  chain.Finish();


Filtering Micro-DSTs

There are certainly times when only a portion of a micro-DST is of interest to a user. One then may wish to create a sub-micro-DST (or nano-DST) containing only that portion which is of interest. This allows greater portability for that portion of the micro-DST, and also improves access speed. To do this, two makers are needed: one to create the new sub-micro-DST, and one to read an already existing micro-DST. The maker for the new data must be instantiated before the maker for the existing data. After both have been instantiated, the maker for the new data must be told that it will be a subset of the old - this is done via the SubDst() member function, and should be implemented in the following style:
  StStrangeMuDstMaker strangeNewDst("strangeNewMuDst");
  strangeNewDst.DoV0();
  strangeNewDst.SetWrite("newEvMuDst.root");

  StStrangeMuDstMaker strangeOldDst("strangeOldMuDst");
  strangeOldDst.SetRead("oldEvMuDst.root");

  strangeNewDst.SubDst(strangeOldDst);
The filenames used above are up to the user, and different combinations of strangeness micro-DST components other than simply the V0 components can be used. Note that DoV0() need not be called for the old micro-DST maker, as the new one knows to at least try to make sure that the reading of the V0's is turned on.

If one is filtering the micro-DST, new or additional cuts are probably also being used. If this is the case, one should consider adding cuts information to the strangeness cuts list stored with the micro-DST. More information on how to do this is available for the cuts.

The selection process for what goes into the sub-micro-DST can be viewed under two paradigms. The first paradigm encompasses a user who wants to take only those vertices which satisfy some criteria, and the second paradigm is for those who want to keep events which satisfy a criteria. Both selection processes are done during the event loop, between the Make() and Clear() member function calls of the chain. The event for the old micro-DST is read in during Make(), selections are made, and the new sub-micro-DST is written during Clear(). Selections are made via the Select*() member functions:

Vertex Selection Paradigm:
In this model, only certain vertices from within an event are selected to go on the sub-micro-DST. The vertices are selected with the following six member functions: In these calls, j specifies the index of the vertex within that event to be selected. For example, one might like to select only those cascades with a transverse momentum above 500 Mev/c. This would be done in a loop as follows:
  // Loop over events
  for( Int_t i=0; i < Nevents; i++ ) {
    if( chain.Make(i) ) break;
    // Loop over Xi vertices
    for( Int_t j=0; j < strangeOldDst.GetNXi(); j++ ) {
      StXiMuDst *xim = strangeDst.GetXi(j);
      if (xim->ptXi() > 0.5) strangeNewDst.SelectXi(j);
    }
    if( i != Nevents) chain.Clear();
  }
  chain.Finish();
The user can also specify that all of a specific vertex type can be selected for a particular event by setting the input parameter j to be a negative number. Thus, SelectV0(-1) would select all V0 candidates in the event, but would select no other vertex types for copying.

Similarly, one can unselect a specific vertex that one has already selected using the Unselect*() member functions. This is particularly useful if there are only a few vertices one would not wish to select, in which case one would use Select*(-1) followed by a few specific uses of Unselect*(j). Again, a negative input unselects the entire event.

One note of caution: Select*(j) does not check to see if a vertex has already been selected, so selecting one twice results in two copies.

Event Selection Paradigm:
If one wishes to select events as a whole for a sub-micro-DST, this can be done with the SelectEvent() member function. This essentially copies the entire event from the old micro-DST to the new sub-micro-DST, but is only true to the extent that what micro-DST components are turned on for the new maker (via the Do*() member functions, see Specifying Content for Micro-DSTs). Here is an example where the user takes all events which have more than 10 Kinks:
  // Loop over events
  for( Int_t i=0; i < Nevents; i++ ) {
    if( chain.Make(i) ) break;
    if (strangeOldDst.GetNKink() > 10) strangeNewDst.SelectEvent();
    if( i != Nevents) chain.Clear();
  }
  chain.Finish();
This is really the equivalent of using the vertex selection paradigm with input parameter j = -1 for all active components.

It is important to understand that an event which has not been selected will still have its event branch data copied to the filtered micro-DST, even if none of the secondary vertex data is copied. UnselectEvent() will restore any event in which any selections have already been made to this same state of copying the event branch only. If you truly want to eliminate an event entirely from the filtered micro-DST, then use AbortEvent(). Here is an example:

    StStrangeEvMuDst* ev = strageOlDst.GetEvent();
    if (fabs(ev->primaryVertexZ())<100.) {
      strangeNewDst.SelectEvent();
    } else {
      strangeNewDst.AbortEvent();
    }
Both paradigms are shown in the following macro, with the vertex selection paradigm scheme run by default, and an event selection example commented out:
$STAR/StRoot/macros/analysis/filterStrangeMuDst.C.


Tools from StStrangeMuDstMaker


Buffer Sizes

The maker allocates a certain amount of buffer memory for each variable stored. Optimally, one would like this buffer to be as large as possible, reducing the time spent for I/O, and consequently reducing the number of keys Root uses in its files (thereby reducing the file size). However, with limited memory available on computers running the job, making this buffer large makes the executable very large, and may lead to significant virtual memory paging (disk swapping), thereby reducing the performance. If other users are running jobs, the amount of available memory may be very limited. Default values of the buffer sizes have been chosen which accomodate the usage of about 200 MB of memory for running with all branches activated. But in some situations, the user may wish to changes these limits. To do this, the following functions are provided: where b is the buffer size in bytes for each of the variables associated with the particular branch. The default values are 1024000 (1 MB) for the V0 and Xi branches, and 128000 (128 KB) for the Kink branch. The event branch currently uses 64 KB per variable. The total memory used is then a function of which branches are activated, as well as how many variables are in each branch. The above functions must be called between the instantiation of the maker, and the call to Init().


Adding Content

You might wish to add more information to your own MuDsts. This is done quite simply using the filtering code. An example is provided here:

Let's say I want to add a branch which contains a little bit of information that describes an event. In my example, I've made a new class called "StA" with two values that describe the event. This is what its header file looks like:

#ifndef StA_hh
#define StA_hh
#include "TObject.h"
class StA : public TObject {
  public:
   StA() {}
   virtual ~StA() {}
   Int_t value1;
   Float_t value2;
   ClassDef(StA,1)
};
#endif
Now I must modify the filterStrangeMuDst.C macro to include the new information. Here are 3 steps:

  1. Load the library for the new class:
      gSystem->Load("myNewLibrary");
    

  2. Assign a new branch (I called mine "ABranch") for the new MuDst after the chain.Init() call, but before the event loop:
       StA::Class()->IgnoreTObjectStreamer(); // Don't store extra stuff
       TTree* tree = strangeNewDst.GetTree();
       TClonesArray* AArray = new TClonesArray("StA",1);
       tree->Branch("ABranch",&AArray,64000,99);
       StA* Ainstance = new((*AArray)[0]) StA();
    
  3. Fill the data in the new branch during the event loop. Here is a simple example where I count Lambda candidates and their mean pt:
         Int_t totalLams = 0;
         Float_t totalLamPt = 0.;
         for( Int_t j = 0; j < strangeOldDst.GetNV0(); j++ ) {
           StV0MuDst *v0m = strangeOldDst.GetV0(j);
           Float_t ml = v0m->massLambda();
           if ((ml > 1.11) && (ml < 1.12)) {
             totalLams++;
             totalLamPt += v0m->ptV0();
           }
         }
         Ainstance->value1 = totalLams;
         Ainstance->value2 = totalLamPt/totalLams;
    
That's all there is to it. Your output MuDst will contain the new branch with the new data. Simply opening the new file in Root will allow you to do things like making histograms with the new data using TTree::Draw(). Reading the actual class back will require adding code like this after the chain.Init() of readStrangeMuDst.C:
   TTree* tree = strangeNewDst.GetTree();
   TClonesArray* AArray = new TClonesArray("StA",0);
   tree->SetBranchStatus("ABranch.*",1);
   tree->SetBranchAddress("ABranch",&AArray);
and using the following code during the event loop:
     StA* Ainstance = (StA*) (AArray->At(0));