StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
rdGeant.C
1 class StChain;
2 StChain *chain=0;
3 
4 int rdGeant( int maxEve=2){
5 
6  TString geantFile = "/star/data26/reco/pp200/pythia_6.203/default/minbias/year2003/hadronic_on/trs_ic/rcf1200_2723_2000evts.geant.root";
7 
8 
9  gROOT->LoadMacro("$STAR/StRoot/StMuDSTMaker/COMMON/macros/loadSharedLibraries.C");
10  loadSharedLibraries();
11  cout << " loading done " << endl;
12 
13 
14  // create chain
15  chain = new StChain("StChain");
16  //chain->SetDebug();
17  //chain->PrintInfo();
18 
19  // set up maker in read mode
20  StIOMaker* IOMk = new StIOMaker("IO","r",geantFile,"bfcTree");
21  IOMk->SetDebug();
22  IOMk->SetIOMode("r");
23  IOMk->SetBranch("*",0,"0");
24  IOMk->SetBranch("geantBranch",0,"r");
25 
26  chain->Init();
27  chain->ls(3);
28 
29  // StMuDebug::setLevel(1); // switch of some debug output
30 
31  int eventCounter=0;
32 
33  printf(" requested maxEve=%d\n",maxEve);
34  //---------------------------------------------------
35  while ( 1) {// loop over events
36  eventCounter++;;
37  if(eventCounter >maxEve) break;
38  chain->Clear();
39  int stat = chain->Make();
40  if(stat) break;
41 
42  // Access to geant-tables .......................
43 
44  St_DataSet* Event = chain->GetDataSet("geant");
45  //Event->ls(3);
46  St_DataSetIter geantDstI(Event);
47 
48  // Event generator info ........................
49  St_g2t_event *Pg2t_event=(St_g2t_event *) geantDstI("g2t_event");
50  //Pg2t_event->Print();
51  g2t_event_st *g2t_event1=Pg2t_event->GetTable();
52  printf("nr=%d %p\n",Pg2t_event->GetNRows(),g2t_event1);
53  int k1= g2t_event1->eg_label;
54  int k2= g2t_event1->n_event;
55  int k3= g2t_event1->subprocess_id;
56 
57  printf("eg_label=%d n_event=%d subprocess_id=%d\n", k1,k2,k3);
58 
59 
60  // This is an example to access the particle table directly.
61  St_particle *particleTabPtr = (St_particle *) geantDstI("particle");
62  particle_st* particleTable = particleTabPtr->GetTable();
63 
64  for (int i=0; i<particleTabPtr->GetNRows();++i) {
65  if(i>10) break;
66  cout << "track " << i << endl;
67  cout << " id = " << particleTable[i].idhep << endl;
68  cout << " px = " << particleTable[i].phep[0] << endl;
69  cout << " py = " << particleTable[i].phep[1] << endl;
70  cout << " pz = " << particleTable[i].phep[2] << endl;
71  cout << " e = " << particleTable[i].phep[3] << endl;
72  cout << " m = " << particleTable[i].phep[4] << endl;
73  cout << " moth1 = " << particleTable[i].jmohep[0] << endl;
74  cout << " moth2 = " << particleTable[i].jmohep[1] << endl;
75  }
76 
77  }
78 
79 
80 
81 }
82 
83 
84 
85 
86 
87 // /star/data31/reco/ppMinBias/FullField/P03if/2002/...
88 // /star/data31/reco/ppMinBias/ReversedFullField/P03if/2002/...
89 
virtual void SetIOMode(Option_t *iomode="w")
number of transactions
Definition: StIOInterFace.h:35
virtual void Clear(Option_t *option="")
User defined functions.
Definition: StChain.cxx:77
virtual void ls(Option_t *option="") const
Definition: TDataSet.cxx:495
virtual Int_t Make()
Definition: StChain.cxx:110
Definition: AgUStep.h:26