G4Ellipsoid.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.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 // $Id: G4Ellipsoid.cc 67011 2013-01-29 16:17:41Z gcosmo $
00027 //
00028 // class G4Ellipsoid
00029 //
00030 // Implementation for G4Ellipsoid class
00031 //
00032 // History:
00033 //
00034 // 10.11.99 G.Horton-Smith  -- first writing, based on G4Sphere class
00035 // 25.02.05 G.Guerrieri -- Modified for future Geant4 release
00036 //
00037 // --------------------------------------------------------------------
00038 
00039 #include "globals.hh"
00040 
00041 #include "G4Ellipsoid.hh"
00042 
00043 #include "G4VoxelLimits.hh"
00044 #include "G4AffineTransform.hh"
00045 #include "G4GeometryTolerance.hh"
00046 
00047 #include "meshdefs.hh"
00048 
00049 #include "Randomize.hh"
00050 
00051 #include "G4VGraphicsScene.hh"
00052 #include "G4Polyhedron.hh"
00053 #include "G4NURBS.hh"
00054 #include "G4NURBSbox.hh"
00055 #include "G4VisExtent.hh"
00056 
00057 using namespace CLHEP;
00058 
00060 //
00061 // constructor - check parameters, convert angles so 0<sphi+dpshi<=2_PI
00062 //             - note if pDPhi>2PI then reset to 2PI
00063 
00064 G4Ellipsoid::G4Ellipsoid(const G4String& pName,
00065                                G4double pxSemiAxis,
00066                                G4double pySemiAxis,
00067                                G4double pzSemiAxis,
00068                                G4double pzBottomCut,
00069                                G4double pzTopCut)
00070   : G4VSolid(pName), fpPolyhedron(0), fCubicVolume(0.), fSurfaceArea(0.),
00071     zBottomCut(0.), zTopCut(0.)
00072 {
00073   // note: for users that want to use the full ellipsoid it is useful
00074   // to include a default for the cuts 
00075 
00076   kRadTolerance = G4GeometryTolerance::GetInstance()->GetRadialTolerance();
00077 
00078   // Check Semi-Axis
00079   if ( (pxSemiAxis<=0.) || (pySemiAxis<=0.) || (pzSemiAxis<=0.) )
00080   {
00081      std::ostringstream message;
00082      message << "Invalid semi-axis - " << GetName();
00083      G4Exception("G4Ellipsoid::G4Ellipsoid()", "GeomSolids0002",
00084                  FatalErrorInArgument, message);
00085   }
00086   SetSemiAxis(pxSemiAxis, pySemiAxis, pzSemiAxis);
00087 
00088   if ( pzBottomCut == 0 && pzTopCut == 0 )
00089   {
00090      SetZCuts(-pzSemiAxis, pzSemiAxis);
00091   }
00092   else if ( (pzBottomCut < pzSemiAxis) && (pzTopCut > -pzSemiAxis)
00093          && (pzBottomCut < pzTopCut) )
00094   {
00095      SetZCuts(pzBottomCut, pzTopCut);
00096   }
00097   else
00098   {
00099      std::ostringstream message;
00100      message << "Invalid z-coordinate for cutting plane - " << GetName();
00101      G4Exception("G4Ellipsoid::G4Ellipsoid()", "GeomSolids0002",
00102                  FatalErrorInArgument, message);
00103   }
00104 }
00105 
00107 //
00108 // Fake default constructor - sets only member data and allocates memory
00109 //                            for usage restricted to object persistency.
00110 //
00111 G4Ellipsoid::G4Ellipsoid( __void__& a )
00112   : G4VSolid(a), fpPolyhedron(0), kRadTolerance(0.), fCubicVolume(0.),
00113     fSurfaceArea(0.), xSemiAxis(0.), ySemiAxis(0.), zSemiAxis(0.),
00114     semiAxisMax(0.), zBottomCut(0.), zTopCut(0.)
00115 {
00116 }
00117 
00119 //
00120 // Destructor
00121 
00122 G4Ellipsoid::~G4Ellipsoid()
00123 {
00124 }
00125 
00127 //
00128 // Copy constructor
00129 
00130 G4Ellipsoid::G4Ellipsoid(const G4Ellipsoid& rhs)
00131   : G4VSolid(rhs),
00132     fpPolyhedron(0), kRadTolerance(rhs.kRadTolerance),
00133     fCubicVolume(rhs.fCubicVolume), fSurfaceArea(rhs.fSurfaceArea),
00134     xSemiAxis(rhs.xSemiAxis), ySemiAxis(rhs.ySemiAxis),
00135     zSemiAxis(rhs.zSemiAxis), semiAxisMax(rhs.semiAxisMax),
00136     zBottomCut(rhs.zBottomCut), zTopCut(rhs.zTopCut)
00137 {
00138 }
00139 
00141 //
00142 // Assignment operator
00143 
00144 G4Ellipsoid& G4Ellipsoid::operator = (const G4Ellipsoid& rhs) 
00145 {
00146    // Check assignment to self
00147    //
00148    if (this == &rhs)  { return *this; }
00149 
00150    // Copy base class data
00151    //
00152    G4VSolid::operator=(rhs);
00153 
00154    // Copy data
00155    //
00156    fpPolyhedron = 0; kRadTolerance = rhs.kRadTolerance;
00157    fCubicVolume = rhs.fCubicVolume; fSurfaceArea = rhs.fSurfaceArea;
00158    xSemiAxis = rhs.xSemiAxis; ySemiAxis = rhs.ySemiAxis;
00159    zSemiAxis = rhs.zSemiAxis; semiAxisMax = rhs.semiAxisMax;
00160    zBottomCut = rhs.zBottomCut; zTopCut = rhs.zTopCut;
00161 
00162    return *this;
00163 }
00164 
00166 //
00167 // Calculate extent under transform and specified limit
00168 
00169 G4bool
00170 G4Ellipsoid::CalculateExtent(const EAxis pAxis,
00171                              const G4VoxelLimits& pVoxelLimit,
00172                              const G4AffineTransform& pTransform,
00173                                    G4double& pMin, G4double& pMax) const
00174 {
00175   if (!pTransform.IsRotated())
00176   {
00177     // Special case handling for unrotated solid ellipsoid
00178     // Compute x/y/z mins and maxs for bounding box respecting limits,
00179     // with early returns if outside limits. Then switch() on pAxis,
00180     // and compute exact x and y limit for x/y case
00181 
00182     G4double xoffset,xMin,xMax;
00183     G4double yoffset,yMin,yMax;
00184     G4double zoffset,zMin,zMax;
00185 
00186     G4double maxDiff,newMin,newMax;
00187     G4double xoff,yoff;
00188 
00189     xoffset=pTransform.NetTranslation().x();
00190     xMin=xoffset - xSemiAxis;
00191     xMax=xoffset + xSemiAxis;
00192     if (pVoxelLimit.IsXLimited())
00193     {
00194       if ( (xMin>pVoxelLimit.GetMaxXExtent()+kCarTolerance)
00195         || (xMax<pVoxelLimit.GetMinXExtent()-kCarTolerance) )
00196       {
00197         return false;
00198       }
00199       else
00200       {
00201         if (xMin<pVoxelLimit.GetMinXExtent())
00202         {
00203           xMin=pVoxelLimit.GetMinXExtent();
00204         }
00205         if (xMax>pVoxelLimit.GetMaxXExtent())
00206         {
00207           xMax=pVoxelLimit.GetMaxXExtent();
00208         }
00209       }
00210     }
00211 
00212     yoffset=pTransform.NetTranslation().y();
00213     yMin=yoffset - ySemiAxis;
00214     yMax=yoffset + ySemiAxis;
00215     if (pVoxelLimit.IsYLimited())
00216     {
00217       if ( (yMin>pVoxelLimit.GetMaxYExtent()+kCarTolerance)
00218         || (yMax<pVoxelLimit.GetMinYExtent()-kCarTolerance) )
00219       {
00220         return false;
00221       }
00222       else
00223       {
00224         if (yMin<pVoxelLimit.GetMinYExtent())
00225         {
00226           yMin=pVoxelLimit.GetMinYExtent();
00227         }
00228         if (yMax>pVoxelLimit.GetMaxYExtent())
00229         {
00230           yMax=pVoxelLimit.GetMaxYExtent();
00231         }
00232       }
00233     }
00234 
00235     zoffset=pTransform.NetTranslation().z();
00236     zMin=zoffset + (-zSemiAxis > zBottomCut ? -zSemiAxis : zBottomCut);
00237     zMax=zoffset + ( zSemiAxis < zTopCut ? zSemiAxis : zTopCut);
00238     if (pVoxelLimit.IsZLimited())
00239     {
00240       if ( (zMin>pVoxelLimit.GetMaxZExtent()+kCarTolerance)
00241         || (zMax<pVoxelLimit.GetMinZExtent()-kCarTolerance) )
00242       {
00243         return false;
00244       }
00245       else
00246       {
00247         if (zMin<pVoxelLimit.GetMinZExtent())
00248         {
00249           zMin=pVoxelLimit.GetMinZExtent();
00250         }
00251         if (zMax>pVoxelLimit.GetMaxZExtent())
00252         {
00253           zMax=pVoxelLimit.GetMaxZExtent();
00254         }
00255       }
00256     }
00257 
00258     // if here, then known to cut bounding box around ellipsoid
00259     //
00260     xoff = (xoffset < xMin) ? (xMin-xoffset)
00261          : (xoffset > xMax) ? (xoffset-xMax) : 0.0;
00262     yoff = (yoffset < yMin) ? (yMin-yoffset)
00263          : (yoffset > yMax) ? (yoffset-yMax) : 0.0;
00264 
00265     // detailed calculations
00266     // NOTE: does not use X or Y offsets to adjust Z range,
00267     // and does not use Z offset to adjust X or Y range,
00268     // which is consistent with G4Sphere::CalculateExtent behavior
00269     //
00270     switch (pAxis)
00271     {
00272       case kXAxis:
00273         if (yoff==0.)
00274         {
00275           // YZ limits cross max/min x => no change
00276           //
00277           pMin=xMin;
00278           pMax=xMax;
00279         }
00280         else
00281         {
00282           // YZ limits don't cross max/min x => compute max delta x,
00283           // hence new mins/maxs
00284           //
00285           maxDiff= 1.0-sqr(yoff/ySemiAxis);
00286           if (maxDiff < 0.0) { return false; }
00287           maxDiff= xSemiAxis * std::sqrt(maxDiff);
00288           newMin=xoffset-maxDiff;
00289           newMax=xoffset+maxDiff;
00290           pMin=(newMin<xMin) ? xMin : newMin;
00291           pMax=(newMax>xMax) ? xMax : newMax;
00292         }
00293         break;
00294       case kYAxis:
00295         if (xoff==0.)
00296         {
00297           // XZ limits cross max/min y => no change
00298           //
00299           pMin=yMin;
00300           pMax=yMax;
00301         }
00302         else
00303         {
00304           // XZ limits don't cross max/min y => compute max delta y,
00305           // hence new mins/maxs
00306           //
00307           maxDiff= 1.0-sqr(xoff/xSemiAxis);
00308           if (maxDiff < 0.0) { return false; }
00309           maxDiff= ySemiAxis * std::sqrt(maxDiff);
00310           newMin=yoffset-maxDiff;
00311           newMax=yoffset+maxDiff;
00312           pMin=(newMin<yMin) ? yMin : newMin;
00313           pMax=(newMax>yMax) ? yMax : newMax;
00314         }
00315         break;
00316       case kZAxis:
00317         pMin=zMin;
00318         pMax=zMax;
00319         break;
00320       default:
00321         break;
00322     }
00323   
00324     pMin-=kCarTolerance;
00325     pMax+=kCarTolerance;
00326     return true;
00327   }
00328   else  // not rotated
00329   {
00330     G4int i,j,noEntries,noBetweenSections;
00331     G4bool existsAfterClip=false;
00332 
00333     // Calculate rotated vertex coordinates
00334 
00335     G4int noPolygonVertices=0;
00336     G4ThreeVectorList* vertices =
00337       CreateRotatedVertices(pTransform,noPolygonVertices);
00338 
00339     pMin=+kInfinity;
00340     pMax=-kInfinity;
00341 
00342     noEntries=vertices->size(); // noPolygonVertices*noPhiCrossSections
00343     noBetweenSections=noEntries-noPolygonVertices;
00344     
00345     G4ThreeVectorList ThetaPolygon;
00346     for (i=0;i<noEntries;i+=noPolygonVertices)
00347     {
00348       for(j=0;j<(noPolygonVertices/2)-1;j++)
00349       {
00350         ThetaPolygon.push_back((*vertices)[i+j]);  
00351         ThetaPolygon.push_back((*vertices)[i+j+1]);  
00352         ThetaPolygon.push_back((*vertices)[i+noPolygonVertices-2-j]);
00353         ThetaPolygon.push_back((*vertices)[i+noPolygonVertices-1-j]);
00354         CalculateClippedPolygonExtent(ThetaPolygon,pVoxelLimit,pAxis,pMin,pMax);
00355         ThetaPolygon.clear();
00356       }
00357     }
00358     for (i=0;i<noBetweenSections;i+=noPolygonVertices)
00359     {
00360       for(j=0;j<noPolygonVertices-1;j++)
00361       {
00362         ThetaPolygon.push_back((*vertices)[i+j]);  
00363         ThetaPolygon.push_back((*vertices)[i+j+1]);  
00364         ThetaPolygon.push_back((*vertices)[i+noPolygonVertices+j+1]);
00365         ThetaPolygon.push_back((*vertices)[i+noPolygonVertices+j]);
00366         CalculateClippedPolygonExtent(ThetaPolygon,pVoxelLimit,pAxis,pMin,pMax);
00367         ThetaPolygon.clear();
00368       }
00369       ThetaPolygon.push_back((*vertices)[i+noPolygonVertices-1]);
00370       ThetaPolygon.push_back((*vertices)[i]);
00371       ThetaPolygon.push_back((*vertices)[i+noPolygonVertices]);
00372       ThetaPolygon.push_back((*vertices)[i+2*noPolygonVertices-1]);
00373       CalculateClippedPolygonExtent(ThetaPolygon,pVoxelLimit,pAxis,pMin,pMax);
00374       ThetaPolygon.clear();
00375     }
00376     if ( (pMin!=kInfinity) || (pMax!=-kInfinity) )
00377     {
00378       existsAfterClip=true;
00379     
00380       // Add 2*tolerance to avoid precision troubles
00381       //
00382       pMin-=kCarTolerance;
00383       pMax+=kCarTolerance;
00384 
00385     }
00386     else
00387     {
00388       // Check for case where completely enveloping clipping volume
00389       // If point inside then we are confident that the solid completely
00390       // envelopes the clipping volume. Hence set min/max extents according
00391       // to clipping volume extents along the specified axis.
00392       //
00393       G4ThreeVector
00394       clipCentre((pVoxelLimit.GetMinXExtent()+pVoxelLimit.GetMaxXExtent())*0.5,
00395                  (pVoxelLimit.GetMinYExtent()+pVoxelLimit.GetMaxYExtent())*0.5,
00396                  (pVoxelLimit.GetMinZExtent()+pVoxelLimit.GetMaxZExtent())*0.5);
00397 
00398       if (Inside(pTransform.Inverse().TransformPoint(clipCentre))!=kOutside)
00399       {
00400         existsAfterClip=true;
00401         pMin=pVoxelLimit.GetMinExtent(pAxis);
00402         pMax=pVoxelLimit.GetMaxExtent(pAxis);
00403       }
00404     }
00405     delete vertices;
00406     return existsAfterClip;
00407   }
00408 }
00409 
00411 //
00412 // Return whether point inside/outside/on surface
00413 // Split into radius, phi, theta checks
00414 // Each check modifies `in', or returns as approprate
00415 
00416 EInside G4Ellipsoid::Inside(const G4ThreeVector& p) const
00417 {
00418   G4double rad2oo,  // outside surface outer tolerance
00419            rad2oi;  // outside surface inner tolerance
00420   EInside in;
00421 
00422   static const G4double halfRadTolerance=kRadTolerance*0.5;
00423 
00424   // check this side of z cut first, because that's fast
00425   //
00426   if (p.z() < zBottomCut-halfRadTolerance) { return in=kOutside; }
00427   if (p.z() > zTopCut+halfRadTolerance)    { return in=kOutside; }
00428 
00429   rad2oo= sqr(p.x()/(xSemiAxis+halfRadTolerance))
00430         + sqr(p.y()/(ySemiAxis+halfRadTolerance))
00431         + sqr(p.z()/(zSemiAxis+halfRadTolerance));
00432 
00433   if (rad2oo > 1.0)  { return in=kOutside; }
00434     
00435   rad2oi= sqr(p.x()*(1.0+halfRadTolerance/xSemiAxis)/xSemiAxis)
00436       + sqr(p.y()*(1.0+halfRadTolerance/ySemiAxis)/ySemiAxis)
00437       + sqr(p.z()*(1.0+halfRadTolerance/zSemiAxis)/zSemiAxis);
00438 
00439   // Check radial surfaces
00440   //  sets `in' (already checked for rad2oo > 1.0)
00441   //
00442   if (rad2oi < 1.0)
00443   {
00444     in = ( (p.z() < zBottomCut+halfRadTolerance)
00445         || (p.z() > zTopCut-halfRadTolerance) ) ? kSurface : kInside;
00446     if ( rad2oi > 1.0-halfRadTolerance )  { in=kSurface; }
00447   }
00448   else 
00449   {
00450     in = kSurface;
00451   }
00452   return in;
00453 
00454 }
00455 
00457 //
00458 // Return unit normal of surface closest to p not protected against p=0
00459 
00460 G4ThreeVector G4Ellipsoid::SurfaceNormal( const G4ThreeVector& p) const
00461 {
00462   G4double distR, distZBottom, distZTop;
00463 
00464   // normal vector with special magnitude:  parallel to normal, units 1/length
00465   // norm*p == 1.0 if on surface, >1.0 if outside, <1.0 if inside
00466   //
00467   G4ThreeVector norm(p.x()/(xSemiAxis*xSemiAxis),
00468                      p.y()/(ySemiAxis*ySemiAxis),
00469                      p.z()/(zSemiAxis*zSemiAxis));
00470   G4double radius = 1.0/norm.mag();
00471 
00472   // approximate distance to curved surface
00473   //
00474   distR = std::fabs( (p*norm - 1.0) * radius ) / 2.0;
00475 
00476   // Distance to z-cut plane
00477   //
00478   distZBottom = std::fabs( p.z() - zBottomCut );
00479   distZTop = std::fabs( p.z() - zTopCut );
00480 
00481   if ( (distZBottom < distR) || (distZTop < distR) )
00482   {
00483     return G4ThreeVector(0.,0.,(distZBottom < distZTop) ? -1.0 : 1.0);
00484   }
00485   return ( norm *= radius );
00486 }
00487 
00489 //
00490 // Calculate distance to shape from outside, along normalised vector
00491 // - return kInfinity if no intersection, or intersection distance <= tolerance
00492 //
00493 
00494 G4double G4Ellipsoid::DistanceToIn( const G4ThreeVector& p,
00495                                     const G4ThreeVector& v  ) const
00496 {
00497   static const G4double halfCarTolerance=kCarTolerance*0.5;
00498   static const G4double halfRadTolerance=kRadTolerance*0.5;
00499 
00500   G4double distMin = std::min(xSemiAxis,ySemiAxis);
00501   const G4double dRmax = 100.*std::min(distMin,zSemiAxis);
00502   distMin= kInfinity;
00503 
00504   // check to see if Z plane is relevant
00505   if (p.z() <= zBottomCut+halfCarTolerance)
00506   {
00507     if (v.z() <= 0.0) { return distMin; }
00508     G4double distZ = (zBottomCut - p.z()) / v.z();
00509 
00510     if ( (distZ > -halfRadTolerance) && (Inside(p+distZ*v) != kOutside) )
00511     {
00512       // early exit since can't intercept curved surface if we reach here
00513       if ( std::fabs(distZ) < halfRadTolerance ) { distZ=0.; }
00514       return distMin= distZ;
00515     }
00516   }
00517   if (p.z() >= zTopCut-halfCarTolerance)
00518   {
00519     if (v.z() >= 0.0) { return distMin;}
00520     G4double distZ = (zTopCut - p.z()) / v.z();
00521     if ( (distZ > -halfRadTolerance) && (Inside(p+distZ*v) != kOutside) )
00522     {
00523       // early exit since can't intercept curved surface if we reach here
00524       if ( std::fabs(distZ) < halfRadTolerance ) { distZ=0.; }
00525       return distMin= distZ;
00526     }
00527   }
00528   // if fZCut1 <= p.z() <= fZCut2, then must hit curved surface
00529 
00530   // now check curved surface intercept
00531   G4double A,B,C;
00532 
00533   A= sqr(v.x()/xSemiAxis) + sqr(v.y()/ySemiAxis) + sqr(v.z()/zSemiAxis);
00534   C= sqr(p.x()/xSemiAxis) + sqr(p.y()/ySemiAxis) + sqr(p.z()/zSemiAxis) - 1.0;
00535   B= 2.0 * ( p.x()*v.x()/(xSemiAxis*xSemiAxis)
00536            + p.y()*v.y()/(ySemiAxis*ySemiAxis)
00537            + p.z()*v.z()/(zSemiAxis*zSemiAxis) );
00538 
00539   C= B*B - 4.0*A*C;
00540   if (C > 0.0)
00541   {    
00542     G4double distR= (-B - std::sqrt(C)) / (2.0*A);
00543     G4double intZ = p.z()+distR*v.z();
00544     if ( (distR > halfRadTolerance)
00545       && (intZ >= zBottomCut-halfRadTolerance)
00546       && (intZ <= zTopCut+halfRadTolerance) )
00547     { 
00548       distMin = distR;
00549     }
00550     else if( (distR >- halfRadTolerance)
00551             && (intZ >= zBottomCut-halfRadTolerance)
00552             && (intZ <= zTopCut+halfRadTolerance) )
00553     {
00554       // p is on the curved surface, DistanceToIn returns 0 or kInfinity:
00555       // DistanceToIn returns 0, if second root is positive (means going inside)
00556       // If second root is negative, DistanceToIn returns kInfinity (outside)
00557       //
00558       distR = (-B + std::sqrt(C) ) / (2.0*A);
00559       if(distR>0.) { distMin=0.; }
00560     }
00561     else
00562     {
00563       distR= (-B + std::sqrt(C)) / (2.0*A);
00564       intZ = p.z()+distR*v.z();
00565       if ( (distR > halfRadTolerance)
00566         && (intZ >= zBottomCut-halfRadTolerance)
00567         && (intZ <= zTopCut+halfRadTolerance) )
00568       {
00569         G4ThreeVector norm=SurfaceNormal(p);
00570         if (norm.dot(v)<0.) { distMin = distR; }
00571       }
00572     }
00573     if ( (distMin!=kInfinity) && (distMin>dRmax) ) 
00574     {                    // Avoid rounding errors due to precision issues on
00575                          // 64 bits systems. Split long distances and recompute
00576       G4double fTerm = distMin-std::fmod(distMin,dRmax);
00577       distMin = fTerm + DistanceToIn(p+fTerm*v,v);
00578     }
00579   }
00580   
00581   if (std::fabs(distMin)<halfRadTolerance) { distMin=0.; }
00582   return distMin;
00583 } 
00584 
00586 //
00587 // Calculate distance (<= actual) to closest surface of shape from outside
00588 // - Return 0 if point inside
00589 
00590 G4double G4Ellipsoid::DistanceToIn(const G4ThreeVector& p) const
00591 {
00592   G4double distR, distZ;
00593 
00594   // normal vector:  parallel to normal, magnitude 1/(characteristic radius)
00595   //
00596   G4ThreeVector norm(p.x()/(xSemiAxis*xSemiAxis),
00597                      p.y()/(ySemiAxis*ySemiAxis),
00598                      p.z()/(zSemiAxis*zSemiAxis));
00599   G4double radius= 1.0/norm.mag();
00600 
00601   // approximate distance to curved surface ( <= actual distance )
00602   //
00603   distR= (p*norm - 1.0) * radius / 2.0;
00604 
00605   // Distance to z-cut plane
00606   //
00607   distZ= zBottomCut - p.z();
00608   if (distZ < 0.0)
00609   {
00610     distZ = p.z() - zTopCut;
00611   }
00612 
00613   // Distance to closest surface from outside
00614   //
00615   if (distZ < 0.0)
00616   {
00617     return (distR < 0.0) ? 0.0 : distR;
00618   }
00619   else if (distR < 0.0)
00620   {
00621     return distZ;
00622   }
00623   else
00624   {
00625     return (distZ < distR) ? distZ : distR;
00626   }
00627 }
00628 
00630 //
00631 // Calculate distance to surface of shape from `inside', allowing for tolerance
00632 
00633 G4double G4Ellipsoid::DistanceToOut(const G4ThreeVector& p,
00634                                     const G4ThreeVector& v,
00635                                     const G4bool calcNorm,
00636                                           G4bool *validNorm,
00637                                           G4ThreeVector *n  ) const
00638 {
00639   G4double distMin;
00640   enum surface_e {kPlaneSurf, kCurvedSurf, kNoSurf} surface;
00641   
00642   distMin= kInfinity;
00643   surface= kNoSurf;
00644 
00645   // check to see if Z plane is relevant
00646   //
00647   if (v.z() < 0.0)
00648   {
00649     G4double distZ = (zBottomCut - p.z()) / v.z();
00650     if (distZ < 0.0)
00651     {
00652       distZ= 0.0;
00653       if (!calcNorm) {return 0.0;}
00654     }
00655     distMin= distZ;
00656     surface= kPlaneSurf;
00657   }
00658   if (v.z() > 0.0)
00659   {
00660     G4double distZ = (zTopCut - p.z()) / v.z();
00661     if (distZ < 0.0)
00662     {
00663       distZ= 0.0;
00664       if (!calcNorm) {return 0.0;}
00665     }
00666     distMin= distZ;
00667     surface= kPlaneSurf;
00668   }
00669 
00670   // normal vector:  parallel to normal, magnitude 1/(characteristic radius)
00671   //
00672   G4ThreeVector nearnorm(p.x()/(xSemiAxis*xSemiAxis),
00673                          p.y()/(ySemiAxis*ySemiAxis),
00674                          p.z()/(zSemiAxis*zSemiAxis));
00675   
00676   // now check curved surface intercept
00677   //
00678   G4double A,B,C;
00679   
00680   A= sqr(v.x()/xSemiAxis) + sqr(v.y()/ySemiAxis) + sqr(v.z()/zSemiAxis);
00681   C= (p * nearnorm) - 1.0;
00682   B= 2.0 * (v * nearnorm);
00683 
00684   C= B*B - 4.0*A*C;
00685   if (C > 0.0)
00686   {
00687     G4double distR= (-B + std::sqrt(C) ) / (2.0*A);
00688     if (distR < 0.0)
00689     {
00690       distR= 0.0;
00691       if (!calcNorm) {return 0.0;}
00692     }
00693     if (distR < distMin)
00694     {
00695       distMin= distR;
00696       surface= kCurvedSurf;
00697     }
00698   }
00699 
00700   // set normal if requested
00701   //
00702   if (calcNorm)
00703   {
00704     if (surface == kNoSurf)
00705     {
00706       *validNorm = false;
00707     }
00708     else
00709     {
00710       *validNorm = true;
00711       switch (surface)
00712       {
00713         case kPlaneSurf:
00714           *n= G4ThreeVector(0.,0.,(v.z() > 0.0 ? 1. : -1.));
00715           break;
00716         case kCurvedSurf:
00717         {
00718           G4ThreeVector pexit= p + distMin*v;
00719           G4ThreeVector truenorm(pexit.x()/(xSemiAxis*xSemiAxis),
00720                                  pexit.y()/(ySemiAxis*ySemiAxis),
00721                                  pexit.z()/(zSemiAxis*zSemiAxis));
00722           truenorm *= 1.0/truenorm.mag();
00723           *n= truenorm;
00724         } break;
00725         default:           // Should never reach this case ...
00726           DumpInfo();
00727           std::ostringstream message;
00728           G4int oldprc = message.precision(16);
00729           message << "Undefined side for valid surface normal to solid."
00730                   << G4endl
00731                   << "Position:"  << G4endl
00732                   << "   p.x() = "   << p.x()/mm << " mm" << G4endl
00733                   << "   p.y() = "   << p.y()/mm << " mm" << G4endl
00734                   << "   p.z() = "   << p.z()/mm << " mm" << G4endl
00735                   << "Direction:" << G4endl << G4endl
00736                   << "   v.x() = "   << v.x() << G4endl
00737                   << "   v.y() = "   << v.y() << G4endl
00738                   << "   v.z() = "   << v.z() << G4endl
00739                   << "Proposed distance :" << G4endl
00740                   << "   distMin = "    << distMin/mm << " mm";
00741           message.precision(oldprc);
00742           G4Exception("G4Ellipsoid::DistanceToOut(p,v,..)",
00743                       "GeomSolids1002", JustWarning, message);
00744           break;
00745       }
00746     }
00747   }
00748    
00749   return distMin;
00750 }
00751 
00753 //
00754 // Calculate distance (<=actual) to closest surface of shape from inside
00755 
00756 G4double G4Ellipsoid::DistanceToOut(const G4ThreeVector& p) const
00757 {
00758   G4double distR, distZ;
00759 
00760 #ifdef G4SPECSDEBUG
00761   if( Inside(p) == kOutside )
00762   {
00763      DumpInfo();
00764      std::ostringstream message;
00765      G4int oldprc = message.precision(16);
00766      message << "Point p is outside !?" << G4endl
00767              << "Position:"  << G4endl
00768              << "   p.x() = "   << p.x()/mm << " mm" << G4endl
00769              << "   p.y() = "   << p.y()/mm << " mm" << G4endl
00770              << "   p.z() = "   << p.z()/mm << " mm";
00771      message.precision(oldprc) ;
00772      G4Exception("G4Ellipsoid::DistanceToOut(p)", "GeomSolids1002",
00773                  JustWarning, message);
00774   }
00775 #endif
00776 
00777   // Normal vector:  parallel to normal, magnitude 1/(characteristic radius)
00778   //
00779   G4ThreeVector norm(p.x()/(xSemiAxis*xSemiAxis),
00780                      p.y()/(ySemiAxis*ySemiAxis),
00781                      p.z()/(zSemiAxis*zSemiAxis));
00782 
00783   // the following is a safe inlined "radius= min(1.0/norm.mag(),p.mag())
00784   //
00785   G4double radius= p.mag();
00786   G4double tmp= norm.mag();
00787   if ( (tmp > 0.0) && (1.0 < radius*tmp) ) {radius = 1.0/tmp;}
00788 
00789   // Approximate distance to curved surface ( <= actual distance )
00790   //
00791   distR = (1.0 - p*norm) * radius / 2.0;
00792     
00793   // Distance to z-cut plane
00794   //
00795   distZ = p.z() - zBottomCut;
00796   if (distZ < 0.0) {distZ= zTopCut - p.z();}
00797 
00798   // Distance to closest surface from inside
00799   //
00800   if ( (distZ < 0.0) || (distR < 0.0) )
00801   {
00802     return 0.0;
00803   }
00804   else
00805   {
00806     return (distZ < distR) ? distZ : distR;
00807   }
00808 }
00809 
00811 //
00812 // Create a List containing the transformed vertices
00813 // Ordering [0-3] -fDz cross section
00814 //          [4-7] +fDz cross section such that [0] is below [4],
00815 //                                             [1] below [5] etc.
00816 // Note:
00817 //  Caller has deletion resposibility
00818 //  Potential improvement: For last slice, use actual ending angle
00819 //                         to avoid rounding error problems.
00820 
00821 G4ThreeVectorList*
00822 G4Ellipsoid::CreateRotatedVertices(const G4AffineTransform& pTransform,
00823                                          G4int& noPolygonVertices) const
00824 {
00825   G4ThreeVectorList *vertices;
00826   G4ThreeVector vertex;
00827   G4double meshAnglePhi, meshRMaxFactor,
00828            crossAnglePhi, coscrossAnglePhi, sincrossAnglePhi, sAnglePhi;
00829   G4double meshTheta, crossTheta, startTheta;
00830   G4double rMaxX, rMaxY, rMaxZ, rMaxMax, rx, ry, rz;
00831   G4int crossSectionPhi, noPhiCrossSections, crossSectionTheta, noThetaSections;
00832 
00833   // Phi cross sections
00834   //
00835   noPhiCrossSections=G4int (twopi/kMeshAngleDefault)+1;  // = 9!
00836     
00837 /*
00838   if (noPhiCrossSections<kMinMeshSections)        // <3
00839   {
00840     noPhiCrossSections=kMinMeshSections;
00841   }
00842   else if (noPhiCrossSections>kMaxMeshSections)   // >37
00843   {
00844     noPhiCrossSections=kMaxMeshSections;
00845   }
00846 */
00847   meshAnglePhi=twopi/(noPhiCrossSections-1);
00848     
00849   // Set start angle such that mesh will be at fRMax
00850   // on the x axis. Will give better extent calculations when not rotated.
00851     
00852   sAnglePhi = -meshAnglePhi*0.5;
00853 
00854   // Theta cross sections
00855     
00856   noThetaSections = G4int(pi/kMeshAngleDefault)+3;  //  = 7!
00857 
00858 /*
00859   if (noThetaSections<kMinMeshSections)       // <3
00860   {
00861     noThetaSections=kMinMeshSections;
00862   }
00863   else if (noThetaSections>kMaxMeshSections)  // >37
00864   {
00865     noThetaSections=kMaxMeshSections;
00866   }
00867 */
00868   meshTheta= pi/(noThetaSections-2);
00869     
00870   // Set start angle such that mesh will be at fRMax
00871   // on the z axis. Will give better extent calculations when not rotated.
00872     
00873   startTheta = -meshTheta*0.5;
00874 
00875   meshRMaxFactor =  1.0/std::cos(0.5*
00876                     std::sqrt(meshAnglePhi*meshAnglePhi+meshTheta*meshTheta));
00877   rMaxMax= (xSemiAxis > ySemiAxis ? xSemiAxis : ySemiAxis);
00878   if (zSemiAxis > rMaxMax) rMaxMax= zSemiAxis;
00879   rMaxX= xSemiAxis + rMaxMax*(meshRMaxFactor-1.0);
00880   rMaxY= ySemiAxis + rMaxMax*(meshRMaxFactor-1.0);
00881   rMaxZ= zSemiAxis + rMaxMax*(meshRMaxFactor-1.0);
00882   G4double* cosCrossTheta = new G4double[noThetaSections];
00883   G4double* sinCrossTheta = new G4double[noThetaSections];    
00884   vertices=new G4ThreeVectorList(noPhiCrossSections*noThetaSections);
00885   if (vertices && cosCrossTheta && sinCrossTheta)
00886   {
00887     for (crossSectionTheta=0; crossSectionTheta<noThetaSections;
00888          crossSectionTheta++)
00889     {
00890       // Compute sine and cosine table (for historical reasons)
00891       //
00892       crossTheta=startTheta+crossSectionTheta*meshTheta;
00893       cosCrossTheta[crossSectionTheta]=std::cos(crossTheta);
00894       sinCrossTheta[crossSectionTheta]=std::sin(crossTheta);
00895     }
00896     for (crossSectionPhi=0; crossSectionPhi<noPhiCrossSections;
00897          crossSectionPhi++)
00898     {
00899       crossAnglePhi=sAnglePhi+crossSectionPhi*meshAnglePhi;
00900       coscrossAnglePhi=std::cos(crossAnglePhi);
00901       sincrossAnglePhi=std::sin(crossAnglePhi);
00902       for (crossSectionTheta=0; crossSectionTheta<noThetaSections;
00903            crossSectionTheta++)
00904       {
00905         // Compute coordinates of cross section at section crossSectionPhi
00906         //
00907         rx= sinCrossTheta[crossSectionTheta]*coscrossAnglePhi*rMaxX;
00908         ry= sinCrossTheta[crossSectionTheta]*sincrossAnglePhi*rMaxY;
00909         rz= cosCrossTheta[crossSectionTheta]*rMaxZ;
00910         if (rz < zBottomCut)
00911           { rz= zBottomCut; }
00912         if (rz > zTopCut)
00913           { rz= zTopCut; }
00914         vertex= G4ThreeVector(rx,ry,rz);
00915         vertices->push_back(pTransform.TransformPoint(vertex));
00916       }    // Theta forward     
00917     }    // Phi
00918     noPolygonVertices = noThetaSections ;
00919   }
00920   else
00921   {
00922     DumpInfo();
00923     G4Exception("G4Ellipsoid::CreateRotatedVertices()",
00924                 "GeomSolids0003", FatalException,
00925                 "Error in allocation of vertices. Out of memory !");
00926   }
00927 
00928   delete[] cosCrossTheta;
00929   delete[] sinCrossTheta;
00930 
00931   return vertices;
00932 }
00933 
00935 //
00936 // G4EntityType
00937 
00938 G4GeometryType G4Ellipsoid::GetEntityType() const
00939 {
00940   return G4String("G4Ellipsoid");
00941 }
00942 
00944 //
00945 // Make a clone of the object
00946 
00947 G4VSolid* G4Ellipsoid::Clone() const
00948 {
00949   return new G4Ellipsoid(*this);
00950 }
00951 
00953 //
00954 // Stream object contents to an output stream
00955 
00956 std::ostream& G4Ellipsoid::StreamInfo( std::ostream& os ) const
00957 {
00958   G4int oldprc = os.precision(16);
00959   os << "-----------------------------------------------------------\n"
00960      << "    *** Dump for solid - " << GetName() << " ***\n"
00961      << "    ===================================================\n"
00962      << " Solid type: G4Ellipsoid\n"
00963      << " Parameters: \n"
00964 
00965      << "    semi-axis x: " << xSemiAxis/mm << " mm \n"
00966      << "    semi-axis y: " << ySemiAxis/mm << " mm \n"
00967      << "    semi-axis z: " << zSemiAxis/mm << " mm \n"
00968      << "    max semi-axis: " << semiAxisMax/mm << " mm \n"
00969      << "    lower cut plane level z: " << zBottomCut/mm << " mm \n"
00970      << "    upper cut plane level z: " << zTopCut/mm << " mm \n"
00971      << "-----------------------------------------------------------\n";
00972   os.precision(oldprc);
00973 
00974   return os;
00975 }
00976 
00978 //
00979 // GetPointOnSurface
00980 
00981 G4ThreeVector G4Ellipsoid::GetPointOnSurface() const
00982 {
00983   G4double aTop, aBottom, aCurved, chose, xRand, yRand, zRand, phi;
00984   G4double cosphi, sinphi, costheta, sintheta, alpha, beta, max1, max2, max3;
00985 
00986   max1  = xSemiAxis > ySemiAxis ? xSemiAxis : ySemiAxis;
00987   max1  = max1 > zSemiAxis ? max1 : zSemiAxis;
00988   if (max1 == xSemiAxis)      { max2 = ySemiAxis; max3 = zSemiAxis; }
00989   else if (max1 == ySemiAxis) { max2 = xSemiAxis; max3 = zSemiAxis; }
00990   else                        { max2 = xSemiAxis; max3 = ySemiAxis; }
00991 
00992   phi   = RandFlat::shoot(0.,twopi);
00993   
00994   cosphi = std::cos(phi);   sinphi = std::sin(phi);
00995   costheta = RandFlat::shoot(zBottomCut,zTopCut)/zSemiAxis;
00996   sintheta = std::sqrt(1.-sqr(costheta));
00997   
00998   alpha = 1.-sqr(max2/max1); beta  = 1.-sqr(max3/max1);
00999   
01000   aTop    = pi*xSemiAxis*ySemiAxis*(1 - sqr(zTopCut/zSemiAxis));
01001   aBottom = pi*xSemiAxis*ySemiAxis*(1 - sqr(zBottomCut/zSemiAxis));
01002   
01003   // approximation
01004   // from:" http://www.citr.auckland.ac.nz/techreports/2004/CITR-TR-139.pdf"
01005   aCurved = 4.*pi*max1*max2*(1.-1./6.*(alpha+beta)-
01006                             1./120.*(3.*sqr(alpha)+2.*alpha*beta+3.*sqr(beta)));
01007 
01008   aCurved *= 0.5*(1.2*zTopCut/zSemiAxis - 1.2*zBottomCut/zSemiAxis);
01009   
01010   if( ( zTopCut >= zSemiAxis && zBottomCut <= -1.*zSemiAxis )
01011    || ( zTopCut == 0 && zBottomCut ==0 ) )
01012   {
01013     aTop = 0; aBottom = 0;
01014   }
01015   
01016   chose = RandFlat::shoot(0.,aTop + aBottom + aCurved); 
01017   
01018   if(chose < aCurved)
01019   { 
01020     xRand = xSemiAxis*sintheta*cosphi;
01021     yRand = ySemiAxis*sintheta*sinphi;
01022     zRand = zSemiAxis*costheta;
01023     return G4ThreeVector (xRand,yRand,zRand); 
01024   }
01025   else if(chose >= aCurved && chose < aCurved + aTop)
01026   {
01027     xRand = RandFlat::shoot(-1.,1.)*xSemiAxis
01028           * std::sqrt(1-sqr(zTopCut/zSemiAxis));
01029     yRand = RandFlat::shoot(-1.,1.)*ySemiAxis
01030           * std::sqrt(1.-sqr(zTopCut/zSemiAxis)-sqr(xRand/xSemiAxis));
01031     zRand = zTopCut;
01032     return G4ThreeVector (xRand,yRand,zRand);
01033   }
01034   else
01035   {
01036     xRand = RandFlat::shoot(-1.,1.)*xSemiAxis
01037           * std::sqrt(1-sqr(zBottomCut/zSemiAxis));
01038     yRand = RandFlat::shoot(-1.,1.)*ySemiAxis
01039           * std::sqrt(1.-sqr(zBottomCut/zSemiAxis)-sqr(xRand/xSemiAxis)); 
01040     zRand = zBottomCut;
01041     return G4ThreeVector (xRand,yRand,zRand);
01042   }
01043 }
01044 
01046 //
01047 // Methods for visualisation
01048 
01049 void G4Ellipsoid::DescribeYourselfTo (G4VGraphicsScene& scene) const
01050 {
01051   scene.AddSolid(*this);
01052 }
01053 
01054 G4VisExtent G4Ellipsoid::GetExtent() const
01055 {
01056   // Define the sides of the box into which the G4Ellipsoid instance would fit.
01057   //
01058   return G4VisExtent (-semiAxisMax, semiAxisMax,
01059                       -semiAxisMax, semiAxisMax,
01060                       -semiAxisMax, semiAxisMax);
01061 }
01062 
01063 G4NURBS* G4Ellipsoid::CreateNURBS () const
01064 {
01065   // Box for now!!!
01066   //
01067   return new G4NURBSbox(semiAxisMax, semiAxisMax, semiAxisMax);
01068 }
01069 
01070 G4Polyhedron* G4Ellipsoid::CreatePolyhedron () const
01071 {
01072   return new G4PolyhedronEllipsoid(xSemiAxis, ySemiAxis, zSemiAxis,
01073                                    zBottomCut, zTopCut);
01074 }
01075 
01076 G4Polyhedron* G4Ellipsoid::GetPolyhedron () const
01077 {
01078   if (!fpPolyhedron ||
01079       fpPolyhedron->GetNumberOfRotationStepsAtTimeOfCreation() !=
01080       fpPolyhedron->GetNumberOfRotationSteps())
01081     {
01082       delete fpPolyhedron;
01083       fpPolyhedron = CreatePolyhedron();
01084     }
01085   return fpPolyhedron;
01086 }

Generated on Mon May 27 17:48:07 2013 for Geant4 by  doxygen 1.4.7