ROOT logo
// 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
      
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);
  Int_t nentries = (Int_t)T->GetEntries();
  
  for (Int_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");
  Int_t nentries = (Int_t)T->GetEntries();
  
  const Int_t LOOPMAX=10;
  JetEvent *events[LOOPMAX];
  Int_t loop;
  for (loop=0;loop<LOOPMAX;loop++) events[loop] = 0;
  for (Int_t ev=0;ev<nev;ev++) {
      if (ev%10 == 0) printf("building pileup: %d\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);
}
 jets.C:1
 jets.C:2
 jets.C:3
 jets.C:4
 jets.C:5
 jets.C:6
 jets.C:7
 jets.C:8
 jets.C:9
 jets.C:10
 jets.C:11
 jets.C:12
 jets.C:13
 jets.C:14
 jets.C:15
 jets.C:16
 jets.C:17
 jets.C:18
 jets.C:19
 jets.C:20
 jets.C:21
 jets.C:22
 jets.C:23
 jets.C:24
 jets.C:25
 jets.C:26
 jets.C:27
 jets.C:28
 jets.C:29
 jets.C:30
 jets.C:31
 jets.C:32
 jets.C:33
 jets.C:34
 jets.C:35
 jets.C:36
 jets.C:37
 jets.C:38
 jets.C:39
 jets.C:40
 jets.C:41
 jets.C:42
 jets.C:43
 jets.C:44
 jets.C:45
 jets.C:46
 jets.C:47
 jets.C:48
 jets.C:49
 jets.C:50
 jets.C:51
 jets.C:52
 jets.C:53
 jets.C:54
 jets.C:55
 jets.C:56
 jets.C:57
 jets.C:58
 jets.C:59
 jets.C:60
 jets.C:61
 jets.C:62
 jets.C:63
 jets.C:64
 jets.C:65
 jets.C:66
 jets.C:67
 jets.C:68
 jets.C:69
 jets.C:70
 jets.C:71
 jets.C:72
 jets.C:73
 jets.C:74