Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
tree2.C
Go to the documentation of this file.
1 /// \file
2 /// \ingroup tutorial_tree
3 /// \notebook
4 /// This example illustrates how to make a Tree from variables or arrays
5 /// in a C struct - without a dictionary, by creating the branches for
6 /// builtin types (int, float, double) and arrays explicitly.
7 /// See tree2a.C for the same example using a class with dictionary
8 /// instead of a C-struct.
9 ///
10 /// In this example, we are mapping a C struct to one of the Geant3
11 /// common blocks /gctrak/. In the real life, this common will be filled
12 /// by Geant3 at each step and only the Tree Fill function should be called.
13 /// The example emulates the Geant3 step routines.
14 ///
15 /// to run the example, do:
16 /// ~~~
17 /// .x tree2.C to execute with the Cling interpreter
18 /// .x tree2.C++ to execute with native compiler
19 /// ~~~
20 /// \macro_code
21 ///
22 /// \author Rene Brun
23 
24 #include "TFile.h"
25 #include "TTree.h"
26 #include "TH2.h"
27 #include "TRandom.h"
28 #include "TCanvas.h"
29 #include "TMath.h"
30 
31 const Int_t MAXMEC = 30;
32 
33 typedef struct {
34  Float_t vect[7];
35  Float_t getot;
36  Float_t gekin;
37  Float_t vout[7];
38  Int_t nmec;
39  Int_t lmec[MAXMEC];
40  Int_t namec[MAXMEC];
41  Int_t nstep;
42  Int_t pid;
43  Float_t destep;
44  Float_t destel;
45  Float_t safety;
46  Float_t sleng;
47  Float_t step;
48  Float_t snext;
49  Float_t sfield;
50  Float_t tofg;
51  Float_t gekrat;
52  Float_t upwght;
53 } Gctrak_t;
54 
55 
56 void helixStep(Float_t step, Float_t *vect, Float_t *vout)
57 {
58  // extrapolate track in constant field
59  Float_t field = 20; //magnetic field in kilogauss
60  enum Evect {kX,kY,kZ,kPX,kPY,kPZ,kPP};
61  vout[kPP] = vect[kPP];
62  Float_t h4 = field*2.99792e-4;
63  Float_t rho = -h4/vect[kPP];
64  Float_t tet = rho*step;
65  Float_t tsint = tet*tet/6;
66  Float_t sintt = 1 - tsint;
67  Float_t sint = tet*sintt;
68  Float_t cos1t = tet/2;
69  Float_t f1 = step*sintt;
70  Float_t f2 = step*cos1t;
71  Float_t f3 = step*tsint*vect[kPZ];
72  Float_t f4 = -tet*cos1t;
73  Float_t f5 = sint;
74  Float_t f6 = tet*cos1t*vect[kPZ];
75  vout[kX] = vect[kX] + (f1*vect[kPX] - f2*vect[kPY]);
76  vout[kY] = vect[kY] + (f1*vect[kPY] + f2*vect[kPX]);
77  vout[kZ] = vect[kZ] + (f1*vect[kPZ] + f3);
78  vout[kPX] = vect[kPX] + (f4*vect[kPX] - f5*vect[kPY]);
79  vout[kPY] = vect[kPY] + (f4*vect[kPY] + f5*vect[kPX]);
80  vout[kPZ] = vect[kPZ] + (f4*vect[kPZ] + f6);
81 }
82 
83 void tree2w()
84 {
85  //create a Tree file tree2.root
86 
87  //create the file, the Tree and a few branches with
88  //a subset of gctrak
89  TFile f("tree2.root","recreate");
90  TTree t2("t2","a Tree with data from a fake Geant3");
91  Gctrak_t gstep;
92  t2.Branch("vect",gstep.vect,"vect[7]/F");
93  t2.Branch("getot",&gstep.getot);
94  t2.Branch("gekin",&gstep.gekin);
95  t2.Branch("nmec",&gstep.nmec);
96  t2.Branch("lmec",gstep.lmec,"lmec[nmec]/I");
97  t2.Branch("destep",&gstep.destep);
98  t2.Branch("pid",&gstep.pid);
99 
100  //Initialize particle parameters at first point
101  Float_t px,py,pz,p,charge=0;
102  Float_t vout[7];
103  Float_t mass = 0.137;
104  Bool_t newParticle = kTRUE;
105  gstep.step = 0.1;
106  gstep.destep = 0;
107  gstep.nmec = 0;
108  gstep.pid = 0;
109 
110  //transport particles
111  for (Int_t i=0;i<10000;i++) {
112  //generate a new particle if necessary
113  if (newParticle) {
114  px = gRandom->Gaus(0,.02);
115  py = gRandom->Gaus(0,.02);
116  pz = gRandom->Gaus(0,.02);
117  p = TMath::Sqrt(px*px+py*py+pz*pz);
118  charge = 1; if (gRandom->Rndm() < 0.5) charge = -1;
119  gstep.pid += 1;
120  gstep.vect[0] = 0;
121  gstep.vect[1] = 0;
122  gstep.vect[2] = 0;
123  gstep.vect[3] = px/p;
124  gstep.vect[4] = py/p;
125  gstep.vect[5] = pz/p;
126  gstep.vect[6] = p*charge;
127  gstep.getot = TMath::Sqrt(p*p + mass*mass);
128  gstep.gekin = gstep.getot - mass;
129  newParticle = kFALSE;
130  }
131 
132  // fill the Tree with current step parameters
133  t2.Fill();
134 
135  //transport particle in magnetic field
136  helixStep(gstep.step, gstep.vect, vout); //make one step
137 
138  //apply energy loss
139  gstep.destep = gstep.step*gRandom->Gaus(0.0002,0.00001);
140  gstep.gekin -= gstep.destep;
141  gstep.getot = gstep.gekin + mass;
142  gstep.vect[6] = charge*TMath::Sqrt(gstep.getot*gstep.getot - mass*mass);
143  gstep.vect[0] = vout[0];
144  gstep.vect[1] = vout[1];
145  gstep.vect[2] = vout[2];
146  gstep.vect[3] = vout[3];
147  gstep.vect[4] = vout[4];
148  gstep.vect[5] = vout[5];
149  gstep.nmec = (Int_t)(5*gRandom->Rndm());
150  for (Int_t l=0;l<gstep.nmec;l++) gstep.lmec[l] = l;
151  if (gstep.gekin < 0.001) newParticle = kTRUE;
152  if (TMath::Abs(gstep.vect[2]) > 30) newParticle = kTRUE;
153  }
154 
155  //save the Tree header. The file will be automatically closed
156  //when going out of the function scope
157  t2.Write();
158 }
159 
160 void tree2r()
161 {
162  //read the Tree generated by tree2w and fill one histogram
163  //we are only interested by the destep branch.
164 
165  //note that we use "new" to create the TFile and TTree objects !
166  //because we want to keep these objects alive when we leave
167  //this function.
168  TFile *f = new TFile("tree2.root");
169  TTree *t2 = (TTree*)f->Get("t2");
170  static Float_t destep;
171  TBranch *b_destep = t2->GetBranch("destep");
172  b_destep->SetAddress(&destep);
173 
174  //create one histogram
175  TH1F *hdestep = new TH1F("hdestep","destep in Mev",100,1e-5,3e-5);
176 
177  //read only the destep branch for all entries
178  Long64_t nentries = t2->GetEntries();
179  for (Long64_t i=0;i<nentries;i++) {
180  b_destep->GetEntry(i);
181  hdestep->Fill(destep);
182  }
183 
184  //we do not close the file.
185  //We want to keep the generated histograms
186  //We fill a 3-d scatter plot with the particle step coordinates
187  TCanvas *c1 = new TCanvas("c1","c1",600,800);
188  c1->SetFillColor(42);
189  c1->Divide(1,2);
190  c1->cd(1);
191  hdestep->SetFillColor(45);
192  hdestep->Fit("gaus");
193  c1->cd(2);
194  gPad->SetFillColor(37);
195  t2->SetMarkerColor(kRed);
196  t2->Draw("vect[0]:vect[1]:vect[2]");
197 
198  // Allow to use the TTree after the end of the function.
199  t2->ResetBranchAddresses();
200 }
201 
202 void tree2() {
203  tree2w();
204  tree2r();
205 }