This example illustrates how to make a Tree from variables or arrays in a C struct - without a dictionary, by creating the branches for builtin types (int, float, double) and arrays explicitly.
See tree2a.C for the same example using a class with dictionary instead of a C-struct.
In this example, we are mapping a C struct to one of the Geant3 common blocks /gctrak/. In the real life, this common will be filled by Geant3 at each step and only the Tree Fill function should be called. The example emulates the Geant3 step routines.
to run the example, do:
.x tree2.C to execute with the Cling interpreter
.x tree2.C++ to execute with native compiler
typedef struct {
} Gctrak_t;
{
enum Evect {kX,kY,kZ,kPX,kPY,kPZ,kPP};
vout[kPP] = vect[kPP];
vout[kX] = vect[kX] + (f1*vect[kPX] - f2*vect[kPY]);
vout[kY] = vect[kY] + (f1*vect[kPY] + f2*vect[kPX]);
vout[kZ] = vect[kZ] + (f1*vect[kPZ] + f3);
vout[kPX] = vect[kPX] + (f4*vect[kPX] - f5*vect[kPY]);
vout[kPY] = vect[kPY] + (f4*vect[kPY] + f5*vect[kPX]);
vout[kPZ] = vect[kPZ] + (f4*vect[kPZ] + f6);
}
void tree2w()
{
TFile f(
"tree2.root",
"recreate");
TTree t2(
"t2",
"a Tree with data from a fake Geant3");
Gctrak_t gstep;
t2.Branch("vect",gstep.vect,"vect[7]/F");
t2.Branch("getot",&gstep.getot);
t2.Branch("gekin",&gstep.gekin);
t2.Branch("nmec",&gstep.nmec);
t2.Branch("lmec",gstep.lmec,"lmec[nmec]/I");
t2.Branch("destep",&gstep.destep);
t2.Branch("pid",&gstep.pid);
gstep.step = 0.1;
gstep.destep = 0;
gstep.nmec = 0;
gstep.pid = 0;
for (
Int_t i=0;i<10000;i++) {
if (newParticle) {
px = gRandom->Gaus(0,.02);
py = gRandom->Gaus(0,.02);
pz = gRandom->Gaus(0,.02);
charge = 1; if (gRandom->Rndm() < 0.5) charge = -1;
gstep.pid += 1;
gstep.vect[0] = 0;
gstep.vect[1] = 0;
gstep.vect[2] = 0;
gstep.vect[3] = px/p;
gstep.vect[4] = py/p;
gstep.vect[5] = pz/p;
gstep.vect[6] = p*charge;
gstep.gekin = gstep.getot - mass;
}
t2.Fill();
helixStep(gstep.step, gstep.vect, vout);
gstep.destep = gstep.step*gRandom->Gaus(0.0002,0.00001);
gstep.gekin -= gstep.destep;
gstep.getot = gstep.gekin + mass;
gstep.vect[6] = charge*
TMath::Sqrt(gstep.getot*gstep.getot - mass*mass);
gstep.vect[0] = vout[0];
gstep.vect[1] = vout[1];
gstep.vect[2] = vout[2];
gstep.vect[3] = vout[3];
gstep.vect[4] = vout[4];
gstep.vect[5] = vout[5];
gstep.nmec = (
Int_t)(5*gRandom->Rndm());
for (
Int_t l=0;l<gstep.nmec;l++) gstep.lmec[l] = l;
if (gstep.gekin < 0.001) newParticle =
kTRUE;
}
}
void tree2r()
{
TH1F *hdestep =
new TH1F(
"hdestep",
"destep in Mev",100,1e-5,3e-5);
}
t2->
Draw(
"vect[0]:vect[1]:vect[2]");
}
void tree2() {
tree2w();
tree2r();
}
- Author
- Rene Brun
Definition in file tree2.C.