]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/TTreeStream.cxx
Adding includes now needed by ROOT
[u/mrichter/AliRoot.git] / STEER / TTreeStream.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 //
19 //  marian.ivanov@cern.ch
20 //
21 //  ------------------------------------------------------------------------------------------------
22 //  TTreeStream
23 //  Standard stream (cout) like input for the tree
24 //  Run and see TTreeStreamer::Test() - to see TTreeStreamer functionality
25 //  ------------------------------------------------------------------------------------------------  
26 //
27 //  -------------------------------------------------------------------------------------------------
28 //  TTreeSRedirector
29 //  Redirect file to  different TTreeStreams  
30 //  Run and see   TTreeSRedirector::Test() as an example of TTreeSRedirectorer functionality 
31 // 
32
33 #include <TClass.h>
34 #include <TFile.h>
35 #include <TObjArray.h>
36 #include <TTree.h>
37 #include <TTreeStream.h>
38
39 ClassImp(TTreeDataElement)
40 ClassImp(TTreeStream)
41 ClassImp(TTreeSRedirector)
42
43
44
45 void TTreeStream::Test()
46 {
47   //
48   // 
49   TFile *ftest = new TFile("teststreamer.root","recreate");
50   if (!ftest) ftest = new TFile("teststreamer.root","new");
51   //
52   //create to streems Tree1 and Tree2
53   TTreeStream stream1("Tree1");
54   TTreeStream stream2("Tree2");
55   //
56   Char_t ch='s';
57   Float_t f=3.;
58   Float_t f2=1;
59   TObject *po  = new TObject;
60   TObject *po2 = new TObject;
61   for (Int_t i=0;i<100000;i++) {
62     f=i*100;
63     po->SetUniqueID(i);
64     po2->SetUniqueID(i*100);
65     ch=i%120;
66     //
67     //    Stream the data
68     //    The data layout of stream is defined during first invocation of streamer.
69     //    Endl is the trigger which define the end of structure.
70     // 
71     //    The name of branch can be specified using strings with = at the the end
72     //    if string is not specified automatic convention is u (sed B0, B1, ...Bn)
73     stream1<<"i="<<i<<"ch="<<ch<<"f="<<f<<"po="<<po<<"\n";
74     f  = 1./(100.1+i);
75     f2 = -f;     
76     //3.) just another example - we can fill the same tree with different objects
77     //
78     stream2<<f<<po<<"\n";
79     stream2<<f2<<po2<<"\n";
80   }
81   //
82   //4.) Close the streeamers (Write the streamed tree's to the file) and close the corresponding file.
83   //
84   stream1.Close();
85   stream2.Close();
86   ftest->Close();
87   delete ftest;
88   //
89   //5.) and now see results  in file tteststreamer.root
90 }
91
92
93
94 void TTreeSRedirector::Test()
95 {
96   //
97   //Example test function to show functionality of TTreeSRedirector
98   //
99   //
100   //1.)create the  redirector associated with file (testredirector.root)
101   //
102   //
103   TTreeSRedirector *pmistream= new TTreeSRedirector("testredirector.root");
104   TTreeSRedirector &mistream = *pmistream;
105   Char_t ch='s';
106   Float_t f=3.;
107   Float_t f2=1;
108   TObject *po  = new TObject;
109   TObject *po2 = new TObject;
110   for (Int_t i=0;i<100000;i++) {
111     f=i*100;
112     po->SetUniqueID(i);
113     po2->SetUniqueID(i*100);
114     ch=i%120;
115     //
116     //2.) create the tree with identifier specified by first argument
117     //                                layout specified by sequence of arguments
118     //                                Tree identifier has to be specified as first argument !!! 
119     //    if the tree and layout was already defined the consistency if layout is checked
120     //                                if the data are consisten fill given tree 
121     //    the name of branch can be specified using strings with = at the the end
122     //    if string is not specified use automatic convention  B0, B1, ...Bn
123     mistream<<"TreeIdentifier"<<"i="<<i<<"ch="<<ch<<"f="<<f<<"po="<<po<<"\n";
124     f  = 1./(100.1+i);
125     f2 = -f; 
126     
127     //3.) just another example - we can fill the same tree with different objects
128     //
129     mistream<<"TreeK"<<f<<po<<"\n";
130     mistream<<"TreeK"<<f2<<po2<<"\n";
131   }
132   //
133   //4.) write the streamed tree's to the file and close the corresponding file in destructor
134   //
135   delete pmistream;
136   //
137   //5.) and now see results in file testredirector.root 
138 }
139
140
141 TTreeSRedirector::TTreeSRedirector(const char *fname) :
142   fFile(new TFile(fname,"recreate")),
143   fDataLayouts(0)
144 {
145   //
146   // Constructor
147   //
148   if (!fFile){
149     fFile = new TFile(fname,"new");
150   }
151 }
152
153 TTreeSRedirector::~TTreeSRedirector()
154 {
155   //
156   // Destructor
157   //
158   Close();       //write the tree to the selected file
159   fFile->Close();
160   delete fFile;
161 }
162 void TTreeSRedirector::StoreObject(TObject* object){
163   //
164   //
165   //
166   TFile * backup = gFile;
167   fFile->cd();
168   object->Write();
169   if (backup) backup->cd();
170 }
171
172
173
174 TTreeStream  & TTreeSRedirector::operator<<(Int_t id)
175 {
176   //
177   // return reference to the data layout with given identifier
178   // if not existing - creates new
179   if (!fDataLayouts) fDataLayouts = new TObjArray(10000);
180   TTreeStream *clayout=0;
181   Int_t entries = fDataLayouts->GetEntriesFast();
182   for (Int_t i=0;i<entries;i++){
183     TTreeStream * layout = (TTreeStream*)fDataLayouts->At(i);
184     if (!layout) continue;
185     if (layout->fId==id) {
186       clayout = layout;
187       break;
188     }
189   }
190   if (!clayout){
191     fFile->cd();
192     char chname[100];
193     sprintf(chname,"Tree%d",id);
194     clayout = new TTreeStream(chname);
195     clayout->fId=id;
196     fDataLayouts->AddAt(clayout,entries);
197   }
198   return *clayout;
199 }
200
201
202 TTreeStream  & TTreeSRedirector::operator<<(const char* name)
203 {
204   //
205   // return reference to the data layout with given identifier
206   // if not existing - creates new
207   if (!fDataLayouts) fDataLayouts = new TObjArray(10000);
208   TTreeStream *clayout=(TTreeStream*)fDataLayouts->FindObject(name);
209   Int_t entries = fDataLayouts->GetEntriesFast();
210
211   if (!clayout){
212     fFile->cd();
213     clayout = new TTreeStream(name);
214     clayout->fId=-1;
215     clayout->SetName(name);
216     fDataLayouts->AddAt(clayout,entries);    
217   }
218   return *clayout;
219 }
220
221
222
223
224 void TTreeSRedirector::Close(){
225   //
226   //
227   TFile * backup = gFile;
228   fFile->cd();
229   if (fDataLayouts){
230     Int_t entries = fDataLayouts->GetEntriesFast();
231     for (Int_t i=0;i<entries;i++){
232       TTreeStream * layout = (TTreeStream*)fDataLayouts->At(i);
233       if (layout){
234         if (layout->fTree) layout->fTree->Write(layout->GetName());
235       }
236     }
237     delete fDataLayouts;
238     fDataLayouts=0;
239   }
240   if (backup) backup->cd();
241 }
242
243
244
245 //-------------------------------------------------------------
246 TTreeDataElement:: TTreeDataElement(Char_t type) :
247   TNamed(),
248   fType(type),
249   fDType(0),
250   fClass(0),
251   fPointer(0)
252 {
253   //
254   //
255   //
256 }
257
258 TTreeDataElement:: TTreeDataElement(TDataType* type) :
259   TNamed(),
260   fType(0),
261   fDType(type),
262   fClass(0),
263   fPointer(0)
264 {
265   //
266   //
267   //
268 }
269
270 TTreeDataElement:: TTreeDataElement(TClass* cl) :
271   TNamed(),
272   fType(0),
273   fDType(0),
274   fClass(cl),
275   fPointer(0)
276 {
277   //
278   //
279   //
280 }
281
282 //-------------------------------------------------------------------
283 TTreeStream::TTreeStream(const char *treename):
284   TNamed(treename,treename),
285   fElements(0),
286   fBranches(0),
287   fTree(new TTree(treename, treename)),
288   fCurrentIndex(0),
289   fId(0),
290   fNextName(),
291   fNextNameCounter(),
292   fStatus(0)
293 {
294   //
295   // Standard ctor
296   //
297 }
298
299 TTreeStream::~TTreeStream()
300 {
301   //
302   // Class dtor
303   //
304   fElements->Delete();
305   fBranches->Clear();
306   delete fElements;
307   delete fBranches;
308 }
309
310 void TTreeStream::Close()
311 {
312   //
313   // Flush data to disk and close
314   //
315   fTree->Write();
316 }
317
318 Int_t TTreeStream::CheckIn(Char_t type, void *pointer)
319 {
320   //
321   // Insert object of given type
322   //
323   if (!fElements) fElements = new TObjArray(1000);
324   TTreeDataElement* element = (TTreeDataElement*)fElements->At(fCurrentIndex);
325   if (!element) {
326     element = new TTreeDataElement(type);
327     //
328     char name[1000];
329     if (fNextName.Length()>0){
330       if (fNextNameCounter==0){
331         sprintf(name,"%s",(const char*)fNextName);
332       }
333       if (fNextNameCounter>0){
334         sprintf(name,"%s%d",(const char*)fNextName,fNextNameCounter);
335       }      
336     }
337     else{
338       sprintf(name,"B%d.",fCurrentIndex);
339     }
340     element->SetName(name);
341     //
342     element->SetPointer(pointer);
343     fElements->AddAt(element,fCurrentIndex);
344     fCurrentIndex++;
345     return 0; //new element added
346   }
347   if (element->GetType()!=type){
348     fStatus++;
349     return 1; //mismatched data element
350   }
351   element->SetPointer(pointer);
352   fCurrentIndex++;
353   return 0;
354 }
355
356 Int_t TTreeStream::CheckIn(TObject *o){
357   //
358   // Insert TObject
359   //
360   if (!o) return 0;
361   if (!fElements) fElements = new TObjArray(1000);
362   TTreeDataElement* element = (TTreeDataElement*)fElements->At(fCurrentIndex);
363   if (!element) {
364     element = new TTreeDataElement(o->IsA());
365     //
366     char name[1000];
367     if (fNextName.Length()>0){
368       if (fNextNameCounter==0){
369         sprintf(name,"%s",(const char*)fNextName);
370       }
371       if (fNextNameCounter>0){
372         sprintf(name,"%s%d",(const char*)fNextName,fNextNameCounter);
373       }      
374     }
375     else{
376       sprintf(name,"B%d",fCurrentIndex);
377     }
378     element->SetName(name);
379
380     element->SetPointer(o);
381     fElements->AddAt(element,fCurrentIndex);
382     fCurrentIndex++;
383     return 0; //new element added
384   }
385   if (element->fClass!=o->IsA()){
386     fStatus++;
387     return 1; //mismatched data element
388   }
389   element->SetPointer(o);
390   fCurrentIndex++;
391   return 0;  
392 }
393
394 void TTreeStream::BuildTree(){
395   //
396   // Build the Tree
397   //
398   if (fTree->GetEntries()>0) return;
399   fTree = new TTree(GetName(),GetName());
400   Int_t entries = fElements->GetEntriesFast();  
401   fBranches = new TObjArray(entries);
402   
403   for (Int_t i=0;i<entries;i++){
404     //
405     TTreeDataElement* element = (TTreeDataElement*)fElements->At(i);
406     char bname1[1000];
407     if (element->GetName()[0]==0){
408       sprintf(bname1,"B%d",i);
409     }
410     else{
411       sprintf(bname1,element->GetName());
412     }
413     if (element->fClass){
414       if (element->fClass->GetBaseClass("TClonesArray")){
415         TBranch * br = fTree->Branch(bname1,element->fClass->GetName(),&(element->fPointer));
416         fBranches->AddAt(br,i);
417       }else
418         {
419           TBranch * br = fTree->Branch(bname1,element->fClass->GetName(),&(element->fPointer));
420           fBranches->AddAt(br,i);
421         }
422     }
423     if (element->GetType()>0){
424       char bname2[1000];
425       sprintf(bname2,"B%d/%c",i,element->GetType());
426       TBranch * br = fTree->Branch(bname1,element->fPointer,bname2);
427       fBranches->AddAt(br,i);
428     }
429   }
430 }
431
432 void TTreeStream::Fill(){
433   //
434   // Fill the tree
435   //
436   if (fTree) { 
437     Int_t entries=fElements->GetEntriesFast();
438     if (entries>fTree->GetNbranches()) BuildTree();
439     for (Int_t i=0;i<entries;i++){    
440       TTreeDataElement* el  = (TTreeDataElement*)fElements->At(i);
441       if (!el) continue;
442       if (!el->GetType()) continue;
443       TBranch      * br  = (TBranch*)fBranches->At(i);
444       if (br &&el){
445         if (el->GetType())  br->SetAddress(el->fPointer);
446       }
447     }
448     if (fStatus==0) fTree->Fill(); //fill only in case of non conflicts
449     fStatus=0;
450   }
451 }
452
453 TTreeStream & TTreeStream::Endl()
454 {
455   //
456   // Perform pseudo endl operation
457   //
458   if (fTree->GetNbranches()==0) BuildTree();
459   Fill();
460   fStatus =0;
461   fCurrentIndex=0;
462   return *this;
463 }
464
465
466 TTreeStream  &TTreeStream::operator<<(Char_t *name)
467 {
468   //
469   // Endl 
470   //
471   if (name[0]=='\n'){
472     return Endl();
473   }
474   //
475   //if tree was already defined ignore
476   if (fTree->GetEntries()>0) return *this;
477   //check branch name if tree was not 
478   //
479   Int_t last=0;
480   for (last=0;;last++){
481     if (name[last]==0) break;    
482   }
483   
484   if (last>0&&name[last-1]=='='){
485     fNextName = name;
486     fNextName[last-1]=0;
487     fNextNameCounter=0;
488   }
489   return *this;
490 }
491