ROOT logo
// Title:   Sample Analysis Macro for looking at data on STAR NTuples
// Author:  Jim Thomas     jhthomas@lbl.gov
// Modified: Mikolaj Krzewicki  mikolaj.krzewicki@cern.ch
// Date:    04-Aug-2010
//

#include "TCanvas.h"
#include "TH1F.h"
#include "TMath.h"
#include "TNtuple.h"
#include "TLeaf.h"

void  readStarEvents()
{
  gSystem->Load("libTree.so");
  gSystem->Load("libVMC.so");
  gSystem->Load("libPhysics.so");
  gSystem->Load("libPWGflowBase");

  Long64_t EventCounter = 0 ;
  Int_t    TrackCounter = 0 ;
  Int_t    PrintInfo    = 1 ;                               // Print event and track data (1/0) (on/off)
  Int_t    PDG_ID           ;                               // for PID test using Particle Data Group PID names

  TH1F* myHistogram = new TH1F("Pt","Transverse Momentum", 100, 0.0, 3.0) ;
  TCanvas* myCanvas = new TCanvas("c1","c1",150,50,500,500) ;
  myCanvas -> cd() ;
  myHistogram -> Draw() ;                // Prepare the histogram and canvas for updates

  AliStarEventReader*  starReader = new AliStarEventReader( "/data/alice3/jthomas/testData/") ;
  AliStarEvent* starEvent = starReader->GetEvent();

  while ( starReader->GetNextEvent() )                                // Get next event
  {
    if ( !starReader->AcceptEvent(starEvent) ) continue;                           // Test if the event is good

    if ( PrintInfo == 1 ) starEvent->Print() ;  // Print basic information for this event

    if ((EventCounter%100) == 0)
    {
      TrackCounter = 0 ;

      //loop over track
      for ( Int_t j = 0 ; j < starEvent->GetNumberOfTracks() ; j++ )
      {
        AliStarTrack* starTrack = starEvent->GetTrack(j);             // Get next track
        if ( !starReader->AcceptTrack(starTrack) ) continue;               // Test if the track is good

        // Do something useful with the track
        if ( TrackCounter<5 && PrintInfo==1 ) starTrack->Print(); ; // Print the track data

        PDG_ID = starReader->ParticleID(starTrack);         // Assign Particle Data Book ID number to track (very simple PID algorithm)
        // Note this is an overly simplified ID algorithm and should really be extended for serious work
        if ( PDG_ID == 211 ) myHistogram->Fill( starTrack->GetPt() ) ;
        TrackCounter ++       ;               // Count the number of accepted and analyzed tracks
      }
      
      //starEvent->Print("all"); //if you want to print info for all tracks
      
      
      if ( PrintInfo == 1 )                         // Talk to the user and decide whether to continue, or stop.
      {
        cout << endl <<  "Enter 0 to quit, 1 to continue, 2 to continue without printing " << endl ;
        Int_t k = 0 ;
        cin >> k ;
        if ( k == 2 ) PrintInfo = 0 ;
        if (k <= 0) return ;
      }
    }
    EventCounter ++ ;
    if ( (EventCounter%10000) == 0 )   cout << EventCounter << endl ;
    if ( (EventCounter%10000) == 0 )   myCanvas->Update() ;
  }

  myCanvas->Update() ;                                      // Final update of histograms
  delete starReader ;
  starReader = NULL ;                              // Prepare to exit
}
 readStarEvents.C:1
 readStarEvents.C:2
 readStarEvents.C:3
 readStarEvents.C:4
 readStarEvents.C:5
 readStarEvents.C:6
 readStarEvents.C:7
 readStarEvents.C:8
 readStarEvents.C:9
 readStarEvents.C:10
 readStarEvents.C:11
 readStarEvents.C:12
 readStarEvents.C:13
 readStarEvents.C:14
 readStarEvents.C:15
 readStarEvents.C:16
 readStarEvents.C:17
 readStarEvents.C:18
 readStarEvents.C:19
 readStarEvents.C:20
 readStarEvents.C:21
 readStarEvents.C:22
 readStarEvents.C:23
 readStarEvents.C:24
 readStarEvents.C:25
 readStarEvents.C:26
 readStarEvents.C:27
 readStarEvents.C:28
 readStarEvents.C:29
 readStarEvents.C:30
 readStarEvents.C:31
 readStarEvents.C:32
 readStarEvents.C:33
 readStarEvents.C:34
 readStarEvents.C:35
 readStarEvents.C:36
 readStarEvents.C:37
 readStarEvents.C:38
 readStarEvents.C:39
 readStarEvents.C:40
 readStarEvents.C:41
 readStarEvents.C:42
 readStarEvents.C:43
 readStarEvents.C:44
 readStarEvents.C:45
 readStarEvents.C:46
 readStarEvents.C:47
 readStarEvents.C:48
 readStarEvents.C:49
 readStarEvents.C:50
 readStarEvents.C:51
 readStarEvents.C:52
 readStarEvents.C:53
 readStarEvents.C:54
 readStarEvents.C:55
 readStarEvents.C:56
 readStarEvents.C:57
 readStarEvents.C:58
 readStarEvents.C:59
 readStarEvents.C:60
 readStarEvents.C:61
 readStarEvents.C:62
 readStarEvents.C:63
 readStarEvents.C:64
 readStarEvents.C:65
 readStarEvents.C:66
 readStarEvents.C:67
 readStarEvents.C:68
 readStarEvents.C:69
 readStarEvents.C:70
 readStarEvents.C:71
 readStarEvents.C:72
 readStarEvents.C:73
 readStarEvents.C:74
 readStarEvents.C:75
 readStarEvents.C:76
 readStarEvents.C:77
 readStarEvents.C:78
 readStarEvents.C:79