Logo ROOT  
Reference Guide
tree2a.C File Reference

Detailed Description

View in nbviewer Open in SWAN This example is the same as tree2.C, but uses a class instead of a C-struct.

In this example, we are mapping a class 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 to execute with native compiler:

.x tree2a.C+

Note that since IO is involved, ACLiC has to be invoked to create the dictionary of class Gctrak.

#include "TROOT.h"
#include "TFile.h"
#include "TTree.h"
#include "TBrowser.h"
#include "TH2.h"
#include "TMath.h"
#include "TRandom.h"
#include "TCanvas.h"
const Int_t MAXMEC = 30;
class Gctrak : public TObject {
public:
Float_t vect[7];
Float_t getot;
Float_t gekin;
Float_t vout[7]; //! not persistent
Int_t nmec;
Int_t *lmec; //[nmec]
Int_t *namec; //[nmec]
Int_t nstep; //! not persistent
Int_t pid;
Float_t destep;
Float_t destel; //! not persistent
Float_t safety; //! not persistent
Float_t sleng; //! not persistent
Float_t step; //! not persistent
Float_t snext; //! not persistent
Float_t sfield; //! not persistent
Float_t tofg; //! not persistent
Float_t gekrat; //! not persistent
Float_t upwght; //! not persistent
Gctrak() {lmec=0; namec=0;}
ClassDef(Gctrak,1)
};
void helixStep(Float_t step, Float_t *vect, Float_t *vout)
{
// extrapolate track in constant field
Float_t field = 20; //magnetic field in kilogauss
enum Evect {kX,kY,kZ,kPX,kPY,kPZ,kPP};
vout[kPP] = vect[kPP];
Float_t h4 = field*2.99792e-4;
Float_t rho = -h4/vect[kPP];
Float_t tet = rho*step;
Float_t tsint = tet*tet/6;
Float_t sintt = 1 - tsint;
Float_t sint = tet*sintt;
Float_t cos1t = tet/2;
Float_t f1 = step*sintt;
Float_t f2 = step*cos1t;
Float_t f3 = step*tsint*vect[kPZ];
Float_t f4 = -tet*cos1t;
Float_t f5 = sint;
Float_t f6 = tet*cos1t*vect[kPZ];
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 tree2aw()
{
//create a Tree file tree2.root
//create the file, the Tree and a few branches with
//a subset of gctrak
TFile f("tree2.root","recreate");
TTree t2("t2","a Tree with data from a fake Geant3");
Gctrak *gstep = new Gctrak;
t2.Branch("track",&gstep,8000,1);
//Initialize particle parameters at first point
Float_t px,py,pz,p,charge=0;
Float_t vout[7];
Float_t mass = 0.137;
Bool_t newParticle = kTRUE;
gstep->lmec = new Int_t[MAXMEC];
gstep->namec = new Int_t[MAXMEC];
gstep->step = 0.1;
gstep->destep = 0;
gstep->nmec = 0;
gstep->pid = 0;
//transport particles
for (Int_t i=0;i<10000;i++) {
//generate a new particle if necessary
if (newParticle) {
px = gRandom->Gaus(0,.02);
py = gRandom->Gaus(0,.02);
pz = gRandom->Gaus(0,.02);
p = TMath::Sqrt(px*px+py*py+pz*pz);
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->getot = TMath::Sqrt(p*p + mass*mass);
gstep->gekin = gstep->getot - mass;
newParticle = kFALSE;
}
// fill the Tree with current step parameters
t2.Fill();
//transport particle in magnetic field
helixStep(gstep->step, gstep->vect, vout); //make one step
//apply energy loss
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;
gstep->namec[l] = l+100;
}
if (gstep->gekin < 0.001) newParticle = kTRUE;
if (TMath::Abs(gstep->vect[2]) > 30) newParticle = kTRUE;
}
//save the Tree header. The file will be automatically closed
//when going out of the function scope
t2.Write();
}
void tree2ar()
{
//read the Tree generated by tree2w and fill one histogram
//we are only interested by the destep branch.
//note that we use "new" to create the TFile and TTree objects !
//because we want to keep these objects alive when we leave
//this function.
TFile *f = new TFile("tree2.root");
TTree *t2 = (TTree*)f->Get("t2");
Gctrak *gstep = 0;
t2->SetBranchAddress("track",&gstep);
TBranch *b_destep = t2->GetBranch("destep");
//create one histogram
TH1F *hdestep = new TH1F("hdestep","destep in Mev",100,1e-5,3e-5);
//read only the destep branch for all entries
for (Long64_t i=0;i<nentries;i++) {
b_destep->GetEntry(i);
hdestep->Fill(gstep->destep);
}
//we do not close the file.
//We want to keep the generated histograms
//We fill a 3-d scatter plot with the particle step coordinates
TCanvas *c1 = new TCanvas("c1","c1",600,800);
c1->SetFillColor(42);
c1->Divide(1,2);
c1->cd(1);
hdestep->SetFillColor(45);
hdestep->Fit("gaus");
c1->cd(2);
gPad->SetFillColor(37);
t2->Draw("vect[0]:vect[1]:vect[2]");
if (gROOT->IsBatch()) return;
// invoke the x3d viewer
gPad->GetViewer3D("x3d");
}
void tree2a() {
tree2aw();
tree2ar();
}
#define f(i)
Definition: RSha256.hxx:104
#define e(i)
Definition: RSha256.hxx:103
int Int_t
Definition: RtypesCore.h:41
const Bool_t kFALSE
Definition: RtypesCore.h:88
bool Bool_t
Definition: RtypesCore.h:59
long long Long64_t
Definition: RtypesCore.h:69
float Float_t
Definition: RtypesCore.h:53
const Bool_t kTRUE
Definition: RtypesCore.h:87
#define ClassDef(name, id)
Definition: Rtypes.h:326
@ kRed
Definition: Rtypes.h:64
int nentries
Definition: THbookFile.cxx:89
#define gROOT
Definition: TROOT.h:415
R__EXTERN TRandom * gRandom
Definition: TRandom.h:62
#define gPad
Definition: TVirtualPad.h:286
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition: TAttFill.h:37
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition: TAttMarker.h:38
A TTree is a list of TBranches.
Definition: TBranch.h:91
virtual Int_t GetEntry(Long64_t entry=0, Int_t getall=0)
Read all leaves of entry and return total number of bytes read.
Definition: TBranch.cxx:1579
The Canvas class.
Definition: TCanvas.h:31
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format.
Definition: TFile.h:48
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:571
virtual TFitResultPtr Fit(const char *formula, Option_t *option="", Option_t *goption="", Double_t xmin=0, Double_t xmax=0)
Fit histogram with function fname.
Definition: TH1.cxx:3808
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition: TH1.cxx:3275
Mother of all ROOT objects.
Definition: TObject.h:37
virtual Double_t Gaus(Double_t mean=0, Double_t sigma=1)
Samples a random number from the standard Normal (Gaussian) Distribution with the given mean and sigm...
Definition: TRandom.cxx:263
virtual Double_t Rndm()
Machine independent random number generator.
Definition: TRandom.cxx:541
A TTree represents a columnar dataset.
Definition: TTree.h:72
virtual TBranch * GetBranch(const char *name)
Return pointer to the branch with the given name in this tree or its friends.
Definition: TTree.cxx:5170
virtual Int_t SetBranchAddress(const char *bname, void *add, TBranch **ptr=0)
Change branch address, dealing with clone trees properly.
Definition: TTree.cxx:8127
virtual Long64_t GetEntries() const
Definition: TTree.h:450
virtual void Draw(Option_t *opt)
Default Draw method for all objects.
Definition: TTree.h:419
return c1
Definition: legend1.C:41
TF1 * f1
Definition: legend1.C:11
Double_t Sqrt(Double_t x)
Definition: TMath.h:681
Short_t Abs(Short_t d)
Definition: TMathBase.h:120
auto * l
Definition: textangle.C:4
#define snext(osub1, osub2)
Definition: triangle.c:1167
Author
Rene Brun

Definition in file tree2a.C.