TGeoChecker.cxx

Go to the documentation of this file.
00001 // @(#)root/geom:$Id: TGeoChecker.cxx 36535 2010-11-08 14:41:54Z agheata $
00002 // Author: Andrei Gheata   01/11/01
00003 // CheckGeometry(), CheckOverlaps() by Mihaela Gheata
00004 
00005 /*************************************************************************
00006  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers.               *
00007  * All rights reserved.                                                  *
00008  *                                                                       *
00009  * For the licensing terms see $ROOTSYS/LICENSE.                         *
00010  * For the list of contributors see $ROOTSYS/README/CREDITS.             *
00011  *************************************************************************/
00012 
00013 //______________________________________________________________________________
00014 //   TGeoChecker - Geometry checking package.
00015 //===============
00016 //
00017 // TGeoChecker class provides several geometry checking methods. There are two
00018 // types of tests that can be performed. One is based on random sampling or
00019 // ray-tracing and provides a visual check on how navigation methods work for
00020 // a given geometry. The second actually checks the validity of the geometry
00021 // definition in terms of overlapping/extruding objects. Both types of checks
00022 // can be done for a given branch (starting with a given volume) as well as for 
00023 // the geometry as a whole.
00024 //
00025 // 1. TGeoChecker::CheckPoint(Double_t x, Double_t y, Double_t z)
00026 //
00027 // This method can be called direcly from the TGeoManager class and print a
00028 // report on how a given point is classified by the modeller: which is the
00029 // full path to the deepest node containing it, and the (under)estimation 
00030 // of the distance to the closest boundary (safety).
00031 //
00032 // 2. TGeoChecker::RandomPoints(Int_t npoints)
00033 //
00034 // Can be called from TGeoVolume class. It first draws the volume and its
00035 // content with the current visualization settings. Then randomly samples points
00036 // in its bounding box, plotting in the geometry display only the points 
00037 // classified as belonging to visible volumes. 
00038 //
00039 // 3. TGeoChecker::RandomRays(Int_t nrays, Double_t startx, starty, startz)
00040 //
00041 // Can be called and acts in the same way as the previous, but instead of points,
00042 // rays having random isotropic directions are generated from the given point.
00043 // A raytracing algorithm propagates all rays untill they exit geometry, plotting
00044 // all segments crossing visible nodes in the same color as these.
00045 //
00046 // 4. TGeoChecker::Test(Int_t npoints)
00047 //
00048 // Implementation of TGeoManager::Test(). Computes the time for the modeller
00049 // to find out "Where am I?" for a given number of random points.
00050 //
00051 // 5. TGeoChecker::LegoPlot(ntheta, themin, themax, nphi, phimin, phimax,...)
00052 // 
00053 // Implementation of TGeoVolume::LegoPlot(). Draws a spherical radiation length 
00054 // lego plot for a given volume, in a given theta/phi range.
00055 //
00056 // 6. TGeoChecker::Weigth(Double_t precision)
00057 //
00058 // Implementation of TGeoVolume::Weigth(). Estimates the total weigth of a given 
00059 // volume by matrial sampling. Accepts as input the desired precision.
00060 //
00061 // Overlap checker
00062 //-----------------
00063 //
00064 //Begin_Html
00065 /*
00066 <img src="gif/t_checker.jpg">
00067 */
00068 //End_Html
00069 
00070 #include "TVirtualPad.h"
00071 #include "TCanvas.h"
00072 #include "TStyle.h"
00073 #include "TFile.h"
00074 #include "TNtuple.h"
00075 #include "TH2.h"
00076 #include "TRandom3.h"
00077 #include "TPolyMarker3D.h"
00078 #include "TPolyLine3D.h"
00079 #include "TStopwatch.h"
00080 
00081 #include "TGeoVoxelFinder.h"
00082 #include "TGeoBBox.h"
00083 #include "TGeoPcon.h"
00084 #include "TGeoManager.h"
00085 #include "TGeoOverlap.h"
00086 #include "TGeoPainter.h"
00087 #include "TGeoChecker.h"
00088 
00089 #include "TBuffer3D.h"
00090 #include "TBuffer3DTypes.h"
00091 #include "TMath.h"
00092 
00093 #include <stdlib.h>
00094 
00095 // statics and globals
00096 
00097 ClassImp(TGeoChecker)
00098 
00099 //______________________________________________________________________________
00100 TGeoChecker::TGeoChecker()
00101             :TObject(),
00102              fGeoManager(NULL),
00103              fVsafe(NULL),
00104              fBuff1(NULL),
00105              fBuff2(NULL),
00106              fFullCheck(kFALSE),
00107              fVal1(NULL),
00108              fVal2(NULL),
00109              fFlags(NULL),
00110              fTimer(NULL),
00111              fSelectedNode(NULL),
00112              fNchecks(0),
00113              fNmeshPoints(1000)
00114 {
00115 // Default constructor
00116 }
00117 
00118 //______________________________________________________________________________
00119 TGeoChecker::TGeoChecker(TGeoManager *geom)
00120             :TObject(),
00121              fGeoManager(geom),
00122              fVsafe(NULL),
00123              fBuff1(NULL),
00124              fBuff2(NULL),
00125              fFullCheck(kFALSE),
00126              fVal1(NULL),
00127              fVal2(NULL),
00128              fFlags(NULL),
00129              fTimer(NULL),
00130              fSelectedNode(NULL),
00131              fNchecks(0),
00132              fNmeshPoints(1000)
00133 {
00134 // Constructor for a given geometry
00135    fBuff1 = new TBuffer3D(TBuffer3DTypes::kGeneric,500,3*500,0,0,0,0);
00136    fBuff2 = new TBuffer3D(TBuffer3DTypes::kGeneric,500,3*500,0,0,0,0);   
00137 }
00138 
00139 //______________________________________________________________________________
00140 TGeoChecker::~TGeoChecker()
00141 {
00142 // Destructor
00143    if (fBuff1) delete fBuff1;
00144    if (fBuff2) delete fBuff2;
00145    if (fTimer) delete fTimer;
00146 }
00147 
00148 //______________________________________________________________________________
00149 void TGeoChecker::OpProgress(const char *opname, Long64_t current, Long64_t size, TStopwatch *watch, Bool_t last, Bool_t refresh)
00150 {
00151 // Print current operation progress.
00152    static Long64_t icount = 0;
00153    static TString oname;
00154    static TString nname;
00155    static Long64_t ocurrent = 0;
00156    static Long64_t osize = 0;
00157    static Int_t oseconds = 0;
00158    static TStopwatch *owatch = 0;
00159    static Bool_t oneoftwo = kFALSE;
00160    static Int_t nrefresh = 0;
00161    const char symbol[4] = {'=','\\','|','/'}; 
00162    char progress[11] = "          ";
00163    Int_t ichar = icount%4;
00164    
00165    if (!refresh) {
00166       nrefresh = 0;
00167       if (!size) return;
00168       owatch = watch;
00169       oname = opname;
00170       ocurrent = TMath::Abs(current);
00171       osize = TMath::Abs(size);
00172       if (ocurrent > osize) ocurrent=osize;
00173    } else {
00174       nrefresh++;
00175       if (!osize) return;
00176    }     
00177    icount++;
00178    Double_t time = 0.;
00179    Int_t hours = 0;
00180    Int_t minutes = 0;
00181    Int_t seconds = 0;
00182    if (owatch && !last) {
00183       owatch->Stop();
00184       time = owatch->RealTime();
00185       hours = (Int_t)(time/3600.);
00186       time -= 3600*hours;
00187       minutes = (Int_t)(time/60.);
00188       time -= 60*minutes;
00189       seconds = (Int_t)time;
00190       if (refresh)  {
00191          if (oseconds==seconds) {
00192             owatch->Continue();
00193             return;
00194          }
00195          oneoftwo = !oneoftwo;   
00196       }
00197       oseconds = seconds;   
00198    }
00199    if (refresh && oneoftwo) {
00200       nname = oname;
00201       if (fNchecks <= 0) fNchecks = nrefresh+1;
00202       Int_t pctdone = (Int_t)(100.*nrefresh/fNchecks);
00203       oname = TString::Format("     == %d%% ==", pctdone);
00204    }         
00205    Double_t percent = 100.0*ocurrent/osize;
00206    Int_t nchar = Int_t(percent/10);
00207    if (nchar>10) nchar=10;
00208    Int_t i;
00209    for (i=0; i<nchar; i++)  progress[i] = '=';
00210    progress[nchar] = symbol[ichar];
00211    for (i=nchar+1; i<10; i++) progress[i] = ' ';
00212    progress[10] = '\0';
00213    oname += "                    ";
00214    oname.Remove(20);
00215    if(size<10000) fprintf(stderr, "%s [%10s] %4lld ", oname.Data(), progress, ocurrent);
00216    else if(size<100000) fprintf(stderr, "%s [%10s] %5lld ",oname.Data(), progress, ocurrent);
00217    else fprintf(stderr, "%s [%10s] %7lld ",oname.Data(), progress, ocurrent);
00218    if (time>0.) fprintf(stderr, "[%6.2f %%]   TIME %.2d:%.2d:%.2d             \r", percent, hours, minutes, seconds);
00219    else fprintf(stderr, "[%6.2f %%]\r", percent);
00220    if (refresh && oneoftwo) oname = nname;
00221    if (owatch) owatch->Continue();
00222    if (last) {
00223       icount = 0;
00224       owatch = 0;
00225       ocurrent = 0;
00226       osize = 0;
00227       oseconds = 0;
00228       oneoftwo = kFALSE;
00229       nrefresh = 0;
00230       fprintf(stderr, "\n");
00231    }   
00232 }   
00233 
00234 //_____________________________________________________________________________
00235 void TGeoChecker::CheckBoundaryErrors(Int_t ntracks, Double_t radius)
00236 {
00237 // Check pushes and pulls needed to cross the next boundary with respect to the
00238 // position given by FindNextBoundary. If radius is not mentioned the full bounding
00239 // box will be sampled.
00240    gRandom = new TRandom3();
00241    TGeoVolume *tvol = fGeoManager->GetTopVolume();
00242    Info("CheckBoundaryErrors", "Top volume is %s",tvol->GetName());  
00243    const TGeoShape *shape = tvol->GetShape();
00244    TGeoBBox *box = (TGeoBBox *)shape;
00245    Double_t dl[3];
00246    Double_t ori[3];
00247    Double_t xyz[3];
00248    Double_t nxyz[3];
00249    Double_t dir[3];
00250    Double_t relp;
00251    Char_t path[1024];
00252    Char_t cdir[10];
00253 
00254    // Tree part
00255    TFile *f=new TFile("geobugs.root","recreate");
00256    TTree *bug=new TTree("bug","Geometrical problems");
00257    bug->Branch("pos",xyz,"xyz[3]/D");
00258    bug->Branch("dir",dir,"dir[3]/D");
00259    bug->Branch("push",&relp,"push/D");
00260    bug->Branch("path",&path,"path/C");
00261    bug->Branch("cdir",&cdir,"cdir/C");
00262   
00263    dl[0] = box->GetDX();
00264    dl[1] = box->GetDY();
00265    dl[2] = box->GetDZ();
00266    ori[0] = (box->GetOrigin())[0];
00267    ori[1] = (box->GetOrigin())[1];
00268    ori[2] = (box->GetOrigin())[2];
00269    if (radius>0)
00270       dl[0] = dl[1] = dl[2] = radius;
00271   
00272    TH1::AddDirectory(kFALSE);
00273    TH1F *hnew = new TH1F("hnew","Precision pushing",30,-20.,10.);
00274    TH1F *hold = new TH1F("hold","Precision pulling", 30,-20.,10.);
00275    TH2F *hplotS = new TH2F("hplotS","Problematic points",100,-dl[0],dl[0],100,-dl[1],dl[1]);
00276    gStyle->SetOptStat(111111);
00277 
00278    TGeoNode *node = 0;
00279    Long_t igen=0;
00280    Long_t itry=0;
00281    Long_t n100 = ntracks/100;
00282    Double_t rad = TMath::Sqrt(dl[0]*dl[0]+dl[1]*dl[1]);
00283    printf("Random box : %f, %f, %f, %f, %f, %f\n", ori[0], ori[1], ori[2], dl[0], dl[1], dl[2]);
00284    printf("Start... %i points\n", ntracks);
00285    if (!fTimer) fTimer = new TStopwatch();
00286    fTimer->Reset();
00287    fTimer->Start();
00288    while (igen<ntracks) {
00289       Double_t phi1  = TMath::TwoPi()*gRandom->Rndm();
00290       Double_t r     = rad*gRandom->Rndm();
00291       xyz[0] = ori[0] + r*TMath::Cos(phi1);
00292       xyz[1] = ori[1] + r*TMath::Sin(phi1);
00293       Double_t z = (1.-2.*gRandom->Rndm());
00294       xyz[2] = ori[2]+dl[2]*z*TMath::Abs(z);
00295       ++itry;
00296       fGeoManager->SetCurrentPoint(xyz);
00297       node = fGeoManager->FindNode();
00298       if (!node || node==fGeoManager->GetTopNode()) continue;
00299       ++igen;
00300       if (n100 && !(igen%n100)) 
00301          OpProgress("Sampling progress:",igen, ntracks, fTimer);
00302       Double_t cost = 1.-2.*gRandom->Rndm();
00303       Double_t sint = TMath::Sqrt((1.+cost)*(1.-cost));
00304       Double_t phi  = TMath::TwoPi()*gRandom->Rndm();
00305       dir[0] = sint * TMath::Cos(phi);
00306       dir[1] = sint * TMath::Sin(phi);
00307       dir[2] = cost;
00308       fGeoManager->SetCurrentDirection(dir);
00309       fGeoManager->FindNextBoundary();
00310       Double_t step = fGeoManager->GetStep();
00311 
00312       relp = 1.e-21;
00313       for(Int_t i=0; i<30; ++i) {
00314          relp *=10.;
00315          for(Int_t j=0; j<3; ++j) nxyz[j]=xyz[j]+step*(1.+relp)*dir[j];
00316          if(!fGeoManager->IsSameLocation(nxyz[0],nxyz[1],nxyz[2])) {
00317             hnew->Fill(i-20.);
00318             if(i>15) {
00319                const Double_t* norm = fGeoManager->FindNormal();
00320                strncpy(path,fGeoManager->GetPath(),1024);
00321                path[1023] = '\0';
00322                Double_t dotp = norm[0]*dir[0]+norm[1]*dir[1]+norm[2]*dir[2];
00323                printf("Forward error i=%d p=%5.4f %5.4f %5.4f s=%5.4f dot=%5.4f path=%s\n",
00324                        i,xyz[0],xyz[1],xyz[2],step,dotp,path);
00325                hplotS->Fill(xyz[0],xyz[1],(Double_t)i);
00326                strncpy(cdir,"Forward",10);
00327                bug->Fill();
00328             }
00329                  break;
00330          }
00331       }
00332       
00333       relp = -1.e-21;
00334       for(Int_t i=0; i<30; ++i) {
00335          relp *=10.;
00336          for(Int_t j=0; j<3; ++j) nxyz[j]=xyz[j]+step*(1.+relp)*dir[j];
00337          if(fGeoManager->IsSameLocation(nxyz[0],nxyz[1],nxyz[2])) {
00338             hold->Fill(i-20.);
00339             if(i>15) {
00340                const Double_t* norm = fGeoManager->FindNormal();
00341                strncpy(path,fGeoManager->GetPath(),1024);
00342                path[1023] = '\0';
00343                Double_t dotp = norm[0]*dir[0]+norm[1]*dir[1]+norm[2]*dir[2];
00344                printf("Backward error i=%d p=%5.4f %5.4f %5.4f s=%5.4f dot=%5.4f path=%s\n",
00345                        i,xyz[0],xyz[1],xyz[2],step,dotp,path);
00346                strncpy(cdir,"Backward",10);
00347                bug->Fill();
00348             }
00349             break;
00350          }
00351       }
00352    }
00353    fTimer->Stop();
00354 
00355    printf("CPU time/point = %5.2emus: Real time/point = %5.2emus\n",
00356                1000000.*fTimer->CpuTime()/itry,1000000.*fTimer->RealTime()/itry);
00357    bug->Write();
00358    delete bug;
00359    bug=0;
00360    delete f;
00361 
00362    CheckBoundaryReference();
00363 
00364    printf("Effic = %3.1f%%\n",(100.*igen)/itry);
00365    TCanvas *c1 = new TCanvas("c1","Results",600,800);
00366    c1->Divide(1,2);
00367    c1->cd(1);
00368    gPad->SetLogy();
00369    hold->Draw();
00370    c1->cd(2);
00371    gPad->SetLogy();
00372    hnew->Draw();
00373    /*TCanvas *c3 = */new TCanvas("c3","Plot",600,600);
00374    hplotS->Draw("cont0");
00375 }   
00376 
00377 //_____________________________________________________________________________
00378 void TGeoChecker::CheckBoundaryReference(Int_t icheck)
00379 {
00380 // Check the boundary errors reference file created by CheckBoundaryErrors method.
00381 // The shape for which the crossing failed is drawn with the starting point in red
00382 // and the extrapolated point to boundary (+/- failing push/pull) in yellow.
00383    Double_t xyz[3];
00384    Double_t nxyz[3];
00385    Double_t dir[3];
00386    Double_t lnext[3];
00387    Double_t push;
00388    Char_t path[1024];
00389    Char_t cdir[10];
00390    // Tree part
00391    TFile *f=new TFile("geobugs.root","read");
00392    TTree *bug=(TTree*)f->Get("bug");
00393    bug->SetBranchAddress("pos",xyz);
00394    bug->SetBranchAddress("dir",dir);
00395    bug->SetBranchAddress("push",&push);
00396    bug->SetBranchAddress("path",&path);
00397    bug->SetBranchAddress("cdir",&cdir);
00398 
00399    Int_t nentries = (Int_t)bug->GetEntries();
00400    printf("nentries %d\n",nentries);
00401    if (icheck<0) {
00402       for (Int_t i=0;i<nentries;i++) {
00403          bug->GetEntry(i);
00404          printf("%-9s error push=%g p=%5.4f %5.4f %5.4f s=%5.4f dot=%5.4f path=%s\n",
00405                      cdir,push,xyz[0],xyz[1],xyz[2],1.,1.,path);
00406       }
00407    } else {
00408       if (icheck>=nentries) return;
00409       Int_t idebug = TGeoManager::GetVerboseLevel();
00410       TGeoManager::SetVerboseLevel(5);
00411       bug->GetEntry(icheck);
00412       printf("%-9s error push=%g p=%5.4f %5.4f %5.4f s=%5.4f dot=%5.4f path=%s\n",
00413                   cdir,push,xyz[0],xyz[1],xyz[2],1.,1.,path);
00414       fGeoManager->SetCurrentPoint(xyz);
00415       fGeoManager->SetCurrentDirection(dir);
00416       fGeoManager->FindNode();
00417       TGeoNode *next = fGeoManager->FindNextBoundary();
00418       Double_t step = fGeoManager->GetStep();
00419       for (Int_t j=0; j<3; j++) nxyz[j]=xyz[j]+step*(1.+0.1*push)*dir[j];
00420       Bool_t change = !fGeoManager->IsSameLocation(nxyz[0],nxyz[1],nxyz[2]);
00421       printf("step=%g in: %s\n", step, fGeoManager->GetPath());
00422       printf("  -> next = %s push=%g  change=%d\n", next->GetName(),push, (UInt_t)change);
00423       next->GetVolume()->InspectShape();
00424       next->GetVolume()->DrawOnly();
00425       if (next != fGeoManager->GetCurrentNode()) {
00426          Int_t index1 = fGeoManager->GetCurrentVolume()->GetIndex(next);
00427          if (index1>=0) fGeoManager->CdDown(index1);
00428       }
00429       TPolyMarker3D *pm = new TPolyMarker3D();
00430       fGeoManager->MasterToLocal(xyz, lnext);
00431       pm->SetNextPoint(lnext[0], lnext[1], lnext[2]);
00432       pm->SetMarkerStyle(2);
00433       pm->SetMarkerSize(0.2);
00434       pm->SetMarkerColor(kRed);
00435       pm->Draw("SAME");
00436       TPolyMarker3D *pm1 = new TPolyMarker3D();
00437       for (Int_t j=0; j<3; j++) nxyz[j]=xyz[j]+step*dir[j];
00438       fGeoManager->MasterToLocal(nxyz, lnext);
00439       pm1->SetNextPoint(lnext[0], lnext[1], lnext[2]);
00440       pm1->SetMarkerStyle(2);
00441       pm1->SetMarkerSize(0.2);
00442       pm1->SetMarkerColor(kYellow);
00443       pm1->Draw("SAME");
00444       TGeoManager::SetVerboseLevel(idebug);
00445    }   
00446    delete bug;
00447    delete f;
00448 }   
00449 
00450 //______________________________________________________________________________
00451 void TGeoChecker::CheckGeometryFull(Bool_t checkoverlaps, Bool_t checkcrossings, Int_t ntracks, const Double_t *vertex)
00452 {
00453 // Geometry checking. Opional overlap checkings (by sampling and by mesh). Optional
00454 // boundary crossing check + timing per volume.
00455 //
00456 // STAGE 1: extensive overlap checking by sampling per volume. Stdout need to be 
00457 //  checked by user to get report, then TGeoVolume::CheckOverlaps(0.01, "s") can 
00458 //  be called for the suspicious volumes.
00459 // STAGE2 : normal overlap checking using the shapes mesh - fills the list of 
00460 //  overlaps.
00461 // STAGE3 : shooting NRAYS rays from VERTEX and counting the total number of 
00462 //  crossings per volume (rays propagated from boundary to boundary until 
00463 //  geometry exit). Timing computed and results stored in a histo.
00464 // STAGE4 : shooting 1 mil. random rays inside EACH volume and calling 
00465 //  FindNextBoundary() + Safety() for each call. The timing is normalized by the 
00466 //  number of crossings computed at stage 2 and presented as percentage. 
00467 //  One can get a picture on which are the most "burned" volumes during 
00468 //  transportation from geometry point of view. Another plot of the timing per 
00469 //  volume vs. number of daughters is produced. 
00470 // All histos are saved in the file statistics.root
00471    Int_t nuid = fGeoManager->GetListOfUVolumes()->GetEntries();
00472    if (!fTimer) fTimer = new TStopwatch();
00473    Int_t i;
00474    Double_t value;
00475    fFlags = new Bool_t[nuid];
00476    memset(fFlags, 0, nuid*sizeof(Bool_t));
00477    TGeoVolume *vol;
00478    TCanvas *c = new TCanvas("overlaps", "Overlaps by sampling", 800,800);
00479    new TRandom3();
00480 
00481 // STAGE 1: Overlap checking by sampling per volume
00482    if (checkoverlaps) {
00483       printf("====================================================================\n");
00484       printf("STAGE 1: Overlap checking by sampling within 10 microns\n");
00485       printf("====================================================================\n");
00486       fGeoManager->CheckOverlaps(0.001, "s"); 
00487 
00488    // STAGE 2: Global overlap/extrusion checking
00489       printf("====================================================================\n");
00490       printf("STAGE 2: Global overlap/extrusion checking within 10 microns\n");
00491       printf("====================================================================\n");
00492       fGeoManager->CheckOverlaps(0.001);
00493    }   
00494    
00495    if (!checkcrossings) {
00496       delete [] fFlags;
00497       fFlags = 0;
00498       delete c;
00499       return;
00500    }   
00501    
00502    fVal1 = new Double_t[nuid];
00503    fVal2 = new Double_t[nuid];
00504    memset(fVal1, 0, nuid*sizeof(Double_t));
00505    memset(fVal2, 0, nuid*sizeof(Double_t));
00506    // STAGE 3: How many crossings per volume in a realistic case ?
00507    // Ignore volumes with no daughters
00508 
00509    // Generate rays from vertex in phi=[0,2*pi] theta=[0,pi]
00510 //   Int_t ntracks = 1000000;
00511    printf("====================================================================\n");
00512    printf("STAGE 3: Propagating %i tracks starting from vertex\n and conting number of boundary crossings...\n", ntracks);
00513    printf("====================================================================\n");
00514    Int_t nbound = 0;
00515    Double_t theta, phi;
00516    Double_t point[3], dir[3];
00517    memset(point, 0, 3*sizeof(Double_t));
00518    if (vertex) memcpy(point, vertex, 3*sizeof(Double_t));
00519    
00520    fTimer->Reset();
00521    fTimer->Start();
00522    for (i=0; i<ntracks; i++) {
00523       phi = 2.*TMath::Pi()*gRandom->Rndm();
00524       theta= TMath::ACos(1.-2.*gRandom->Rndm());
00525       dir[0]=TMath::Sin(theta)*TMath::Cos(phi);
00526       dir[1]=TMath::Sin(theta)*TMath::Sin(phi);
00527       dir[2]=TMath::Cos(theta);
00528       if ((i%100)==0) OpProgress("Transporting tracks",i, ntracks, fTimer); 
00529 //      if ((i%1000)==0) printf("... remaining tracks %i\n", ntracks-i);
00530       nbound += PropagateInGeom(point,dir);
00531    }
00532    fTimer->Stop();
00533    Double_t time1 = fTimer->CpuTime() *1.E6;
00534    Double_t time2 = time1/ntracks;
00535    Double_t time3 = time1/nbound;
00536    OpProgress("Transporting tracks",ntracks, ntracks, fTimer, kTRUE); 
00537    printf("Time for crossing %i boundaries: %g [ms]\n", nbound, time1);
00538    printf("Time per track for full geometry traversal: %g [ms], per crossing: %g [ms]\n", time2, time3);
00539 
00540 // STAGE 4: How much time per volume:
00541    
00542    printf("====================================================================\n");
00543    printf("STAGE 4: How much navigation time per volume per next+safety call\n");
00544    printf("====================================================================\n");
00545    TGeoIterator next(fGeoManager->GetTopVolume());
00546    TGeoNode*current;
00547    TString path;
00548    vol = fGeoManager->GetTopVolume();
00549    memset(fFlags, 0, nuid*sizeof(Bool_t));
00550    TStopwatch timer;
00551    timer.Start();
00552    i = 0;
00553    char volname[30];
00554    strncpy(volname, vol->GetName(),15); 
00555    volname[15] = '\0';
00556    OpProgress(volname,i++, nuid, &timer); 
00557    Score(vol, 1, TimingPerVolume(vol)); 
00558    while ((current=next())) {
00559       vol = current->GetVolume();
00560       Int_t uid = vol->GetNumber();
00561       if (fFlags[uid]) continue;
00562       fFlags[uid] = kTRUE;
00563       next.GetPath(path);
00564       fGeoManager->cd(path.Data());
00565       strncpy(volname, vol->GetName(),15); 
00566       volname[15] = '\0';
00567       OpProgress(volname,i++, nuid, &timer); 
00568       Score(vol,1,TimingPerVolume(vol));
00569    }   
00570    OpProgress("STAGE 4 completed",i, nuid, &timer, kTRUE); 
00571    // Draw some histos
00572    Double_t time_tot_pertrack = 0.;
00573    TCanvas *c1 = new TCanvas("c2","ncrossings",10,10,900,500);
00574    c1->SetGrid();
00575    c1->SetTopMargin(0.15);
00576    TFile *f = new TFile("statistics.root", "RECREATE");
00577    TH1F *h = new TH1F("h","number of boundary crossings per volume",3,0,3);
00578    h->SetStats(0);
00579    h->SetFillColor(38);
00580    h->SetBit(TH1::kCanRebin);
00581    
00582    memset(fFlags, 0, nuid*sizeof(Bool_t));
00583    for (i=0; i<nuid; i++) {
00584       vol = fGeoManager->GetVolume(i);
00585       if (!vol->GetNdaughters()) continue;
00586       time_tot_pertrack += fVal1[i]*fVal2[i];
00587       h->Fill(vol->GetName(), (Int_t)fVal1[i]);
00588    }
00589    time_tot_pertrack /= ntracks;
00590    h->LabelsDeflate();
00591    h->LabelsOption(">","X");
00592    h->Draw();   
00593 
00594 
00595    TCanvas *c2 = new TCanvas("c3","time spent per volume in navigation",10,10,900,500);
00596    c2->SetGrid();
00597    c2->SetTopMargin(0.15);
00598    TH2F *h2 = new TH2F("h2", "time per FNB call vs. ndaughters", 100, 0,100,100,0,15);
00599    h2->SetStats(0);
00600    h2->SetMarkerStyle(2);
00601    TH1F *h1 = new TH1F("h1","percent of time spent per volume",3,0,3);
00602    h1->SetStats(0);
00603    h1->SetFillColor(38);
00604    h1->SetBit(TH1::kCanRebin);
00605    for (i=0; i<nuid; i++) {
00606       vol = fGeoManager->GetVolume(i);
00607       if (!vol->GetNdaughters()) continue;
00608       value = fVal1[i]*fVal2[i]/ntracks/time_tot_pertrack;
00609       h1->Fill(vol->GetName(), value);
00610       h2->Fill(vol->GetNdaughters(), fVal2[i]);
00611    }     
00612    h1->LabelsDeflate();
00613    h1->LabelsOption(">","X");
00614    h1->Draw();
00615    TCanvas *c3 = new TCanvas("c4","timing vs. ndaughters",10,10,900,500);
00616    c3->SetGrid();
00617    c3->SetTopMargin(0.15);
00618    h2->Draw();   
00619    f->Write();
00620    delete [] fFlags;
00621    fFlags = 0;
00622    delete [] fVal1;
00623    fVal1 = 0;
00624    delete [] fVal2;
00625    fVal2 = 0;
00626    delete fTimer;
00627    fTimer = 0;
00628    delete c;
00629 }
00630 
00631 //______________________________________________________________________________
00632 Int_t TGeoChecker::PropagateInGeom(Double_t *start, Double_t *dir)
00633 {
00634 // Propagate from START along DIR from boundary to boundary until exiting 
00635 // geometry. Fill array of hits.
00636    fGeoManager->InitTrack(start, dir);
00637    TGeoNode *current = 0;
00638    Int_t nzero = 0;
00639    Int_t nhits = 0;
00640    while (!fGeoManager->IsOutside()) {
00641       if (nzero>3) {
00642       // Problems in trying to cross a boundary
00643          printf("Error in trying to cross boundary of %s\n", current->GetName());
00644          return nhits;
00645       }
00646       current = fGeoManager->FindNextBoundaryAndStep(TGeoShape::Big(), kFALSE);
00647       if (!current || fGeoManager->IsOutside()) return nhits;
00648       Double_t step = fGeoManager->GetStep();
00649       if (step<2.*TGeoShape::Tolerance()) {
00650          nzero++;
00651          continue;
00652       } 
00653       else nzero = 0;
00654       // Generate the hit
00655       nhits++;
00656       TGeoVolume *vol = current->GetVolume();
00657       Score(vol,0,1.);
00658       Int_t iup = 1;
00659       TGeoNode *mother = fGeoManager->GetMother(iup++);
00660       while (mother && mother->GetVolume()->IsAssembly()) {
00661          Score(mother->GetVolume(), 0, 1.);
00662          mother = fGeoManager->GetMother(iup++);
00663       }   
00664    }
00665    return nhits;
00666 }      
00667 
00668 //______________________________________________________________________________
00669 void TGeoChecker::Score(TGeoVolume *vol, Int_t ifield, Double_t value)
00670 {
00671 // Score a hit for VOL
00672    Int_t uid = vol->GetNumber();
00673    switch (ifield) {
00674       case 0:
00675          fVal1[uid] += value;
00676          break;
00677       case 1:
00678          fVal2[uid] += value;
00679    }
00680 }   
00681 
00682 //______________________________________________________________________________
00683 void TGeoChecker::SetNmeshPoints(Int_t npoints) {
00684 // Set number of points to be generated on the shape outline when checking for overlaps.
00685    fNmeshPoints = npoints;
00686    if (npoints<1000) {
00687       Error("SetNmeshPoints", "Cannot allow less than 1000 points for checking - set to 1000");
00688       fNmeshPoints = 1000;
00689    }
00690 }   
00691 
00692 //______________________________________________________________________________
00693 Double_t TGeoChecker::TimingPerVolume(TGeoVolume *vol)
00694 {
00695 // Compute timing per "FindNextBoundary" + "Safety" call. Volume must be
00696 // in the current path.
00697    fTimer->Reset();
00698    const TGeoShape *shape = vol->GetShape();
00699    TGeoBBox *box = (TGeoBBox *)shape;
00700    Double_t dx = box->GetDX();
00701    Double_t dy = box->GetDY();
00702    Double_t dz = box->GetDZ();
00703    Double_t ox = (box->GetOrigin())[0];
00704    Double_t oy = (box->GetOrigin())[1];
00705    Double_t oz = (box->GetOrigin())[2];
00706    Double_t point[3], dir[3], lpt[3], ldir[3];
00707    Double_t pstep = 0.;
00708    pstep = TMath::Max(pstep,dz);
00709    Double_t theta, phi;
00710    Int_t idaughter;
00711    fTimer->Start();
00712    Double_t dist;
00713    Bool_t inside;
00714    for (Int_t i=0; i<1000000; i++) {
00715       lpt[0] = ox-dx+2*dx*gRandom->Rndm();
00716       lpt[1] = oy-dy+2*dy*gRandom->Rndm();
00717       lpt[2] = oz-dz+2*dz*gRandom->Rndm();
00718       fGeoManager->GetCurrentMatrix()->LocalToMaster(lpt,point);
00719       fGeoManager->SetCurrentPoint(point[0],point[1],point[2]);
00720       phi = 2*TMath::Pi()*gRandom->Rndm();
00721       theta= TMath::ACos(1.-2.*gRandom->Rndm());
00722       ldir[0]=TMath::Sin(theta)*TMath::Cos(phi);
00723       ldir[1]=TMath::Sin(theta)*TMath::Sin(phi);
00724       ldir[2]=TMath::Cos(theta);
00725       fGeoManager->GetCurrentMatrix()->LocalToMasterVect(ldir,dir);
00726       fGeoManager->SetCurrentDirection(dir);
00727       fGeoManager->SetStep(pstep);
00728       fGeoManager->ResetState();
00729       inside = kTRUE;
00730       dist = TGeoShape::Big();
00731       if (!vol->IsAssembly()) {
00732          inside = vol->Contains(lpt);
00733          if (!inside) {
00734             dist = vol->GetShape()->DistFromOutside(lpt,ldir,3,pstep); 
00735 //            if (dist>=pstep) continue;
00736          } else {   
00737             vol->GetShape()->DistFromInside(lpt,ldir,3,pstep);
00738          }   
00739             
00740          if (!vol->GetNdaughters()) vol->GetShape()->Safety(lpt, inside);
00741       }   
00742       if (vol->GetNdaughters()) {
00743          fGeoManager->Safety();
00744          fGeoManager->FindNextDaughterBoundary(point,dir,idaughter,kFALSE);
00745       }   
00746    }
00747    fTimer->Stop();
00748    Double_t time_per_track = fTimer->CpuTime();
00749    Int_t uid = vol->GetNumber();
00750    Int_t ncrossings = (Int_t)fVal1[uid];
00751    if (!vol->GetNdaughters())
00752       printf("Time for volume %s (shape=%s): %g [ms] ndaughters=%d ncross=%d\n", vol->GetName(), vol->GetShape()->GetName(), time_per_track, vol->GetNdaughters(), ncrossings);
00753    else
00754       printf("Time for volume %s (assemb=%d): %g [ms] ndaughters=%d ncross=%d\n", vol->GetName(), vol->IsAssembly(), time_per_track, vol->GetNdaughters(), ncrossings);
00755    return time_per_track;
00756 }
00757 
00758 //______________________________________________________________________________
00759 void TGeoChecker::CheckGeometry(Int_t nrays, Double_t startx, Double_t starty, Double_t startz) const
00760 {
00761 // Shoot nrays with random directions from starting point (startx, starty, startz)
00762 // in the reference frame of this volume. Track each ray until exiting geometry, then
00763 // shoot backwards from exiting point and compare boundary crossing points.
00764    Int_t i, j;
00765    Double_t start[3], end[3];
00766    Double_t dir[3];
00767    Double_t dummy[3];
00768    Double_t eps = 0.;
00769    Double_t *array1 = new Double_t[3*1000];
00770    Double_t *array2 = new Double_t[3*1000];
00771    TObjArray *pma = new TObjArray();
00772    TPolyMarker3D *pm;
00773    pm = new TPolyMarker3D();
00774    pm->SetMarkerColor(2); // error > eps
00775    pm->SetMarkerStyle(8);
00776    pm->SetMarkerSize(0.4);
00777    pma->AddAt(pm, 0);
00778    pm = new TPolyMarker3D();
00779    pm->SetMarkerColor(4); // point not found
00780    pm->SetMarkerStyle(8);
00781    pm->SetMarkerSize(0.4);
00782    pma->AddAt(pm, 1);
00783    pm = new TPolyMarker3D();
00784    pm->SetMarkerColor(6); // extra point back
00785    pm->SetMarkerStyle(8);
00786    pm->SetMarkerSize(0.4);
00787    pma->AddAt(pm, 2);
00788    Int_t nelem1, nelem2;
00789    Int_t dim1=1000, dim2=1000;
00790    if ((startx==0) && (starty==0) && (startz==0)) eps=1E-3;
00791    start[0] = startx+eps;
00792    start[1] = starty+eps;
00793    start[2] = startz+eps;
00794    Int_t n10=nrays/10;
00795    Double_t theta,phi;
00796    Double_t dw, dwmin, dx, dy, dz;
00797    Int_t ist1, ist2, ifound;
00798    for (i=0; i<nrays; i++) {
00799       if (n10) {
00800          if ((i%n10) == 0) printf("%i percent\n", Int_t(100*i/nrays));
00801       }
00802       phi = 2*TMath::Pi()*gRandom->Rndm();
00803       theta= TMath::ACos(1.-2.*gRandom->Rndm());
00804       dir[0]=TMath::Sin(theta)*TMath::Cos(phi);
00805       dir[1]=TMath::Sin(theta)*TMath::Sin(phi);
00806       dir[2]=TMath::Cos(theta);
00807       // shoot direct ray
00808       nelem1=nelem2=0;
00809 //      printf("DIRECT %i\n", i);
00810       array1 = ShootRay(&start[0], dir[0], dir[1], dir[2], array1, nelem1, dim1);
00811       if (!nelem1) continue;
00812 //      for (j=0; j<nelem1; j++) printf("%i : %f %f %f\n", j, array1[3*j], array1[3*j+1], array1[3*j+2]);
00813       memcpy(&end[0], &array1[3*(nelem1-1)], 3*sizeof(Double_t));
00814       // shoot ray backwards
00815 //      printf("BACK %i\n", i);
00816       array2 = ShootRay(&end[0], -dir[0], -dir[1], -dir[2], array2, nelem2, dim2, &start[0]);
00817       if (!nelem2) {
00818          printf("#### NOTHING BACK ###########################\n");
00819          for (j=0; j<nelem1; j++) {
00820             pm = (TPolyMarker3D*)pma->At(0);
00821             pm->SetNextPoint(array1[3*j], array1[3*j+1], array1[3*j+2]);
00822          }
00823          continue;
00824       }             
00825 //      printf("BACKWARDS\n");
00826       Int_t k=nelem2>>1;
00827       for (j=0; j<k; j++) {
00828          memcpy(&dummy[0], &array2[3*j], 3*sizeof(Double_t));
00829          memcpy(&array2[3*j], &array2[3*(nelem2-1-j)], 3*sizeof(Double_t));
00830          memcpy(&array2[3*(nelem2-1-j)], &dummy[0], 3*sizeof(Double_t));
00831       }
00832 //      for (j=0; j<nelem2; j++) printf("%i : %f ,%f ,%f   \n", j, array2[3*j], array2[3*j+1], array2[3*j+2]);         
00833       if (nelem1!=nelem2) printf("### DIFFERENT SIZES : nelem1=%i nelem2=%i ##########\n", nelem1, nelem2);
00834       ist1 = ist2 = 0;
00835       // check first match
00836       
00837       dx = array1[3*ist1]-array2[3*ist2];
00838       dy = array1[3*ist1+1]-array2[3*ist2+1];
00839       dz = array1[3*ist1+2]-array2[3*ist2+2];
00840       dw = dx*dir[0]+dy*dir[1]+dz*dir[2];
00841       fGeoManager->SetCurrentPoint(&array1[3*ist1]);
00842       fGeoManager->FindNode();
00843 //      printf("%i : %s (%g, %g, %g)\n", ist1, fGeoManager->GetPath(), 
00844 //             array1[3*ist1], array1[3*ist1+1], array1[3*ist1+2]);
00845       if (TMath::Abs(dw)<1E-4) {
00846 //         printf("   matching %i (%g, %g, %g)\n", ist2, array2[3*ist2], array2[3*ist2+1], array2[3*ist2+2]);
00847          ist2++;
00848       } else {
00849          printf("### NOT MATCHING %i f:(%f, %f, %f) b:(%f %f %f) DCLOSE=%f\n", ist2, array1[3*ist1], array1[3*ist1+1], array1[3*ist1+2], array2[3*ist2], array2[3*ist2+1], array2[3*ist2+2],dw);
00850          pm = (TPolyMarker3D*)pma->At(0);
00851          pm->SetNextPoint(array2[3*ist2], array2[3*ist2+1], array2[3*ist2+2]);
00852          if (dw<0) {
00853             // first boundary missed on way back
00854          } else {
00855             // first boundary different on way back
00856             ist2++;
00857          }
00858       }      
00859       
00860       while ((ist1<nelem1-1) && (ist2<nelem2)) {
00861          fGeoManager->SetCurrentPoint(&array1[3*ist1+3]);
00862          fGeoManager->FindNode();
00863 //         printf("%i : %s (%g, %g, %g)\n", ist1+1, fGeoManager->GetPath(), 
00864 //                array1[3*ist1+3], array1[3*ist1+4], array1[3*ist1+5]);
00865          
00866          dx = array1[3*ist1+3]-array1[3*ist1];
00867          dy = array1[3*ist1+4]-array1[3*ist1+1];
00868          dz = array1[3*ist1+5]-array1[3*ist1+2];
00869          // distance to next point
00870          dwmin = dx+dir[0]+dy*dir[1]+dz*dir[2];
00871          while (ist2<nelem2) {
00872             ifound = 0;
00873             dx = array2[3*ist2]-array1[3*ist1];
00874             dy = array2[3*ist2+1]-array1[3*ist1+1];
00875             dz = array2[3*ist2+2]-array1[3*ist1+2];
00876             dw = dx+dir[0]+dy*dir[1]+dz*dir[2];
00877             if (TMath::Abs(dw-dwmin)<1E-4) {
00878                ist1++;
00879                ist2++;
00880                break;
00881             }   
00882             if (dw<dwmin) {
00883             // point found on way back. Check if close enough to ist1+1
00884                ifound++;
00885                dw = dwmin-dw;
00886                if (dw<1E-4) {
00887                // point is matching ist1+1
00888 //                  printf("   matching %i (%g, %g, %g) DCLOSE=%g\n", ist2, array2[3*ist2], array2[3*ist2+1], array2[3*ist2+2], dw);
00889                   ist2++;
00890                   ist1++;
00891                   break;
00892                } else {
00893                // extra boundary found on way back   
00894                   fGeoManager->SetCurrentPoint(&array2[3*ist2]);
00895                   fGeoManager->FindNode();
00896                   pm = (TPolyMarker3D*)pma->At(2);
00897                   pm->SetNextPoint(array2[3*ist2], array2[3*ist2+1], array2[3*ist2+2]);
00898                   printf("### EXTRA BOUNDARY %i :  %s found at DCLOSE=%f\n", ist2, fGeoManager->GetPath(), dw);
00899                   ist2++;
00900                   continue;
00901                }
00902             } else {
00903                if (!ifound) {
00904                   // point ist1+1 not found on way back
00905                   fGeoManager->SetCurrentPoint(&array1[3*ist1+3]);
00906                   fGeoManager->FindNode();
00907                   pm = (TPolyMarker3D*)pma->At(1);
00908                   pm->SetNextPoint(array2[3*ist1+3], array2[3*ist1+4], array2[3*ist1+5]);
00909                   printf("### BOUNDARY MISSED BACK #########################\n");
00910                   ist1++;
00911                   break;
00912                } else {
00913                   ist1++;
00914                   break;
00915                }
00916             }
00917          }               
00918       }                          
00919    }   
00920    pm = (TPolyMarker3D*)pma->At(0);
00921    pm->Draw("SAME");
00922    pm = (TPolyMarker3D*)pma->At(1);
00923    pm->Draw("SAME");
00924    pm = (TPolyMarker3D*)pma->At(2);
00925    pm->Draw("SAME");
00926    if (gPad) {
00927       gPad->Modified();
00928       gPad->Update();
00929    }      
00930    delete [] array1;
00931    delete [] array2;
00932 }
00933 
00934 //______________________________________________________________________________
00935 void TGeoChecker::CleanPoints(Double_t *points, Int_t &numPoints) const
00936 {
00937 // Clean-up the mesh of pcon/pgon from useless points
00938    Int_t ipoint = 0;
00939    Int_t j, k=0;
00940    Double_t rsq;
00941    for (Int_t i=0; i<numPoints; i++) {
00942       j = 3*i;
00943       rsq = points[j]*points[j]+points[j+1]*points[j+1];
00944       if (rsq < 1.e-10) continue;
00945       points[k] = points[j];
00946       points[k+1] = points[j+1];
00947       points[k+2] = points[j+2];
00948       ipoint++;
00949       k = 3*ipoint;
00950    }
00951    numPoints = ipoint;
00952 }      
00953 
00954 //______________________________________________________________________________
00955 TGeoOverlap *TGeoChecker::MakeCheckOverlap(const char *name, TGeoVolume *vol1, TGeoVolume *vol2, TGeoMatrix *mat1, TGeoMatrix *mat2, Bool_t isovlp, Double_t ovlp)
00956 {
00957 // Check if the 2 non-assembly volume candidates overlap/extrude. Returns overlap object.
00958    TGeoOverlap *nodeovlp = 0;
00959    Int_t numPoints1 = fBuff1->NbPnts();
00960    Int_t numSegs1   = fBuff1->NbSegs();
00961    Int_t numPols1   = fBuff1->NbPols();
00962    Int_t numPoints2 = fBuff2->NbPnts();
00963    Int_t numSegs2   = fBuff2->NbSegs();
00964    Int_t numPols2   = fBuff2->NbPols();
00965    Int_t ip;
00966    Bool_t extrude, isextrusion, isoverlapping;
00967    Double_t *points1 = fBuff1->fPnts; 
00968    Double_t *points2 = fBuff2->fPnts;
00969    Double_t local[3], local1[3];
00970    Double_t point[3];
00971    Double_t safety = TGeoShape::Big();
00972    Double_t tolerance = TGeoShape::Tolerance();
00973    if (vol1->IsAssembly() || vol2->IsAssembly()) return nodeovlp;
00974    TGeoShape *shape1 = vol1->GetShape();
00975    TGeoShape *shape2 = vol2->GetShape();
00976    OpProgress("refresh", 0,0,NULL,kFALSE,kTRUE);   
00977    shape1->GetMeshNumbers(numPoints1, numSegs1, numPols1);
00978    if (fBuff1->fID != (TObject*)shape1) {
00979       // Fill first buffer.
00980       fBuff1->SetRawSizes(TMath::Max(numPoints1,fNmeshPoints), 3*TMath::Max(numPoints1,fNmeshPoints), 0, 0, 0, 0);
00981       points1 = fBuff1->fPnts;
00982       if (shape1->GetPointsOnSegments(fNmeshPoints, points1)) {
00983          numPoints1 = fNmeshPoints;
00984       } else {   
00985          shape1->SetPoints(points1);
00986       }   
00987 //      if (shape1->InheritsFrom(TGeoPcon::Class())) {
00988 //         CleanPoints(points1, numPoints1);
00989 //         fBuff1->SetRawSizes(numPoints1, 3*numPoints1,0, 0, 0, 0);
00990 //      }   
00991       fBuff1->fID = shape1;
00992    }   
00993    shape2->GetMeshNumbers(numPoints2, numSegs2, numPols2);
00994    if (fBuff2->fID != (TObject*)shape2) {
00995       // Fill second buffer.
00996       fBuff2->SetRawSizes(TMath::Max(numPoints2,fNmeshPoints), 3*TMath::Max(numPoints2,fNmeshPoints), 0, 0, 0,0);
00997       points2 = fBuff2->fPnts;
00998       if (shape2->GetPointsOnSegments(fNmeshPoints, points2)) {
00999          numPoints2 = fNmeshPoints;
01000       } else {   
01001          shape2->SetPoints(points2);
01002       }   
01003 //      if (shape2->InheritsFrom(TGeoPcon::Class())) {
01004 //         CleanPoints(points2, numPoints2);
01005 //         fBuff1->SetRawSizes(numPoints2, 3*numPoints2,0,0,0,0);
01006 //      }   
01007       fBuff2->fID = shape2;
01008    }   
01009       
01010    if (!isovlp) {
01011    // Extrusion case. Test vol2 extrude vol1.
01012       isextrusion=kFALSE;
01013       // loop all points of the daughter
01014       for (ip=0; ip<numPoints2; ip++) {
01015          memcpy(local, &points2[3*ip], 3*sizeof(Double_t));
01016          if (TMath::Abs(local[0])<tolerance && TMath::Abs(local[1])<tolerance) continue;
01017          mat2->LocalToMaster(local, point);
01018          mat1->MasterToLocal(point, local);
01019          extrude = !shape1->Contains(local);
01020          if (extrude) {
01021             safety = shape1->Safety(local, kFALSE);
01022             if (safety<ovlp) extrude=kFALSE;
01023          }    
01024          if (extrude) {
01025             if (!isextrusion) {
01026                isextrusion = kTRUE;
01027                nodeovlp = new TGeoOverlap(name, vol1, vol2, mat1,mat2,kFALSE,safety);
01028                nodeovlp->SetNextPoint(point[0],point[1],point[2]);
01029                fGeoManager->AddOverlap(nodeovlp);
01030             } else {
01031                if (safety>nodeovlp->GetOverlap()) nodeovlp->SetOverlap(safety);
01032                nodeovlp->SetNextPoint(point[0],point[1],point[2]);
01033            }   
01034          }
01035       }
01036       // loop all points of the mother
01037       for (ip=0; ip<numPoints1; ip++) {
01038          memcpy(local, &points1[3*ip], 3*sizeof(Double_t));
01039          if (local[0]<1e-10 && local[1]<1e-10) continue;
01040          mat1->LocalToMaster(local, point);
01041          mat2->MasterToLocal(point, local1);
01042          extrude = shape2->Contains(local1);
01043          if (extrude) {
01044             // skip points on mother mesh that have no neghbourhood ouside mother
01045             safety = shape1->Safety(local,kTRUE);
01046             if (safety>1E-6) {
01047                extrude = kFALSE;
01048             } else {   
01049                safety = shape2->Safety(local1,kTRUE);
01050                if (safety<ovlp) extrude=kFALSE;
01051             }   
01052          }   
01053          if (extrude) {
01054             if (!isextrusion) {
01055                isextrusion = kTRUE;
01056                nodeovlp = new TGeoOverlap(name, vol1,vol2,mat1,mat2,kFALSE,safety);
01057                nodeovlp->SetNextPoint(point[0],point[1],point[2]);
01058                fGeoManager->AddOverlap(nodeovlp);
01059             } else {
01060                if (safety>nodeovlp->GetOverlap()) nodeovlp->SetOverlap(safety);
01061                nodeovlp->SetNextPoint(point[0],point[1],point[2]);
01062             }   
01063          }
01064       }
01065       return nodeovlp;
01066    }            
01067    // Check overlap
01068    Bool_t overlap;
01069    overlap = kFALSE;
01070    isoverlapping = kFALSE;
01071    // loop all points of first candidate
01072    for (ip=0; ip<numPoints1; ip++) {
01073       memcpy(local, &points1[3*ip], 3*sizeof(Double_t));
01074       if (local[0]<1e-10 && local[1]<1e-10) continue;
01075       mat1->LocalToMaster(local, point);
01076       mat2->MasterToLocal(point, local); // now point in local reference of second
01077       overlap = shape2->Contains(local);
01078       if (overlap) {
01079          safety = shape2->Safety(local, kTRUE);
01080          if (safety<ovlp) overlap=kFALSE;
01081       }    
01082       if (overlap) {
01083          if (!isoverlapping) {
01084             isoverlapping = kTRUE;
01085             nodeovlp = new TGeoOverlap(name,vol1,vol2,mat1,mat2,kTRUE,safety);
01086             nodeovlp->SetNextPoint(point[0],point[1],point[2]);
01087             fGeoManager->AddOverlap(nodeovlp);
01088          } else {
01089             if (safety>nodeovlp->GetOverlap()) nodeovlp->SetOverlap(safety);
01090             nodeovlp->SetNextPoint(point[0],point[1],point[2]);
01091          }     
01092       }
01093    }
01094    // loop all points of second candidate
01095    for (ip=0; ip<numPoints2; ip++) {
01096       memcpy(local, &points2[3*ip], 3*sizeof(Double_t));
01097       if (local[0]<1e-10 && local[1]<1e-10) continue;
01098       mat2->LocalToMaster(local, point);
01099       mat1->MasterToLocal(point, local); // now point in local reference of first
01100       overlap = shape1->Contains(local);
01101       if (overlap) {
01102          safety = shape1->Safety(local, kTRUE);
01103          if (safety<ovlp) overlap=kFALSE;
01104       }    
01105       if (overlap) {
01106          if (!isoverlapping) {
01107             isoverlapping = kTRUE;
01108             nodeovlp = new TGeoOverlap(name,vol1,vol2,mat1,mat2,kTRUE,safety);
01109             nodeovlp->SetNextPoint(point[0],point[1],point[2]);
01110             fGeoManager->AddOverlap(nodeovlp);
01111          } else {
01112             if (safety>nodeovlp->GetOverlap()) nodeovlp->SetOverlap(safety);
01113             nodeovlp->SetNextPoint(point[0],point[1],point[2]);
01114          }     
01115       }
01116    }
01117    return nodeovlp;  
01118 }
01119 
01120 //______________________________________________________________________________
01121 void TGeoChecker::CheckOverlapsBySampling(TGeoVolume *vol, Double_t /* ovlp */, Int_t npoints) const
01122 {
01123 // Check illegal overlaps for volume VOL within a limit OVLP by sampling npoints
01124 // inside the volume shape.
01125    Int_t nd = vol->GetNdaughters();
01126    if (nd<2) return;
01127    TGeoVoxelFinder *voxels = vol->GetVoxels();
01128    if (!voxels) return;
01129    if (voxels->NeedRebuild()) {
01130       voxels->Voxelize();
01131       vol->FindOverlaps();
01132    }   
01133    TGeoBBox *box = (TGeoBBox*)vol->GetShape();
01134    TGeoShape *shape;
01135    TGeoNode *node;
01136    Double_t dx = box->GetDX();
01137    Double_t dy = box->GetDY();
01138    Double_t dz = box->GetDZ();
01139    Double_t pt[3];
01140    Double_t local[3];
01141    Int_t *check_list = 0;
01142    Int_t ncheck = 0;
01143    const Double_t *orig = box->GetOrigin();
01144    Int_t ipoint = 0;
01145    Int_t itry = 0;
01146    Int_t iovlp = 0;
01147    Int_t id=0, id0=0, id1=0;
01148    Bool_t in, incrt;
01149    Double_t safe;
01150    TString name1 = "";
01151    TString name2 = "";
01152    TGeoOverlap **flags = 0;
01153    TGeoNode *node1, *node2;
01154    Int_t novlps = 0;
01155    TGeoHMatrix mat1, mat2;
01156    if (!gRandom) new TRandom3();
01157    while (ipoint < npoints) {
01158    // Shoot randomly in the bounding box.
01159       pt[0] = orig[0] - dx + 2.*dx*gRandom->Rndm();
01160       pt[1] = orig[1] - dy + 2.*dy*gRandom->Rndm();
01161       pt[2] = orig[2] - dz + 2.*dz*gRandom->Rndm();
01162       if (!vol->Contains(pt)) {
01163          itry++;
01164          if (itry>10000 && !ipoint) {
01165             Error("CheckOverlapsBySampling", "No point inside volume!!! - aborting");
01166             break;
01167          }   
01168          continue;
01169       }   
01170       // Check if the point is inside one or more daughters
01171       in = kFALSE;
01172       ipoint++;
01173       check_list = voxels->GetCheckList(pt, ncheck);
01174       if (!check_list || ncheck<2) continue;
01175       for (id=0; id<ncheck; id++) {
01176          id0 = check_list[id];
01177          node  = vol->GetNode(id0);
01178          // Ignore MANY's
01179          if (node->IsOverlapping()) continue;
01180          node->GetMatrix()->MasterToLocal(pt,local);
01181          shape = node->GetVolume()->GetShape();
01182          incrt = shape->Contains(local);
01183          if (!incrt) continue;
01184          if (!in) {
01185             in = kTRUE;
01186             id1 = id0;
01187             continue;
01188          }
01189          // The point is inside 2 or more daughters, check safety
01190          safe = shape->Safety(local, kTRUE);
01191 //         if (safe < ovlp) continue;
01192          // We really have found an overlap -> store the point in a container
01193          iovlp++;
01194          if (!novlps) {
01195             flags = new TGeoOverlap*[nd*nd];
01196             memset(flags, 0, nd*nd*sizeof(TGeoOverlap*));
01197          }
01198          TGeoOverlap *nodeovlp = flags[nd*id1+id0];
01199          if (!nodeovlp) {
01200             novlps++;
01201             node1 = vol->GetNode(id1);
01202             name1 = node1->GetName();
01203             mat1 = node1->GetMatrix();
01204             Int_t cindex = node1->GetVolume()->GetCurrentNodeIndex();
01205             while (cindex >= 0) {
01206                node1 = node1->GetVolume()->GetNode(cindex);
01207                name1 += TString::Format("/%s", node1->GetName());
01208                mat1.Multiply(node1->GetMatrix());
01209                cindex = node1->GetVolume()->GetCurrentNodeIndex();
01210             }   
01211             node2 = vol->GetNode(id0);
01212             name2 = node2->GetName();
01213             mat2 = node2->GetMatrix();
01214             cindex = node2->GetVolume()->GetCurrentNodeIndex();
01215             while (cindex >= 0) {
01216                node2 = node2->GetVolume()->GetNode(cindex);
01217                name2 += TString::Format("/%s", node2->GetName());
01218                mat2.Multiply(node2->GetMatrix());
01219                cindex = node2->GetVolume()->GetCurrentNodeIndex();
01220             }   
01221             nodeovlp = new TGeoOverlap(TString::Format("Volume %s: node %s overlapping %s", 
01222                vol->GetName(), name1.Data(), name2.Data()), node1->GetVolume(),node2->GetVolume(),
01223                &mat1,&mat2, kTRUE, safe);
01224             flags[nd*id1+id0] = nodeovlp;
01225             fGeoManager->AddOverlap(nodeovlp);
01226          } 
01227          // Max 100 points per marker
01228          if (nodeovlp->GetPolyMarker()->GetN()<100) nodeovlp->SetNextPoint(pt[0],pt[1],pt[2]);
01229          if (nodeovlp->GetOverlap()<safe) nodeovlp->SetOverlap(safe);
01230       }
01231    }
01232 
01233    if (flags) delete [] flags;
01234    if (!novlps) return;
01235    Double_t capacity = vol->GetShape()->Capacity();
01236    capacity *= Double_t(iovlp)/Double_t(npoints);
01237    Double_t err = 1./TMath::Sqrt(Double_t(iovlp));
01238    Info("CheckOverlapsBySampling", "#Found %d overlaps adding-up to %g +/- %g [cm3] for daughters of %s",
01239          novlps, capacity, err*capacity, vol->GetName());
01240 }   
01241 
01242 //______________________________________________________________________________
01243 Int_t TGeoChecker::NChecksPerVolume(TGeoVolume *vol)
01244 {
01245 // Compute number of overlaps combinations to check per volume
01246    if (vol->GetFinder()) return 0;
01247    UInt_t nd = vol->GetNdaughters();
01248    if (!nd) return 0;
01249    Bool_t is_assembly = vol->IsAssembly();
01250    TGeoIterator next1(vol);
01251    TGeoIterator next2(vol);
01252    Int_t nchecks = 0;
01253    TGeoNode *node;
01254    UInt_t id;
01255    if (!is_assembly) {      
01256       while ((node=next1())) {
01257          if (node->IsOverlapping()) {
01258             next1.Skip();
01259             continue;
01260          }   
01261          if (!node->GetVolume()->IsAssembly()) {
01262             nchecks++;
01263             next1.Skip();
01264          }
01265       }            
01266    }
01267    // now check if the daughters overlap with each other
01268    if (nd<2) return nchecks;
01269    TGeoVoxelFinder *vox = vol->GetVoxels();
01270    if (!vox) return nchecks;
01271 
01272 
01273    TGeoNode *node1, *node01, *node02;
01274    Int_t novlp;
01275    Int_t *ovlps;
01276    Int_t ko;
01277    UInt_t io;
01278    for (id=0; id<nd; id++) {  
01279       node01 = vol->GetNode(id);
01280       if (node01->IsOverlapping()) continue;
01281       vox->FindOverlaps(id);
01282       ovlps = node01->GetOverlaps(novlp);
01283       if (!ovlps) continue;
01284       for (ko=0; ko<novlp; ko++) { // loop all possible overlapping candidates
01285          io = ovlps[ko];           // (node1, shaped, matrix1, points, fBuff1)
01286          if (io<=id) continue;
01287          node02 = vol->GetNode(io);
01288          if (node02->IsOverlapping()) continue;        
01289         
01290          // We have to check node against node1, but they may be assemblies
01291          if (node01->GetVolume()->IsAssembly()) {
01292             next1.Reset(node01->GetVolume());
01293             while ((node=next1())) {
01294                if (!node->GetVolume()->IsAssembly()) {
01295                   if (node02->GetVolume()->IsAssembly()) {
01296                      next2.Reset(node02->GetVolume());
01297                      while ((node1=next2())) {
01298                         if (!node1->GetVolume()->IsAssembly()) {
01299                            nchecks++;
01300                            next2.Skip();
01301                         }
01302                      }
01303                   } else {
01304                      nchecks++;
01305                   }
01306                   next1.Skip();
01307                }
01308             }
01309          } else {
01310             // node not assembly         
01311             if (node02->GetVolume()->IsAssembly()) { 
01312                next2.Reset(node02->GetVolume());
01313                while ((node1=next2())) {
01314                   if (!node1->GetVolume()->IsAssembly()) {
01315                      nchecks++;
01316                      next2.Skip();
01317                   }
01318                }
01319             } else {
01320                // node1 also not assembly
01321                nchecks++;
01322             }
01323          }                         
01324       }                
01325       node01->SetOverlaps(0,0);
01326    }   
01327    return nchecks;
01328 }      
01329 
01330 //______________________________________________________________________________
01331 void TGeoChecker::CheckOverlaps(const TGeoVolume *vol, Double_t ovlp, Option_t *option)
01332 {
01333 // Check illegal overlaps for volume VOL within a limit OVLP.
01334    if (vol->GetFinder()) return;
01335    UInt_t nd = vol->GetNdaughters();
01336    if (!nd) return;
01337    TGeoShape::SetTransform(gGeoIdentity);
01338    fNchecks = NChecksPerVolume((TGeoVolume*)vol);
01339    Bool_t sampling = kFALSE;
01340    TString opt(option);
01341    opt.ToLower();
01342    if (opt.Contains("s")) sampling = kTRUE;
01343    if (opt.Contains("f")) fFullCheck = kTRUE;
01344    else                   fFullCheck = kFALSE;
01345    if (sampling) {
01346       opt = opt.Strip(TString::kLeading, 's');
01347       Int_t npoints = atoi(opt.Data());
01348       if (!npoints) npoints = 1000000;
01349       CheckOverlapsBySampling((TGeoVolume*)vol, ovlp, npoints);
01350       return;
01351    }   
01352    Bool_t is_assembly = vol->IsAssembly();
01353    TGeoIterator next1((TGeoVolume*)vol);
01354    TGeoIterator next2((TGeoVolume*)vol);
01355    TString path;
01356    // first, test if daughters extrude their container
01357    TGeoNode * node, *nodecheck;
01358    TGeoChecker *checker = (TGeoChecker*)this;
01359 
01360 //   TGeoOverlap *nodeovlp = 0;
01361    UInt_t id;
01362    Int_t level;
01363 // Check extrusion only for daughters of a non-assembly volume
01364    if (!is_assembly) {      
01365       while ((node=next1())) {
01366          if (node->IsOverlapping()) {
01367             next1.Skip();
01368             continue;
01369          }   
01370          if (!node->GetVolume()->IsAssembly()) {
01371             if (fSelectedNode) {
01372             // We have to check only overlaps of the selected node (or real daughters if an assembly)
01373                if ((fSelectedNode != node) && (!fSelectedNode->GetVolume()->IsAssembly())) {
01374                   next1.Skip();
01375                   continue;
01376                }
01377                if (node != fSelectedNode) {   
01378                   level = next1.GetLevel();
01379                   while ((nodecheck=next1.GetNode(level--))) {
01380                      if (nodecheck == fSelectedNode) break;
01381                   }
01382                   if (!nodecheck) {   
01383                      next1.Skip();
01384                      continue;
01385                   }   
01386                }   
01387             }
01388             next1.GetPath(path);
01389             checker->MakeCheckOverlap(TString::Format("%s extruded by: %s", vol->GetName(),path.Data()),
01390                                  (TGeoVolume*)vol,node->GetVolume(),gGeoIdentity,(TGeoMatrix*)next1.GetCurrentMatrix(),kFALSE,ovlp);
01391             next1.Skip();
01392          }
01393       }            
01394    }
01395 
01396    // now check if the daughters overlap with each other
01397    if (nd<2) return;
01398    TGeoVoxelFinder *vox = vol->GetVoxels();
01399    if (!vox) {
01400       Warning("CheckOverlaps", "Volume %s with %i daughters but not voxelized", vol->GetName(),nd);
01401       return;
01402    } 
01403    if (vox->NeedRebuild()) {
01404       vox->Voxelize();
01405       vol->FindOverlaps();
01406    }     
01407    TGeoNode *node1, *node01, *node02;
01408    TGeoHMatrix hmat1, hmat2;
01409    TString path1;
01410    Int_t novlp;
01411    Int_t *ovlps;
01412    Int_t ko;
01413    UInt_t io;
01414    for (id=0; id<nd; id++) {  
01415       node01 = vol->GetNode(id);
01416       if (node01->IsOverlapping()) continue;
01417       vox->FindOverlaps(id);
01418       ovlps = node01->GetOverlaps(novlp);
01419       if (!ovlps) continue;
01420       next1.SetTopName(node01->GetName());
01421       path = node01->GetName();
01422       for (ko=0; ko<novlp; ko++) { // loop all possible overlapping candidates
01423          io = ovlps[ko];           // (node1, shaped, matrix1, points, fBuff1)
01424          if (io<=id) continue;
01425          node02 = vol->GetNode(io);
01426          if (node02->IsOverlapping()) continue;        
01427          // Try to fasten-up things...
01428 //         if (!TGeoBBox::AreOverlapping((TGeoBBox*)node01->GetVolume()->GetShape(), node01->GetMatrix(),
01429 //                                       (TGeoBBox*)node02->GetVolume()->GetShape(), node02->GetMatrix())) continue;
01430          next2.SetTopName(node02->GetName());
01431          path1 = node02->GetName();
01432         
01433          // We have to check node against node1, but they may be assemblies
01434          if (node01->GetVolume()->IsAssembly()) {
01435             next1.Reset(node01->GetVolume());
01436             while ((node=next1())) {
01437                if (!node->GetVolume()->IsAssembly()) {
01438                   next1.GetPath(path);
01439                   hmat1 = node01->GetMatrix();
01440                   hmat1 *= *next1.GetCurrentMatrix();
01441                   if (node02->GetVolume()->IsAssembly()) {
01442                      next2.Reset(node02->GetVolume());
01443                      while ((node1=next2())) {
01444                         if (!node1->GetVolume()->IsAssembly()) {
01445                            if (fSelectedNode) {
01446                            // We have to check only overlaps of the selected node (or real daughters if an assembly)
01447                               if ((fSelectedNode != node) && (fSelectedNode != node1) && (!fSelectedNode->GetVolume()->IsAssembly())) {
01448                                  next2.Skip();
01449                                  continue;
01450                               }   
01451                               if ((fSelectedNode != node1) && (fSelectedNode != node)) {
01452                                  level = next2.GetLevel();
01453                                  while ((nodecheck=next2.GetNode(level--))) {
01454                                     if (nodecheck == fSelectedNode) break;
01455                                  }
01456                                  if (node02 == fSelectedNode) nodecheck = node02;
01457                                  if (!nodecheck) {   
01458                                     level = next1.GetLevel();
01459                                     while ((nodecheck=next1.GetNode(level--))) {
01460                                        if (nodecheck == fSelectedNode) break;
01461                                     }   
01462                                  }
01463                                  if (node01 == fSelectedNode) nodecheck = node01;
01464                                  if (!nodecheck) {   
01465                                     next2.Skip();
01466                                     continue;
01467                                  }   
01468                               }   
01469                            }
01470                            next2.GetPath(path1);
01471                            hmat2 = node02->GetMatrix();
01472                            hmat2 *= *next2.GetCurrentMatrix();
01473                            checker->MakeCheckOverlap(TString::Format("%s/%s overlapping %s/%s", vol->GetName(),path.Data(),vol->GetName(),path1.Data()),
01474                                               node->GetVolume(),node1->GetVolume(),&hmat1,&hmat2,kTRUE,ovlp);  
01475                            next2.Skip();
01476                         }
01477                      }
01478                   } else {
01479                      if (fSelectedNode) {
01480                      // We have to check only overlaps of the selected node (or real daughters if an assembly)
01481                         if ((fSelectedNode != node) && (fSelectedNode != node02) && (!fSelectedNode->GetVolume()->IsAssembly())) {
01482                            next1.Skip();
01483                            continue;
01484                         }   
01485                         if ((fSelectedNode != node) && (fSelectedNode != node02)) {
01486                            level = next1.GetLevel();
01487                            while ((nodecheck=next1.GetNode(level--))) {
01488                               if (nodecheck == fSelectedNode) break;
01489                            }
01490                            if (node01 == fSelectedNode) nodecheck = node01;
01491                            if (!nodecheck) {   
01492                               next1.Skip();
01493                               continue;
01494                            }   
01495                         }   
01496                      }
01497                      checker->MakeCheckOverlap(TString::Format("%s/%s overlapping %s/%s", vol->GetName(),path.Data(),vol->GetName(),path1.Data()),
01498                                         node->GetVolume(),node02->GetVolume(),&hmat1,node02->GetMatrix(),kTRUE,ovlp);  
01499                   }
01500                   next1.Skip();
01501                }
01502             }
01503          } else {
01504             // node not assembly         
01505             if (node02->GetVolume()->IsAssembly()) { 
01506                next2.Reset(node02->GetVolume());
01507                while ((node1=next2())) {
01508                   if (!node1->GetVolume()->IsAssembly()) {
01509                      if (fSelectedNode) {
01510                      // We have to check only overlaps of the selected node (or real daughters if an assembly)
01511                         if ((fSelectedNode != node1) && (fSelectedNode != node01) && (!fSelectedNode->GetVolume()->IsAssembly())) {
01512                            next2.Skip();
01513                            continue;
01514                         }   
01515                         if ((fSelectedNode != node1) && (fSelectedNode != node01)) {
01516                            level = next2.GetLevel();
01517                            while ((nodecheck=next2.GetNode(level--))) {
01518                               if (nodecheck == fSelectedNode) break;
01519                            }
01520                            if (node02 == fSelectedNode) nodecheck = node02;
01521                            if (!nodecheck) {   
01522                               next2.Skip();
01523                               continue;
01524                            }   
01525                         }   
01526                      }
01527                      next2.GetPath(path1);
01528                      hmat2 = node02->GetMatrix();
01529                      hmat2 *= *next2.GetCurrentMatrix();
01530                      checker->MakeCheckOverlap(TString::Format("%s/%s overlapping %s/%s", vol->GetName(),path.Data(),vol->GetName(),path1.Data()),
01531                                         node01->GetVolume(),node1->GetVolume(),node01->GetMatrix(),&hmat2,kTRUE,ovlp);  
01532                      next2.Skip();
01533                   }
01534                }
01535             } else {
01536                // node1 also not assembly
01537                if (fSelectedNode && (fSelectedNode != node01) && (fSelectedNode != node02)) continue;
01538                checker->MakeCheckOverlap(TString::Format("%s/%s overlapping %s/%s", vol->GetName(),path.Data(),vol->GetName(),path1.Data()),
01539                                   node01->GetVolume(),node02->GetVolume(),node01->GetMatrix(),node02->GetMatrix(),kTRUE,ovlp);  
01540             }
01541          }                         
01542       }                
01543       node01->SetOverlaps(0,0);
01544    }
01545 }
01546 
01547 //______________________________________________________________________________
01548 void TGeoChecker::PrintOverlaps() const
01549 {
01550 // Print the current list of overlaps held by the manager class.
01551    TIter next(fGeoManager->GetListOfOverlaps());
01552    TGeoOverlap *ov;
01553    printf("=== Overlaps for %s ===\n", fGeoManager->GetName());
01554    while ((ov=(TGeoOverlap*)next())) ov->PrintInfo();
01555 }
01556 
01557 //______________________________________________________________________________
01558 void TGeoChecker::CheckPoint(Double_t x, Double_t y, Double_t z, Option_t *)
01559 {
01560 //--- Draw point (x,y,z) over the picture of the daughers of the volume containing this point.
01561 //   Generates a report regarding the path to the node containing this point and the distance to
01562 //   the closest boundary.
01563 
01564    Double_t point[3];
01565    Double_t local[3];
01566    point[0] = x;
01567    point[1] = y;
01568    point[2] = z;
01569    TGeoVolume *vol = fGeoManager->GetTopVolume();
01570    if (fVsafe) {
01571       TGeoNode *old = fVsafe->GetNode("SAFETY_1");
01572       if (old) fVsafe->GetNodes()->RemoveAt(vol->GetNdaughters()-1);
01573    }   
01574 //   if (vol != fGeoManager->GetMasterVolume()) fGeoManager->RestoreMasterVolume();
01575    TGeoNode *node = fGeoManager->FindNode(point[0], point[1], point[2]);
01576    fGeoManager->MasterToLocal(point, local);
01577    // get current node
01578    printf("===  Check current point : (%g, %g, %g) ===\n", point[0], point[1], point[2]);
01579    printf("  - path : %s\n", fGeoManager->GetPath());
01580    // get corresponding volume
01581    if (node) vol = node->GetVolume();
01582    // compute safety distance (distance to boundary ignored)
01583    Double_t close = fGeoManager->Safety();
01584    printf("Safety radius : %f\n", close);
01585    if (close>1E-4) {
01586       TGeoVolume *sph = fGeoManager->MakeSphere("SAFETY", vol->GetMedium(), 0, close, 0,180,0,360);
01587       sph->SetLineColor(2);
01588       sph->SetLineStyle(3);
01589       vol->AddNode(sph,1,new TGeoTranslation(local[0], local[1], local[2]));
01590       fVsafe = vol;
01591    }
01592    TPolyMarker3D *pm = new TPolyMarker3D();
01593    pm->SetMarkerColor(2);
01594    pm->SetMarkerStyle(8);
01595    pm->SetMarkerSize(0.5);
01596    pm->SetNextPoint(local[0], local[1], local[2]);
01597    if (vol->GetNdaughters()<2) fGeoManager->SetTopVisible();
01598    else fGeoManager->SetTopVisible(kFALSE);
01599    fGeoManager->SetVisLevel(1);
01600    if (!vol->IsVisible()) vol->SetVisibility(kTRUE);
01601    vol->Draw();
01602    pm->Draw("SAME");
01603    gPad->Modified();
01604    gPad->Update();
01605 }  
01606 
01607 //______________________________________________________________________________
01608 TH2F *TGeoChecker::LegoPlot(Int_t ntheta, Double_t themin, Double_t themax,
01609                             Int_t nphi,   Double_t phimin, Double_t phimax,
01610                             Double_t /*rmin*/, Double_t /*rmax*/, Option_t *option)
01611 {
01612 // Generate a lego plot fot the top volume, according to option.
01613    TH2F *hist = new TH2F("lego", option, nphi, phimin, phimax, ntheta, themin, themax);
01614    
01615    Double_t degrad = TMath::Pi()/180.;
01616    Double_t theta, phi, step, matprop, x;
01617    Double_t start[3];
01618    Double_t dir[3];
01619    TGeoNode *startnode, *endnode;
01620    Int_t i;  // loop index for phi
01621    Int_t j;  // loop index for theta
01622    Int_t ntot = ntheta * nphi;
01623    Int_t n10 = ntot/10;
01624    Int_t igen = 0, iloop=0;
01625    printf("=== Lego plot sph. => nrays=%i\n", ntot);
01626    for (i=1; i<=nphi; i++) {
01627       for (j=1; j<=ntheta; j++) {
01628          igen++;
01629          if (n10) {
01630             if ((igen%n10) == 0) printf("%i percent\n", Int_t(100*igen/ntot));
01631          }  
01632          x = 0;
01633          theta = hist->GetYaxis()->GetBinCenter(j);
01634          phi   = hist->GetXaxis()->GetBinCenter(i)+1E-3;
01635          start[0] = start[1] = start[2] = 1E-3;
01636          dir[0]=TMath::Sin(theta*degrad)*TMath::Cos(phi*degrad);
01637          dir[1]=TMath::Sin(theta*degrad)*TMath::Sin(phi*degrad);
01638          dir[2]=TMath::Cos(theta*degrad);
01639          fGeoManager->InitTrack(&start[0], &dir[0]);
01640          startnode = fGeoManager->GetCurrentNode();
01641          if (fGeoManager->IsOutside()) startnode=0;
01642          if (startnode) {
01643             matprop = startnode->GetVolume()->GetMaterial()->GetRadLen();
01644          } else {
01645             matprop = 0.;
01646          }      
01647          fGeoManager->FindNextBoundary();
01648 //         fGeoManager->IsStepEntering();
01649          // find where we end-up
01650          endnode = fGeoManager->Step();
01651          step = fGeoManager->GetStep();
01652          while (step<1E10) {
01653             // now see if we can make an other step
01654             iloop=0;
01655             while (!fGeoManager->IsEntering()) {
01656                iloop++;
01657                fGeoManager->SetStep(1E-3);
01658                step += 1E-3;
01659                endnode = fGeoManager->Step();
01660             }
01661             if (iloop>1000) printf("%i steps\n", iloop);   
01662             if (matprop>0) {
01663                x += step/matprop;
01664             }   
01665             if (endnode==0 && step>1E10) break;
01666             // generate an extra step to cross boundary
01667             startnode = endnode;    
01668             if (startnode) {
01669                matprop = startnode->GetVolume()->GetMaterial()->GetRadLen();
01670             } else {
01671                matprop = 0.;
01672             }      
01673             
01674             fGeoManager->FindNextBoundary();
01675             endnode = fGeoManager->Step();
01676             step = fGeoManager->GetStep();
01677          }
01678          hist->Fill(phi, theta, x); 
01679       }
01680    }
01681    return hist;           
01682 }
01683 
01684 //______________________________________________________________________________
01685 void TGeoChecker::RandomPoints(TGeoVolume *vol, Int_t npoints, Option_t *option)
01686 {
01687 // Draw random points in the bounding box of a volume.
01688    if (!vol) return;
01689    gRandom = new TRandom3();
01690    vol->VisibleDaughters(kTRUE);
01691    vol->Draw();
01692    TString opt = option;
01693    opt.ToLower();
01694    TObjArray *pm = new TObjArray(128);
01695    TPolyMarker3D *marker = 0;
01696    const TGeoShape *shape = vol->GetShape();
01697    TGeoBBox *box = (TGeoBBox *)shape;
01698    Double_t dx = box->GetDX();
01699    Double_t dy = box->GetDY();
01700    Double_t dz = box->GetDZ();
01701    Double_t ox = (box->GetOrigin())[0];
01702    Double_t oy = (box->GetOrigin())[1];
01703    Double_t oz = (box->GetOrigin())[2];
01704    Double_t *xyz = new Double_t[3];
01705    printf("Random box : %f, %f, %f\n", dx, dy, dz);
01706    TGeoNode *node = 0;
01707    printf("Start... %i points\n", npoints);
01708    Int_t i=0;
01709    Int_t igen=0;
01710    Int_t ic = 0;
01711    Int_t n10 = npoints/10;
01712    Double_t ratio=0;
01713    while (igen<npoints) {
01714       xyz[0] = ox-dx+2*dx*gRandom->Rndm();
01715       xyz[1] = oy-dy+2*dy*gRandom->Rndm();
01716       xyz[2] = oz-dz+2*dz*gRandom->Rndm();
01717       fGeoManager->SetCurrentPoint(xyz);
01718       igen++;
01719       if (n10) {
01720          if ((igen%n10) == 0) printf("%i percent\n", Int_t(100*igen/npoints));
01721       }  
01722       node = fGeoManager->FindNode();
01723       if (!node) continue;
01724       if (!node->IsOnScreen()) continue;
01725       // draw only points in overlapping/non-overlapping volumes
01726       if (opt.Contains("many") && !node->IsOverlapping()) continue;
01727       if (opt.Contains("only") && node->IsOverlapping()) continue;
01728       ic = node->GetColour();
01729       if ((ic<0) || (ic>=128)) ic = 1;
01730       marker = (TPolyMarker3D*)pm->At(ic);
01731       if (!marker) {
01732          marker = new TPolyMarker3D();
01733          marker->SetMarkerColor(ic);
01734 //         marker->SetMarkerStyle(8);
01735 //         marker->SetMarkerSize(0.4);
01736          pm->AddAt(marker, ic);
01737       }
01738       marker->SetNextPoint(xyz[0], xyz[1], xyz[2]);
01739       i++;
01740    }
01741    printf("Number of visible points : %i\n", i);
01742    ratio = (Double_t)i/(Double_t)igen;
01743    printf("efficiency : %g\n", ratio);
01744    for (Int_t m=0; m<128; m++) {
01745       marker = (TPolyMarker3D*)pm->At(m);
01746       if (marker) marker->Draw("SAME");
01747    }
01748    fGeoManager->GetTopVolume()->VisibleDaughters(kFALSE);
01749    printf("---Daughters of %s made invisible.\n", fGeoManager->GetTopVolume()->GetName());
01750    printf("---Make them visible with : gGeoManager->GetTopVolume()->VisibleDaughters();\n");
01751    delete pm;
01752    delete [] xyz;
01753 }   
01754 
01755 //______________________________________________________________________________
01756 void TGeoChecker::RandomRays(Int_t nrays, Double_t startx, Double_t starty, Double_t startz)
01757 {
01758 // Randomly shoot nrays from point (startx,starty,startz) and plot intersections 
01759 // with surfaces for current top node.
01760    TObjArray *pm = new TObjArray(128);
01761    TPolyLine3D *line = 0;
01762    TPolyLine3D *normline = 0;
01763    gRandom = new TRandom3();
01764    TGeoVolume *vol=fGeoManager->GetTopVolume();
01765 //   vol->VisibleDaughters(kTRUE);
01766 
01767    Double_t start[3];
01768    Double_t dir[3];
01769    Int_t istep= 0;
01770    const Double_t *point = fGeoManager->GetCurrentPoint();
01771    vol->Draw();
01772    printf("Start... %i rays\n", nrays);
01773    TGeoNode *startnode, *endnode;
01774    Bool_t vis1,vis2;
01775    Int_t i=0;
01776    Int_t ipoint;
01777    Int_t itot=0;
01778    Int_t n10=nrays/10;
01779    Double_t theta,phi, step, normlen;
01780    const Double_t *normal;
01781    Double_t dx = ((TGeoBBox*)vol->GetShape())->GetDX();
01782    Double_t dy = ((TGeoBBox*)vol->GetShape())->GetDY();
01783    Double_t dz = ((TGeoBBox*)vol->GetShape())->GetDZ();
01784    normlen = TMath::Max(dx,dy);
01785    normlen = TMath::Max(normlen,dz);
01786    normlen *= 0.1;
01787    while (itot<nrays) {
01788       itot++;
01789       ipoint = 0;
01790       if (n10) {
01791          if ((itot%n10) == 0) printf("%i percent\n", Int_t(100*itot/nrays));
01792       }
01793       start[0] = startx;
01794       start[1] = starty;
01795       start[2] = startz;
01796       phi = 2*TMath::Pi()*gRandom->Rndm();
01797       theta= TMath::ACos(1.-2.*gRandom->Rndm());
01798       dir[0]=TMath::Sin(theta)*TMath::Cos(phi);
01799       dir[1]=TMath::Sin(theta)*TMath::Sin(phi);
01800       dir[2]=TMath::Cos(theta);
01801       startnode = fGeoManager->InitTrack(start[0],start[1],start[2], dir[0],dir[1],dir[2]);
01802       line = 0;
01803       if (fGeoManager->IsOutside()) startnode=0;
01804       vis1 = (startnode)?(startnode->IsOnScreen()):kFALSE;
01805       if (vis1) {
01806          line = new TPolyLine3D(2);
01807          line->SetLineColor(startnode->GetVolume()->GetLineColor());
01808          line->SetPoint(ipoint++, startx, starty, startz);
01809          i++;
01810          pm->Add(line);
01811       }
01812       // find where we end-up
01813       fGeoManager->FindNextBoundaryAndStep();
01814       step = fGeoManager->GetStep();
01815       endnode = fGeoManager->GetCurrentNode();
01816       normal = fGeoManager->FindNormalFast();
01817       vis2 = (endnode)?(endnode->IsOnScreen()):kFALSE;
01818       while (endnode) {
01819          istep = 0;
01820          vis2 = endnode->IsOnScreen();
01821          if (ipoint>0) {
01822          // old visible node had an entry point -> finish segment
01823             line->SetPoint(ipoint, point[0], point[1], point[2]);
01824             if (!vis2) {
01825                normline = new TPolyLine3D(2);
01826                normline->SetLineColor(kBlue);
01827                normline->SetLineWidth(2);
01828                normline->SetPoint(0, point[0], point[1], point[2]);
01829                normline->SetPoint(1, point[0]+normal[0]*normlen, 
01830                                      point[1]+normal[1]*normlen, 
01831                                      point[2]+normal[2]*normlen);
01832                pm->Add(normline);
01833             }   
01834             ipoint = 0;
01835             line   = 0;
01836          }
01837          if (vis2) {
01838             // create new segment
01839             line = new TPolyLine3D(2);   
01840             line->SetLineColor(endnode->GetVolume()->GetLineColor());
01841             line->SetPoint(ipoint++, point[0], point[1], point[2]);
01842             i++;
01843             normline = new TPolyLine3D(2);
01844             normline->SetLineColor(kBlue);
01845             normline->SetLineWidth(2);
01846             normline->SetPoint(0, point[0], point[1], point[2]);
01847             normline->SetPoint(1, point[0]+normal[0]*normlen, 
01848                                   point[1]+normal[1]*normlen, 
01849                                   point[2]+normal[2]*normlen);
01850             pm->Add(line);
01851             pm->Add(normline);
01852          } 
01853          // generate an extra step to cross boundary
01854          startnode = endnode;    
01855          fGeoManager->FindNextBoundary();
01856          step = fGeoManager->GetStep();
01857          endnode = fGeoManager->Step();
01858       }      
01859    }   
01860    // draw all segments
01861    for (Int_t m=0; m<pm->GetEntriesFast(); m++) {
01862       line = (TPolyLine3D*)pm->At(m);
01863       if (line) line->Draw("SAME");
01864    }
01865    printf("number of segments : %i\n", i);
01866 //   fGeoManager->GetTopVolume()->VisibleDaughters(kFALSE);
01867 //   printf("---Daughters of %s made invisible.\n", fGeoManager->GetTopVolume()->GetName());
01868 //   printf("---Make them visible with : gGeoManager->GetTopVolume()->VisibleDaughters();\n");
01869    delete pm;
01870 }
01871 
01872 //______________________________________________________________________________
01873 TGeoNode *TGeoChecker::SamplePoints(Int_t npoints, Double_t &dist, Double_t epsil,
01874                                     const char* g3path)
01875 {
01876 // shoot npoints randomly in a box of 1E-5 arround current point.
01877 // return minimum distance to points outside
01878    // make sure that path to current node is updated
01879    // get the response of tgeo
01880    TGeoNode *node = fGeoManager->FindNode();
01881    TGeoNode *nodegeo = 0;
01882    TGeoNode *nodeg3 = 0;
01883    TGeoNode *solg3 = 0;
01884    if (!node) {dist=-1; return 0;}
01885    gRandom = new TRandom3();
01886    Bool_t hasg3 = kFALSE;
01887    if (strlen(g3path)) hasg3 = kTRUE;
01888    TString geopath = fGeoManager->GetPath();
01889    dist = 1E10;
01890    TString common = "";
01891    // cd to common path
01892    Double_t point[3];
01893    Double_t closest[3];
01894    TGeoNode *node1 = 0;
01895    TGeoNode *node_close = 0;
01896    dist = 1E10;
01897    Double_t dist1 = 0;
01898    // initialize size of random box to epsil
01899    Double_t eps[3];
01900    eps[0] = epsil; eps[1]=epsil; eps[2]=epsil;
01901    const Double_t *pointg = fGeoManager->GetCurrentPoint();
01902    if (hasg3) {
01903       TString spath = geopath;
01904       TString name = "";
01905       Int_t index=0;
01906       while (index>=0) {
01907          index = spath.Index("/", index+1);
01908          if (index>0) {
01909             name = spath(0, index);
01910             if (strstr(g3path, name.Data())) {
01911                common = name;
01912                continue;
01913             } else break;
01914          }
01915       }
01916       // if g3 response was given, cd to common path
01917       if (strlen(common.Data())) {
01918          while (strcmp(fGeoManager->GetPath(), common.Data()) && fGeoManager->GetLevel()) {
01919             nodegeo = fGeoManager->GetCurrentNode();
01920             fGeoManager->CdUp();
01921          }
01922          fGeoManager->cd(g3path);
01923          solg3 = fGeoManager->GetCurrentNode();
01924          while (strcmp(fGeoManager->GetPath(), common.Data()) && fGeoManager->GetLevel()) {
01925             nodeg3 = fGeoManager->GetCurrentNode();
01926             fGeoManager->CdUp();
01927          }
01928          if (!nodegeo) return 0;
01929          if (!nodeg3) return 0;
01930          fGeoManager->cd(common.Data());
01931          fGeoManager->MasterToLocal(fGeoManager->GetCurrentPoint(), &point[0]);
01932          Double_t xyz[3], local[3];
01933          for (Int_t i=0; i<npoints; i++) {
01934             xyz[0] = point[0] - eps[0] + 2*eps[0]*gRandom->Rndm();
01935             xyz[1] = point[1] - eps[1] + 2*eps[1]*gRandom->Rndm();
01936             xyz[2] = point[2] - eps[2] + 2*eps[2]*gRandom->Rndm();
01937             nodeg3->MasterToLocal(&xyz[0], &local[0]);
01938             if (!nodeg3->GetVolume()->Contains(&local[0])) continue;
01939             dist1 = TMath::Sqrt((xyz[0]-point[0])*(xyz[0]-point[0])+
01940                    (xyz[1]-point[1])*(xyz[1]-point[1])+(xyz[2]-point[2])*(xyz[2]-point[2]));
01941             if (dist1<dist) {
01942             // save node and closest point
01943                dist = dist1;
01944                node_close = solg3;
01945                // make the random box smaller
01946                eps[0] = TMath::Abs(point[0]-pointg[0]);
01947                eps[1] = TMath::Abs(point[1]-pointg[1]);
01948                eps[2] = TMath::Abs(point[2]-pointg[2]);
01949             }
01950          }
01951       }
01952       if (!node_close) dist = -1;
01953       return node_close;
01954    }
01955 
01956 //   gRandom = new TRandom3();
01957    // save current point
01958    memcpy(&point[0], pointg, 3*sizeof(Double_t));
01959    for (Int_t i=0; i<npoints; i++) {
01960       // generate a random point in MARS
01961       fGeoManager->SetCurrentPoint(point[0] - eps[0] + 2*eps[0]*gRandom->Rndm(),
01962                                    point[1] - eps[1] + 2*eps[1]*gRandom->Rndm(),
01963                                    point[2] - eps[2] + 2*eps[2]*gRandom->Rndm());
01964       // check if new node is different from the old one
01965       if (node1!=node) {
01966          dist1 = TMath::Sqrt((point[0]-pointg[0])*(point[0]-pointg[0])+
01967                  (point[1]-pointg[1])*(point[1]-pointg[1])+(point[2]-pointg[2])*(point[2]-pointg[2]));
01968          if (dist1<dist) {
01969             dist = dist1;
01970             node_close = node1;
01971             memcpy(&closest[0], pointg, 3*sizeof(Double_t));
01972             // make the random box smaller
01973             eps[0] = TMath::Abs(point[0]-pointg[0]);
01974             eps[1] = TMath::Abs(point[1]-pointg[1]);
01975             eps[2] = TMath::Abs(point[2]-pointg[2]);
01976          }
01977       }
01978    }
01979    // restore the original point and path
01980    fGeoManager->FindNode(point[0], point[1], point[2]);  // really needed ?
01981    if (!node_close) dist=-1;
01982    return node_close;
01983 }
01984 
01985 //______________________________________________________________________________
01986 Double_t *TGeoChecker::ShootRay(Double_t *start, Double_t dirx, Double_t diry, Double_t dirz, Double_t *array, Int_t &nelem, Int_t &dim, Double_t *endpoint) const
01987 {
01988 // Shoot one ray from start point with direction (dirx,diry,dirz). Fills input array
01989 // with points just after boundary crossings.
01990 //   Int_t array_dimension = 3*dim;
01991    nelem = 0;
01992    Int_t istep = 0;
01993    if (!dim) {
01994       printf("empty input array\n");
01995       return array;
01996    }   
01997 //   fGeoManager->CdTop();
01998    const Double_t *point = fGeoManager->GetCurrentPoint();
01999    TGeoNode *endnode;
02000    Bool_t is_entering;
02001    Double_t step, forward;
02002    Double_t dir[3];
02003    dir[0] = dirx;
02004    dir[1] = diry;
02005    dir[2] = dirz;
02006    fGeoManager->InitTrack(start, &dir[0]);
02007    fGeoManager->GetCurrentNode();
02008 //   printf("Start : (%f,%f,%f)\n", point[0], point[1], point[2]);
02009    fGeoManager->FindNextBoundary();
02010    step = fGeoManager->GetStep();
02011 //   printf("---next : at step=%f\n", step);
02012    if (step>1E10) return array;
02013    endnode = fGeoManager->Step();
02014    is_entering = fGeoManager->IsEntering();
02015    while (step<1E10) {
02016       if (endpoint) {
02017          forward = dirx*(endpoint[0]-point[0])+diry*(endpoint[1]-point[1])+dirz*(endpoint[2]-point[2]);
02018          if (forward<1E-3) {
02019 //            printf("exit : Passed start point. nelem=%i\n", nelem); 
02020             return array;
02021          }
02022       }
02023       if (is_entering) {
02024          if (nelem>=dim) {
02025             Double_t *temparray = new Double_t[3*(dim+20)];
02026             memcpy(temparray, array, 3*dim*sizeof(Double_t));
02027             delete [] array;
02028             array = temparray;
02029                   dim += 20;
02030          }
02031          memcpy(&array[3*nelem], point, 3*sizeof(Double_t)); 
02032 //         printf("%i (%f, %f, %f) step=%f\n", nelem, point[0], point[1], point[2], step);
02033          nelem++; 
02034       } else {
02035          if (endnode==0 && step>1E10) {
02036 //            printf("exit : NULL endnode. nelem=%i\n", nelem); 
02037             return array;
02038          }    
02039          if (!fGeoManager->IsEntering()) {
02040 //            if (startnode) printf("stepping %f from (%f, %f, %f) inside %s...\n", step,point[0], point[1], point[2], startnode->GetName());
02041 //            else printf("stepping %f from (%f, %f, %f) OUTSIDE...\n", step,point[0], point[1], point[2]);
02042             istep = 0;
02043          }    
02044          while (!fGeoManager->IsEntering()) {
02045             istep++;
02046             if (istep>1E3) {
02047 //               Error("ShootRay", "more than 1000 steps. Step was %f", step);
02048                nelem = 0;
02049                return array;
02050             }   
02051             fGeoManager->SetStep(1E-5);
02052             endnode = fGeoManager->Step();
02053          }
02054          if (istep>0) printf("%i steps\n", istep);   
02055          if (nelem>=dim) {
02056             Double_t *temparray = new Double_t[3*(dim+20)];
02057             memcpy(temparray, array, 3*dim*sizeof(Double_t));
02058             delete [] array;
02059             array = temparray;
02060                   dim += 20;
02061          }
02062          memcpy(&array[3*nelem], point, 3*sizeof(Double_t)); 
02063 //         if (endnode) printf("%i (%f, %f, %f) step=%f\n", nelem, point[0], point[1], point[2], step);
02064          nelem++;   
02065          is_entering = kTRUE;
02066       }
02067       fGeoManager->FindNextBoundary();
02068       step = fGeoManager->GetStep();
02069 //      printf("---next at step=%f\n", step);
02070       endnode = fGeoManager->Step();
02071       is_entering = fGeoManager->IsEntering();
02072    }
02073    return array;
02074 //   printf("exit : INFINITE step. nelem=%i\n", nelem);
02075 }
02076 
02077 //______________________________________________________________________________
02078 void TGeoChecker::Test(Int_t npoints, Option_t *option)
02079 {
02080    // Check time of finding "Where am I" for n points.
02081    gRandom= new TRandom3();
02082    Bool_t recheck = !strcmp(option, "RECHECK");
02083    if (recheck) printf("RECHECK\n");
02084    const TGeoShape *shape = fGeoManager->GetTopVolume()->GetShape();
02085    Double_t dx = ((TGeoBBox*)shape)->GetDX();
02086    Double_t dy = ((TGeoBBox*)shape)->GetDY();
02087    Double_t dz = ((TGeoBBox*)shape)->GetDZ();
02088    Double_t ox = (((TGeoBBox*)shape)->GetOrigin())[0];
02089    Double_t oy = (((TGeoBBox*)shape)->GetOrigin())[1];
02090    Double_t oz = (((TGeoBBox*)shape)->GetOrigin())[2];
02091    Double_t *xyz = new Double_t[3*npoints];
02092    TStopwatch *timer = new TStopwatch();
02093    printf("Random box : %f, %f, %f\n", dx, dy, dz);
02094    timer->Start(kFALSE);
02095    Int_t i;
02096    for (i=0; i<npoints; i++) {
02097       xyz[3*i] = ox-dx+2*dx*gRandom->Rndm();
02098       xyz[3*i+1] = oy-dy+2*dy*gRandom->Rndm();
02099       xyz[3*i+2] = oz-dz+2*dz*gRandom->Rndm();
02100    }
02101    timer->Stop();
02102    printf("Generation time :\n");
02103    timer->Print();
02104    timer->Reset();
02105    TGeoNode *node, *node1;
02106    printf("Start... %i points\n", npoints);
02107    timer->Start(kFALSE);
02108    for (i=0; i<npoints; i++) {
02109       fGeoManager->SetCurrentPoint(xyz+3*i);
02110       if (recheck) fGeoManager->CdTop();
02111       node = fGeoManager->FindNode();
02112       if (recheck) {
02113          node1 = fGeoManager->FindNode();
02114          if (node1 != node) {
02115             printf("Difference for x=%g y=%g z=%g\n", xyz[3*i], xyz[3*i+1], xyz[3*i+2]);
02116             printf(" from top : %s\n", node->GetName());
02117             printf(" redo     : %s\n", fGeoManager->GetPath());
02118          }
02119       }
02120    }
02121    timer->Stop();
02122    timer->Print();
02123    delete [] xyz;
02124    delete timer;
02125 }
02126 
02127 //______________________________________________________________________________
02128 void TGeoChecker::TestOverlaps(const char* path)
02129 {
02130 //--- Geometry overlap checker based on sampling. 
02131    if (fGeoManager->GetTopVolume()!=fGeoManager->GetMasterVolume()) fGeoManager->RestoreMasterVolume();
02132    printf("Checking overlaps for path :\n");
02133    if (!fGeoManager->cd(path)) return;
02134    TGeoNode *checked = fGeoManager->GetCurrentNode();
02135    checked->InspectNode();
02136    // shoot 1E4 points in the shape of the current volume
02137    gRandom= new TRandom3();
02138    Int_t npoints = 1000000;
02139    Double_t big = 1E6;
02140    Double_t xmin = big;
02141    Double_t xmax = -big;
02142    Double_t ymin = big;
02143    Double_t ymax = -big;
02144    Double_t zmin = big;
02145    Double_t zmax = -big;
02146    TObjArray *pm = new TObjArray(128);
02147    TPolyMarker3D *marker = 0;
02148    TPolyMarker3D *markthis = new TPolyMarker3D();
02149    markthis->SetMarkerColor(5);
02150    TNtuple *ntpl = new TNtuple("ntpl","random points","x:y:z");
02151    TGeoShape *shape = fGeoManager->GetCurrentNode()->GetVolume()->GetShape();
02152    Double_t *point = new Double_t[3];
02153    Double_t dx = ((TGeoBBox*)shape)->GetDX();
02154    Double_t dy = ((TGeoBBox*)shape)->GetDY();
02155    Double_t dz = ((TGeoBBox*)shape)->GetDZ();
02156    Double_t ox = (((TGeoBBox*)shape)->GetOrigin())[0];
02157    Double_t oy = (((TGeoBBox*)shape)->GetOrigin())[1];
02158    Double_t oz = (((TGeoBBox*)shape)->GetOrigin())[2];
02159    Double_t *xyz = new Double_t[3*npoints];
02160    Int_t i=0;
02161    printf("Generating %i points inside %s\n", npoints, fGeoManager->GetPath());
02162    while (i<npoints) {
02163       point[0] = ox-dx+2*dx*gRandom->Rndm();
02164       point[1] = oy-dy+2*dy*gRandom->Rndm();
02165       point[2] = oz-dz+2*dz*gRandom->Rndm();
02166       if (!shape->Contains(point)) continue;
02167       // convert each point to MARS
02168 //      printf("local  %9.3f %9.3f %9.3f\n", point[0], point[1], point[2]);
02169       fGeoManager->LocalToMaster(point, &xyz[3*i]);
02170 //      printf("master %9.3f %9.3f %9.3f\n", xyz[3*i], xyz[3*i+1], xyz[3*i+2]);
02171       xmin = TMath::Min(xmin, xyz[3*i]);
02172       xmax = TMath::Max(xmax, xyz[3*i]);
02173       ymin = TMath::Min(ymin, xyz[3*i+1]);
02174       ymax = TMath::Max(ymax, xyz[3*i+1]);
02175       zmin = TMath::Min(zmin, xyz[3*i+2]);
02176       zmax = TMath::Max(zmax, xyz[3*i+2]);
02177       i++;
02178    }
02179    delete [] point;
02180    ntpl->Fill(xmin,ymin,zmin);
02181    ntpl->Fill(xmax,ymin,zmin);
02182    ntpl->Fill(xmin,ymax,zmin);
02183    ntpl->Fill(xmax,ymax,zmin);
02184    ntpl->Fill(xmin,ymin,zmax);
02185    ntpl->Fill(xmax,ymin,zmax);
02186    ntpl->Fill(xmin,ymax,zmax);
02187    ntpl->Fill(xmax,ymax,zmax);
02188    ntpl->Draw("z:y:x");
02189 
02190    // shoot the poins in the geometry
02191    TGeoNode *node;
02192    TString cpath;
02193    Int_t ic=0;
02194    TObjArray *overlaps = new TObjArray();
02195    printf("using FindNode...\n");
02196    for (Int_t j=0; j<npoints; j++) {
02197       // always start from top level (testing only)
02198       fGeoManager->CdTop();
02199       fGeoManager->SetCurrentPoint(&xyz[3*j]);
02200       node = fGeoManager->FindNode();
02201       cpath = fGeoManager->GetPath();
02202       if (cpath.Contains(path)) {
02203          markthis->SetNextPoint(xyz[3*j], xyz[3*j+1], xyz[3*j+2]);
02204          continue;
02205       }
02206       // current point is found in an overlapping node
02207       if (!node) ic=128;
02208       else ic = node->GetVolume()->GetLineColor();
02209       if (ic >= 128) ic = 0;
02210       marker = (TPolyMarker3D*)pm->At(ic);
02211       if (!marker) {
02212          marker = new TPolyMarker3D();
02213          marker->SetMarkerColor(ic);
02214          pm->AddAt(marker, ic);
02215       }
02216       // draw the overlapping point
02217       marker->SetNextPoint(xyz[3*j], xyz[3*j+1], xyz[3*j+2]);
02218       if (node) {
02219          if (overlaps->IndexOf(node) < 0) overlaps->Add(node);
02220       }
02221    }
02222    // draw all overlapping points
02223 //   for (Int_t m=0; m<128; m++) {
02224 //      marker = (TPolyMarker3D*)pm->At(m);
02225 //      if (marker) marker->Draw("SAME");
02226 //   }
02227    markthis->Draw("SAME");
02228    if (gPad) gPad->Update();
02229    // display overlaps
02230    if (overlaps->GetEntriesFast()) {
02231       printf("list of overlapping nodes :\n");
02232       for (i=0; i<overlaps->GetEntriesFast(); i++) {
02233          node = (TGeoNode*)overlaps->At(i);
02234          if (node->IsOverlapping()) printf("%s  MANY\n", node->GetName());
02235          else printf("%s  ONLY\n", node->GetName());
02236       }
02237    } else printf("No overlaps\n");
02238    delete ntpl;
02239    delete pm;
02240    delete [] xyz;
02241    delete overlaps;
02242 }
02243 
02244 //______________________________________________________________________________
02245 Double_t TGeoChecker::Weight(Double_t precision, Option_t *option)
02246 {
02247 // Estimate weight of top level volume with a precision SIGMA(W)/W
02248 // better than PRECISION. Option can be "v" - verbose (default).
02249    TList *matlist = fGeoManager->GetListOfMaterials();
02250    Int_t nmat = matlist->GetSize();
02251    if (!nmat) return 0;
02252    Int_t *nin = new Int_t[nmat];
02253    memset(nin, 0, nmat*sizeof(Int_t));
02254    gRandom = new TRandom3();
02255    TString opt = option;
02256    opt.ToLower();
02257    Bool_t isverbose = opt.Contains("v");
02258    TGeoBBox *box = (TGeoBBox *)fGeoManager->GetTopVolume()->GetShape();
02259    Double_t dx = box->GetDX();
02260    Double_t dy = box->GetDY();
02261    Double_t dz = box->GetDZ();
02262    Double_t ox = (box->GetOrigin())[0];
02263    Double_t oy = (box->GetOrigin())[1];
02264    Double_t oz = (box->GetOrigin())[2];
02265    Double_t x,y,z;
02266    TGeoNode *node;
02267    TGeoMaterial *mat;
02268    Double_t vbox = 0.000008*dx*dy*dz; // m3
02269    Bool_t end = kFALSE;
02270    Double_t weight=0, sigma, eps, dens;
02271    Double_t eps0=1.;
02272    Int_t indmat;
02273    Int_t igen=0;
02274    Int_t iin = 0;
02275    while (!end) {
02276       x = ox-dx+2*dx*gRandom->Rndm();
02277       y = oy-dy+2*dy*gRandom->Rndm();
02278       z = oz-dz+2*dz*gRandom->Rndm();
02279       node = fGeoManager->FindNode(x,y,z);
02280       igen++;
02281       if (!node) continue;
02282       mat = node->GetVolume()->GetMedium()->GetMaterial();
02283       indmat = mat->GetIndex();
02284       if (indmat<0) continue;
02285       nin[indmat]++;
02286       iin++;
02287       if ((iin%100000)==0 || igen>1E8) {
02288          weight = 0;
02289          sigma = 0;
02290          for (indmat=0; indmat<nmat; indmat++) {
02291             mat = (TGeoMaterial*)matlist->At(indmat);
02292             dens = mat->GetDensity(); //  [g/cm3]
02293             if (dens<1E-2) dens=0;
02294             dens *= 1000.;            // [kg/m3]
02295             weight += dens*Double_t(nin[indmat]);
02296             sigma  += dens*dens*nin[indmat];
02297          }
02298                sigma = TMath::Sqrt(sigma);
02299                eps = sigma/weight;
02300                weight *= vbox/Double_t(igen);
02301                sigma *= vbox/Double_t(igen);
02302                if (eps<precision || igen>1E8) {
02303                   if (isverbose) {
02304                      printf("=== Weight of %s : %g +/- %g [kg]\n", 
02305                             fGeoManager->GetTopVolume()->GetName(), weight, sigma);
02306                   }
02307                   end = kTRUE;                      
02308                } else {
02309                   if (isverbose && eps<0.5*eps0) {
02310                      printf("%8dK: %14.7g kg  %g %%\n", 
02311                        igen/1000, weight, eps*100);
02312                eps0 = eps;
02313             } 
02314          }
02315       }
02316    }      
02317    delete [] nin;
02318    return weight;
02319 }
02320 //______________________________________________________________________________
02321 Double_t TGeoChecker::CheckVoxels(TGeoVolume *vol, TGeoVoxelFinder *voxels, Double_t *xyz, Int_t npoints)
02322 {
02323 // count voxel timing
02324    TStopwatch timer;
02325    Double_t time;
02326    TGeoShape *shape = vol->GetShape();
02327    TGeoNode *node;
02328    TGeoMatrix *matrix;
02329    Double_t *point;
02330    Double_t local[3];
02331    Int_t *checklist;
02332    Int_t ncheck;
02333 
02334    timer.Start();
02335    for (Int_t i=0; i<npoints; i++) {
02336       point = xyz + 3*i;
02337       if (!shape->Contains(point)) continue;
02338       checklist = voxels->GetCheckList(point, ncheck);
02339       if (!checklist) continue;
02340       if (!ncheck) continue;
02341       for (Int_t id=0; id<ncheck; id++) {
02342          node = vol->GetNode(checklist[id]);
02343          matrix = node->GetMatrix();
02344          matrix->MasterToLocal(point, &local[0]);
02345          if (node->GetVolume()->GetShape()->Contains(&local[0])) break;
02346       }   
02347    }
02348    time = timer.CpuTime();
02349    return time;
02350 }   
02351 
02352 //______________________________________________________________________________
02353 Bool_t TGeoChecker::TestVoxels(TGeoVolume * /*vol*/, Int_t /*npoints*/)
02354 {
02355 // Returns optimal voxelization type for volume vol.
02356 //   kFALSE - cartesian
02357 //   kTRUE  - cylindrical
02358    return kFALSE;
02359 }

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