hvector.C

Go to the documentation of this file.
00001 //
00002 // This tutorials demonstrate how to store and restore simple vectors
00003 // in a TTree
00004 //
00005 
00006 #include <vector>
00007 
00008 #include "TFile.h"
00009 #include "TTree.h"
00010 #include "TCanvas.h"
00011 #include "TFrame.h"
00012 #include "TH1F.h"
00013 #include "TBenchmark.h"
00014 #include "TRandom.h"
00015 #include "TSystem.h"
00016 
00017 #ifdef __MAKECINT__
00018 #pragma link C++ class vector<float>+;
00019 #endif
00020 
00021 void write() 
00022 {
00023   
00024    TFile *f = TFile::Open("hvector.root","RECREATE");
00025    
00026    if (!f) { return; }
00027 
00028    // Create one histograms
00029    TH1F *hpx = new TH1F("hpx","This is the px distribution",100,-4,4);
00030    hpx->SetFillColor(48);
00031 
00032    std::vector<float> vpx;
00033    std::vector<float> vpy;
00034    std::vector<float> vpz;
00035    std::vector<float> vrand;
00036 
00037    // Create a TTree
00038    TTree *t = new TTree("tvec","Tree with vectors");
00039    t->Branch("vpx",&vpx);
00040    t->Branch("vpy",&vpy);
00041    t->Branch("vpz",&vpz);
00042    t->Branch("vrand",&vrand);
00043 
00044 
00045   // Create a new canvas.
00046    TCanvas *c1 = new TCanvas("c1","Dynamic Filling Example",200,10,700,500);
00047    c1->SetFillColor(42);
00048    c1->GetFrame()->SetFillColor(21);
00049    c1->GetFrame()->SetBorderSize(6);
00050    c1->GetFrame()->SetBorderMode(-1);
00051 
00052    gRandom->SetSeed();
00053    const Int_t kUPDATE = 1000;
00054    for (Int_t i = 0; i < 25000; i++) {
00055       Int_t npx = (Int_t)(gRandom->Rndm(1)*15);
00056 
00057       vpx.clear();
00058       vpy.clear();
00059       vpz.clear();
00060       vrand.clear();
00061 
00062       for (Int_t j = 0; j < npx; ++j) {
00063 
00064          Float_t px,py,pz;
00065          gRandom->Rannor(px,py);
00066          pz = px*px + py*py;
00067          Float_t random = gRandom->Rndm(1);
00068  
00069          hpx->Fill(px);
00070           
00071          vpx.push_back(px);
00072          vpy.push_back(py);
00073          vpz.push_back(pz);
00074          vrand.push_back(random);
00075 
00076       }
00077       if (i && (i%kUPDATE) == 0) {
00078          if (i == kUPDATE) hpx->Draw();
00079          c1->Modified();
00080          c1->Update();
00081          if (gSystem->ProcessEvents())
00082             break;
00083       }
00084       t->Fill();
00085    }
00086    f->Write();
00087    
00088    delete f;
00089 }
00090 
00091 
00092 void read() 
00093 {
00094   
00095    TFile *f = TFile::Open("hvector.root","READ");
00096    
00097    if (!f) { return; }
00098 
00099    TTree *t; f->GetObject("tvec",t);
00100 
00101    std::vector<float> *vpx = 0;
00102 
00103   // Create a new canvas.
00104    TCanvas *c1 = new TCanvas("c1","Dynamic Filling Example",200,10,700,500);
00105    c1->SetFillColor(42);
00106    c1->GetFrame()->SetFillColor(21);
00107    c1->GetFrame()->SetBorderSize(6);
00108    c1->GetFrame()->SetBorderMode(-1);
00109 
00110    const Int_t kUPDATE = 1000;
00111    
00112    TBranch *bvpx = 0;
00113    t->SetBranchAddress("vpx",&vpx,&bvpx);
00114    
00115 
00116    // Create one histograms
00117    TH1F *h = new TH1F("h","This is the px distribution",100,-4,4);
00118    h->SetFillColor(48);
00119 
00120    for (Int_t i = 0; i < 25000; i++) {
00121       
00122       Long64_t tentry = t->LoadTree(i);
00123       bvpx->GetEntry(tentry);
00124       
00125       for (UInt_t j = 0; j < vpx->size(); ++j) {
00126  
00127          h->Fill(vpx->at(j));          
00128 
00129       }
00130       if (i && (i%kUPDATE) == 0) {
00131          if (i == kUPDATE) h->Draw();
00132          c1->Modified();
00133          c1->Update();
00134          if (gSystem->ProcessEvents())
00135             break;
00136       }
00137    }
00138 
00139    // Since we passed the address of a local variable we need
00140    // to remove it.
00141    t->ResetBranchAddresses();
00142 }
00143 
00144 
00145 void hvector() 
00146 {
00147    gBenchmark->Start("hvector");
00148 
00149    write();
00150    read();
00151    
00152    gBenchmark->Show("hvector");   
00153 }

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