MainEvent.cxx

Go to the documentation of this file.
00001 // @(#)root/test:$Id: MainEvent.cxx 32212 2010-02-04 15:32:12Z rdm $
00002 // Author: Rene Brun   19/01/97
00003 
00004 ////////////////////////////////////////////////////////////////////////
00005 //
00006 //             A simple example with a ROOT tree
00007 //             =================================
00008 //
00009 //  This program creates :
00010 //    - a ROOT file
00011 //    - a tree
00012 //  Additional arguments can be passed to the program to control the flow
00013 //  of execution. (see comments describing the arguments in the code).
00014 //      Event  nevent comp split fill
00015 //  All arguments are optional. Default is:
00016 //      Event  400      1    1     1
00017 //
00018 //  In this example, the tree consists of one single "super branch"
00019 //  The statement ***tree->Branch("event", &event, 64000,split);*** below
00020 //  will parse the structure described in Event.h and will make
00021 //  a new branch for each data member of the class if split is set to 1.
00022 //    - 9 branches corresponding to the basic types fType, fNtrack,fNseg,
00023 //           fNvertex,fFlag,fTemperature,fMeasures,fMatrix,fClosesDistance.
00024 //    - 3 branches corresponding to the members of the subobject EventHeader.
00025 //    - one branch for each data member of the class Track of TClonesArray.
00026 //    - one branch for the TRefArray of high Pt tracks
00027 //    - one branch for the TRefArray of muon tracks
00028 //    - one branch for the reference pointer to the last track
00029 //    - one branch for the object fH (histogram of class TH1F).
00030 //
00031 //  if split = 0 only one single branch is created and the complete event
00032 //  is serialized in one single buffer.
00033 //  if split = -2 the event is split using the old TBranchObject mechanism
00034 //  if split = -1 the event is streamed using the old TBranchObject mechanism
00035 //  if split > 0  the event is split using the new TBranchElement mechanism.
00036 //
00037 //  if comp = 0 no compression at all.
00038 //  if comp = 1 event is compressed.
00039 //  if comp = 2 same as 1. In addition branches with floats in the TClonesArray
00040 //                         are also compressed.
00041 //  The 4th argument fill can be set to 0 if one wants to time
00042 //     the percentage of time spent in creating the event structure and
00043 //     not write the event in the file.
00044 //  In this example, one loops over nevent events.
00045 //  The branch "event" is created at the first event.
00046 //  The branch address is set for all other events.
00047 //  For each event, the event header is filled and ntrack tracks
00048 //  are generated and added to the TClonesArray list.
00049 //  For each event the event histogram is saved as well as the list
00050 //  of all tracks.
00051 //
00052 //  The two TRefArray contain only references to the original tracks owned by
00053 //  the TClonesArray fTracks.
00054 //
00055 //  The number of events can be given as the first argument to the program.
00056 //  By default 400 events are generated.
00057 //  The compression option can be activated/deactivated via the second argument.
00058 //
00059 //   ---Running/Linking instructions----
00060 //  This program consists of the following files and procedures.
00061 //    - Event.h event class description
00062 //    - Event.C event class implementation
00063 //    - MainEvent.C the main program to demo this class might be used (this file)
00064 //    - EventCint.C  the CINT dictionary for the event and Track classes
00065 //        this file is automatically generated by rootcint (see Makefile),
00066 //        when the class definition in Event.h is modified.
00067 //
00068 //   ---Analyzing the Event.root file with the interactive root
00069 //        example of a simple session
00070 //   Root > TFile f("Event.root")
00071 //   Root > T.Draw("fNtrack")   //histogram the number of tracks per event
00072 //   Root > T.Draw("fPx")       //histogram fPx for all tracks in all events
00073 //   Root > T.Draw("fXfirst:fYfirst","fNtrack>600")
00074 //                              //scatter-plot for x versus y of first point of each track
00075 //   Root > T.Draw("fH.GetRMS()")  //histogram of the RMS of the event histogram
00076 //
00077 //   Look also in the same directory at the following macros:
00078 //     - eventa.C  an example how to read the tree
00079 //     - eventb.C  how to read events conditionally
00080 //
00081 ////////////////////////////////////////////////////////////////////////
00082 
00083 #include <stdlib.h>
00084 
00085 #include "Riostream.h"
00086 #include "TROOT.h"
00087 #include "TFile.h"
00088 #include "TNetFile.h"
00089 #include "TRandom.h"
00090 #include "TTree.h"
00091 #include "TBranch.h"
00092 #include "TClonesArray.h"
00093 #include "TStopwatch.h"
00094 
00095 #include "Event.h"
00096 
00097 
00098 //______________________________________________________________________________
00099 int main(int argc, char **argv)
00100 {
00101    Int_t nevent = 400;     // by default create 400 events
00102    Int_t comp   = 1;       // by default file is compressed
00103    Int_t split  = 1;       // by default, split Event in sub branches
00104    Int_t write  = 1;       // by default the tree is filled
00105    Int_t hfill  = 0;       // by default histograms are not filled
00106    Int_t read   = 0;
00107    Int_t arg4   = 1;
00108    Int_t arg5   = 600;     //default number of tracks per event
00109    Int_t netf   = 0;
00110    Int_t punzip = 0;
00111 
00112    if (argc > 1)  nevent = atoi(argv[1]);
00113    if (argc > 2)  comp   = atoi(argv[2]);
00114    if (argc > 3)  split  = atoi(argv[3]);
00115    if (argc > 4)  arg4   = atoi(argv[4]);
00116    if (argc > 5)  arg5   = atoi(argv[5]);
00117    if (arg4 ==  0) { write = 0; hfill = 0; read = 1;}
00118    if (arg4 ==  1) { write = 1; hfill = 0;}
00119    if (arg4 ==  2) { write = 0; hfill = 0;}
00120    if (arg4 == 10) { write = 0; hfill = 1;}
00121    if (arg4 == 11) { write = 1; hfill = 1;}
00122    if (arg4 == 20) { write = 0; read  = 1;}  //read sequential
00123    if (arg4 == 21) { write = 0; read  = 1;  punzip = 1;}  //read sequential + parallel unzipping
00124    if (arg4 == 25) { write = 0; read  = 2;}  //read random
00125    if (arg4 >= 30) { netf  = 1; }            //use TNetFile
00126    if (arg4 == 30) { write = 0; read  = 1;}  //netfile + read sequential
00127    if (arg4 == 35) { write = 0; read  = 2;}  //netfile + read random
00128    if (arg4 == 36) { write = 1; }            //netfile + write sequential
00129    Int_t branchStyle = 1; //new style by default
00130    if (split < 0) {branchStyle = 0; split = -1-split;}
00131 
00132    TFile *hfile;
00133    TTree *tree;
00134    Event *event = 0;
00135 
00136    // Fill event, header and tracks with some random numbers
00137    //   Create a timer object to benchmark this loop
00138    TStopwatch timer;
00139    timer.Start();
00140    Long64_t nb = 0;
00141    Int_t ev;
00142    Int_t bufsize;
00143    Double_t told = 0;
00144    Double_t tnew = 0;
00145    Int_t printev = 100;
00146    if (arg5 < 100) printev = 1000;
00147    if (arg5 < 10)  printev = 10000;
00148 
00149 //         Read case
00150    if (read) {
00151       if (netf) {
00152          hfile = new TNetFile("root://localhost/root/test/EventNet.root");
00153       } else
00154          hfile = new TFile("Event.root");
00155       tree = (TTree*)hfile->Get("T");
00156       TBranch *branch = tree->GetBranch("event");
00157       branch->SetAddress(&event);
00158       Int_t nentries = (Int_t)tree->GetEntries();
00159       nevent = TMath::Min(nevent,nentries);
00160       if (read == 1) {  //read sequential
00161          //by setting the read cache to -1 we set it to the AutoFlush value when writing
00162          Int_t cachesize = -1; 
00163          if (punzip) tree->SetParallelUnzip();
00164          tree->SetCacheSize(cachesize);
00165          tree->SetCacheLearnEntries(1); //one entry is sufficient to learn
00166          tree->SetCacheEntryRange(0,nevent);
00167          for (ev = 0; ev < nevent; ev++) {
00168             tree->LoadTree(ev);  //this call is required when using the cache
00169             if (ev%printev == 0) {
00170                tnew = timer.RealTime();
00171                printf("event:%d, rtime=%f s\n",ev,tnew-told);
00172                told=tnew;
00173                timer.Continue();
00174             }
00175             nb += tree->GetEntry(ev);        //read complete event in memory
00176          }
00177       } else {    //read random
00178          Int_t evrandom;
00179          for (ev = 0; ev < nevent; ev++) {
00180             if (ev%printev == 0) cout<<"event="<<ev<<endl;
00181             evrandom = Int_t(nevent*gRandom->Rndm(1));
00182             nb += tree->GetEntry(evrandom);  //read complete event in memory
00183          }
00184       }
00185    } else {
00186 //         Write case
00187       // Create a new ROOT binary machine independent file.
00188       // Note that this file may contain any kind of ROOT objects, histograms,
00189       // pictures, graphics objects, detector geometries, tracks, events, etc..
00190       // This file is now becoming the current directory.
00191       if (netf) {
00192          hfile = new TNetFile("root://localhost/root/test/EventNet.root","RECREATE","TTree benchmark ROOT file");
00193       } else
00194          hfile = new TFile("Event.root","RECREATE","TTree benchmark ROOT file");
00195       hfile->SetCompressionLevel(comp);
00196 
00197      // Create histogram to show write_time in function of time
00198      Float_t curtime = -0.5;
00199      Int_t ntime = nevent/printev;
00200      TH1F *htime = new TH1F("htime","Real-Time to write versus time",ntime,0,ntime);
00201      HistogramManager *hm = 0;
00202      if (hfill) {
00203         TDirectory *hdir = new TDirectory("histograms", "all histograms");
00204         hm = new HistogramManager(hdir);
00205      }
00206 
00207      // Create a ROOT Tree and one superbranch
00208       tree = new TTree("T","An example of a ROOT tree");
00209       tree->SetAutoSave(1000000000); // autosave when 1 Gbyte written
00210       tree->SetCacheSize(10000000);  //set a 10 MBytes cache (useless when writing local files)
00211       bufsize = 64000;
00212       if (split)  bufsize /= 4;
00213       event = new Event();
00214       TTree::SetBranchStyle(branchStyle);
00215       TBranch *branch = tree->Branch("event", &event, bufsize,split);
00216       branch->SetAutoDelete(kFALSE);
00217       if(split >= 0 && branchStyle) tree->BranchRef();
00218       Float_t ptmin = 1;
00219 
00220       for (ev = 0; ev < nevent; ev++) {
00221          if (ev%printev == 0) {
00222             tnew = timer.RealTime();
00223             printf("event:%d, rtime=%f s\n",ev,tnew-told);
00224             htime->Fill(curtime,tnew-told);
00225             curtime += 1;
00226             told=tnew;
00227             timer.Continue();
00228          }
00229 
00230          event->Build(ev, arg5, ptmin);
00231 
00232          if (write) nb += tree->Fill();  //fill the tree
00233 
00234          if (hm) hm->Hfill(event);      //fill histograms
00235       }
00236       if (write) {
00237          hfile = tree->GetCurrentFile(); //just in case we switched to a new file
00238          hfile->Write();
00239          tree->Print();
00240       }
00241    }
00242 
00243    //  Stop timer and print results
00244    timer.Stop();
00245    Float_t mbytes = 0.000001*nb;
00246    Double_t rtime = timer.RealTime();
00247    Double_t ctime = timer.CpuTime();
00248 
00249 
00250    printf("\n%d events and %lld bytes processed.\n",nevent,nb);
00251    printf("RealTime=%f seconds, CpuTime=%f seconds\n",rtime,ctime);
00252    if (read) {
00253       tree->PrintCacheStats();
00254       printf("You read %f Mbytes/Realtime seconds\n",mbytes/rtime);
00255       printf("You read %f Mbytes/Cputime seconds\n",mbytes/ctime);
00256    } else {
00257       printf("compression level=%d, split=%d, arg4=%d\n",comp,split,arg4);
00258       printf("You write %f Mbytes/Realtime seconds\n",mbytes/rtime);
00259       printf("You write %f Mbytes/Cputime seconds\n",mbytes/ctime);
00260       //printf("file compression factor = %f\n",hfile.GetCompressionFactor());
00261    }
00262    hfile->Close();
00263    return 0;
00264 }

Generated on Tue Jul 5 15:15:03 2011 for ROOT_528-00b_version by  doxygen 1.5.1