tree2.C

Go to the documentation of this file.
00001 #include "TFile.h"
00002 #include "TTree.h"
00003 #include "TBrowser.h"
00004 #include "TH2.h"
00005 #include "TRandom.h"
00006 #include "TCanvas.h"
00007 
00008 // This example illustrates how to make a Tree from variables or arrays
00009 // in a C struct. Use of C structs is strongly discouraged and one should
00010 // use classes instead. However support for C structs is important for 
00011 // legacy applications written in C or Fortran.
00012 //    see tree2a.C for the same example using a class instead of a C-struct.
00013 //
00014 // In this example, we are mapping a C struct to one of the Geant3
00015 // common blocks /gctrak/. In the real life, this common will be filled
00016 // by Geant3 at each step and only the Tree Fill function should be called.
00017 // The example emulates the Geant3 step routines.
00018 //
00019 // to run the example, do:
00020 // .x tree2.C   to execute with the CINT interpreter
00021 // .x tree2.C++ to execute with native compiler
00022 //
00023 //  Author: Rene Brun
00024 
00025 const Int_t MAXMEC = 30;
00026 //      PARAMETER (MAXMEC=30) 
00027 //      COMMON/GCTRAK/VECT(7),GETOT,GEKIN,VOUT(7),NMEC,LMEC(MAXMEC) 
00028 //     + ,NAMEC(MAXMEC),NSTEP ,PID,DESTEP,DESTEL,SAFETY,SLENG 
00029 //     + ,STEP  ,SNEXT ,SFIELD,TOFG  ,GEKRAT,UPWGHT
00030 typedef struct { 
00031   Float_t  vect[7]; 
00032   Float_t  getot; 
00033   Float_t  gekin; 
00034   Float_t  vout[7]; 
00035   Int_t    nmec; 
00036   Int_t    lmec[MAXMEC]; 
00037   Int_t    namec[MAXMEC]; 
00038   Int_t    nstep; 
00039   Int_t    pid; 
00040   Float_t  destep; 
00041   Float_t  destel; 
00042   Float_t  safety; 
00043   Float_t  sleng; 
00044   Float_t  step; 
00045   Float_t  snext; 
00046   Float_t  sfield; 
00047   Float_t  tofg; 
00048   Float_t  gekrat; 
00049   Float_t  upwght; 
00050 } Gctrak_t; 
00051 
00052 
00053 void helixStep(Float_t step, Float_t *vect, Float_t *vout)
00054 {
00055   // extrapolate track in constant field
00056    Float_t field = 20;      //magnetic field in kilogauss
00057    enum Evect {kX,kY,kZ,kPX,kPY,kPZ,kPP};
00058    vout[kPP] = vect[kPP];
00059    Float_t h4    = field*2.99792e-4;
00060    Float_t rho   = -h4/vect[kPP];
00061    Float_t tet   = rho*step;
00062    Float_t tsint = tet*tet/6;
00063    Float_t sintt = 1 - tsint;
00064    Float_t sint  = tet*sintt;
00065    Float_t cos1t = tet/2;
00066    Float_t f1 = step*sintt;
00067    Float_t f2 = step*cos1t;
00068    Float_t f3 = step*tsint*vect[kPZ];
00069    Float_t f4 = -tet*cos1t;
00070    Float_t f5 = sint;
00071    Float_t f6 = tet*cos1t*vect[kPZ];
00072    vout[kX]   = vect[kX]  + (f1*vect[kPX] - f2*vect[kPY]);
00073    vout[kY]   = vect[kY]  + (f1*vect[kPY] + f2*vect[kPX]);
00074    vout[kZ]   = vect[kZ]  + (f1*vect[kPZ] + f3);
00075    vout[kPX]  = vect[kPX] + (f4*vect[kPX] - f5*vect[kPY]);
00076    vout[kPY]  = vect[kPY] + (f4*vect[kPY] + f5*vect[kPX]);
00077    vout[kPZ]  = vect[kPZ] + (f4*vect[kPZ] + f6);   
00078 }
00079 
00080 void tree2w()
00081 {
00082    //create a Tree file tree2.root
00083    
00084    //create the file, the Tree and a few branches with 
00085    //a subset of gctrak
00086    TFile f("tree2.root","recreate");
00087    TTree t2("t2","a Tree with data from a fake Geant3");
00088    Gctrak_t gstep;
00089    t2.Branch("vect",gstep.vect,"vect[7]/F");
00090    t2.Branch("getot",&gstep.getot,"getot/F");
00091    t2.Branch("gekin",&gstep.gekin,"gekin/F");
00092    t2.Branch("nmec",&gstep.nmec,"nmec/I");
00093    t2.Branch("lmec",gstep.lmec,"lmec[nmec]/I");
00094    t2.Branch("destep",&gstep.destep,"destep/F");
00095    t2.Branch("pid",&gstep.pid,"pid/I");
00096    
00097    //Initialize particle parameters at first point
00098    Float_t px,py,pz,p,charge=0;
00099    Float_t vout[7];
00100    Float_t mass  = 0.137;
00101    Bool_t newParticle = kTRUE;
00102    gstep.step    = 0.1;
00103    gstep.destep  = 0;
00104    gstep.nmec    = 0;
00105    gstep.pid     = 0;
00106       
00107    //transport particles 
00108    for (Int_t i=0;i<10000;i++) {
00109       //generate a new particle if necessary
00110       if (newParticle) {
00111          px = gRandom->Gaus(0,.02);
00112          py = gRandom->Gaus(0,.02);
00113          pz = gRandom->Gaus(0,.02);
00114          p  = TMath::Sqrt(px*px+py*py+pz*pz);
00115          charge = 1; if (gRandom->Rndm() < 0.5) charge = -1;
00116          gstep.pid    += 1;
00117          gstep.vect[0] = 0;
00118          gstep.vect[1] = 0;
00119          gstep.vect[2] = 0;
00120          gstep.vect[3] = px/p;
00121          gstep.vect[4] = py/p;
00122          gstep.vect[5] = pz/p;
00123          gstep.vect[6] = p*charge;
00124          gstep.getot   = TMath::Sqrt(p*p + mass*mass);
00125          gstep.gekin   = gstep.getot - mass;
00126          newParticle = kFALSE;
00127       }
00128       
00129       // fill the Tree with current step parameters
00130       t2.Fill();
00131       
00132       //transport particle in magnetic field
00133       helixStep(gstep.step, gstep.vect, vout); //make one step
00134       
00135       //apply energy loss
00136       gstep.destep = gstep.step*gRandom->Gaus(0.0002,0.00001);
00137       gstep.gekin -= gstep.destep;
00138       gstep.getot   = gstep.gekin + mass;
00139       gstep.vect[6] = charge*TMath::Sqrt(gstep.getot*gstep.getot - mass*mass);
00140       gstep.vect[0] = vout[0];
00141       gstep.vect[1] = vout[1];
00142       gstep.vect[2] = vout[2];
00143       gstep.vect[3] = vout[3];
00144       gstep.vect[4] = vout[4];
00145       gstep.vect[5] = vout[5];
00146       gstep.nmec    = (Int_t)(5*gRandom->Rndm());
00147       for (Int_t l=0;l<gstep.nmec;l++) gstep.lmec[l] = l;
00148       if (gstep.gekin < 0.001)            newParticle = kTRUE;
00149       if (TMath::Abs(gstep.vect[2]) > 30) newParticle = kTRUE;
00150    }
00151   
00152    //save the Tree header. The file will be automatically closed
00153    //when going out of the function scope
00154    t2.Write();
00155 }
00156 
00157 void tree2r()
00158 {
00159    //read the Tree generated by tree2w and fill one histogram
00160    //we are only interested by the destep branch.
00161      
00162    //note that we use "new" to create the TFile and TTree objects !
00163    //because we want to keep these objects alive when we leave 
00164    //this function.
00165    TFile *f = new TFile("tree2.root");
00166    TTree *t2 = (TTree*)f->Get("t2");
00167    static Float_t destep;
00168    TBranch *b_destep = t2->GetBranch("destep");
00169    b_destep->SetAddress(&destep);
00170    
00171    //create one histogram
00172    TH1F *hdestep   = new TH1F("hdestep","destep in Mev",100,1e-5,3e-5);
00173    
00174    //read only the destep branch for all entries
00175    Int_t nentries = (Int_t)t2->GetEntries();
00176    for (Int_t i=0;i<nentries;i++) {
00177       b_destep->GetEntry(i); 
00178       hdestep->Fill(destep);
00179    }
00180   
00181    //we do not close the file. 
00182    //We want to keep the generated histograms
00183    //We fill a 3-d scatter plot with the particle step coordinates
00184    TCanvas *c1 = new TCanvas("c1","c1",600,800);
00185    c1->SetFillColor(42);
00186    c1->Divide(1,2);
00187    c1->cd(1);
00188    hdestep->SetFillColor(45);
00189    hdestep->Fit("gaus");
00190    c1->cd(2);
00191    gPad->SetFillColor(37);
00192    t2->SetMarkerColor(kRed);
00193    t2->Draw("vect[0]:vect[1]:vect[2]");
00194    if (gROOT->IsBatch()) return;
00195    
00196    // invoke the x3d viewer
00197    gPad->GetViewer3D("x3d");
00198 }   
00199 
00200 void tree2() {
00201    tree2w();
00202    tree2r();
00203 }

Generated on Tue Jul 5 15:45:11 2011 for ROOT_528-00b_version by  doxygen 1.5.1