KaliVeda  1.13/01
Heavy-Ion Analysis Toolkit
KVEventListMaker.cpp
Go to the documentation of this file.
1 /*
2 $Id: KVEventListMaker.cpp,v 1.3 2009/01/30 09:27:06 ebonnet Exp $
3 $Revision: 1.3 $
4 $Date: 2009/01/30 09:27:06 $
5 */
6 
7 //Created by KVClassFactory on Thu Mar 20 11:58:48 2008
8 //Author: Eric Bonnet
9 
10 #include "KVEventListMaker.h"
11 #include "TFile.h"
12 #include "TTree.h"
13 #include "TEventList.h"
14 #include "KVList.h"
15 #include "KVString.h"
16 
18 
19 
20 
22 void KVEventListMaker::Process()
23 {
24 
25  if (!IsReady()) return;
26 
27 //open file where tree is stored
28  TFile file(GetFileName().Data(), "update");
29  KVList kevtlist;
30 
31  TTree* tt = NULL;
32 
33 //just count the number of branches
34  KVString lname = GetBranchName();
35  Int_t nbre = 0;
36  lname.Begin(" ");
37  while (!lname.End()) {
38  KVString stamp = lname.Next();
39  nbre += 1;
40  }
41  if (nbre == 0) return;
42 
43  std::vector<Int_t> variable(nbre);
44  KVString evtname;
45  Bool_t ok = kFALSE;
46 //check if the tree are there
47  if ((tt = (TTree*)file.Get(GetTreeName().Data()))) {
48  Int_t nentries = tt->GetEntries();
49  printf("nbre d entree %d\n", nentries);
50  TBranch* bb = NULL;
51 
52  nbre = 0;
53  lname.Begin(" ");
54  //loop on branche names validity, if there is one wrong name
55  //the program will exit
56  while (!lname.End()) {
57  KVString stamp = lname.Next();
58  if ((bb = tt->GetBranch(stamp.Data()))) {
59  tt->SetBranchAddress(stamp.Data(), &variable[nbre]);
60  nbre += 1;
61  ok = kTRUE;
62  }
63  }
64 
65  if (!ok) return;
66 
67  TEventList* el = NULL;
68  //loop on tree entries
69  for (Int_t kk = 0; kk < nentries; kk += 1) {
70  if (kk % 10000 == 0) printf("%d evts lus\n", kk);
71 
72  tt->GetEntry(kk);
73  evtname = "";
74  nbre = 0;
75  lname.Begin(" ");
76  //compute the name of the TEventList
77  //using branch name and branch value
78  while (!lname.End()) {
79  KVString stamp = lname.Next();
80  if ((bb = tt->GetBranch(stamp.Data()))) {
81  KVString val;
82  val.Form("%s_%d_", stamp.Data(), variable[nbre]);
83  evtname += val;
84  nbre += 1;
85  }
86  }
87 
88  if (!(el = kevtlist.get_object<TEventList>(evtname.Data()))) {
89  printf("creation de %s \n", evtname.Data());
90  kevtlist.Add(new TEventList(evtname.Data()));
91  el = (TEventList*)kevtlist.Last();
92  }
93  el->Enter(kk);
94  }
95 
96  //write TEventList created
97  if (ktag_tree) {
98  for (Int_t nn = 0; nn < kevtlist.GetEntries(); nn += 1) {
99  KVString tampname;
100  tampname.Form("%s_%s", GetTreeName().Data(), kevtlist.At(nn)->GetName());
101  ((TEventList*) kevtlist.At(nn))->SetName(tampname.Data());
102  kevtlist.At(nn)->Write();
103  }
104  }
105  else {
106  kevtlist.Write();
107  }
108  }
109  else printf("%s n existe pas\n", GetTreeName().Data());
110 }
111 
112 
int Int_t
ClassImp(KVPartitionList) void KVPartitionList
Initialisation.
const Bool_t kFALSE
bool Bool_t
const Bool_t kTRUE
int nentries
Compute TEventList for TTree.
Extended TList class which owns its objects by default.
Definition: KVList.h:27
T * get_object(const TString &name) const
virtual TObject * Last() const
virtual TObject * At(Int_t idx) const
virtual void Add(TObject *obj)
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition: KVString.h:72
void Begin(TString delim) const
Definition: KVString.cpp:565
Bool_t End() const
Definition: KVString.cpp:634
KVString Next(Bool_t strip_whitespace=kFALSE) const
Definition: KVString.cpp:695
virtual Int_t GetEntries() const
virtual Int_t Write(const char *name=0, Int_t option=0, Int_t bufsize=0)
virtual void Enter(Long64_t entry)
virtual Int_t Write(const char *name=0, Int_t option=0, Int_t bufsize=0)
virtual const char * GetName() const
const char * Data() const
void Form(const char *fmt,...)
gr SetName("gr")
llvm::StringRef GetFileName(const clang::Decl &decl, const cling::Interpreter &interp)
auto * tt