// Program to check a TGeo geometry // The first time you run this program, the geometry files will be taken // from http://root.cern.ch/files // // How the program works // If the file _ref.root does not exist, it is generated. The file // contains a TTree with Npoints (default=100000) obtained with the following // algorithm: // -a point is generated with a uniform distribution x,y,z in the master volume // -a direction theta, phi is generated uniformly in -2pi_ref.root (generated typically with a previous version // of the TGeo classes), the Npoints in the Tree are used to perform the // same operation with the new version. // In case of a disagreement, an error message is reported. // // The ReadRef case is also used as a benchmark // The ROOTMARKS reported are relative to a Linux/P IV 2.8 GHz gcc3.2.3 machine // normalized at 800 ROOTMARKS when running with CINT. // // To run this script, do // stressGeometry // or stressGeometry * // or stressGeometry alice // or from the ROOT command line // root > .L stressGeometry.cxx or .L stressGeometry.cxx+ // root > stressGeometry(exp_name); // where exp_name is the geometry file name without .root // OR simply: stressGeometry(); to run tests for a set of geometries // // Authors: Rene Brun, Andrei Gheata, 22 march 2005 #include "TStopwatch.h" #include "TGeoManager.h" #include "TGeoMatrix.h" #include "TGeoNode.h" #include "TGeoMedium.h" #include "TGeoMaterial.h" #include "TGeoBBox.h" #include "TROOT.h" #include "TFile.h" #include "TTree.h" #include "TRandom3.h" #include "TVectorD.h" #include "TCanvas.h" #include "TError.h" #include "TApplication.h" #include "TMath.h" #include "TSystem.h" #include "TVirtualGeoConverter.h" // Total and reference times Double_t tpstot = 0; Double_t tpsref = 112.1; //time including the generation of the ref files Bool_t testfailed = kFALSE; #ifndef __CINT__ void stressGeometry(const char*, Bool_t, Bool_t); int main(int argc, char **argv) { gROOT->SetBatch(); TApplication theApp("App", &argc, argv); Bool_t vecgeom = kFALSE; TString geom = "*"; if (argc > 1) geom = argv[1]; geom.ToLower(); if (geom == "all") geom = "*"; printf("geom: %s\n", geom.Data()); if (argc > 1) { for (Int_t iarg=1; iargConvertGeometry(); fname = TString::Format("files/%s_ref_%d.root", exps[i],versions[i]); if (gen_ref || !TFile::Open(Form("http://root.cern.ch/files/%s_ref_%d.root",exps[i],versions[i]),"CACHEREAD")) { if (!gen_ref) fprintf(stderr,"File: %s does not exist, generating it\n", fname.Data()); else fprintf(stderr,"Generating reference file %s\n", fname.Data()); WriteRef(i); } ReadRef(i); } if (all && tpstot>0) { Float_t rootmarks = 800*tpsref/tpstot; Bool_t UNIX = strcmp(gSystem->GetName(), "Unix") == 0; if (UNIX) { TString sp = gSystem->GetFromPipe("uname -a"); sp.Resize(60); printf("* SYS: %s\n",sp.Data()); if (strstr(gSystem->GetBuildNode(),"Linux")) { sp = gSystem->GetFromPipe("lsb_release -d -s"); printf("* SYS: %s\n",sp.Data()); } if (strstr(gSystem->GetBuildNode(),"Darwin")) { sp = gSystem->GetFromPipe("sw_vers -productVersion"); sp += " Mac OS X "; printf("* SYS: %s\n",sp.Data()); } } else { const char *os = gSystem->Getenv("OS"); if (!os) fprintf(stderr,"* SYS: Windows 95\n"); else fprintf(stderr,"* SYS: %s %s \n",os,gSystem->Getenv("PROCESSOR_IDENTIFIER")); } fprintf(stderr,"******************************************************************\n"); if (testfailed) fprintf(stderr,"* stressGeometry found bad points ............. FAILED\n"); else fprintf(stderr,"* stressGeometry .................................. OK\n"); fprintf(stderr,"******************************************************************\n"); fprintf(stderr,"* CPU time in ReadRef = %6.2f seconds\n",tpstot); fprintf(stderr,"* ROOTMARKS =%6.1f * Root%-8s %d/%d\n",rootmarks,gROOT->GetVersion(),gROOT->GetVersionDate(),gROOT->GetVersionTime()); } fprintf(stderr,"******************************************************************\n"); } void ReadRef(Int_t kexp) { TString fname; TFile *f = 0; //use ref_[version[i]] files if (!gen_ref) fname = TString::Format("http://root.cern.ch/files/%s_ref_%d.root", exps[kexp],versions[kexp]); else fname.Format("files/%s_ref_%d.root", exps[kexp],versions[kexp]); f = TFile::Open(fname,"CACHEREAD"); if (!f) { fprintf(stderr,"Reference file %s not found ! Skipping.\n", fname.Data()); return; } // fprintf(stderr,"Reference file %s found\n", fname.Data()); fname = TString::Format("%s_diff.root", exps[kexp]); TFile fdiff(fname,"RECREATE"); TTree *TD = new TTree("TD","TGeo stress diff"); TD->Branch("p",&p.x,"x/D:y/D:z/D:theta/D:phi/D:rad[4]/F"); TTree *T = (TTree*)f->Get("T"); T->SetBranchAddress("p",&p.x); Long64_t nentries = T->GetEntries(); TVectorD *vref = (TVectorD *)T->GetUserInfo()->At(0); if (!vref) { fprintf(stderr," ERROR: User info not found, regenerate reference file\n"); return; } TVectorD vect(4); TVectorD vect_ref = *vref; Int_t nbound; Float_t length, safe, rad; Float_t diff; Float_t diffmax = 0.01; // percent of rad! Int_t nbad = 0; vect(0) = 0;//gGeoManager->Weight(0.01, "va"); TStopwatch sw; for (Long64_t i=0;iGetEntry(i); nbound = 0; length = 0.; safe = 0.; rad = 0.; FindRad(p.x,p.y,p.z, p.theta, p.phi, nbound, length, safe, rad); vect(1) += Double_t(nbound); vect(2) += length; vect(3) += rad; diff = 0; diff += TMath::Abs(length-p.length); diff += TMath::Abs(safe-p.safe); diff += TMath::Abs(rad-p.rad); if (((p.rad>0) && (TMath::Abs(rad-p.rad)/p.rad)>diffmax) || TMath::Abs(nbound-p.nbound)>100) { nbad++; if (nbad < 10) { fprintf(stderr," ==>Point %lld differs with diff = %g, x=%g, y=%g, z=%g\n",i,diff,p.x,p.y,p.z); fprintf(stderr," p.nbound=%d, p.length=%g, p.safe=%g, p.rad=%g\n", p.nbound,p.length,p.safe,p.rad); fprintf(stderr," nbound=%d, length=%g, safe=%g, rad=%g\n", nbound,length,safe,rad); } TD->Fill(); p.nbound = nbound; p.length = length; p.safe = safe; p.rad = rad; TD->Fill(); } } diff = 0.; //for (Int_t j=1; j<4; j++) diff += TMath::Abs(vect_ref(j)-vect(j)); diff += TMath::Abs(vect_ref(3)-vect(3))/vect_ref(3); if (diff > diffmax) { // fprintf(stderr,"Total weight=%g ref=%g\n", vect(0), vect_ref(0)); fprintf(stderr,"Total nbound=%g ref=%g\n", vect(1), vect_ref(1)); fprintf(stderr,"Total length=%g ref=%g\n", vect(2), vect_ref(2)); fprintf(stderr,"Total rad=%g ref=%g\n", vect(3), vect_ref(3)); nbad++; } if (nbad) { testfailed = kTRUE; TD->AutoSave(); TD->Print(); } delete TD; delete f; Double_t cp = sw.CpuTime(); tpstot += cp; if (nbad > 0) fprintf(stderr,"* stress %-15s found %5d bad points ............. failed\n",exps[kexp],nbad); else fprintf(stderr,"* stress %-15s: time/ref = %6.2f/%6.2f............ OK\n",exps[kexp],cp,cp_brun[kexp]); } void WriteRef(Int_t kexp) { TRandom3 r; // Double_t theta, phi; Double_t point[3]; TVectorD vect(4); TGeoShape *top = gGeoManager->GetMasterVolume()->GetShape(); // TGeoBBox *box = (TGeoBBox*)top; Double_t xmax = boxes[kexp][0]; //box->GetDX(); // 300; Double_t ymax = boxes[kexp][1]; //box->GetDY(); // 300; Double_t zmax = boxes[kexp][2]; //box->GetDZ(); // 500; TString fname(TString::Format("files/%s_ref_%d.root", exps[kexp], versions[kexp])); TFile f(fname,"recreate"); TTree *T = new TTree("T","TGeo stress"); T->Branch("p",&p.x,"x/D:y/D:z/D:theta/D:phi/D:nbound/I:length/F:safe/F:rad/F"); T->GetUserInfo()->Add(&vect); Long64_t Npoints = 10000; Long64_t i = 0; vect(0) = 0; //gGeoManager->Weight(0.01, "va"); while (iContains(point)) { p.phi = 2*TMath::Pi()*r.Rndm(); p.theta = TMath::ACos(1.-2.*r.Rndm()); FindRad(p.x,p.y,p.z, p.theta, p.phi, p.nbound, p.length, p.safe, p.rad); vect(1) += Double_t(p.nbound); vect(2) += p.length; vect(3) += p.rad; T->Fill(); i++; } } T->AutoSave(); T->GetUserInfo()->Remove(&vect); // T->Print(); delete T; } 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) { Double_t xp = TMath::Sin(theta)*TMath::Cos(phi); Double_t yp = TMath::Sin(theta)*TMath::Sin(phi); Double_t zp = TMath::Cos(theta); Double_t snext; TString path; Double_t pt[3]; Double_t loc[3]; Double_t epsil = 1.E-2; Double_t lastrad = 0.; Int_t ismall = 0; nbound = 0; length = 0.; safe = 0.; rad = 0.; TGeoMedium *med; TGeoShape *shape; TGeoNode *lastnode; gGeoManager->InitTrack(x,y,z,xp,yp,zp); if (verbose) { fprintf(stderr,"Track: (%15.10f,%15.10f,%15.10f,%15.10f,%15.10f,%15.10f)\n", x,y,z,xp,yp,zp); path = gGeoManager->GetPath(); } TGeoNode *nextnode = gGeoManager->GetCurrentNode(); safe = gGeoManager->Safety(); while (nextnode) { med = nextnode->GetVolume()->GetMedium(); shape = nextnode->GetVolume()->GetShape(); lastnode = nextnode; nextnode = gGeoManager->FindNextBoundaryAndStep(); snext = gGeoManager->GetStep(); if (snext<1.e-8) { ismall++; if ((ismall<3) && (lastnode != nextnode)) { // First try to cross a very thin layer length += snext; nextnode = gGeoManager->FindNextBoundaryAndStep(); snext = gGeoManager->GetStep(); if (snext<1.E-8) continue; // We managed to cross the layer ismall = 0; } else { // Relocate point if (ismall > 3) { fprintf(stderr,"ERROR: Small steps in: %s shape=%s\n",gGeoManager->GetPath(), shape->ClassName()); return; } memcpy(pt,gGeoManager->GetCurrentPoint(),3*sizeof(Double_t)); const Double_t *dir = gGeoManager->GetCurrentDirection(); for (Int_t i=0;i<3;i++) pt[i] += epsil*dir[i]; snext = epsil; length += snext; rad += lastrad*snext; gGeoManager->CdTop(); nextnode = gGeoManager->FindNode(pt[0],pt[1],pt[2]); if (gGeoManager->IsOutside()) return; TGeoMatrix *mat = gGeoManager->GetCurrentMatrix(); mat->MasterToLocal(pt,loc); if (!gGeoManager->GetCurrentVolume()->Contains(loc)) { // fprintf(stderr,"Woops - out\n"); gGeoManager->CdUp(); nextnode = gGeoManager->GetCurrentNode(); } continue; } } else { ismall = 0; } nbound++; length += snext; if (med) { Double_t radlen = med->GetMaterial()->GetRadLen(); if (radlen>1.e-5 && radlen<1.e10) { lastrad = med->GetMaterial()->GetDensity()/radlen; rad += lastrad*snext; } else { lastrad = 0.; } if (verbose) { fprintf(stderr," STEP #%d: %s\n",nbound, path.Data()); fprintf(stderr," step=%g length=%g rad=%g %s\n", snext,length, med->GetMaterial()->GetDensity()*snext/med->GetMaterial()->GetRadLen(),med->GetName()); path = gGeoManager->GetPath(); } } } } void InspectDiff(const char* exp="alice",Long64_t ientry=-1) { Int_t nbound = 0; Float_t length = 0.; Float_t safe = 0.; Float_t rad = 0.; TString fname(TString::Format("%s.root",exp)); if (gSystem->AccessPathName(fname)) { TGeoManager::Import(Form("http://root.cern.ch/files/%s",fname.Data())); } else { TGeoManager::Import(fname); } fname = TString::Format("%s_diff.root",exp); TFile f(fname); if (f.IsZombie()) return; TTree *TD = (TTree*)f.Get("TD"); TD->SetBranchAddress("p",&p.x); Long64_t nentries = TD->GetEntries(); nentries = nentries>>1; if (ientry>=0 && ientryGetEntry(2*ientry); fprintf(stderr," NEW: nbound=%d length=%g safe=%g rad=%g\n", p.nbound,p.length,p.safe,p.rad); TD->GetEntry(2*ientry+1); fprintf(stderr," OLD: nbound=%d length=%g safe=%g rad=%g\n", p.nbound,p.length,p.safe,p.rad); FindRad(p.x,p.y,p.z, p.theta, p.phi, nbound,length,safe,rad, kTRUE); return; } for (Long64_t i=0;iGetEntry(2*i); fprintf(stderr," NEW: nbound=%d length=%g safe=%g rad=%g\n", p.nbound,p.length,p.safe,p.rad); TD->GetEntry(2*i+1); fprintf(stderr," OLD: nbound=%d length=%g safe=%g rad=%g\n", p.nbound,p.length,p.safe,p.rad); FindRad(p.x,p.y,p.z, p.theta, p.phi, nbound,length,safe,rad, kTRUE); } } void InspectRef(const char *exp, Int_t vers) { // Inspect current reference. TString fname(TString::Format("%s_ref_%d.root", exp, vers)); if (gSystem->AccessPathName(fname)) { fprintf(stderr,"ERROR: file %s does not exist\n", fname.Data()); return; } TFile f(fname); if (f.IsZombie()) return; TTree *T = (TTree*)f.Get("T"); Long64_t nentries = T->GetEntries(); fname.Format("Stress test for %s geometry", exp); TCanvas *c = new TCanvas("stress", fname,700,800); c->Divide(2,2,0.005,0.005); c->cd(1); gPad->SetLogy(); T->Draw("p.nbound","","", nentries, 0); c->cd(2); gPad->SetLogy(); T->Draw("p.length","","", nentries, 0); c->cd(3); gPad->SetLogy(); T->Draw("p.safe","","", nentries, 0); c->cd(4); gPad->SetLogy(); T->Draw("p.rad","","", nentries, 0); c->cd(0); c->SetFillColor(kYellow); TVectorD *vref = (TVectorD *)T->GetUserInfo()->At(0); TVectorD vect = *vref; fprintf(stderr,"=====================================\n"); // fprintf(stderr,"Total weight: %g [kg]\n", vect(0)); fprintf(stderr,"Total nbound: %g boundaries crossed\n", vect(1)); fprintf(stderr,"Total length: %g [m]\n", 0.01*vect(2)); fprintf(stderr,"Total nradlen: %f\n", vect(3)); fprintf(stderr,"=====================================\n"); }