tree2a.C

Go to the documentation of this file.
00001 #include "TROOT.h"
00002 #include "TFile.h"
00003 #include "TTree.h"
00004 #include "TBrowser.h"
00005 #include "TH2.h"
00006 #include "TMath.h"
00007 #include "TRandom.h"
00008 #include "TCanvas.h"
00009 
00010 //+ This example is the same as tree2.C, but uses a class instead of a C-struct.
00011 // In this example, we are mapping a class to one of the Geant3
00012 // common blocks /gctrak/. In the real life, this common will be filled
00013 // by Geant3 at each step and only the Tree Fill function should be called.
00014 // The example emulates the Geant3 step routines.
00015 //
00016 // to run the example, do:
00017 // .x tree2a.C+ to execute with native compiler
00018 //
00019 //  Note that this example cannot be run under CINT (ie .x tree2a.c)
00020 //  because an interpreted class cannot derive from a compiled class.
00021 //
00022 //  Author: Rene Brun
00023 
00024 const Int_t MAXMEC = 30;
00025 
00026 class Gctrak : public TObject {
00027 public:
00028   Float_t  vect[7]; 
00029   Float_t  getot; 
00030   Float_t  gekin; 
00031   Float_t  vout[7];   //! not persistent
00032   Int_t    nmec; 
00033   Int_t   *lmec;      //[nmec]
00034   Int_t   *namec;     //[nmec]
00035   Int_t    nstep;     //! not persistent
00036   Int_t    pid; 
00037   Float_t  destep; 
00038   Float_t  destel;    //! not persistent
00039   Float_t  safety;    //! not persistent
00040   Float_t  sleng;     //! not persistent
00041   Float_t  step;      //! not persistent
00042   Float_t  snext;     //! not persistent
00043   Float_t  sfield;    //! not persistent 
00044   Float_t  tofg;      //! not persistent
00045   Float_t  gekrat;    //! not persistent 
00046   Float_t  upwght;    //! not persistent 
00047   
00048   Gctrak() {lmec=0; namec=0;}
00049   
00050   ClassDef(Gctrak,1)
00051 }; 
00052 
00053 
00054 void helixStep(Float_t step, Float_t *vect, Float_t *vout)
00055 {
00056   // extrapolate track in constant field
00057    Float_t field = 20;      //magnetic field in kilogauss
00058    enum Evect {kX,kY,kZ,kPX,kPY,kPZ,kPP};
00059    vout[kPP] = vect[kPP];
00060    Float_t h4    = field*2.99792e-4;
00061    Float_t rho   = -h4/vect[kPP];
00062    Float_t tet   = rho*step;
00063    Float_t tsint = tet*tet/6;
00064    Float_t sintt = 1 - tsint;
00065    Float_t sint  = tet*sintt;
00066    Float_t cos1t = tet/2;
00067    Float_t f1 = step*sintt;
00068    Float_t f2 = step*cos1t;
00069    Float_t f3 = step*tsint*vect[kPZ];
00070    Float_t f4 = -tet*cos1t;
00071    Float_t f5 = sint;
00072    Float_t f6 = tet*cos1t*vect[kPZ];
00073    vout[kX]   = vect[kX]  + (f1*vect[kPX] - f2*vect[kPY]);
00074    vout[kY]   = vect[kY]  + (f1*vect[kPY] + f2*vect[kPX]);
00075    vout[kZ]   = vect[kZ]  + (f1*vect[kPZ] + f3);
00076    vout[kPX]  = vect[kPX] + (f4*vect[kPX] - f5*vect[kPY]);
00077    vout[kPY]  = vect[kPY] + (f4*vect[kPY] + f5*vect[kPX]);
00078    vout[kPZ]  = vect[kPZ] + (f4*vect[kPZ] + f6);   
00079 }
00080 
00081 void tree2aw()
00082 {
00083    //create a Tree file tree2.root
00084    
00085    //create the file, the Tree and a few branches with 
00086    //a subset of gctrak
00087    TFile f("tree2.root","recreate");
00088    TTree t2("t2","a Tree with data from a fake Geant3");
00089    Gctrak *gstep = new Gctrak;
00090    t2.Branch("track",&gstep,8000,1);
00091    
00092    //Initialize particle parameters at first point
00093    Float_t px,py,pz,p,charge=0;
00094    Float_t vout[7];
00095    Float_t mass  = 0.137;
00096    Bool_t newParticle = kTRUE;
00097    gstep->lmec    = new Int_t[MAXMEC];
00098    gstep->namec   = new Int_t[MAXMEC];
00099    gstep->step    = 0.1;
00100    gstep->destep  = 0;
00101    gstep->nmec    = 0;
00102    gstep->pid     = 0;
00103       
00104    //transport particles 
00105    for (Int_t i=0;i<10000;i++) {
00106       //generate a new particle if necessary
00107       if (newParticle) {
00108          px = gRandom->Gaus(0,.02);
00109          py = gRandom->Gaus(0,.02);
00110          pz = gRandom->Gaus(0,.02);
00111          p  = TMath::Sqrt(px*px+py*py+pz*pz);
00112          charge = 1; if (gRandom->Rndm() < 0.5) charge = -1;
00113          gstep->pid    += 1;
00114          gstep->vect[0] = 0;
00115          gstep->vect[1] = 0;
00116          gstep->vect[2] = 0;
00117          gstep->vect[3] = px/p;
00118          gstep->vect[4] = py/p;
00119          gstep->vect[5] = pz/p;
00120          gstep->vect[6] = p*charge;
00121          gstep->getot   = TMath::Sqrt(p*p + mass*mass);
00122          gstep->gekin   = gstep->getot - mass;
00123          newParticle = kFALSE;
00124       }
00125       
00126       // fill the Tree with current step parameters
00127       t2.Fill();
00128       
00129       //transport particle in magnetic field
00130       helixStep(gstep->step, gstep->vect, vout); //make one step
00131       
00132       //apply energy loss
00133       gstep->destep = gstep->step*gRandom->Gaus(0.0002,0.00001);
00134       gstep->gekin -= gstep->destep;
00135       gstep->getot   = gstep->gekin + mass;
00136       gstep->vect[6] = charge*TMath::Sqrt(gstep->getot*gstep->getot - mass*mass);
00137       gstep->vect[0] = vout[0];
00138       gstep->vect[1] = vout[1];
00139       gstep->vect[2] = vout[2];
00140       gstep->vect[3] = vout[3];
00141       gstep->vect[4] = vout[4];
00142       gstep->vect[5] = vout[5];
00143       gstep->nmec    = (Int_t)(5*gRandom->Rndm());
00144       for (Int_t l=0;l<gstep->nmec;l++) {
00145          gstep->lmec[l] = l;
00146          gstep->namec[l] = l+100;
00147       }
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 tree2ar()
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    Gctrak *gstep = 0;
00168    t2->SetBranchAddress("track",&gstep);
00169    TBranch *b_destep = t2->GetBranch("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(gstep->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 tree2a() {
00201    tree2aw();
00202    tree2ar();
00203 }

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