StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Example_tags.C
1 void Example_tags(TString topDir = "/star/rcf/GC/daq/tags")
2 {
4 // //
5 // Example_tags.C //
6 // //
7 // shows how to use the STAR tags files //
8 // Input: top level directory //
9 // //
10 // what it does: //
11 // 1. creates TChain from all tags files down from the topDir //
12 // 2. loops over all events in the chain //
13 // //
14 // owner: Alexandre V. Vaniachine <AVVaniachine@lbl.gov> //
16 
17  gSystem->Load("libTable");
18  gSystem->Load("St_base");
19  // start benchmarks
20  gBenchmark = new TBenchmark();
21  gBenchmark->Start("total");
22 
23  // set loop optimization level
24  gROOT->ProcessLine(".O4");
25  // gather all files from the same top directory into one chain
26  // topDir must end with "/"
27  topDir +='/';
28  St_FileSet dirs(topDir);
29  St_DataSetIter next(&dirs,0);
30  St_DataSet *set = 0;
31  TChain chain("Tag");
32  while ( (set = next()) ) {
33  if (strcmp(set->GetTitle(),"file") ||
34  !(strstr(set->GetName(),".tags.root"))) continue;
35  chain.Add(gSystem->ConcatFileName(topDir,set->Path()));
36  }
37  UInt_t nEvents = chain->GetEntries();
38  cout<<"chained "<<nEvents<<" events "<<endl;
39 
40  TObjArray *files = chain.GetListOfFiles();
41  UInt_t nFiles = files->GetEntriesFast();
42  cout << "chained " << nFiles << " files from " << topDir << endl;
43 
44  TObjArray *leaves = chain.GetListOfLeaves();
45  Int_t nleaves = leaves->GetEntriesFast();
46 
47  TString tableName = " ";
48  TObjArray *tagTable = new TObjArray;
49 
50  Int_t tableCount = 0;
51  Int_t *tableIndex = new Int_t[nleaves];
52  Int_t tagCount = 0;
53 
54  // decode tag table names
55  for (Int_t l=0;l<nleaves;l++) {
56  TLeaf *leaf = (TLeaf*)leaves->UncheckedAt(l);
57  tagCount+=leaf->GetNdata();
58  TBranch *branch = leaf->GetBranch();
59  // new tag table name
60  if ( strstr(branch->GetName(), tableName.Data()) == 0 ) {
61  tableName = branch->GetName();
62  // the tableName is encoded in the branch Name before the "."
63  tableName.Resize(tableName->Last('.'));
64  tagTable->AddLast(new TObjString(tableName.Data()));
65  tableCount++;
66  }
67  tableIndex[l]=tableCount-1;
68  }
69  cout << " tot num tables, tags = " << tableCount << " "
70  << tagCount << endl << endl;
71 
72  //EXAMPLE 1: how to print out names of all tags and values for first event
73  for (l=0;l<nleaves;l++) {
74  leaf = (TLeaf*)leaves->UncheckedAt(l);
75  branch = leaf->GetBranch();
76  branch->GetEntry();
77  // tag comment is in the title
78  TString Title = leaf->GetTitle();
79  Int_t dim = leaf->GetNdata();
80  if (dim==1) {
81  Title.ReplaceAll('['," '");
82  Title.ReplaceAll(']',"'");
83  }
84  cout << "\n Table: ";
85  cout.width(10);
86  cout << ((TObjString*)tagTable->UncheckedAt(tableIndex[l]))->GetString()
87  <<" -- has tag: " << Title << endl;
88  for (Int_t i=0;i<dim;i++) {
89  cout <<" "<< leaf->GetName();
90  if (dim>1) cout << '['<<i<<']';
91  cout << " = " << leaf->GetValue(i) << endl;
92  }
93  }
94 
95  // EXAMPLE 2: how to make a plot
96  c1 = new TCanvas("c1","Beam-Gas Rejection",600,1000);
97  gStyle->SetMarkerStyle(8);
98  chain->Draw("n_trk_tpc[0]:n_trk_tpc[1]");
99 
100  // EXAMPLE 3: how to make a selection (write selected event numbers on the plot)
101  Int_t ncoll=0;
102  char aevent[10];
103  TText t(0,0,"a");
104  t.SetTextFont(52);
105  t.SetTextSize(0.02);
106  Float_t cut = 0.35;
107  cout <<"\n Events with ntrk>400 and |asim|<"<<cut<<endl;
108  //loop over all events: READ ONLY n_trk_tpc AND mEventNumber BRANCHES!
109  gBenchmark->Start("loop");
110  for (Int_t i=0;i<nFiles;i++) {
111  chain.LoadTree(*(chain.GetTreeOffset()+i));
112  TTree *T = chain.GetTree();
113  //must renew leaf pointer for each tree
114  TLeaf *ntrk = T->GetLeaf("n_trk_tpc");
115  TLeaf *run = T->GetLeaf("mRunNumber");
116  TLeaf *event = T->GetLeaf("mEventNumber");
117  for (Int_t j=0; j<T->GetEntries(); j++){
118  ntrk->GetBranch()->GetEntry(j);
119  event->GetBranch()->GetEntry(j);
120  run->GetBranch()->GetEntry(j);
121  Int_t Nm=ntrk->GetValue(0);
122  Int_t Np=ntrk->GetValue(1);
123  Int_t Ntrk = Np+Nm;
124  // avoid division by 0
125  Float_t asim = Np-Nm;
126  if (Ntrk>0) asim /= Ntrk;
127  if (-cut < asim&&asim < cut && Ntrk>400) {
128  cout<<" Run "<<(UInt_t)run->GetValue()
129  <<", Event "<<event->GetValue() <<endl;
130  ncoll++;
131  sprintf(aevent,"%d",event->GetValue());
132  t.DrawText(Np+10,Nm+10,aevent);
133  }
134  }
135  }
136  gBenchmark->Stop("loop");
137  t.SetTextSize(0.05);
138  t.DrawText(50,2550,"ntrk>400 and |(Np-Nm)/(Np+Nm)| < 0.35 ");
139  t.DrawText(500,-300,"Ntrk with tanl<0 ");
140  cout << " Selected " << ncoll << " collision candidates out of "
141  << nEvents << " events" << endl;
142  // stop timer and print benchmarks
143  gBenchmark->Print("loop");
144  gBenchmark->Stop("total");
145  gBenchmark->Print("total");
146 }
virtual const char * GetName() const
special overload
Definition: StMaker.cxx:237
virtual TString Path() const
return the full path of this data set
Definition: TDataSet.cxx:626