G4TessellatedSolid.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration and of QinetiQ Ltd,   *
00020 // * subject to DEFCON 705 IPR conditions.                            *
00021 // * By using,  copying,  modifying or  distributing the software (or *
00022 // * any work based  on the software)  you  agree  to acknowledge its *
00023 // * use  in  resulting  scientific  publications,  and indicate your *
00024 // * acceptance of all terms of the Geant4 Software license.          *
00025 // ********************************************************************
00026 //
00027 // $Id: G4TessellatedSolid.cc 69790 2013-05-15 12:39:10Z gcosmo $
00028 //
00029 // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
00030 //
00031 // CHANGE HISTORY
00032 // --------------
00033 // 12 October 2012,   M Gayer, CERN, complete rewrite reducing memory
00034 //                    requirements more than 50% and speedup by a factor of
00035 //                    tens or more depending on the number of facets, thanks
00036 //                    to voxelization of surface and improvements. 
00037 //                    Speedup factor of thousands for solids with number of
00038 //                    facets in hundreds of thousands.
00039 //
00040 // 22 August 2011,    I Hrivnacova, Orsay, fix in DistanceToOut(p) and
00041 //                    DistanceToIn(p) to exactly compute distance from facet
00042 //                    avoiding use of 'outgoing' flag shortcut variant.
00043 //
00044 // 04 August 2011,    T Nikitina, CERN, added SetReferences() to
00045 //                    CreatePolyhedron() for Visualization of Boolean Operations  
00046 //
00047 // 12 April 2010,     P R Truscott, QinetiQ, bug fixes to treat optical
00048 //                    photon transport, in particular internal reflection
00049 //                    at surface.
00050 //
00051 // 14 November 2007,  P R Truscott, QinetiQ & Stan Seibert, U Texas
00052 //                    Bug fixes to CalculateExtent
00053 //
00054 // 17 September 2007, P R Truscott, QinetiQ Ltd & Richard Holmberg
00055 //                    Updated extensively prior to this date to deal with
00056 //                    concaved tessellated surfaces, based on the algorithm
00057 //                    of Richard Holmberg.  This had been slightly modified
00058 //                    to determine with inside the geometry by projecting
00059 //                    random rays from the point provided.  Now random rays
00060 //                    are predefined rather than making use of random
00061 //                    number generator at run-time.
00062 //
00063 // 22 November 2005, F Lei
00064 //  - Changed ::DescribeYourselfTo(), line 464
00065 //  - added GetPolyHedron()
00066 // 
00067 // 31 October 2004, P R Truscott, QinetiQ Ltd, UK
00068 //  - Created.
00069 //
00070 // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
00071 
00072 #include <iostream>
00073 #include <stack>
00074 #include <iostream>
00075 #include <iomanip>
00076 #include <fstream>
00077 #include <algorithm>
00078 #include <list>
00079 
00080 #include "G4TessellatedSolid.hh"
00081 
00082 #include "geomdefs.hh"
00083 #include "Randomize.hh"
00084 #include "G4SystemOfUnits.hh"
00085 #include "G4PhysicalConstants.hh"
00086 #include "G4GeometryTolerance.hh"
00087 #include "G4VFacet.hh"
00088 #include "G4VoxelLimits.hh"
00089 #include "G4AffineTransform.hh"
00090 
00091 #include "G4PolyhedronArbitrary.hh"
00092 #include "G4VGraphicsScene.hh"
00093 #include "G4VisExtent.hh"
00094 
00095 using namespace std;
00096 
00098 //
00099 // Standard contructor has blank name and defines no fFacets.
00100 //
00101 G4TessellatedSolid::G4TessellatedSolid () : G4VSolid("dummy")
00102 {
00103   Initialize();
00104 }
00105 
00107 //
00108 // Alternative constructor. Simple define name and geometry type - no fFacets
00109 // to detine.
00110 //
00111 G4TessellatedSolid::G4TessellatedSolid (const G4String &name)
00112   : G4VSolid(name)
00113 {
00114   Initialize();
00115 }
00116 
00118 //
00119 // Fake default constructor - sets only member data and allocates memory
00120 //                            for usage restricted to object persistency.
00121 //
00122 G4TessellatedSolid::G4TessellatedSolid( __void__& a) : G4VSolid(a)
00123 {
00124   Initialize();
00125   fMinExtent.set(0,0,0);
00126   fMaxExtent.set(0,0,0);
00127 }
00128 
00129 
00131 G4TessellatedSolid::~G4TessellatedSolid ()
00132 {
00133   DeleteObjects ();
00134 }
00135 
00137 //
00138 // Copy constructor.
00139 //
00140 G4TessellatedSolid::G4TessellatedSolid (const G4TessellatedSolid &ts)
00141   : G4VSolid(ts), fpPolyhedron(0)
00142 {
00143   Initialize();
00144 
00145   CopyObjects(ts);
00146 }
00147 
00149 //
00150 // Assignment operator.
00151 //
00152 G4TessellatedSolid&
00153 G4TessellatedSolid::operator= (const G4TessellatedSolid &ts)
00154 {
00155   if (&ts == this) return *this;
00156 
00157   // Copy base class data
00158   G4VSolid::operator=(ts);
00159 
00160   DeleteObjects ();
00161 
00162   Initialize();
00163 
00164   CopyObjects (ts);
00165 
00166   return *this;
00167 }
00168 
00170 //
00171 void G4TessellatedSolid::Initialize()
00172 {
00173   kCarToleranceHalf = 0.5*kCarTolerance;
00174 
00175   fpPolyhedron = 0; fCubicVolume = 0.; fSurfaceArea = 0.;
00176 
00177   fGeometryType = "G4TessellatedSolid";
00178   fSolidClosed  = false;
00179 
00180   fMinExtent.set(kInfinity,kInfinity,kInfinity);
00181   fMaxExtent.set(-kInfinity,-kInfinity,-kInfinity);
00182 
00183   SetRandomVectors();
00184 }
00185 
00187 //
00188 void G4TessellatedSolid::DeleteObjects ()
00189 {
00190   G4int size = fFacets.size();
00191   for (G4int i = 0; i < size; ++i)  { delete fFacets[i]; }
00192   fFacets.clear();
00193 }
00194 
00196 //
00197 void G4TessellatedSolid::CopyObjects (const G4TessellatedSolid &ts)
00198 {
00199   G4ThreeVector reductionRatio;
00200   G4int fmaxVoxels = fVoxels.GetMaxVoxels(reductionRatio);
00201   if (fmaxVoxels < 0)
00202     fVoxels.SetMaxVoxels(reductionRatio);
00203   else
00204     fVoxels.SetMaxVoxels(fmaxVoxels);
00205 
00206   G4int n = ts.GetNumberOfFacets();
00207   for (G4int i = 0; i < n; ++i)
00208   {
00209     G4VFacet *facetClone = (ts.GetFacet(i))->GetClone();
00210     AddFacet(facetClone);
00211   }
00212   if (ts.GetSolidClosed()) SetSolidClosed(true);
00213 }
00214 
00216 //
00217 // Add a facet to the facet list.
00218 // Note that you can add, but you cannot delete.
00219 //
00220 G4bool G4TessellatedSolid::AddFacet (G4VFacet *aFacet)
00221 {
00222   // Add the facet to the vector.
00223   //
00224   if (fSolidClosed)
00225   {
00226     G4Exception("G4TessellatedSolid::AddFacet()", "GeomSolids1002",
00227                 JustWarning, "Attempt to add facets when solid is closed.");
00228     return false;
00229   }
00230   else if (aFacet->IsDefined())
00231   {
00232     set<G4VertexInfo,G4VertexComparator>::iterator begin
00233       = fFacetList.begin(), end = fFacetList.end(), pos, it;
00234     G4ThreeVector p = aFacet->GetCircumcentre();
00235     G4VertexInfo value;
00236     value.id = fFacetList.size();
00237     value.mag2 = p.x() + p.y() + p.z();
00238 
00239     G4bool found = false;
00240     if (!OutsideOfExtent(p, kCarTolerance))
00241     {
00242       G4double kCarTolerance3 = 3 * kCarTolerance;
00243       pos = fFacetList.lower_bound(value);
00244 
00245       it = pos;
00246       while (!found && it != end)
00247       {
00248         G4int id = (*it).id;
00249         G4VFacet *facet = fFacets[id];
00250         G4ThreeVector q = facet->GetCircumcentre();
00251         if ((found = (facet == aFacet))) break;
00252         G4double dif = q.x() + q.y() + q.z() - value.mag2;
00253         if (dif > kCarTolerance3) break;
00254         it++;
00255       }
00256 
00257       if (fFacets.size() > 1)
00258       {
00259         it = pos;
00260         while (!found && it != begin)
00261         {
00262           --it;
00263           G4int id = (*it).id;
00264           G4VFacet *facet = fFacets[id];          
00265           G4ThreeVector q = facet->GetCircumcentre();
00266           found = (facet == aFacet);
00267           if (found) break;
00268           G4double dif = value.mag2 - (q.x() + q.y() + q.z());
00269           if (dif > kCarTolerance3) break;
00270         }
00271       }
00272     }
00273 
00274     if (!found)
00275     {
00276       fFacets.push_back(aFacet);
00277       fFacetList.insert(value);
00278     }
00279 
00280     return true;
00281   }
00282   else
00283   {
00284     G4Exception("G4TessellatedSolid::AddFacet()", "GeomSolids1002",
00285                 JustWarning, "Attempt to add facet not properly defined.");    
00286     aFacet->StreamInfo(G4cout);
00287     return false;
00288   }
00289 }
00290 
00292 //
00293 G4int G4TessellatedSolid::SetAllUsingStack(const std::vector<G4int> &voxel,
00294                                            const std::vector<G4int> &max,
00295                                            G4bool status, G4SurfBits &checked)
00296 {
00297   vector<G4int> xyz = voxel;
00298   stack<vector<G4int> > pos;
00299   pos.push(xyz);
00300   G4int filled = 0;
00301   G4int cc = 0, nz = 0;
00302 
00303   vector<G4int> candidates;
00304 
00305   while (!pos.empty())
00306   {
00307     xyz = pos.top();
00308     pos.pop();
00309     G4int index = fVoxels.GetVoxelsIndex(xyz);
00310     if (!checked[index])
00311     {
00312       checked.SetBitNumber(index, true);
00313       cc++;
00314 
00315       if (fVoxels.IsEmpty(index))
00316       {
00317         filled++;
00318 
00319         fInsides.SetBitNumber(index, status);
00320 
00321         for (G4int i = 0; i <= 2; ++i)
00322         {
00323           if (xyz[i] < max[i] - 1)
00324           {
00325             xyz[i]++;
00326             pos.push(xyz);
00327             xyz[i]--;
00328           }
00329 
00330           if (xyz[i] > 0)
00331           {
00332             xyz[i]--;
00333             pos.push(xyz);
00334             xyz[i]++;
00335           }
00336         }
00337       }
00338       else
00339       {
00340         nz++;
00341       }
00342     }
00343   }
00344   return filled;
00345 }
00346 
00348 //
00349 void G4TessellatedSolid::PrecalculateInsides()
00350 {
00351   vector<G4int> voxel(3), maxVoxels(3);
00352   for (G4int i = 0; i <= 2; ++i) maxVoxels[i] = fVoxels.GetBoundary(i).size();
00353   G4int size = maxVoxels[0] * maxVoxels[1] * maxVoxels[2];
00354 
00355   G4SurfBits checked(size-1);
00356   fInsides.Clear();
00357   fInsides.ResetBitNumber(size-1);
00358 
00359   G4ThreeVector point;
00360   for (voxel[2] = 0; voxel[2] < maxVoxels[2] - 1; ++voxel[2])
00361   {
00362     for (voxel[1] = 0; voxel[1] < maxVoxels[1] - 1; ++voxel[1])
00363     {
00364       for (voxel[0] = 0; voxel[0] < maxVoxels[0] - 1; ++voxel[0])
00365       {
00366         G4int index = fVoxels.GetVoxelsIndex(voxel);
00367         if (!checked[index] && fVoxels.IsEmpty(index))
00368         {
00369           for (G4int i = 0; i <= 2; ++i) point[i] = fVoxels.GetBoundary(i)[voxel[i]];
00370           G4bool inside = (G4bool) (InsideNoVoxels(point) == kInside);
00371           SetAllUsingStack(voxel, maxVoxels, inside, checked);
00372         }
00373         else checked.SetBitNumber(index);
00374       }
00375     }
00376   }
00377 }
00378 
00380 //
00381 void G4TessellatedSolid::Voxelize ()
00382 {
00383 #ifdef G4SPECSDEBUG
00384   G4cout << "Voxelizing..." << G4endl;
00385 #endif
00386   fVoxels.Voxelize(fFacets);
00387 
00388   if (fVoxels.Empty().GetNbits())
00389   {
00390 #ifdef G4SPECSDEBUG
00391     G4cout << "Precalculating Insides..." << G4endl;
00392 #endif
00393     PrecalculateInsides();
00394   }
00395 }
00396 
00398 //
00399 // Compute extremeFacets, i.e. find those facets that have surface
00400 // planes that bound the volume.
00401 // Note that this is going to reject concaved surfaces as being extreme. Also
00402 // note that if the vertex is on the facet, displacement is zero, so IsInside
00403 // returns true. So will this work??  Need non-equality
00404 // "G4bool inside = displacement < 0.0;"
00405 // or
00406 // "G4bool inside = displacement <= -0.5*kCarTolerance" 
00407 // (Notes from PT 13/08/2007).
00408 //
00409 void G4TessellatedSolid::SetExtremeFacets()
00410 {
00411   G4int size = fFacets.size();
00412   for (G4int j = 0; j < size; ++j)
00413   {
00414     G4VFacet &facet = *fFacets[j];
00415 
00416     G4bool isExtreme = true;
00417     G4int vsize = fVertexList.size();
00418     for (G4int i=0; i < vsize; ++i)
00419     {
00420       if (!facet.IsInside(fVertexList[i]))
00421       {
00422         isExtreme = false;
00423         break;
00424       }
00425     }
00426     if (isExtreme) fExtremeFacets.insert(&facet);
00427   }
00428 }
00429 
00431 //
00432 void G4TessellatedSolid::CreateVertexList()
00433 {
00434   // The algorithm:
00435   // we will have additional vertexListSorted, where all the items will be
00436   // sorted by magnitude of vertice vector.
00437   // New candidate for fVertexList - we will determine the position fo first
00438   // item which would be within its magnitude - 0.5*kCarTolerance.
00439   // We will go trough until we will reach > +0.5 kCarTolerance.
00440   // Comparison (q-p).mag() < 0.5*kCarTolerance will be made.
00441   // They can be just stored in std::vector, with custom insertion based
00442   // on binary search.
00443 
00444   set<G4VertexInfo,G4VertexComparator> vertexListSorted;
00445   set<G4VertexInfo,G4VertexComparator>::iterator begin
00446      = vertexListSorted.begin(), end = vertexListSorted.end(), pos, it;
00447   G4ThreeVector p;
00448   G4VertexInfo value;
00449 
00450   fVertexList.clear();
00451   G4int size = fFacets.size();
00452 
00453   G4double kCarTolerance24 = kCarTolerance * kCarTolerance / 4.0;
00454   G4double kCarTolerance3 = 3 * kCarTolerance;
00455   vector<G4int> newIndex(100);
00456   
00457   for (G4int k = 0; k < size; ++k)
00458   {
00459     G4VFacet &facet = *fFacets[k];
00460     G4int max = facet.GetNumberOfVertices();
00461 
00462     for (G4int i = 0; i < max; ++i)
00463     {
00464       p = facet.GetVertex(i);
00465       value.id = fVertexList.size();
00466       value.mag2 = p.x() + p.y() + p.z();
00467 
00468       G4bool found = false;
00469       G4int id = 0;
00470       if (!OutsideOfExtent(p, kCarTolerance))
00471       {
00472         pos = vertexListSorted.lower_bound(value);
00473         it = pos;
00474         while (it != end)
00475         {
00476           id = (*it).id;
00477           G4ThreeVector q = fVertexList[id];
00478           G4double dif = (q-p).mag2();
00479           found = (dif < kCarTolerance24);
00480           if (found) break;
00481           dif = q.x() + q.y() + q.z() - value.mag2;
00482           if (dif > kCarTolerance3) break;
00483           it++;
00484         }
00485 
00486         if (!found && (fVertexList.size() > 1))
00487         {
00488           it = pos;
00489           while (it != begin)
00490           {
00491             --it;
00492             id = (*it).id;
00493             G4ThreeVector q = fVertexList[id];
00494             G4double dif = (q-p).mag2();
00495             found = (dif < kCarTolerance24);
00496             if (found) break;
00497         dif = value.mag2 - (q.x() + q.y() + q.z());
00498             if (dif > kCarTolerance3) break;
00499           }
00500         }
00501       }
00502 
00503       if (!found)
00504       {
00505 #ifdef G4SPECSDEBUG     
00506         G4cout << p.x() << ":" << p.y() << ":" << p.z() << G4endl;
00507         G4cout << "Adding new vertex #" << i << " of facet " << k
00508                << " id " << value.id << G4endl;
00509         G4cout << "===" << G4endl;      
00510 #endif
00511         fVertexList.push_back(p);
00512         vertexListSorted.insert(value);
00513         begin = vertexListSorted.begin();
00514         end = vertexListSorted.end();
00515         newIndex[i] = value.id;
00516         //
00517         // Now update the maximum x, y and z limits of the volume.
00518         //
00519         if (value.id == 0) fMinExtent = fMaxExtent = p; 
00520         else
00521         {
00522           if (p.x() > fMaxExtent.x()) fMaxExtent.setX(p.x());
00523           else if (p.x() < fMinExtent.x()) fMinExtent.setX(p.x());
00524           if (p.y() > fMaxExtent.y()) fMaxExtent.setY(p.y());
00525           else if (p.y() < fMinExtent.y()) fMinExtent.setY(p.y());
00526           if (p.z() > fMaxExtent.z()) fMaxExtent.setZ(p.z());
00527           else if (p.z() < fMinExtent.z()) fMinExtent.setZ(p.z());
00528         }
00529       }
00530       else
00531       {
00532 #ifdef G4SPECSDEBUG     
00533         G4cout << p.x() << ":" << p.y() << ":" << p.z() << G4endl;
00534         G4cout << "Vertex #" << i << " of facet " << k
00535                << " found, redirecting to " << id << G4endl;
00536         G4cout << "===" << G4endl;
00537 #endif
00538         newIndex[i] = id;
00539       }
00540     }
00541     // only now it is possible to change vertices pointer
00542     //
00543     facet.SetVertices(&fVertexList);
00544     for (G4int i = 0; i < max; i++)
00545         facet.SetVertexIndex(i,newIndex[i]);
00546   }
00547   vector<G4ThreeVector>(fVertexList).swap(fVertexList);
00548   
00549 #ifdef G4SPECSDEBUG
00550   G4double previousValue = 0;
00551   for (set<G4VertexInfo,G4VertexComparator>::iterator res=
00552        vertexListSorted.begin(); res!=vertexListSorted.end(); ++res)
00553   {
00554     G4int id = (*res).id;
00555     G4ThreeVector vec = fVertexList[id];
00556     G4double mvalue = vec.x() + vec.y() + vec.z();
00557     if (previousValue && (previousValue - 1e-9 > mvalue))
00558       G4cout << "Error in CreateVertexList: previousValue " << previousValue 
00559              <<  " is smaller than mvalue " << mvalue << G4endl;
00560     previousValue = mvalue;
00561   }
00562 #endif
00563 }
00564 
00566 //
00567 void G4TessellatedSolid::DisplayAllocatedMemory()
00568 {
00569   G4int without = AllocatedMemoryWithoutVoxels();
00570   G4int with = AllocatedMemory();
00571   G4double ratio = (G4double) with / without;
00572   G4cout << "G4TessellatedSolid - Allocated memory without voxel overhead "
00573          << without << "; with " << with << "; ratio: " << ratio << G4endl;      
00574 }
00575 
00577 //
00578 void G4TessellatedSolid::SetSolidClosed (const G4bool t)
00579 {
00580   if (t)
00581   {
00582 #ifdef G4SPECSDEBUG    
00583     G4cout << "Creating vertex list..." << G4endl;
00584 #endif
00585     CreateVertexList();
00586 
00587 #ifdef G4SPECSDEBUG    
00588     G4cout << "Setting extreme facets..." << G4endl;
00589 #endif
00590     SetExtremeFacets();
00591     
00592 #ifdef G4SPECSDEBUG    
00593     G4cout << "Voxelizing..." << G4endl;
00594 #endif
00595     Voxelize();
00596 
00597 #ifdef G4SPECSDEBUG
00598     DisplayAllocatedMemory();
00599 #endif
00600 
00601   }  
00602   fSolidClosed = t;
00603 }
00604 
00606 //
00607 // GetSolidClosed
00608 //
00609 // Used to determine whether the solid is closed to adding further fFacets.
00610 //
00611 G4bool G4TessellatedSolid::GetSolidClosed () const
00612 {
00613   return fSolidClosed;
00614 }
00615 
00617 //
00618 // operator+=
00619 //
00620 // This operator allows the user to add two tessellated solids together, so
00621 // that the solid on the left then includes all of the facets in the solid
00622 // on the right.  Note that copies of the facets are generated, rather than
00623 // using the original facet set of the solid on the right.
00624 //
00625 G4TessellatedSolid &
00626 G4TessellatedSolid::operator+=(const G4TessellatedSolid &right)
00627 {
00628   G4int size = right.GetNumberOfFacets();
00629   for (G4int i = 0; i < size; ++i)
00630     AddFacet(right.GetFacet(i)->GetClone());
00631 
00632   return *this;
00633 }
00634 
00636 //
00637 // GetNumberOfFacets
00638 //
00639 G4int G4TessellatedSolid::GetNumberOfFacets () const
00640 {
00641   return fFacets.size();
00642 }
00643 
00645 //
00646 EInside G4TessellatedSolid::InsideVoxels(const G4ThreeVector &p) const
00647 {
00648   //
00649   // First the simple test - check if we're outside of the X-Y-Z extremes
00650   // of the tessellated solid.
00651   //
00652   if (OutsideOfExtent(p, kCarTolerance))
00653     return kOutside;
00654 
00655   vector<G4int> startingVoxel(3);
00656   fVoxels.GetVoxel(startingVoxel, p);
00657 
00658   const G4double dirTolerance = 1.0E-14;
00659 
00660   const vector<G4int> &startingCandidates =
00661     fVoxels.GetCandidates(startingVoxel);
00662   G4int limit = startingCandidates.size();
00663   if (limit == 0 && fInsides.GetNbits())
00664   {
00665     G4int index = fVoxels.GetPointIndex(p);
00666     EInside location = fInsides[index] ? kInside : kOutside;
00667     return location;
00668   }
00669 
00670   G4double minDist = kInfinity;
00671 
00672   for(G4int i = 0; i < limit; ++i)
00673   {
00674     G4int candidate = startingCandidates[i];
00675     G4VFacet &facet = *fFacets[candidate];
00676     G4double dist = facet.Distance(p,minDist);
00677     if (dist < minDist) minDist = dist;
00678     if (dist <= kCarToleranceHalf)
00679       return kSurface;
00680   }
00681 
00682   // The following is something of an adaptation of the method implemented by
00683   // Rickard Holmberg augmented with information from Schneider & Eberly,
00684   // "Geometric Tools for Computer Graphics," pp700-701, 2003. In essence,
00685   // we're trying to determine whether we're inside the volume by projecting
00686   // a few rays and determining if the first surface crossed is has a normal
00687   // vector between 0 to pi/2 (out-going) or pi/2 to pi (in-going).
00688   // We should also avoid rays which are nearly within the plane of the
00689   // tessellated surface, and therefore produce rays randomly.
00690   // For the moment, this is a bit over-engineered (belt-braces-and-ducttape).
00691   //
00692   G4double distOut          = kInfinity;
00693   G4double distIn           = kInfinity;
00694   G4double distO            = 0.0;
00695   G4double distI            = 0.0;
00696   G4double distFromSurfaceO = 0.0;
00697   G4double distFromSurfaceI = 0.0;
00698   G4ThreeVector normalO, normalI;
00699   G4bool crossingO          = false;
00700   G4bool crossingI          = false;
00701   EInside location          = kOutside;
00702   G4int sm                  = 0;
00703 
00704   G4bool nearParallel = false;
00705   do
00706   {
00707     // We loop until we find direction where the vector is not nearly parallel
00708     // to the surface of any facet since this causes ambiguities.  The usual
00709     // case is that the angles should be sufficiently different, but there
00710     // are 20 random directions to select from - hopefully sufficient.
00711     //
00712     distOut = distIn = kInfinity;
00713     const G4ThreeVector &v = fRandir[sm];
00714     sm++;
00715     //
00716     // This code could be voxelized by the same algorithm, which is used for
00717     // DistanceToOut(). We will traverse through fVoxels. we will call
00718     // intersect only for those, which would be candidates and was not
00719     // checked before.
00720     //
00721     G4ThreeVector currentPoint = p;
00722     G4ThreeVector direction = v.unit();
00723     // G4SurfBits exclusion(fVoxels.GetBitsPerSlice());
00724     vector<G4int> curVoxel(3);
00725     curVoxel = startingVoxel;
00726     G4double shiftBonus = kCarTolerance;
00727 
00728     G4bool crossed = false;
00729     G4bool started = true;
00730 
00731     do
00732     {
00733       const vector<G4int> &candidates =
00734         started ? startingCandidates : fVoxels.GetCandidates(curVoxel);
00735       started = false;
00736       if (G4int candidatesCount = candidates.size())
00737       {  
00738         for (G4int i = 0 ; i < candidatesCount; ++i)
00739         {
00740           G4int candidate = candidates[i];
00741           // bits.SetBitNumber(candidate);
00742           G4VFacet &facet = *fFacets[candidate];
00743 
00744           crossingO = facet.Intersect(p,v,true,distO,distFromSurfaceO,normalO);
00745           crossingI = facet.Intersect(p,v,false,distI,distFromSurfaceI,normalI);
00746 
00747           if (crossingO || crossingI)
00748           {
00749             crossed = true;
00750 
00751             nearParallel = (crossingO
00752                      && std::fabs(normalO.dot(v))<dirTolerance)
00753                      || (crossingI && std::fabs(normalI.dot(v))<dirTolerance);
00754             if (!nearParallel)
00755             {
00756               if (crossingO && distO > 0.0 && distO < distOut) 
00757                 distOut = distO;
00758               if (crossingI && distI > 0.0 && distI < distIn)  
00759                 distIn  = distI;
00760             }
00761             else break;
00762           }
00763         }
00764         if (nearParallel) break;
00765       }
00766       else
00767       {
00768         if (!crossed)
00769         {
00770           G4int index = fVoxels.GetVoxelsIndex(curVoxel);
00771           G4bool inside = fInsides[index];
00772           location = inside ? kInside : kOutside;
00773           return location;
00774         }
00775       }
00776 
00777       G4double shift=fVoxels.DistanceToNext(currentPoint, direction, curVoxel);
00778       if (shift == kInfinity) break;
00779 
00780       currentPoint += direction * (shift + shiftBonus);
00781     }
00782     while (fVoxels.UpdateCurrentVoxel(currentPoint, direction, curVoxel));
00783 
00784   }
00785   while (nearParallel && sm!=fMaxTries);
00786   //
00787   // Here we loop through the facets to find out if there is an intersection
00788   // between the ray and that facet.  The test if performed separately whether
00789   // the ray is entering the facet or exiting.
00790   //
00791 #ifdef G4VERBOSE
00792   if (sm == fMaxTries)
00793   {
00794     //
00795     // We've run out of random vector directions. If nTries is set sufficiently
00796     // low (nTries <= 0.5*maxTries) then this would indicate that there is
00797     // something wrong with geometry.
00798     //
00799     std::ostringstream message;
00800     G4int oldprc = message.precision(16);
00801     message << "Cannot determine whether point is inside or outside volume!"
00802       << G4endl
00803       << "Solid name       = " << GetName()  << G4endl
00804       << "Geometry Type    = " << fGeometryType  << G4endl
00805       << "Number of facets = " << fFacets.size() << G4endl
00806       << "Position:"  << G4endl << G4endl
00807       << "p.x() = "   << p.x()/mm << " mm" << G4endl
00808       << "p.y() = "   << p.y()/mm << " mm" << G4endl
00809       << "p.z() = "   << p.z()/mm << " mm";
00810     message.precision(oldprc);
00811     G4Exception("G4TessellatedSolid::Inside()",
00812                 "GeomSolids1002", JustWarning, message);
00813   }
00814 #endif
00815 
00816   // In the next if-then-elseif G4String the logic is as follows:
00817   // (1) You don't hit anything so cannot be inside volume, provided volume
00818   //     constructed correctly!
00819   // (2) Distance to inside (ie. nearest facet such that you enter facet) is
00820   //     shorter than distance to outside (nearest facet such that you exit
00821   //     facet) - on condition of safety distance - therefore we're outside.
00822   // (3) Distance to outside is shorter than distance to inside therefore
00823   //     we're inside.
00824   //
00825   if (distIn == kInfinity && distOut == kInfinity)
00826     location = kOutside;
00827   else if (distIn <= distOut - kCarToleranceHalf)
00828     location = kOutside;
00829   else if (distOut <= distIn - kCarToleranceHalf)
00830     location = kInside;
00831 
00832   return location;
00833 }
00834  
00836 //
00837 EInside G4TessellatedSolid::InsideNoVoxels (const G4ThreeVector &p) const
00838 {
00839   //
00840   // First the simple test - check if we're outside of the X-Y-Z extremes
00841   // of the tessellated solid.
00842   //
00843   if (OutsideOfExtent(p, kCarTolerance))
00844     return kOutside;
00845 
00846   const G4double dirTolerance = 1.0E-14;
00847 
00848   G4double minDist = kInfinity;
00849   //
00850   // Check if we are close to a surface
00851   //
00852   G4int size = fFacets.size();
00853   for (G4int i = 0; i < size; ++i)
00854   {
00855     G4VFacet &facet = *fFacets[i];
00856     G4double dist = facet.Distance(p,minDist);
00857     if (dist < minDist) minDist = dist;
00858     if (dist <= kCarToleranceHalf)
00859     {
00860       return kSurface;
00861     }
00862   }
00863   //
00864   // The following is something of an adaptation of the method implemented by
00865   // Rickard Holmberg augmented with information from Schneider & Eberly,
00866   // "Geometric Tools for Computer Graphics," pp700-701, 2003. In essence, we're
00867   // trying to determine whether we're inside the volume by projecting a few
00868   // rays and determining if the first surface crossed is has a normal vector
00869   // between 0 to pi/2 (out-going) or pi/2 to pi (in-going). We should also
00870   // avoid rays which are nearly within the plane of the tessellated surface,
00871   // and therefore produce rays randomly. For the moment, this is a bit
00872   // over-engineered (belt-braces-and-ducttape).
00873   //
00874 #if G4SPECSDEBUG
00875   G4int nTry                = 7;
00876 #else
00877   G4int nTry                = 3;
00878 #endif
00879   G4double distOut          = kInfinity;
00880   G4double distIn           = kInfinity;
00881   G4double distO            = 0.0;
00882   G4double distI            = 0.0;
00883   G4double distFromSurfaceO = 0.0;
00884   G4double distFromSurfaceI = 0.0;
00885   G4ThreeVector normalO(0.0,0.0,0.0);
00886   G4ThreeVector normalI(0.0,0.0,0.0);
00887   G4bool crossingO          = false;
00888   G4bool crossingI          = false;
00889   EInside location          = kOutside;
00890   EInside locationprime     = kOutside;
00891   G4int sm = 0;
00892 
00893   for (G4int i=0; i<nTry; ++i)
00894   {
00895     G4bool nearParallel = false;
00896     do
00897     {
00898       //
00899       // We loop until we find direction where the vector is not nearly parallel
00900       // to the surface of any facet since this causes ambiguities.  The usual
00901       // case is that the angles should be sufficiently different, but there
00902       // are 20 random directions to select from - hopefully sufficient.
00903       //
00904       distOut =  distIn = kInfinity;
00905       G4ThreeVector v = fRandir[sm];
00906       sm++;
00907       vector<G4VFacet*>::const_iterator f = fFacets.begin();
00908 
00909       do
00910       {
00911         //
00912         // Here we loop through the facets to find out if there is an
00913         // intersection between the ray and that facet. The test if performed
00914         // separately whether the ray is entering the facet or exiting.
00915         //
00916         crossingO = ((*f)->Intersect(p,v,true,distO,distFromSurfaceO,normalO));
00917         crossingI = ((*f)->Intersect(p,v,false,distI,distFromSurfaceI,normalI));
00918         if (crossingO || crossingI)
00919         {
00920           nearParallel = (crossingO && std::fabs(normalO.dot(v))<dirTolerance)
00921                       || (crossingI && std::fabs(normalI.dot(v))<dirTolerance);
00922           if (!nearParallel)
00923           {
00924             if (crossingO && distO > 0.0 && distO < distOut) distOut = distO;
00925             if (crossingI && distI > 0.0 && distI < distIn)  distIn  = distI;
00926           }
00927         }
00928       } while (!nearParallel && ++f!=fFacets.end());
00929     } while (nearParallel && sm!=fMaxTries);
00930 
00931 #ifdef G4VERBOSE
00932     if (sm == fMaxTries)
00933     {
00934       //
00935       // We've run out of random vector directions. If nTries is set
00936       // sufficiently low (nTries <= 0.5*maxTries) then this would indicate
00937       // that there is something wrong with geometry.
00938       //
00939       std::ostringstream message;
00940       G4int oldprc = message.precision(16);
00941       message << "Cannot determine whether point is inside or outside volume!"
00942         << G4endl
00943         << "Solid name       = " << GetName()  << G4endl
00944         << "Geometry Type    = " << fGeometryType  << G4endl
00945         << "Number of facets = " << fFacets.size() << G4endl
00946         << "Position:"  << G4endl << G4endl
00947         << "p.x() = "   << p.x()/mm << " mm" << G4endl
00948         << "p.y() = "   << p.y()/mm << " mm" << G4endl
00949         << "p.z() = "   << p.z()/mm << " mm";
00950       message.precision(oldprc);
00951       G4Exception("G4TessellatedSolid::Inside()",
00952         "GeomSolids1002", JustWarning, message);
00953     }
00954 #endif
00955     //
00956     // In the next if-then-elseif G4String the logic is as follows:
00957     // (1) You don't hit anything so cannot be inside volume, provided volume
00958     //     constructed correctly!
00959     // (2) Distance to inside (ie. nearest facet such that you enter facet) is
00960     //     shorter than distance to outside (nearest facet such that you exit
00961     //     facet) - on condition of safety distance - therefore we're outside.
00962     // (3) Distance to outside is shorter than distance to inside therefore
00963     // we're inside.
00964     //
00965     if (distIn == kInfinity && distOut == kInfinity)
00966       locationprime = kOutside;
00967     else if (distIn <= distOut - kCarToleranceHalf)
00968       locationprime = kOutside;
00969     else if (distOut <= distIn - kCarToleranceHalf)
00970       locationprime = kInside;
00971 
00972     if (i == 0) location = locationprime;
00973   }
00974 
00975   return location;
00976 }
00977 
00979 //
00980 // Return the outwards pointing unit normal of the shape for the
00981 // surface closest to the point at offset p.
00982 //
00983 G4bool G4TessellatedSolid::Normal (const G4ThreeVector &p,
00984                                          G4ThreeVector &aNormal) const
00985 {
00986   G4double minDist;
00987   G4VFacet *facet = 0;
00988 
00989   if (fVoxels.GetCountOfVoxels() > 1)
00990   {
00991     vector<G4int> curVoxel(3);
00992     fVoxels.GetVoxel(curVoxel, p);
00993     const vector<G4int> &candidates = fVoxels.GetCandidates(curVoxel);
00994     // fVoxels.GetCandidatesVoxelArray(p, candidates, 0);
00995 
00996     if (G4int limit = candidates.size())
00997     {
00998       minDist = kInfinity;
00999       for(G4int i = 0 ; i < limit ; ++i)
01000       {      
01001         G4int candidate = candidates[i];
01002         G4VFacet &fct = *fFacets[candidate];
01003         G4double dist = fct.Distance(p,minDist);
01004         if (dist < minDist) minDist = dist;
01005         if (dist <= kCarToleranceHalf)
01006         {
01007           aNormal = fct.GetSurfaceNormal();
01008           return true;
01009         }
01010       }
01011     }
01012     minDist = MinDistanceFacet(p, true, facet);
01013   }
01014   else
01015   {
01016     minDist = kInfinity;
01017     G4int size = fFacets.size();
01018     for (G4int i = 0; i < size; ++i)
01019     {
01020       G4VFacet &f = *fFacets[i];
01021       G4double dist = f.Distance(p, minDist);
01022       if (dist < minDist)
01023       {
01024         minDist  = dist;
01025         facet = &f;
01026       }
01027     }
01028   }
01029 
01030   if (minDist != kInfinity)
01031   {
01032     if (facet)  { aNormal = facet->GetSurfaceNormal(); }
01033     return minDist <= kCarToleranceHalf;
01034   }
01035   else
01036   {
01037 #ifdef G4VERBOSE
01038     std::ostringstream message;
01039     message << "Point p is not on surface !?" << G4endl
01040       << "          No facets found for point: " << p << " !" << G4endl
01041       << "          Returning approximated value for normal.";
01042 
01043     G4Exception("G4TessellatedSolid::SurfaceNormal(p)",
01044                 "GeomSolids1002", JustWarning, message );
01045 #endif
01046     aNormal = (p.z() > 0 ? G4ThreeVector(0,0,1) : G4ThreeVector(0,0,-1));
01047     return false;
01048   }
01049 }
01050 
01052 //
01053 // G4double DistanceToIn(const G4ThreeVector& p, const G4ThreeVector& v)
01054 //
01055 // Return the distance along the normalised vector v to the shape,
01056 // from the point at offset p. If there is no intersection, return
01057 // kInfinity. The first intersection resulting from 'leaving' a
01058 // surface/volume is discarded. Hence, this is tolerant of points on
01059 // surface of shape.
01060 //
01061 G4double
01062 G4TessellatedSolid::DistanceToInNoVoxels (const G4ThreeVector &p,
01063                                           const G4ThreeVector &v,
01064                                                 G4double /*aPstep*/) const
01065 {
01066   G4double minDist         = kInfinity;
01067   G4double dist            = 0.0;
01068   G4double distFromSurface = 0.0;
01069   G4ThreeVector normal;
01070 
01071 #if G4SPECSDEBUG
01072   if (Inside(p) == kInside )
01073   {
01074     std::ostringstream message;
01075     G4int oldprc = message.precision(16) ;
01076     message << "Point p is already inside!?" << G4endl
01077       << "Position:"  << G4endl << G4endl
01078       << "   p.x() = "   << p.x()/mm << " mm" << G4endl
01079       << "   p.y() = "   << p.y()/mm << " mm" << G4endl
01080       << "   p.z() = "   << p.z()/mm << " mm" << G4endl
01081       << "DistanceToOut(p) == " << DistanceToOut(p);
01082     message.precision(oldprc) ;
01083     G4Exception("G4TriangularFacet::DistanceToIn(p,v)",
01084                 "GeomSolids1002", JustWarning, message);
01085   }
01086 #endif
01087 
01088   G4int size = fFacets.size();
01089   for (G4int i = 0; i < size; ++i)
01090   {
01091     G4VFacet &facet = *fFacets[i];
01092     if (facet.Intersect(p,v,false,dist,distFromSurface,normal))
01093     {
01094       //
01095       // set minDist to the new distance to current facet if distFromSurface
01096       // is in positive direction and point is not at surface. If the point is
01097       // within 0.5*kCarTolerance of the surface, then force distance to be
01098       // zero and leave member function immediately (for efficiency), as
01099       // proposed by & credit to Akira Okumura.
01100       //
01101       if (distFromSurface > kCarToleranceHalf && dist >= 0.0 && dist < minDist)
01102       {
01103         minDist  = dist;
01104       }
01105       else 
01106         if (-kCarToleranceHalf <= dist && dist <= kCarToleranceHalf)
01107           return 0.0;
01108     }
01109   }
01110   return minDist;
01111 }
01112 
01114 //
01115 G4double
01116 G4TessellatedSolid::DistanceToOutNoVoxels (const G4ThreeVector &p,
01117                                            const G4ThreeVector &v,
01118                                                  G4ThreeVector &aNormalVector,
01119                                                  G4bool &aConvex,
01120                                                  G4double /*aPstep*/) const
01121 {
01122   G4double minDist         = kInfinity;
01123   G4double dist            = 0.0;
01124   G4double distFromSurface = 0.0;
01125   G4ThreeVector normal, minNormal;
01126 
01127 #if G4SPECSDEBUG
01128   if ( Inside(p) == kOutside )
01129   {
01130     std::ostringstream message;
01131     G4int oldprc = message.precision(16) ;
01132     message << "Point p is already outside!?" << G4endl
01133       << "Position:"  << G4endl << G4endl
01134       << "   p.x() = "   << p.x()/mm << " mm" << G4endl
01135       << "   p.y() = "   << p.y()/mm << " mm" << G4endl
01136       << "   p.z() = "   << p.z()/mm << " mm" << G4endl
01137       << "DistanceToIn(p) == " << DistanceToIn(p);
01138     message.precision(oldprc) ;
01139     G4Exception("G4TriangularFacet::DistanceToOut(p)",
01140                 "GeomSolids1002", JustWarning, message);
01141   }
01142 #endif
01143 
01144   G4bool isExtreme = false;
01145   G4int size = fFacets.size();
01146   for (G4int i = 0; i < size; ++i)
01147   {
01148     G4VFacet &facet = *fFacets[i];
01149     if (facet.Intersect(p,v,true,dist,distFromSurface,normal))
01150     {
01151       if (distFromSurface > 0.0 && distFromSurface <= kCarToleranceHalf &&
01152         facet.Distance(p,kCarTolerance) <= kCarToleranceHalf)
01153       {
01154         // We are on a surface. Return zero.
01155         aConvex = (fExtremeFacets.find(&facet) != fExtremeFacets.end());
01156         // Normal(p, aNormalVector);
01157         // aNormalVector = facet.GetSurfaceNormal();
01158         aNormalVector = normal;
01159         return 0.0;
01160       }
01161       if (dist >= 0.0 && dist < minDist)
01162       {
01163         minDist   = dist;
01164         minNormal = normal;
01165         isExtreme = (fExtremeFacets.find(&facet) != fExtremeFacets.end());
01166       }
01167     }
01168   }
01169   if (minDist < kInfinity)
01170   {
01171     aNormalVector = minNormal;
01172     aConvex = isExtreme;
01173     return minDist;
01174   }
01175   else
01176   {
01177     // No intersection found
01178     aConvex = false;
01179     Normal(p, aNormalVector);
01180     return 0.0;
01181   }
01182 }
01183 
01185 //
01186 void G4TessellatedSolid::
01187 DistanceToOutCandidates(const std::vector<G4int> &candidates,
01188                         const G4ThreeVector &aPoint,
01189                         const G4ThreeVector &direction,
01190                               G4double &minDist, G4ThreeVector &minNormal,
01191                               G4int &minCandidate ) const
01192 {
01193   G4int candidatesCount = candidates.size();
01194   G4double dist            = 0.0;
01195   G4double distFromSurface = 0.0;
01196   G4ThreeVector normal;
01197 
01198   for (G4int i = 0 ; i < candidatesCount; ++i)
01199   {
01200     G4int candidate = candidates[i];
01201     G4VFacet &facet = *fFacets[candidate];
01202     if (facet.Intersect(aPoint,direction,true,dist,distFromSurface,normal))
01203     {
01204       if (distFromSurface > 0.0 && distFromSurface <= kCarToleranceHalf
01205        && facet.Distance(aPoint,kCarTolerance) <= kCarToleranceHalf)
01206       {
01207         // We are on a surface
01208         //
01209         minDist = 0.0;
01210         minNormal = normal;
01211         minCandidate = candidate;
01212         break;
01213       }
01214       if (dist >= 0.0 && dist < minDist)
01215       {
01216         minDist = dist;
01217         minNormal = normal;
01218         minCandidate = candidate;
01219       }
01220     }
01221   }
01222 }
01223 
01225 //
01226 G4double
01227 G4TessellatedSolid::DistanceToOutCore(const G4ThreeVector &aPoint,
01228                                       const G4ThreeVector &aDirection,
01229                                             G4ThreeVector &aNormalVector,
01230                                             G4bool &aConvex,
01231                                             G4double aPstep) const
01232 {
01233   G4double minDistance;
01234 
01235   if (fVoxels.GetCountOfVoxels() > 1)
01236   {
01237     minDistance = kInfinity;
01238 
01239     G4ThreeVector currentPoint = aPoint;
01240     G4ThreeVector direction = aDirection.unit();
01241     G4double totalShift = 0;
01242     vector<G4int> curVoxel(3);
01243     if (!fVoxels.Contains(aPoint)) return 0;
01244 
01245     fVoxels.GetVoxel(curVoxel, currentPoint);
01246 
01247     G4double shiftBonus = kCarTolerance;
01248 
01249     const vector<G4int> *old = 0;
01250 
01251     G4int minCandidate = -1;
01252     do
01253     {
01254       const vector<G4int> &candidates = fVoxels.GetCandidates(curVoxel);
01255       if (old == &candidates)
01256         old++;
01257       if (old != &candidates && candidates.size())
01258       {
01259         DistanceToOutCandidates(candidates, aPoint, direction, minDistance,
01260                                 aNormalVector, minCandidate); 
01261         if (minDistance <= totalShift) break; 
01262       }
01263 
01264       G4double shift=fVoxels.DistanceToNext(currentPoint, direction, curVoxel);
01265       if (shift == kInfinity) break;
01266 
01267       totalShift += shift;
01268       if (minDistance <= totalShift) break;
01269 
01270       currentPoint += direction * (shift + shiftBonus);
01271 
01272       old = &candidates;
01273     }
01274     while (fVoxels.UpdateCurrentVoxel(currentPoint, direction, curVoxel));
01275 
01276     if (minCandidate < 0)
01277     {
01278       // No intersection found
01279       minDistance = 0;
01280       aConvex = false;
01281       Normal(aPoint, aNormalVector);
01282     }
01283     else
01284     {
01285       aConvex = (fExtremeFacets.find(fFacets[minCandidate])
01286               != fExtremeFacets.end());
01287     }
01288   }
01289   else
01290   {
01291     minDistance = DistanceToOutNoVoxels(aPoint, aDirection, aNormalVector,
01292                                         aConvex, aPstep);
01293   }
01294   return minDistance;
01295 }
01296 
01298 //
01299 G4double G4TessellatedSolid::
01300 DistanceToInCandidates(const std::vector<G4int> &candidates,
01301                        const G4ThreeVector &aPoint,
01302                        const G4ThreeVector &direction) const
01303 {
01304   G4int candidatesCount = candidates.size();
01305   G4double dist            = 0.0;
01306   G4double distFromSurface = 0.0;
01307   G4ThreeVector normal;
01308 
01309   G4double minDistance = kInfinity;   
01310   for (G4int i = 0 ; i < candidatesCount; ++i)
01311   {
01312     G4int candidate = candidates[i];
01313     G4VFacet &facet = *fFacets[candidate];
01314     if (facet.Intersect(aPoint,direction,false,dist,distFromSurface,normal))
01315     {
01316       //
01317       // Set minDist to the new distance to current facet if distFromSurface is
01318       // in positive direction and point is not at surface. If the point is
01319       // within 0.5*kCarTolerance of the surface, then force distance to be
01320       // zero and leave member function immediately (for efficiency), as
01321       // proposed by & credit to Akira Okumura.
01322       //
01323       if ( (distFromSurface > kCarToleranceHalf)
01324         && (dist >= 0.0) && (dist < minDistance))
01325       {
01326         minDistance  = dist;
01327       }
01328       else if (-kCarToleranceHalf <= dist && dist <= kCarToleranceHalf)
01329       {
01330         return 0.0;
01331       }
01332     }
01333   }
01334   return minDistance;
01335 }
01336 
01338 //
01339 G4double
01340 G4TessellatedSolid::DistanceToInCore(const G4ThreeVector &aPoint,
01341                                      const G4ThreeVector &aDirection,
01342                                            G4double aPstep) const
01343 {
01344   G4double minDistance;
01345 
01346   if (fVoxels.GetCountOfVoxels() > 1)
01347   {
01348     minDistance = kInfinity;
01349     G4ThreeVector currentPoint = aPoint;
01350     G4ThreeVector direction = aDirection.unit();
01351     G4double shift = fVoxels.DistanceToFirst(currentPoint, direction);
01352     if (shift == kInfinity) return shift;
01353     G4double shiftBonus = kCarTolerance;
01354     if (shift) 
01355       currentPoint += direction * (shift + shiftBonus);
01356     // if (!fVoxels.Contains(currentPoint))  return minDistance;
01357     G4double totalShift = shift;
01358 
01359     // G4SurfBits exclusion; // (1/*fVoxels.GetBitsPerSlice()*/);
01360     vector<G4int> curVoxel(3);
01361 
01362     fVoxels.GetVoxel(curVoxel, currentPoint);
01363     do
01364     {
01365       const vector<G4int> &candidates = fVoxels.GetCandidates(curVoxel);
01366       if (candidates.size())
01367       {
01368         G4double distance=DistanceToInCandidates(candidates, aPoint, direction);
01369         if (minDistance > distance) minDistance = distance;
01370         if (distance < totalShift) break;
01371       }
01372 
01373       shift = fVoxels.DistanceToNext(currentPoint, direction, curVoxel);
01374       if (shift == kInfinity /*|| shift == 0*/) break;
01375 
01376       totalShift += shift;
01377       if (minDistance < totalShift) break;
01378 
01379       currentPoint += direction * (shift + shiftBonus);
01380     }
01381     while (fVoxels.UpdateCurrentVoxel(currentPoint, direction, curVoxel));
01382   }
01383   else
01384   {
01385     minDistance = DistanceToInNoVoxels(aPoint, aDirection, aPstep);
01386   }
01387 
01388   return minDistance;
01389 }
01390 
01392 //
01393 G4bool
01394 G4TessellatedSolid::CompareSortedVoxel(const std::pair<G4int, G4double> &l,
01395                                        const std::pair<G4int, G4double> &r)
01396 {
01397   return l.second < r.second;
01398 }
01399 
01401 //
01402 G4double
01403 G4TessellatedSolid::MinDistanceFacet(const G4ThreeVector &p,
01404                                            G4bool simple,
01405                                            G4VFacet * &minFacet) const
01406 {
01407   G4double minDist = kInfinity;
01408 
01409   G4int size = fVoxels.GetVoxelBoxesSize();
01410   vector<pair<G4int, G4double> > voxelsSorted(size);
01411   
01412   pair<G4int, G4double> info;
01413 
01414   for (G4int i = 0; i < size; ++i)
01415   {
01416     const G4VoxelBox &voxelBox = fVoxels.GetVoxelBox(i);
01417 
01418     G4ThreeVector pointShifted = p - voxelBox.pos;
01419     G4double safety = fVoxels.MinDistanceToBox(pointShifted, voxelBox.hlen);
01420     info.first = i;
01421     info.second = safety;
01422     voxelsSorted[i] = info;
01423   }
01424 
01425   std::sort(voxelsSorted.begin(), voxelsSorted.end(),
01426             &G4TessellatedSolid::CompareSortedVoxel);
01427 
01428   for (G4int i = 0; i < size; ++i)
01429   {
01430     const pair<G4int,G4double> &inf = voxelsSorted[i];
01431     G4double dist = inf.second;
01432     if (dist > minDist) break;
01433 
01434     const vector<G4int> &candidates = fVoxels.GetVoxelBoxCandidates(inf.first);
01435     G4int csize = candidates.size();
01436     for (G4int j = 0; j < csize; ++j)
01437     {
01438       G4int candidate = candidates[j];
01439       G4VFacet &facet = *fFacets[candidate];
01440       dist = simple ? facet.Distance(p,minDist)
01441                     : facet.Distance(p,minDist,false);
01442       if (dist < minDist)
01443       {
01444         minDist  = dist;
01445         minFacet = &facet;
01446       }
01447     }
01448   }
01449   return minDist;
01450 }
01451 
01453 //
01454 G4double G4TessellatedSolid::SafetyFromOutside (const G4ThreeVector &p,
01455                                                       G4bool aAccurate) const
01456 {
01457 #if G4SPECSDEBUG
01458   if ( Inside(p) == kInside )
01459   {
01460     std::ostringstream message;
01461     G4int oldprc = message.precision(16) ;
01462     message << "Point p is already inside!?" << G4endl
01463       << "Position:"  << G4endl << G4endl
01464       << "p.x() = "   << p.x()/mm << " mm" << G4endl
01465       << "p.y() = "   << p.y()/mm << " mm" << G4endl
01466       << "p.z() = "   << p.z()/mm << " mm" << G4endl
01467       << "DistanceToOut(p) == " << DistanceToOut(p);
01468     message.precision(oldprc) ;
01469     G4Exception("G4TriangularFacet::DistanceToIn(p)",
01470                 "GeomSolids1002", JustWarning, message);
01471   }
01472 #endif
01473 
01474   G4double minDist;
01475 
01476   if (fVoxels.GetCountOfVoxels() > 1)
01477   {
01478     if (!aAccurate)
01479       return fVoxels.DistanceToBoundingBox(p);
01480 
01481     if (!OutsideOfExtent(p, kCarTolerance))
01482     {
01483       vector<G4int> startingVoxel(3);
01484       fVoxels.GetVoxel(startingVoxel, p);
01485       const vector<G4int> &candidates = fVoxels.GetCandidates(startingVoxel);
01486       if (candidates.size() == 0 && fInsides.GetNbits())
01487       {
01488         G4int index = fVoxels.GetPointIndex(p);
01489         if (fInsides[index]) return 0.;
01490       }
01491     }
01492 
01493     G4VFacet *facet;
01494     minDist = MinDistanceFacet(p, true, facet);
01495   }
01496   else
01497   {
01498     minDist = kInfinity;
01499     G4int size = fFacets.size();
01500     for (G4int i = 0; i < size; ++i)
01501     {
01502       G4VFacet &facet = *fFacets[i];
01503       G4double dist = facet.Distance(p,minDist);
01504       if (dist < minDist) minDist  = dist;
01505     }
01506   }
01507   return minDist;
01508 }
01509 
01511 //
01512 G4double
01513 G4TessellatedSolid::SafetyFromInside (const G4ThreeVector &p, G4bool) const
01514 {  
01515 #if G4SPECSDEBUG
01516   if ( Inside(p) == kOutside )
01517   {
01518     std::ostringstream message;
01519     G4int oldprc = message.precision(16) ;
01520     message << "Point p is already outside!?" << G4endl
01521       << "Position:"  << G4endl << G4endl
01522       << "p.x() = "   << p.x()/mm << " mm" << G4endl
01523       << "p.y() = "   << p.y()/mm << " mm" << G4endl
01524       << "p.z() = "   << p.z()/mm << " mm" << G4endl
01525       << "DistanceToIn(p) == " << DistanceToIn(p);
01526     message.precision(oldprc) ;
01527     G4Exception("G4TriangularFacet::DistanceToOut(p)",
01528                 "GeomSolids1002", JustWarning, message);
01529   }
01530 #endif
01531 
01532   G4double minDist;
01533 
01534   if (OutsideOfExtent(p, kCarTolerance)) return 0.0;
01535 
01536   if (fVoxels.GetCountOfVoxels() > 1)
01537   {
01538     G4VFacet *facet;
01539     minDist = MinDistanceFacet(p, true, facet);
01540   }
01541   else
01542   {
01543     minDist = kInfinity;
01544     G4double dist = 0.0;
01545     G4int size = fFacets.size();
01546     for (G4int i = 0; i < size; ++i)
01547     {
01548       G4VFacet &facet = *fFacets[i];
01549       dist = facet.Distance(p,minDist);
01550       if (dist < minDist) minDist  = dist;
01551     }
01552   }
01553   return minDist;
01554 }
01555 
01557 //
01558 // G4GeometryType GetEntityType() const;
01559 //
01560 // Provide identification of the class of an object
01561 //
01562 G4GeometryType G4TessellatedSolid::GetEntityType () const
01563 {
01564   return fGeometryType;
01565 }
01566 
01568 //
01569 std::ostream &G4TessellatedSolid::StreamInfo(std::ostream &os) const
01570 {
01571   os << G4endl;
01572   os << "Geometry Type    = " << fGeometryType  << G4endl;
01573   os << "Number of facets = " << fFacets.size() << G4endl;
01574 
01575   G4int size = fFacets.size();
01576   for (G4int i = 0; i < size; ++i)
01577   {
01578     os << "FACET #          = " << i + 1 << G4endl;
01579     G4VFacet &facet = *fFacets[i];
01580     facet.StreamInfo(os);
01581   }
01582   os << G4endl;
01583 
01584   return os;
01585 }
01586 
01588 //
01589 // Make a clone of the object
01590 //
01591 G4VSolid* G4TessellatedSolid::Clone() const
01592 {
01593   return new G4TessellatedSolid(*this);
01594 }
01595 
01597 //
01598 // EInside G4TessellatedSolid::Inside (const G4ThreeVector &p) const
01599 //
01600 // This method must return:
01601 //    * kOutside if the point at offset p is outside the shape
01602 //      boundaries plus kCarTolerance/2,
01603 //    * kSurface if the point is <= kCarTolerance/2 from a surface, or
01604 //    * kInside otherwise.
01605 //
01606 EInside G4TessellatedSolid::Inside (const G4ThreeVector &aPoint) const
01607 {
01608   EInside location;
01609 
01610   if (fVoxels.GetCountOfVoxels() > 1)
01611   {
01612     location = InsideVoxels(aPoint);
01613   }
01614   else
01615   {
01616     location = InsideNoVoxels(aPoint);
01617   }
01618   return location;
01619 }
01620 
01622 //
01623 G4ThreeVector G4TessellatedSolid::SurfaceNormal(const G4ThreeVector& p) const
01624 {
01625   G4ThreeVector n;
01626   Normal(p, n);
01627   return n;
01628 }
01629 
01631 //
01632 // G4double DistanceToIn(const G4ThreeVector& p)
01633 //
01634 // Calculate distance to nearest surface of shape from an outside point p. The
01635 // distance can be an underestimate.
01636 //
01637 G4double G4TessellatedSolid::DistanceToIn(const G4ThreeVector& p) const
01638 {
01639   return SafetyFromOutside(p,false);
01640 }
01641 
01643 //
01644 G4double G4TessellatedSolid::DistanceToIn(const G4ThreeVector& p,
01645                                           const G4ThreeVector& v)const
01646 {
01647   return DistanceToInCore(p,v,kInfinity);
01648 }
01649 
01651 //
01652 // G4double DistanceToOut(const G4ThreeVector& p)
01653 //
01654 // Calculate distance to nearest surface of shape from an inside
01655 // point. The distance can be an underestimate.
01656 //
01657 G4double G4TessellatedSolid::DistanceToOut(const G4ThreeVector& p) const
01658 {
01659   return SafetyFromInside(p,false);
01660 }
01661 
01663 //
01664 // G4double DistanceToOut(const G4ThreeVector& p, const G4ThreeVector& v,
01665 //                        const G4bool calcNorm=false,
01666 //                        G4bool *validNorm=0, G4ThreeVector *n=0);
01667 //
01668 // Return distance along the normalised vector v to the shape, from a
01669 // point at an offset p inside or on the surface of the
01670 // shape. Intersections with surfaces, when the point is not greater
01671 // than kCarTolerance/2 from a surface, must be ignored.
01672 //     If calcNorm is true, then it must also set validNorm to either
01673 //     * true, if the solid lies entirely behind or on the exiting
01674 //        surface. Then it must set n to the outwards normal vector
01675 //        (the Magnitude of the vector is not defined).
01676 //     * false, if the solid does not lie entirely behind or on the
01677 //       exiting surface.
01678 // If calcNorm is false, then validNorm and n are unused.
01679 //
01680 G4double G4TessellatedSolid::DistanceToOut(const G4ThreeVector& p,
01681                                            const G4ThreeVector& v,
01682                                            const G4bool calcNorm,
01683                                                  G4bool *validNorm,
01684                                                  G4ThreeVector *norm) const
01685 {
01686   G4ThreeVector n;
01687   G4bool valid;
01688 
01689   G4double dist = DistanceToOutCore(p, v, n, valid);
01690   if (calcNorm)
01691   {
01692     *norm = n;
01693     *validNorm = valid;
01694   }
01695   return dist;
01696 }
01697 
01699 //
01700 void G4TessellatedSolid::DescribeYourselfTo (G4VGraphicsScene& scene) const
01701 {
01702   scene.AddSolid (*this);
01703 }
01704 
01706 //
01707 G4Polyhedron *G4TessellatedSolid::CreatePolyhedron () const
01708 {
01709   G4int nVertices = fVertexList.size();
01710   G4int nFacets   = fFacets.size();
01711   G4PolyhedronArbitrary *polyhedron =
01712     new G4PolyhedronArbitrary (nVertices, nFacets);
01713   for (G4ThreeVectorList::const_iterator v= fVertexList.begin();
01714                                          v!=fVertexList.end(); ++v)
01715   {
01716     polyhedron->AddVertex(*v);
01717   }
01718 
01719   G4int size = fFacets.size();
01720   for (G4int i = 0; i < size; ++i)
01721   {
01722     G4VFacet &facet = *fFacets[i];
01723     G4int v[4];
01724     G4int n = facet.GetNumberOfVertices();
01725     if (n > 4) n = 4;
01726     else if (n == 3) v[3] = 0;
01727     for (G4int j=0; j<n; ++j)
01728     {
01729       G4int k = facet.GetVertexIndex(j);
01730       v[j] = k+1;
01731     }
01732     polyhedron->AddFacet(v[0],v[1],v[2],v[3]);
01733   }
01734   polyhedron->SetReferences();  
01735 
01736   return (G4Polyhedron*) polyhedron;
01737 }
01738 
01740 //
01741 G4NURBS *G4TessellatedSolid::CreateNURBS () const
01742 {
01743   return 0;
01744 }
01745 
01747 //
01748 // GetPolyhedron
01749 //
01750 G4Polyhedron* G4TessellatedSolid::GetPolyhedron () const
01751 {
01752   if (!fpPolyhedron ||
01753     fpPolyhedron->GetNumberOfRotationStepsAtTimeOfCreation() !=
01754     fpPolyhedron->GetNumberOfRotationSteps())
01755   {
01756     delete fpPolyhedron;
01757     fpPolyhedron = CreatePolyhedron();
01758   }
01759   return fpPolyhedron;
01760 }
01761 
01763 //
01764 // CalculateExtent
01765 //
01766 // Based on correction provided by Stan Seibert, University of Texas.
01767 //
01768 G4bool
01769 G4TessellatedSolid::CalculateExtent(const EAxis pAxis,
01770                                     const G4VoxelLimits& pVoxelLimit,
01771                                     const G4AffineTransform& pTransform,
01772                                           G4double& pMin, G4double& pMax) const
01773 {
01774   G4ThreeVectorList transVertexList(fVertexList);
01775   G4int size = fVertexList.size();
01776 
01777   // Put solid into transformed frame
01778   for (G4int i=0; i < size; ++i)
01779   {
01780     pTransform.ApplyPointTransform(transVertexList[i]);
01781   }
01782 
01783   // Find min and max extent in each dimension
01784   G4ThreeVector minExtent(kInfinity, kInfinity, kInfinity);
01785   G4ThreeVector maxExtent(-kInfinity, -kInfinity, -kInfinity);
01786 
01787   size = transVertexList.size();
01788   for (G4int i=0; i< size; ++i)
01789   {
01790     for (G4int axis=G4ThreeVector::X; axis < G4ThreeVector::SIZE; ++axis)
01791     {
01792       G4double coordinate = transVertexList[i][axis];
01793       if (coordinate < minExtent[axis])
01794       { minExtent[axis] = coordinate; }
01795       if (coordinate > maxExtent[axis])
01796       { maxExtent[axis] = coordinate; }
01797     }
01798   }
01799 
01800   // Check for containment and clamp to voxel boundaries
01801   for (G4int axis=G4ThreeVector::X; axis < G4ThreeVector::SIZE; ++axis)
01802   {
01803     EAxis geomAxis = kXAxis; // U geom classes use different index type
01804     switch(axis)
01805     {
01806     case G4ThreeVector::X: geomAxis = kXAxis; break;
01807     case G4ThreeVector::Y: geomAxis = kYAxis; break;
01808     case G4ThreeVector::Z: geomAxis = kZAxis; break;
01809     }
01810     G4bool isLimited = pVoxelLimit.IsLimited(geomAxis);
01811     G4double voxelMinExtent = pVoxelLimit.GetMinExtent(geomAxis);
01812     G4double voxelMaxExtent = pVoxelLimit.GetMaxExtent(geomAxis);
01813 
01814     if (isLimited)
01815     {
01816       if ( minExtent[axis] > voxelMaxExtent+kCarTolerance ||
01817         maxExtent[axis] < voxelMinExtent-kCarTolerance    )
01818       {
01819         return false ;
01820       }
01821       else
01822       {
01823         if (minExtent[axis] < voxelMinExtent)
01824         {
01825           minExtent[axis] = voxelMinExtent ;
01826         }
01827         if (maxExtent[axis] > voxelMaxExtent)
01828         {
01829           maxExtent[axis] = voxelMaxExtent;
01830         }
01831       }
01832     }
01833   }
01834 
01835   // Convert pAxis into G4ThreeVector index
01836   G4int vecAxis=0;
01837   switch(pAxis)
01838   {
01839   case kXAxis: vecAxis = G4ThreeVector::X; break;
01840   case kYAxis: vecAxis = G4ThreeVector::Y; break;
01841   case kZAxis: vecAxis = G4ThreeVector::Z; break;
01842   default: break;
01843   } 
01844 
01845   pMin = minExtent[vecAxis] - kCarTolerance;
01846   pMax = maxExtent[vecAxis] + kCarTolerance;
01847 
01848   return true;
01849 }
01850 
01852 //
01853 G4double G4TessellatedSolid::GetMinXExtent () const
01854 {
01855   return fMinExtent.x();
01856 }
01857 
01859 //
01860 G4double G4TessellatedSolid::GetMaxXExtent () const
01861 {
01862   return fMaxExtent.x();
01863 }
01864 
01866 //
01867 G4double G4TessellatedSolid::GetMinYExtent () const
01868 {
01869   return fMinExtent.y();
01870 }
01871 
01873 //
01874 G4double G4TessellatedSolid::GetMaxYExtent () const
01875 {
01876   return fMaxExtent.y();
01877 }
01878 
01880 //
01881 G4double G4TessellatedSolid::GetMinZExtent () const
01882 {
01883   return fMinExtent.z();
01884 }
01885 
01887 //
01888 G4double G4TessellatedSolid::GetMaxZExtent () const
01889 {
01890   return fMaxExtent.z();
01891 }
01892 
01894 //
01895 G4VisExtent G4TessellatedSolid::GetExtent () const
01896 {
01897   return G4VisExtent (fMinExtent.x(), fMaxExtent.x(), fMinExtent.y(), fMaxExtent.y(), fMinExtent.z(), fMaxExtent.z());
01898 }
01899 
01901 //
01902 G4double G4TessellatedSolid::GetCubicVolume ()
01903 {
01904   if(fCubicVolume != 0.) {;}
01905   else   { fCubicVolume = G4VSolid::GetCubicVolume(); }
01906   return fCubicVolume;
01907 }
01908 
01910 //
01911 G4double G4TessellatedSolid::GetSurfaceArea ()
01912 {
01913   if (fSurfaceArea != 0.) return fSurfaceArea;
01914 
01915   G4int size = fFacets.size();
01916   for (G4int i = 0; i < size; ++i)
01917   {
01918     G4VFacet &facet = *fFacets[i];
01919     fSurfaceArea += facet.GetArea();
01920   }
01921   return fSurfaceArea;
01922 }
01923 
01925 //
01926 G4ThreeVector G4TessellatedSolid::GetPointOnSurface() const
01927 {
01928   // Select randomly a facet and return a random point on it
01929 
01930   G4int i = (G4int) G4RandFlat::shoot(0., fFacets.size());
01931   return fFacets[i]->GetPointOnFace();
01932 }
01933 
01935 //
01936 // SetRandomVectorSet
01937 //
01938 // This is a set of predefined random vectors (if that isn't a contradition
01939 // in terms!) used to generate rays from a user-defined point.  The member
01940 // function Inside uses these to determine whether the point is inside or
01941 // outside of the tessellated solid.  All vectors should be unit vectors.
01942 //
01943 void G4TessellatedSolid::SetRandomVectors ()
01944 {
01945   fRandir.resize(20);
01946   fRandir[0]  =
01947     G4ThreeVector(-0.9577428892113370, 0.2732676269591740, 0.0897405271949221);
01948   fRandir[1]  =
01949     G4ThreeVector(-0.8331264504940770,-0.5162067214954600,-0.1985722492445700);
01950   fRandir[2]  =
01951     G4ThreeVector(-0.1516671651108820, 0.9666292616127460, 0.2064580868390110);
01952   fRandir[3]  =
01953     G4ThreeVector( 0.6570250350323190,-0.6944539025883300, 0.2933460081893360);
01954   fRandir[4]  =
01955     G4ThreeVector(-0.4820456281280320,-0.6331060000098690,-0.6056474264406270);
01956   fRandir[5]  =
01957     G4ThreeVector( 0.7629032554236800 , 0.1016854697539910,-0.6384658864065180);
01958   fRandir[6]  =
01959     G4ThreeVector( 0.7689540409061150, 0.5034929891988220, 0.3939600142169160);
01960   fRandir[7]  =
01961     G4ThreeVector( 0.5765188359255740, 0.5997271636278330,-0.5549354566343150);
01962   fRandir[8]  =
01963     G4ThreeVector( 0.6660632777862070,-0.6362809868288380, 0.3892379937580790);
01964   fRandir[9]  =
01965     G4ThreeVector( 0.3824415020414780, 0.6541792713761380,-0.6525243125110690);
01966   fRandir[10] =
01967     G4ThreeVector(-0.5107726564526760, 0.6020905056811610, 0.6136760679616570);
01968   fRandir[11] =
01969     G4ThreeVector( 0.7459135439578050, 0.6618796061649330, 0.0743530220183488);
01970   fRandir[12] =
01971     G4ThreeVector( 0.1536405855311580, 0.8117477913978260,-0.5634359711967240);
01972   fRandir[13] =
01973     G4ThreeVector( 0.0744395301705579,-0.8707110101772920,-0.4861286795736560);
01974   fRandir[14] =
01975     G4ThreeVector(-0.1665874645185400, 0.6018553940549240,-0.7810369397872780);
01976   fRandir[15] =
01977     G4ThreeVector( 0.7766902003633100, 0.6014617505959970,-0.1870724331097450);
01978   fRandir[16] =
01979     G4ThreeVector(-0.8710128685847430,-0.1434320216603030,-0.4698551243971010);
01980   fRandir[17] =
01981     G4ThreeVector( 0.8901082092766820,-0.4388411398893870, 0.1229871120030100);
01982   fRandir[18] =
01983     G4ThreeVector(-0.6430417431544370,-0.3295938228697690, 0.6912779675984150);
01984   fRandir[19] =
01985     G4ThreeVector( 0.6331124368380410, 0.6306211461665000, 0.4488714875425340);
01986 
01987   fMaxTries = 20;
01988 }
01989 
01991 //
01992 G4int G4TessellatedSolid::AllocatedMemoryWithoutVoxels()
01993 {
01994   G4int base = sizeof(*this);
01995   base += fVertexList.capacity() * sizeof(G4ThreeVector);
01996   base += fRandir.capacity() * sizeof(G4ThreeVector);
01997 
01998   G4int limit = fFacets.size();
01999   for (G4int i = 0; i < limit; i++)
02000   {
02001     G4VFacet &facet = *fFacets[i];
02002     base += facet.AllocatedMemory();
02003   }
02004 
02005   std::set<G4VFacet *>::const_iterator beg, end, it;
02006   beg = fExtremeFacets.begin();
02007   end = fExtremeFacets.end();
02008   for (it = beg; it != end; it++)
02009   {
02010     G4VFacet &facet = *(*it);
02011     base += facet.AllocatedMemory();
02012   }
02013   return base;
02014 }
02015 
02017 //
02018 G4int G4TessellatedSolid::AllocatedMemory()
02019 {
02020   G4int size = AllocatedMemoryWithoutVoxels();
02021   G4int sizeInsides = fInsides.GetNbytes();
02022   G4int sizeVoxels = fVoxels.AllocatedMemory();
02023   size += sizeInsides + sizeVoxels;
02024   return size;
02025 }

Generated on Mon May 27 17:49:57 2013 for Geant4 by  doxygen 1.4.7