1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
19 // marian.ivanov@cern.ch
21 // ------------------------------------------------------------------------------------------------
23 // Standard stream (cout) like input for the tree
24 // Run and see TTreeStreamer::Test() - to see TTreeStreamer functionality
25 // ------------------------------------------------------------------------------------------------
27 // -------------------------------------------------------------------------------------------------
29 // Redirect file to different TTreeStreams
30 // Run and see TTreeSRedirector::Test() as an example of TTreeSRedirectorer functionality
35 #include <TDirectory.h>
36 #include <TObjArray.h>
38 #include <TTreeStream.h>
40 ClassImp(TTreeDataElement)
42 ClassImp(TTreeSRedirector)
46 void TTreeStream::Test()
50 TFile *ftest = new TFile("teststreamer.root","recreate");
51 if (!ftest) ftest = new TFile("teststreamer.root","new");
53 //create to streems Tree1 and Tree2
54 TTreeStream stream1("Tree1");
55 TTreeStream stream2("Tree2");
60 TObject *po = new TObject;
61 TObject *po2 = new TObject;
62 for (Int_t i=0;i<100000;i++) {
65 po2->SetUniqueID(i*100);
69 // The data layout of stream is defined during first invocation of streamer.
70 // Endl is the trigger which define the end of structure.
72 // The name of branch can be specified using strings with = at the the end
73 // if string is not specified automatic convention is u (sed B0, B1, ...Bn)
74 stream1<<"i="<<i<<"ch="<<ch<<"f="<<f<<"po="<<po<<"\n";
77 //3.) just another example - we can fill the same tree with different objects
80 stream2<<f2<<po2<<"\n";
83 //4.) Close the streeamers (Write the streamed tree's to the file) and close the corresponding file.
90 //5.) and now see results in file tteststreamer.root
95 void TTreeSRedirector::Test()
98 //Example test function to show functionality of TTreeSRedirector
101 //1.)create the redirector associated with file (testredirector.root)
104 TTreeSRedirector *pmistream= new TTreeSRedirector("testredirector.root");
105 TTreeSRedirector &mistream = *pmistream;
109 TObject *po = new TObject;
110 TObject *po2 = new TObject;
111 for (Int_t i=0;i<100000;i++) {
114 po2->SetUniqueID(i*100);
117 //2.) create the tree with identifier specified by first argument
118 // layout specified by sequence of arguments
119 // Tree identifier has to be specified as first argument !!!
120 // if the tree and layout was already defined the consistency if layout is checked
121 // if the data are consisten fill given tree
122 // the name of branch can be specified using strings with = at the the end
123 // if string is not specified use automatic convention B0, B1, ...Bn
124 mistream<<"TreeIdentifier"<<"i="<<i<<"ch="<<ch<<"f="<<f<<"po="<<po<<"\n";
128 //3.) just another example - we can fill the same tree with different objects
130 mistream<<"TreeK"<<f<<po<<"\n";
131 mistream<<"TreeK"<<f2<<po2<<"\n";
134 //4.) write the streamed tree's to the file and close the corresponding file in destructor
138 //5.) and now see results in file testredirector.root
142 TTreeSRedirector::TTreeSRedirector(const char *fname) :
143 fFile(new TFile(fname,"recreate")),
150 fFile = new TFile(fname,"new");
154 TTreeSRedirector::~TTreeSRedirector()
159 Close(); //write the tree to the selected file
163 void TTreeSRedirector::StoreObject(TObject* object){
167 TDirectory * backup = gDirectory;
170 if (backup) backup->cd();
175 TTreeStream & TTreeSRedirector::operator<<(Int_t id)
178 // return reference to the data layout with given identifier
179 // if not existing - creates new
180 if (!fDataLayouts) fDataLayouts = new TObjArray(10000);
181 TTreeStream *clayout=0;
182 Int_t entries = fDataLayouts->GetEntriesFast();
183 for (Int_t i=0;i<entries;i++){
184 TTreeStream * layout = (TTreeStream*)fDataLayouts->At(i);
185 if (!layout) continue;
186 if (layout->fId==id) {
192 TDirectory * backup = gDirectory;
195 sprintf(chname,"Tree%d",id);
196 clayout = new TTreeStream(chname);
198 fDataLayouts->AddAt(clayout,entries);
199 if (backup) backup->cd();
205 TTreeStream & TTreeSRedirector::operator<<(const char* name)
208 // return reference to the data layout with given identifier
209 // if not existing - creates new
210 if (!fDataLayouts) fDataLayouts = new TObjArray(10000);
211 TTreeStream *clayout=(TTreeStream*)fDataLayouts->FindObject(name);
212 Int_t entries = fDataLayouts->GetEntriesFast();
215 TDirectory * backup = gDirectory;
217 clayout = new TTreeStream(name);
219 clayout->SetName(name);
220 fDataLayouts->AddAt(clayout,entries);
221 if (backup) backup->cd();
229 void TTreeSRedirector::Close(){
232 TDirectory * backup = gDirectory;
235 Int_t entries = fDataLayouts->GetEntriesFast();
236 for (Int_t i=0;i<entries;i++){
237 TTreeStream * layout = (TTreeStream*)fDataLayouts->At(i);
239 if (layout->fTree) layout->fTree->Write(layout->GetName());
245 if (backup) backup->cd();
250 //-------------------------------------------------------------
251 TTreeDataElement:: TTreeDataElement(Char_t type) :
263 TTreeDataElement:: TTreeDataElement(TDataType* type) :
275 TTreeDataElement:: TTreeDataElement(TClass* cl) :
287 //-------------------------------------------------------------------
288 TTreeStream::TTreeStream(const char *treename):
289 TNamed(treename,treename),
292 fTree(new TTree(treename, treename)),
304 TTreeStream::~TTreeStream()
315 void TTreeStream::Close()
318 // Flush data to disk and close
323 Int_t TTreeStream::CheckIn(Char_t type, void *pointer)
326 // Insert object of given type
328 if (!fElements) fElements = new TObjArray(10000);
329 if (fElements->GetSize()<=fCurrentIndex) fElements->Expand(fCurrentIndex*2);
330 TTreeDataElement* element = (TTreeDataElement*)fElements->At(fCurrentIndex);
332 element = new TTreeDataElement(type);
335 if (fNextName.Length()>0){
336 if (fNextNameCounter==0){
337 sprintf(name,"%s",(const char*)fNextName);
339 if (fNextNameCounter>0){
340 sprintf(name,"%s%d",(const char*)fNextName,fNextNameCounter);
344 sprintf(name,"B%d.",fCurrentIndex);
346 element->SetName(name);
348 element->SetPointer(pointer);
349 fElements->AddAt(element,fCurrentIndex);
351 return 0; //new element added
353 if (element->GetType()!=type){
355 return 1; //mismatched data element
357 element->SetPointer(pointer);
362 Int_t TTreeStream::CheckIn(TObject *o){
367 if (!fElements) fElements = new TObjArray(1000);
368 TTreeDataElement* element = (TTreeDataElement*)fElements->At(fCurrentIndex);
370 element = new TTreeDataElement(o->IsA());
373 if (fNextName.Length()>0){
374 if (fNextNameCounter==0){
375 sprintf(name,"%s",(const char*)fNextName);
377 if (fNextNameCounter>0){
378 sprintf(name,"%s%d",(const char*)fNextName,fNextNameCounter);
382 sprintf(name,"B%d",fCurrentIndex);
384 element->SetName(name);
386 element->SetPointer(o);
387 fElements->AddAt(element,fCurrentIndex);
389 return 0; //new element added
391 if (element->fClass!=o->IsA()){
393 return 1; //mismatched data element
395 element->SetPointer(o);
400 void TTreeStream::BuildTree(){
404 if (fTree && fTree->GetEntries()>0) return;
405 if (!fTree) fTree = new TTree(GetName(),GetName());
406 Int_t entries = fElements->GetEntriesFast();
407 fBranches = new TObjArray(entries);
409 for (Int_t i=0;i<entries;i++){
411 TTreeDataElement* element = (TTreeDataElement*)fElements->At(i);
413 if (element->GetName()[0]==0){
414 sprintf(bname1,"B%d",i);
417 sprintf(bname1,element->GetName());
419 if (element->fClass){
420 if (element->fClass->GetBaseClass("TClonesArray")){
421 TBranch * br = fTree->Branch(bname1,element->fClass->GetName(),&(element->fPointer));
422 fBranches->AddAt(br,i);
425 TBranch * br = fTree->Branch(bname1,element->fClass->GetName(),&(element->fPointer));
426 fBranches->AddAt(br,i);
429 if (element->GetType()>0){
431 sprintf(bname2,"B%d/%c",i,element->GetType());
432 TBranch * br = fTree->Branch(bname1,element->fPointer,bname2);
433 fBranches->AddAt(br,i);
438 void TTreeStream::Fill(){
443 Int_t entries=fElements->GetEntriesFast();
444 if (entries>fTree->GetNbranches()) BuildTree();
445 for (Int_t i=0;i<entries;i++){
446 TTreeDataElement* el = (TTreeDataElement*)fElements->At(i);
448 if (!el->GetType()) continue;
449 TBranch * br = (TBranch*)fBranches->At(i);
451 if (el->GetType()) br->SetAddress(el->fPointer);
454 if (fStatus==0) fTree->Fill(); //fill only in case of non conflicts
459 TTreeStream & TTreeStream::Endl()
462 // Perform pseudo endl operation
464 if (fTree->GetNbranches()==0) BuildTree();
472 TTreeStream &TTreeStream::operator<<(const Char_t *name)
481 //if tree was already defined ignore
482 if (fTree->GetEntries()>0) return *this;
483 //check branch name if tree was not
486 for (last=0;;last++){
487 if (name[last]==0) break;
490 if (last>0&&name[last-1]=='='){