StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
exampleStMuIOMaker.C
1 
2 #include "TH1.h"
3 #include "TChain.h"
4 #include "TSystem.h"
5 #include "TFile.h"
6 #include <iostream>
7 
8 class StMuIOMaker;
9 StMuIOMaker* maker;
10 
11 TH1D globalPt("globalPt","globalPt",100,0.,3.);
12 TH1D primaryPt("primaryPt","primaryPt",100,0.,3.);
13 TH1D l3Pt("l3Pt","l3Pt",100,0.,3.);
14 TH1D refMult("refMult","refMult",100,0.,100.);
15 
16 TH2D etaLength("etaLength","etaLength",20,-1,+1,200,0.,200.);
17 
18 
19 void exampleStMuIOMaker(const char* file="/star/data14/reco/productionCentral/ReversedFullField/P02ge/2001/324/st_physics_2324001_raw_0006.MuDst.root") {
20  gROOT->LoadMacro("$STAR/StRoot/StMuDSTMaker/COMMON/macros/loadSharedLibraries.C");
21  loadSharedLibraries();
22 
24 
25  int counter=0;
26  int iret=0;
27 
28  StMuIOMaker* maker = new StMuIOMaker("MuDst");
29  maker->SetFileAndOpen(file);
30 // maker->SetFileName("st_physics_2313018_raw_0029.MuDst.root");
31  TMemStat memStat("exampleMuIO");
32 
33  StMuDst* mu;
34 // evtid is an unique id of an event in a run which may involve several files!
35 // it is different from sequential index in a run which starts from "0"!
36 // for ( int evtid=80673; evtid<80694; evtid++) {
37 // for ( int evtid=6300; evtid<6321; evtid++) {
38 // iret = maker->Make( StUKey(2271008,evtid) );
39 // iret = maker->Make( StUKey(2313018,evtid) );
40 // alterative
41 for ( int i=0; i<20; i++) {
42 // iret = maker->Make(); // read an event in natural sequential
43  iret = maker->Make(i); // read an event with seqential index
44  if ( iret == kStOK) {
45  StMuDst* mu = maker->muDst();
46  if ( !mu ) continue;
47 // if(i%10 != 0) continue;
48 
49 // take a look at branches of tracks
50  int n;
51  n= mu->globalTracks()->GetEntries();
52  for (int l=0; l<n; l++) globalPt.Fill(mu->globalTracks(l)->pt());
53  n= mu->primaryTracks()->GetEntries();
54  for (int l=0; l<n; l++) primaryPt.Fill(mu->primaryTracks(l)->pt());
55  n= mu->l3Tracks()->GetEntries();
56  for (int l=0; l<n; l++) {l3Pt.Fill( mu->l3Tracks(l)->pt() );}
57 
58 // take a look at event branch
59  StMuEvent* muEvent = mu->event();
60  int referenceMultiplicity = muEvent->refMult();
61  refMult.Fill(referenceMultiplicity);
62 
63  cout << "#" << i << " index=" << maker->currentIndex()
64 // cout << "eventId =" << evtid << " index=" << maker->currentIndex()
65  << " refmult= "<< referenceMultiplicity
66  << " used= "<< memStat.Used()
67  << " size= "<< memStat.ProgSize() << endl;
68  counter++;
69  }
70  }
71  cout << " # of events:" << counter << endl;
72 
73  globalPt.Draw();
74  primaryPt.Draw("same");
75  l3Pt.Draw("same");
76 
77  TFile f("testMuIO.root","RECREATE");
78  globalPt.Write();
79  primaryPt.Write();
80  l3Pt.Write();
81  refMult.Write();
82  f.Close();
83 }
static TObjArray * globalTracks()
returns pointer to the global tracks list
Definition: StMuDst.h:303
StMuDst * muDst()
Definition: StMuDstMaker.h:425
static void setLevel(unsigned int level)
sets the debug level
Definition: StMuDebug.h:74
static TClonesArray * l3Tracks()
returns pointer to the l3Tracks list
Definition: StMuDst.h:307
StMuIOMaker(const char *name="", const char *ioFile="")
Default constructor.
Definition: StMuIOMaker.cxx:28
static TObjArray * primaryTracks()
returns pointer to a list of tracks belonging to the selected primary vertex
Definition: StMuDst.h:301
unsigned short refMult(int vtx_id=-1)
Reference multiplicity of charged particles as defined in StEventUtilities/StuRefMult.hh for vertex vtx_id (-1 is default index from StMuDst)
Definition: StMuEvent.cxx:195
static StMuEvent * event()
returns pointer to current StMuEvent (class holding the event wise information, e.g. event number, run number)
Definition: StMuDst.h:320
Definition: Stypes.h:40