File: jets.C

package info (click to toggle)
root-system 5.34.00-2
  • links: PTS, VCS
  • area: main
  • in suites: wheezy
  • size: 190,532 kB
  • sloc: cpp: 1,401,904; ansic: 217,182; xml: 25,899; sh: 19,215; fortran: 12,570; python: 7,311; makefile: 7,216; ruby: 553; csh: 317; objc: 88; perl: 85; sql: 14; tcl: 4
file content (84 lines) | stat: -rw-r--r-- 2,234 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
// script illustrating the use of a Tree using the JetEvent class.
// The JetEvent class has several collections (TClonesArray)
// and other collections (TRefArray) referencing objects
// in the TClonesArrays.
// The JetEvent class is in $ROOTSYS/tutorials/JetEvent.h,cxx
// to execute the script, do
// .x jets.C
//
//Author: Rene Brun

#include "TFile.h"
#include "TTree.h"
#include "TRandom.h"
#include "TROOT.h"
#include "TSystem.h"
#ifdef __CINT__
#else
#include "JetEvent.h"
#endif
#include "Riostream.h"
    
void write(Int_t nev=100) {
   //write nev Jet events
   TFile f("JetEvent.root","recreate");
   TTree *T = new TTree("T","Event example with Jets");
   JetEvent *event = new JetEvent;
   T->Branch("event","JetEvent",&event,8000,2);
   
   for (Int_t ev=0;ev<nev;ev++) {
      event->Build();
      T->Fill();
   }
   
   T->Print();
   T->Write();
}

void read() {
  //read the JetEvent file
  TFile f("JetEvent.root");
  TTree *T = (TTree*)f.Get("T");
  JetEvent *event = 0;
  T->SetBranchAddress("event", &event);
  Long64_t nentries = T->GetEntries();
  
  for (Long64_t ev=0;ev<nentries;ev++) {
      T->GetEntry(ev);
      if (ev) continue; //dump first event only
      cout << " Event: "<< ev
           << "  Jets: " << event->GetNjet()
	   << "  Tracks: " << event->GetNtrack()
	   << "  Hits A: " << event->GetNhitA()
	   << "  Hits B: " << event->GetNhitB() << endl;
  }
}

void pileup(Int_t nev=200) {
  //make nev pilepup events, each build with LOOPMAX events selected
  //randomly among the nentries
  TFile f("JetEvent.root");
  TTree *T = (TTree*)f.Get("T");
  // Long64_t nentries = T->GetEntries();
  
  const Int_t LOOPMAX=10;
  JetEvent *events[LOOPMAX];
  Int_t loop;
  for (loop=0;loop<LOOPMAX;loop++) events[loop] = 0;
  for (Long64_t ev=0;ev<nev;ev++) {
      if (ev%10 == 0) printf("building pileup: %lld\n",ev);
      for (loop=0;loop<LOOPMAX;loop++) {
         Int_t rev = gRandom->Uniform(LOOPMAX);
         T->SetBranchAddress("event", &events[loop]);
         T->GetEntry(rev);
      }
  }
}

void jets(Int_t nev=100, Int_t npileup=200) {
   gSystem->Load("libPhysics");
   gROOT->ProcessLine(".L $ROOTSYS/tutorials/tree/JetEvent.cxx+");
   write(nev);
   read();
   pileup(npileup);
}