stressGeometry.cxx

Go to the documentation of this file.
00001 // Program to check a TGeo geometry
00002 // The first time you run this program, the geometry files will be taken
00003 // from http://root.cern.ch/files
00004 //   
00005 //    How the program works
00006 // If the file <geom_name>_ref.root does not exist, it is generated. The file
00007 // contains a TTree with Npoints (default=100000) obtained with the following
00008 // algorithm:
00009 //   -a point is generated with a uniform distribution x,y,z in the master volume
00010 //   -a direction theta, phi is generated uniformly in -2pi<phi<2pi and 0<theta<pi
00011 //   -gGeoManager finds the geometry path for the point
00012 //   -the number of boundaries (nbound), total length (length), safety distance
00013 // from the starting point (safe) and number of radiation lengths (rad) from x,y,z 
00014 // is calculated to the exit of the detector. The total number of crossings, detector
00015 // weight and total number of radiation lengths for all tracks are stored as user info in the tree.
00016 //
00017 //  Using the file <geom_name>_ref.root (generated typically with a previous version
00018 //  of the TGeo classes), the Npoints in the Tree are used to perform the
00019 //  same operation with the new version.
00020 //  In case of a disagreement, an error message is reported.
00021 //
00022 //  The ReadRef case is also used as a benchmark
00023 //  The ROOTMARKS reported are relative to a Linux/P IV 2.8 GHz gcc3.2.3 machine
00024 //  normalized at 800 ROOTMARKS when running with CINT.
00025 //
00026 // To run this script, do
00027 //   stressGeometry
00028 // or  stressGeometry *
00029 // or  stressGeometry alice
00030 // or from the ROOT command line
00031 // root > .L stressGeometry.cxx  or .L stressGeometry.cxx+
00032 // root > stressGeometry(exp_name); // where exp_name is the geometry file name without .root
00033 // OR simply: stressGeometry(); to run tests for a set of geometries
00034 //
00035 // Authors: Rene Brun, Andrei Gheata, 22 march 2005
00036 
00037 #include "TStopwatch.h"
00038 #include "TGeoManager.h"
00039 #include "TGeoMatrix.h"
00040 #include "TGeoNode.h"
00041 #include "TGeoMedium.h"
00042 #include "TGeoMaterial.h"
00043 #include "TGeoBBox.h"
00044 #include "TROOT.h"
00045 #include "TFile.h"
00046 #include "TTree.h"
00047 #include "TRandom3.h"
00048 #include "TVectorD.h"
00049 #include "TCanvas.h"
00050 #include "TError.h"
00051 #include "TApplication.h"
00052 #include "TMath.h"
00053 #include "TSystem.h"
00054 
00055 #ifndef __CINT__
00056 void stressGeometry(const char*, Bool_t);
00057 
00058 int main(int argc, char **argv)
00059 {
00060    TApplication theApp("App", &argc, argv);
00061    if (argc > 1) stressGeometry(argv[1],kFALSE);
00062    else          stressGeometry("*",kFALSE);
00063    return 0;
00064 }
00065 
00066 #endif
00067    
00068 // data structure for one point
00069 typedef struct {
00070    Double_t x,y,z,theta,phi;      // Initial track position and direction
00071    Int_t    nbound;               // Number of boundaries crossed until exit
00072    Float_t  length;               // Total length up to exit
00073    Float_t  safe;                 // Safety distance for the initial location
00074    Float_t  rad;                  // Number of radiation lengths up to exit
00075 } p_t;
00076 p_t p;
00077   
00078 const Int_t NG = 33;
00079 const char *exps[NG] = {"aleph",  
00080                         "barres",
00081                         "felix",
00082                         "phenix",
00083                         "chambers",
00084                         "p326",
00085                         "bes",
00086                         "dubna",
00087                         "ganil",
00088                         "e907",
00089                         "phobos2",
00090                         "hermes",
00091                         "na35",
00092                         "na47",
00093                         "na49",
00094                         "wa91",
00095                         "sdc",
00096                         "integral",
00097                         "ams", 
00098                         "brahms",
00099                         "gem",
00100                         "tesla",
00101                         "btev",
00102                         "cdf",  
00103                         "hades2", 
00104                         "lhcbfull",
00105                         "star", 
00106                         "sld",   
00107                         "cms",   
00108                         "alice2",
00109                         "babar2", 
00110                         "belle",
00111                         "atlas" 
00112 };
00113 // The timings below are on my machine PIV 3GHz
00114 const Double_t cp_brun[NG] = {1.9,  //aleph
00115                               0.1,  //barres
00116                               0.12, //felix
00117                               0.62, //phenix
00118                               0.1,  //chambers
00119                               0.19, //p326
00120                               1.2,  //bes
00121                               0.12, //dubna
00122                               0.11, //ganil
00123                               0.17, //e907
00124                               0.22, //phobos2
00125                               0.24, //hermes
00126                               0.14, //na35
00127                               0.21, //na47
00128                               0.23, //na49
00129                               0.16, //wa91
00130                               0.17, //sdc
00131                               0.63, //integral
00132                               0.9,  //ams
00133                               1.1,  //brahms
00134                               1.8,  //gem
00135                               1.5,  //tesla
00136                               1.6,  //btev
00137                               2.2,  //cdf
00138                               1.2,  //hades2
00139                               1.6,  //lhcbfull
00140                               2.7,  //star
00141                               3.3,  //sld
00142                               7.5,  //cms
00143                               8.0,  //alice2
00144                              19.6,  //babar2
00145                              24.1,  //belle
00146                              26.7}; //atlas
00147 // Bounding boxes for experiments
00148 Double_t boxes[NG][3] = {{600,600,500},     // aleph
00149                          {100,100,220},     // barres
00150                          {200,200,12000},   // felix
00151                          {750,750,1000},    // phenix
00152                          {500,500,500},     // chambers
00153                          {201,201,26000},   // p326
00154                          {400,400,240},     // bes
00155                          {500,500,2000},    // dubna
00156                          {500,500,500},     // ganil
00157                          {250,250,2000},    // e907
00158                          {400,40,520},      // phobos2
00159                          {250,250,770},     // hermes
00160                          {310,160,1500},    // na35
00161                          {750,500,3000},    // na47
00162                          {600,200,2000},    // na49
00163                          {175,325,680},     // wa91
00164                          {1400,1400,2100},  // sdc
00165                          {100,100,200},     // integral
00166                          {200,200,200},     // ams
00167                          {50,50,50},        // brahms
00168                          {2000,2000,5000},  // gem
00169                          {1500,1500,1500},  // tesla
00170                          {600,475,1270},    // btev
00171                          {500,500,500},     // cdf
00172                          {250,250,200},     // hades2
00173                          {6700,5000,19000}, // lhcbfull
00174                          {350,350,350},     // star
00175                          {500,500,500},     // sld
00176                          {800,800,1000},    // cms
00177                          {400,400,400},     // alice2
00178                          {300,300,400},     // babar2
00179                          {440,440,538},     // belle
00180                          {1000,1000,1500}   // atlas
00181 };                     
00182 // Total and reference times
00183 Double_t tpstot = 0;
00184 Double_t tpsref = 112.1; //time including the generation of the ref files
00185 Bool_t testfailed = kFALSE;
00186                          
00187 Int_t iexp[NG];
00188 Bool_t gen_ref=kFALSE;
00189 void FindRad(Double_t x, Double_t y, Double_t z,Double_t theta, Double_t phi, Int_t &nbound, Float_t &length, Float_t &safe, Float_t &rad, Bool_t verbose=kFALSE);
00190 void ReadRef(Int_t kexp);
00191 void WriteRef(Int_t kexp);
00192 void InspectRef(const char *exp="alice");
00193 
00194 void stressGeometry(const char *exp="*", Bool_t generate_ref=kFALSE) {
00195    gen_ref = generate_ref;
00196    gErrorIgnoreLevel = 10;
00197    
00198    fprintf(stderr,"******************************************************************\n");
00199    fprintf(stderr,"* STRESS GEOMETRY\n");
00200    TString opt = exp;
00201    opt.ToLower();
00202    Bool_t all = kFALSE;
00203    if (opt.Contains("*")) all = kTRUE;
00204    Int_t i;
00205    for (i=0; i<NG; i++) {
00206       if (all) {
00207          iexp[i] = 1;
00208          continue;
00209       }
00210       if (opt.Contains(exps[i])) iexp[i] = 1;
00211       else                       iexp[i] = 0;
00212    }       
00213    TFile::SetCacheFileDir(".");
00214    TString fname;
00215    for (i=0; i<NG; i++) {
00216       if (!iexp[i]) continue;
00217       fname = TString::Format("%s.root", exps[i]);
00218       if (gGeoManager) {
00219          delete gGeoManager;
00220          gGeoManager = 0;
00221       }   
00222       TGeoManager::Import(Form("http://root.cern.ch/files/%s",fname.Data()));
00223          
00224       fname = TString::Format("files/%s_ref_3.root", exps[i]);
00225       
00226       if (gen_ref || !TFile::Open(Form("http://root.cern.ch/files/%s_ref_3.root",exps[i]),"CACHEREAD")) {
00227          if (!gen_ref) fprintf(stderr,"File: %s does not exist, generating it\n", fname.Data());
00228          else               fprintf(stderr,"Generating reference file %s\n", fname.Data());
00229          WriteRef(i);
00230       }
00231    
00232       ReadRef(i);
00233    }   
00234    if (all && tpstot>0) {
00235       Float_t rootmarks = 800*tpsref/tpstot;
00236       Bool_t UNIX = strcmp(gSystem->GetName(), "Unix") == 0;
00237       if (UNIX) {
00238          TString sp = gSystem->GetFromPipe("uname -a");
00239          sp.Resize(60);
00240          printf("*  SYS: %s\n",sp.Data());
00241          if (strstr(gSystem->GetBuildNode(),"Linux")) {
00242             sp = gSystem->GetFromPipe("lsb_release -d -s");
00243             printf("*  SYS: %s\n",sp.Data());
00244          }
00245          if (strstr(gSystem->GetBuildNode(),"Darwin")) {
00246             sp  = gSystem->GetFromPipe("sw_vers -productVersion");
00247             sp += " Mac OS X ";
00248             printf("*  SYS: %s\n",sp.Data());
00249          }
00250       } else {
00251          const char *os = gSystem->Getenv("OS");
00252          if (!os) fprintf(stderr,"*  SYS: Windows 95\n");
00253          else     fprintf(stderr,"*  SYS: %s %s \n",os,gSystem->Getenv("PROCESSOR_IDENTIFIER"));
00254       }
00255       fprintf(stderr,"******************************************************************\n");
00256       if (testfailed) fprintf(stderr,"*  stressGeometry found bad points ............. FAILED\n");
00257       else          fprintf(stderr,"*  stressGeometry .................................. OK\n");
00258       fprintf(stderr,"******************************************************************\n");
00259       fprintf(stderr,"*  CPU time in ReadRef = %6.2f seconds\n",tpstot);
00260       fprintf(stderr,"*  ROOTMARKS =%6.1f   *  Root%-8s  %d/%d\n",rootmarks,gROOT->GetVersion(),gROOT->GetVersionDate(),gROOT->GetVersionTime());
00261    }
00262    fprintf(stderr,"******************************************************************\n");
00263 }
00264 
00265 void ReadRef(Int_t kexp) {
00266    TStopwatch sw;
00267    TString fname;
00268    TFile *f = 0;
00269    //use ref_3 files from version 5.23/01
00270    if (!gen_ref)
00271       fname = TString::Format("http://root.cern.ch/files/%s_ref_3.root", exps[kexp]);
00272    else
00273       fname.Format("files/%s_ref_3.root", exps[kexp]);
00274    
00275    f = TFile::Open(fname,"CACHEREAD");
00276    if (!f) {
00277       fprintf(stderr,"Reference file %s not found ! Skipping.\n", fname.Data());
00278       return;
00279    }   
00280    fprintf(stderr,"Reference file %s found\n", fname.Data());
00281    fname = TString::Format("%s_diff.root", exps[kexp]);
00282    TFile fdiff(fname,"RECREATE");
00283    TTree *TD = new TTree("TD","TGeo stress diff");
00284    TD->Branch("p",&p.x,"x/D:y/D:z/D:theta/D:phi/D:rad[4]/F");
00285    TTree *T = (TTree*)f->Get("T");
00286    T->SetBranchAddress("p",&p.x);
00287    Long64_t nentries = T->GetEntries();
00288    TVectorD *vref = (TVectorD *)T->GetUserInfo()->At(0);
00289    if (!vref) {
00290       fprintf(stderr," ERROR: User info not found, regenerate reference file\n");
00291       return;
00292    }   
00293    TVectorD vect(4);
00294    TVectorD vect_ref = *vref;
00295    Int_t nbound;
00296    Float_t length, safe, rad;
00297    Float_t diff;
00298    Float_t diffmax = 0.01;  // percent of rad!
00299    Int_t nbad = 0;
00300    vect(0) = 0;//gGeoManager->Weight(0.01, "va");
00301    for (Long64_t i=0;i<nentries;i++) {
00302       T->GetEntry(i);
00303       nbound = 0;
00304       length = 0.;
00305       safe = 0.;
00306       rad = 0.;
00307       FindRad(p.x,p.y,p.z, p.theta, p.phi, nbound, length, safe, rad);
00308       vect(1) += Double_t(nbound);
00309       vect(2) += length;
00310       vect(3) += rad;
00311       diff = 0;
00312       diff += TMath::Abs(length-p.length);
00313       diff += TMath::Abs(safe-p.safe);
00314       diff += TMath::Abs(rad-p.rad);
00315       if (((p.rad>0) && (TMath::Abs(rad-p.rad)/p.rad)>diffmax) || 
00316            TMath::Abs(nbound-p.nbound)>100) {
00317          nbad++;
00318          if (nbad < 10) {
00319             fprintf(stderr," ==>Point %lld differs with diff = %g, x=%g, y=%g, z=%g\n",i,diff,p.x,p.y,p.z);
00320             fprintf(stderr,"    p.nbound=%d, p.length=%g, p.safe=%g, p.rad=%g\n",
00321                         p.nbound,p.length,p.safe,p.rad);
00322             fprintf(stderr,"      nbound=%d,   length=%g,   safe=%g,   rad=%g\n",
00323                         nbound,length,safe,rad);
00324          }
00325          TD->Fill();
00326          p.nbound = nbound;
00327          p.length = length;
00328          p.safe   = safe;
00329          p.rad    = rad;
00330          TD->Fill();
00331       }    
00332    }
00333    diff = 0.;
00334    //for (Int_t j=1; j<4; j++) diff += TMath::Abs(vect_ref(j)-vect(j));
00335    diff += TMath::Abs(vect_ref(3)-vect(3))/vect_ref(3);
00336    if (diff > diffmax) {
00337 //      fprintf(stderr,"Total weight=%g   ref=%g\n", vect(0), vect_ref(0));
00338       fprintf(stderr,"Total nbound=%g   ref=%g\n", vect(1), vect_ref(1));
00339       fprintf(stderr,"Total length=%g   ref=%g\n", vect(2), vect_ref(2));
00340       fprintf(stderr,"Total    rad=%g   ref=%g\n", vect(3), vect_ref(3));
00341       nbad++;  
00342    }   
00343       
00344    if (nbad) {
00345       testfailed = kTRUE;
00346       TD->AutoSave();
00347       TD->Print();
00348    }   
00349    delete TD;
00350    delete f;
00351    
00352    Double_t cp = sw.CpuTime();
00353    tpstot += cp;
00354    if (nbad > 0) fprintf(stderr,"*     stress %-15s  found %5d bad points ............. failed\n",exps[kexp],nbad);
00355    else          fprintf(stderr,"*     stress %-15s: time/ref = %6.2f/%6.2f............ OK\n",exps[kexp],cp,cp_brun[kexp]);
00356 }
00357 
00358 void WriteRef(Int_t kexp) {
00359    TRandom3 r;
00360 //   Double_t theta, phi;
00361    Double_t point[3];
00362    TVectorD vect(4);
00363    TGeoShape *top = gGeoManager->GetMasterVolume()->GetShape();
00364 //   TGeoBBox *box = (TGeoBBox*)top;
00365    Double_t xmax = boxes[kexp][0]; //box->GetDX(); // 300;
00366    Double_t ymax = boxes[kexp][1]; //box->GetDY(); // 300;
00367    Double_t zmax = boxes[kexp][2]; //box->GetDZ(); // 500;
00368    TString fname(TString::Format("files/%s_ref_3.root", exps[kexp]));
00369    TFile f(fname,"recreate");
00370    TTree *T = new TTree("T","TGeo stress");
00371    T->Branch("p",&p.x,"x/D:y/D:z/D:theta/D:phi/D:nbound/I:length/F:safe/F:rad/F");
00372    T->GetUserInfo()->Add(&vect);
00373    Long64_t Npoints = 10000;
00374    Long64_t i = 0;
00375    vect(0) = 0; //gGeoManager->Weight(0.01, "va");
00376    while (i<Npoints) {
00377       p.x  = r.Uniform(-xmax,xmax);
00378       p.y  = r.Uniform(-ymax,ymax);
00379       p.z  = r.Uniform(-zmax,zmax);
00380       point[0] = p.x;
00381       point[1] = p.y;
00382       point[2] = p.z;
00383       if (top->Contains(point)) {
00384          p.phi   =  2*TMath::Pi()*r.Rndm();
00385          p.theta = TMath::ACos(1.-2.*r.Rndm());
00386          FindRad(p.x,p.y,p.z, p.theta, p.phi, p.nbound, p.length, p.safe, p.rad);
00387          vect(1) += Double_t(p.nbound);
00388          vect(2) += p.length;
00389          vect(3) += p.rad;
00390          T->Fill();
00391          i++;
00392       }
00393    }   
00394    T->AutoSave();
00395    T->GetUserInfo()->Remove(&vect);
00396 //   T->Print();
00397    delete T;
00398 }
00399 
00400 void FindRad(Double_t x, Double_t y, Double_t z,Double_t theta, Double_t phi, Int_t &nbound, Float_t &length, Float_t &safe, Float_t &rad, Bool_t verbose) {
00401    Double_t xp  = TMath::Sin(theta)*TMath::Cos(phi);
00402    Double_t yp  = TMath::Sin(theta)*TMath::Sin(phi);
00403    Double_t zp  = TMath::Cos(theta);
00404    Double_t snext;
00405    TString path;
00406    Double_t pt[3];
00407    Double_t loc[3];
00408    Double_t epsil = 1.E-2;
00409    Double_t lastrad = 0.;
00410    Int_t ismall = 0;
00411    nbound = 0;
00412    length = 0.;
00413    safe   = 0.;
00414    rad    = 0.;
00415    TGeoMedium *med;
00416    TGeoShape *shape;
00417    TGeoNode *lastnode;
00418    gGeoManager->InitTrack(x,y,z,xp,yp,zp);
00419    if (verbose) {
00420       fprintf(stderr,"Track: (%15.10f,%15.10f,%15.10f,%15.10f,%15.10f,%15.10f)\n",
00421                        x,y,z,xp,yp,zp);
00422       path = gGeoManager->GetPath();
00423    }                    
00424    TGeoNode *nextnode = gGeoManager->GetCurrentNode();
00425    safe = gGeoManager->Safety();
00426    while (nextnode) {
00427       med = 0;
00428       if (nextnode) med = nextnode->GetVolume()->GetMedium();
00429       else return;      
00430       shape = nextnode->GetVolume()->GetShape();
00431       lastnode = nextnode;
00432       nextnode = gGeoManager->FindNextBoundaryAndStep();
00433       snext  = gGeoManager->GetStep();
00434       if (snext<1.e-8) {
00435          ismall++;
00436          if ((ismall<3) && (lastnode != nextnode)) {
00437             // First try to cross a very thin layer
00438             length += snext;
00439             nextnode = gGeoManager->FindNextBoundaryAndStep();
00440             snext  = gGeoManager->GetStep();
00441             if (snext<1.E-8) continue;
00442             // We managed to cross the layer
00443             ismall = 0;
00444          } else {  
00445             // Relocate point
00446             if (ismall > 3) {
00447                fprintf(stderr,"ERROR: Small steps in: %s shape=%s\n",gGeoManager->GetPath(), shape->ClassName());
00448                return;
00449             }   
00450             memcpy(pt,gGeoManager->GetCurrentPoint(),3*sizeof(Double_t));
00451             const Double_t *dir = gGeoManager->GetCurrentDirection();
00452             for (Int_t i=0;i<3;i++) pt[i] += epsil*dir[i];
00453             snext = epsil;
00454             length += snext;
00455             rad += lastrad*snext;
00456             gGeoManager->CdTop();
00457             nextnode = gGeoManager->FindNode(pt[0],pt[1],pt[2]);
00458             if (gGeoManager->IsOutside()) return;
00459             TGeoMatrix *mat = gGeoManager->GetCurrentMatrix();
00460             mat->MasterToLocal(pt,loc);
00461             if (!gGeoManager->GetCurrentVolume()->Contains(loc)) {
00462 //            fprintf(stderr,"Woops - out\n");
00463                gGeoManager->CdUp();
00464                nextnode = gGeoManager->GetCurrentNode();
00465             }   
00466             continue;
00467          }   
00468       } else {
00469          ismall = 0;
00470       }      
00471       nbound++;
00472       length += snext;
00473       if (med) {
00474          Double_t radlen = med->GetMaterial()->GetRadLen();
00475          if (radlen>1.e-5 && radlen<1.e10) {
00476             lastrad = med->GetMaterial()->GetDensity()/radlen;
00477             rad += lastrad*snext;
00478          } else {
00479             lastrad = 0.;
00480          }      
00481          if (verbose) {
00482             fprintf(stderr," STEP #%d: %s\n",nbound, path.Data());
00483             fprintf(stderr,"    step=%g  length=%g  rad=%g %s\n", snext,length,
00484                    med->GetMaterial()->GetDensity()*snext/med->GetMaterial()->GetRadLen(),med->GetName());
00485             path =  gGeoManager->GetPath();
00486          }   
00487       }
00488    }   
00489 }
00490   
00491 void InspectDiff(const char* exp="alice",Long64_t ientry=-1) {
00492    Int_t nbound = 0;   
00493    Float_t length = 0.;
00494    Float_t safe   = 0.;
00495    Float_t rad    = 0.;
00496    TString fname(TString::Format("%s.root",exp));
00497    if (gSystem->AccessPathName(fname)) {
00498       TGeoManager::Import(Form("http://root.cern.ch/files/%s",fname.Data()));
00499    } else {
00500       TGeoManager::Import(fname);
00501    }
00502    fname = TString::Format("%s_diff.root",exp);   
00503    TFile f(fname);
00504    if (f.IsZombie()) return;
00505    TTree *TD = (TTree*)f.Get("TD");
00506    TD->SetBranchAddress("p",&p.x);
00507    Long64_t nentries = TD->GetEntries();
00508    nentries = nentries>>1;
00509    if (ientry>=0 && ientry<nentries) {
00510       fprintf(stderr,"DIFFERENCE #%lld\n", ientry);
00511       TD->GetEntry(2*ientry);
00512       fprintf(stderr,"   NEW: nbound=%d  length=%g  safe=%g  rad=%g\n", p.nbound,p.length,p.safe,p.rad);
00513       TD->GetEntry(2*ientry+1);
00514       fprintf(stderr,"   OLD: nbound=%d  length=%g  safe=%g  rad=%g\n", p.nbound,p.length,p.safe,p.rad);
00515       FindRad(p.x,p.y,p.z, p.theta, p.phi, nbound,length,safe,rad, kTRUE);
00516       return;
00517    }   
00518    for (Long64_t i=0;i<nentries;i++) {
00519       fprintf(stderr,"DIFFERENCE #%lld\n", i);
00520       TD->GetEntry(2*i);
00521       fprintf(stderr,"   NEW: nbound=%d  length=%g  safe=%g rad=%g\n", p.nbound,p.length,p.safe,p.rad);
00522       TD->GetEntry(2*i+1);
00523       fprintf(stderr,"   OLD: nbound=%d  length=%g  safe=%g rad=%g\n", p.nbound,p.length,p.safe,p.rad);
00524       FindRad(p.x,p.y,p.z, p.theta, p.phi, nbound,length,safe,rad, kTRUE);
00525    }
00526 }   
00527 
00528 void InspectRef(const char *exp) {
00529 // Inspect current reference.
00530    TString fname(TString::Format("%s_ref_3.root", exp));
00531    if (gSystem->AccessPathName(fname)) {
00532       fprintf(stderr,"ERROR: file %s does not exist\n", fname.Data());
00533       return;
00534    }
00535    TFile f(fname);
00536    if (f.IsZombie()) return;
00537    TTree *T = (TTree*)f.Get("T");
00538    Long64_t nentries = T->GetEntries();
00539    fname.Format("Stress test for %s geometry", exp);
00540    TCanvas *c = new TCanvas("stress", fname,700,800);
00541    c->Divide(2,2,0.005,0.005);
00542    c->cd(1);
00543    gPad->SetLogy();
00544    T->Draw("p.nbound","","", nentries, 0);
00545    c->cd(2);
00546    gPad->SetLogy();
00547    T->Draw("p.length","","", nentries, 0);
00548    c->cd(3);
00549    gPad->SetLogy();
00550    T->Draw("p.safe","","", nentries, 0);
00551    c->cd(4);
00552    gPad->SetLogy();
00553    T->Draw("p.rad","","", nentries, 0);
00554    c->cd(0);
00555    c->SetFillColor(kYellow);
00556    TVectorD *vref = (TVectorD *)T->GetUserInfo()->At(0);
00557    TVectorD vect = *vref;
00558    fprintf(stderr,"=====================================\n");
00559 //   fprintf(stderr,"Total weight:  %g [kg]\n", vect(0));
00560    fprintf(stderr,"Total nbound:  %g boundaries crossed\n", vect(1));
00561    fprintf(stderr,"Total length:  %g [m]\n", 0.01*vect(2));
00562    fprintf(stderr,"Total nradlen: %f\n", vect(3));   
00563    fprintf(stderr,"=====================================\n");
00564 }

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