G4Tet.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  intellectual property  of the *
00019 // * Vanderbilt University Free Electron Laser Center                 *
00020 // * Vanderbilt University, Nashville, TN, USA                        *
00021 // * Development supported by:                                        *
00022 // * United States MFEL program  under grant FA9550-04-1-0045         *
00023 // * and NASA under contract number NNG04CT05P                        *
00024 // * Written by Marcus H. Mendenhall and Robert A. Weller.            *
00025 // *                                                                  *
00026 // * Contributed to the Geant4 Core, January, 2005.                   *
00027 // *                                                                  *
00028 // ********************************************************************
00029 //
00030 // $Id: G4Tet.cc 69790 2013-05-15 12:39:10Z gcosmo $
00031 //
00032 // class G4Tet
00033 //
00034 // Implementation for G4Tet class
00035 //
00036 // History:
00037 //
00038 //  20040903 - Marcus Mendenhall, created G4Tet
00039 //  20041101 - Marcus Mendenhall, optimized constant dot products with
00040 //             fCdotNijk values
00041 //  20041101 - MHM removed tracking error by clipping DistanceToOut to 0
00042 //             for surface cases
00043 //  20041101 - MHM many speed optimizations in if statements
00044 //  20041101 - MHM changed vdotn comparisons to 1e-12 instead of 0.0 to
00045 //             avoid nearly-parallel problems
00046 //  20041102 - MHM Added extra distance into solid to DistanceToIn(p,v)
00047 //             hit testing
00048 //  20041102 - MHM added ability to check for degeneracy without throwing
00049 //             G4Exception
00050 //  20041103 - MHM removed many unused variables from class
00051 //  20040803 - Dionysios Anninos, added GetPointOnSurface() method
00052 //  20061112 - MHM added code for G4VSolid GetSurfaceArea()
00053 //  20100920 - Gabriele Cosmo added copy-ctor and operator=()
00054 //
00055 // --------------------------------------------------------------------
00056 
00057 #include "G4Tet.hh"
00058 
00059 const char G4Tet::CVSVers[]="$Id: G4Tet.cc 69790 2013-05-15 12:39:10Z gcosmo $";
00060 
00061 #include "G4VoxelLimits.hh"
00062 #include "G4AffineTransform.hh"
00063 
00064 #include "G4VPVParameterisation.hh"
00065 
00066 #include "Randomize.hh"
00067 
00068 #include "G4VGraphicsScene.hh"
00069 #include "G4Polyhedron.hh"
00070 #include "G4NURBS.hh"
00071 #include "G4NURBSbox.hh"
00072 #include "G4VisExtent.hh"
00073 
00074 #include "G4ThreeVector.hh"
00075 
00076 #include <cmath>
00077 
00078 using namespace CLHEP;
00079 
00081 //
00082 // Constructor - create a tetrahedron
00083 // This class is implemented separately from general polyhedra,
00084 // because the simplex geometry can be computed very quickly,
00085 // which may become important in situations imported from mesh generators,
00086 // in which a very large number of G4Tets are created.
00087 // A Tet has all of its geometrical information precomputed
00088 
00089 G4Tet::G4Tet(const G4String& pName,
00090                    G4ThreeVector anchor,
00091                    G4ThreeVector p2,
00092                    G4ThreeVector p3,
00093                    G4ThreeVector p4, G4bool *degeneracyFlag)
00094   : G4VSolid(pName), fpPolyhedron(0), warningFlag(0)
00095 {
00096   // fV<x><y> is vector from vertex <y> to vertex <x>
00097   //
00098   G4ThreeVector fV21=p2-anchor;
00099   G4ThreeVector fV31=p3-anchor;
00100   G4ThreeVector fV41=p4-anchor;
00101 
00102   // make sure this is a correctly oriented set of points for the tetrahedron
00103   //
00104   G4double signed_vol=fV21.cross(fV31).dot(fV41);
00105   if(signed_vol<0.0)
00106   {
00107     G4ThreeVector temp(p4);
00108     p4=p3;
00109     p3=temp;
00110     temp=fV41;
00111     fV41=fV31;
00112     fV31=temp; 
00113   }
00114   fCubicVolume = std::fabs(signed_vol) / 6.;
00115 
00116   G4ThreeVector fV24=p2-p4;
00117   G4ThreeVector fV43=p4-p3;
00118   G4ThreeVector fV32=p3-p2;
00119 
00120   fXMin=std::min(std::min(std::min(anchor.x(), p2.x()),p3.x()),p4.x());
00121   fXMax=std::max(std::max(std::max(anchor.x(), p2.x()),p3.x()),p4.x());
00122   fYMin=std::min(std::min(std::min(anchor.y(), p2.y()),p3.y()),p4.y());
00123   fYMax=std::max(std::max(std::max(anchor.y(), p2.y()),p3.y()),p4.y());
00124   fZMin=std::min(std::min(std::min(anchor.z(), p2.z()),p3.z()),p4.z());
00125   fZMax=std::max(std::max(std::max(anchor.z(), p2.z()),p3.z()),p4.z());
00126 
00127   fDx=(fXMax-fXMin)*0.5; fDy=(fYMax-fYMin)*0.5; fDz=(fZMax-fZMin)*0.5;
00128 
00129   fMiddle=G4ThreeVector(fXMax+fXMin, fYMax+fYMin, fZMax+fZMin)*0.5;
00130   fMaxSize=std::max(std::max(std::max((anchor-fMiddle).mag(),
00131                                       (p2-fMiddle).mag()),
00132                              (p3-fMiddle).mag()),
00133                     (p4-fMiddle).mag());
00134 
00135   G4bool degenerate=std::fabs(signed_vol) < 1e-9*fMaxSize*fMaxSize*fMaxSize;
00136 
00137   if(degeneracyFlag) *degeneracyFlag=degenerate;
00138   else if (degenerate)
00139   {
00140     G4Exception("G4Tet::G4Tet()", "GeomSolids0002", FatalException,
00141                 "Degenerate tetrahedron not allowed.");
00142   }
00143 
00144   fTol=1e-9*(std::fabs(fXMin)+std::fabs(fXMax)+std::fabs(fYMin)
00145             +std::fabs(fYMax)+std::fabs(fZMin)+std::fabs(fZMax));
00146   //fTol=kCarTolerance;
00147 
00148   fAnchor=anchor;
00149   fP2=p2;
00150   fP3=p3;
00151   fP4=p4;
00152 
00153   G4ThreeVector fCenter123=(anchor+p2+p3)*(1.0/3.0); // face center
00154   G4ThreeVector fCenter134=(anchor+p4+p3)*(1.0/3.0);
00155   G4ThreeVector fCenter142=(anchor+p4+p2)*(1.0/3.0);
00156   G4ThreeVector fCenter234=(p2+p3+p4)*(1.0/3.0);
00157 
00158   // compute area of each triangular face by cross product
00159   // and sum for total surface area
00160 
00161   G4ThreeVector normal123=fV31.cross(fV21);
00162   G4ThreeVector normal134=fV41.cross(fV31);
00163   G4ThreeVector normal142=fV21.cross(fV41);
00164   G4ThreeVector normal234=fV32.cross(fV43);
00165 
00166   fSurfaceArea=(
00167       normal123.mag()+
00168       normal134.mag()+
00169       normal142.mag()+
00170       normal234.mag()
00171   )/2.0;
00172 
00173   fNormal123=normal123.unit();
00174   fNormal134=normal134.unit();
00175   fNormal142=normal142.unit();
00176   fNormal234=normal234.unit();
00177 
00178   fCdotN123=fCenter123.dot(fNormal123);
00179   fCdotN134=fCenter134.dot(fNormal134);
00180   fCdotN142=fCenter142.dot(fNormal142);
00181   fCdotN234=fCenter234.dot(fNormal234);
00182 }
00183 
00185 //
00186 // Fake default constructor - sets only member data and allocates memory
00187 //                            for usage restricted to object persistency.
00188 //
00189 G4Tet::G4Tet( __void__& a )
00190   : G4VSolid(a), fCubicVolume(0.), fSurfaceArea(0.), fpPolyhedron(0),
00191     fAnchor(0,0,0), fP2(0,0,0), fP3(0,0,0), fP4(0,0,0), fMiddle(0,0,0),
00192     fNormal123(0,0,0), fNormal142(0,0,0), fNormal134(0,0,0),
00193     fNormal234(0,0,0), warningFlag(0),
00194     fCdotN123(0.), fCdotN142(0.), fCdotN134(0.), fCdotN234(0.),
00195     fXMin(0.), fXMax(0.), fYMin(0.), fYMax(0.), fZMin(0.), fZMax(0.),
00196     fDx(0.), fDy(0.), fDz(0.), fTol(0.), fMaxSize(0.)
00197 {
00198 }
00199 
00201 //
00202 // Destructor
00203 
00204 G4Tet::~G4Tet()
00205 {
00206   delete fpPolyhedron;
00207 }
00208 
00210 //
00211 // Copy constructor
00212 
00213 G4Tet::G4Tet(const G4Tet& rhs)
00214   : G4VSolid(rhs),
00215     fCubicVolume(rhs.fCubicVolume), fSurfaceArea(rhs.fSurfaceArea),
00216     fpPolyhedron(0), fAnchor(rhs.fAnchor),
00217     fP2(rhs.fP2), fP3(rhs.fP3), fP4(rhs.fP4), fMiddle(rhs.fMiddle),
00218     fNormal123(rhs.fNormal123), fNormal142(rhs.fNormal142),
00219     fNormal134(rhs.fNormal134), fNormal234(rhs.fNormal234),
00220     warningFlag(rhs.warningFlag), fCdotN123(rhs.fCdotN123),
00221     fCdotN142(rhs.fCdotN142), fCdotN134(rhs.fCdotN134),
00222     fCdotN234(rhs.fCdotN234), fXMin(rhs.fXMin), fXMax(rhs.fXMax),
00223     fYMin(rhs.fYMin), fYMax(rhs.fYMax), fZMin(rhs.fZMin), fZMax(rhs.fZMax),
00224     fDx(rhs.fDx), fDy(rhs.fDy), fDz(rhs.fDz), fTol(rhs.fTol),
00225     fMaxSize(rhs.fMaxSize)
00226 {
00227 }
00228 
00229 
00231 //
00232 // Assignment operator
00233 
00234 G4Tet& G4Tet::operator = (const G4Tet& rhs) 
00235 {
00236    // Check assignment to self
00237    //
00238    if (this == &rhs)  { return *this; }
00239 
00240    // Copy base class data
00241    //
00242    G4VSolid::operator=(rhs);
00243 
00244    // Copy data
00245    //
00246    fCubicVolume = rhs.fCubicVolume; fSurfaceArea = rhs.fSurfaceArea;
00247    fpPolyhedron = 0; fAnchor = rhs.fAnchor;
00248    fP2 = rhs.fP2; fP3 = rhs.fP3; fP4 = rhs.fP4; fMiddle = rhs.fMiddle;
00249    fNormal123 = rhs.fNormal123; fNormal142 = rhs.fNormal142;
00250    fNormal134 = rhs.fNormal134; fNormal234 = rhs.fNormal234;
00251    warningFlag = rhs.warningFlag; fCdotN123 = rhs.fCdotN123;
00252    fCdotN142 = rhs.fCdotN142; fCdotN134 = rhs.fCdotN134;
00253    fCdotN234 = rhs.fCdotN234; fXMin = rhs.fXMin; fXMax = rhs.fXMax;
00254    fYMin = rhs.fYMin; fYMax = rhs.fYMax; fZMin = rhs.fZMin; fZMax = rhs.fZMax;
00255    fDx = rhs.fDx; fDy = rhs.fDy; fDz = rhs.fDz; fTol = rhs.fTol;
00256    fMaxSize = rhs.fMaxSize;
00257 
00258    return *this;
00259 }
00260 
00262 //
00263 // CheckDegeneracy
00264 
00265 G4bool G4Tet::CheckDegeneracy( G4ThreeVector anchor,
00266                                G4ThreeVector p2,
00267                                G4ThreeVector p3,
00268                                G4ThreeVector p4 )
00269 {
00270   G4bool result;
00271   G4Tet *object=new G4Tet("temp",anchor,p2,p3,p4,&result);
00272   delete object;
00273   return result;
00274 }
00275 
00277 //
00278 // Dispatch to parameterisation for replication mechanism dimension
00279 // computation & modification.
00280 
00281 void G4Tet::ComputeDimensions(G4VPVParameterisation* ,
00282                               const G4int ,
00283                               const G4VPhysicalVolume* )
00284 {
00285 }
00286 
00288 //
00289 // Calculate extent under transform and specified limit
00290 
00291 G4bool G4Tet::CalculateExtent(const EAxis pAxis,
00292                               const G4VoxelLimits& pVoxelLimit,
00293                               const G4AffineTransform& pTransform,
00294                                     G4double& pMin, G4double& pMax) const
00295 {
00296   G4double xMin,xMax;
00297   G4double yMin,yMax;
00298   G4double zMin,zMax;
00299 
00300   if (pTransform.IsRotated())
00301   {
00302     G4ThreeVector pp0=pTransform.TransformPoint(fAnchor);
00303     G4ThreeVector pp1=pTransform.TransformPoint(fP2);
00304     G4ThreeVector pp2=pTransform.TransformPoint(fP3);
00305     G4ThreeVector pp3=pTransform.TransformPoint(fP4);
00306 
00307     xMin    = std::min(std::min(std::min(pp0.x(), pp1.x()),pp2.x()),pp3.x());
00308     xMax    = std::max(std::max(std::max(pp0.x(), pp1.x()),pp2.x()),pp3.x());
00309     yMin    = std::min(std::min(std::min(pp0.y(), pp1.y()),pp2.y()),pp3.y());
00310     yMax    = std::max(std::max(std::max(pp0.y(), pp1.y()),pp2.y()),pp3.y());
00311     zMin    = std::min(std::min(std::min(pp0.z(), pp1.z()),pp2.z()),pp3.z());
00312     zMax    = std::max(std::max(std::max(pp0.z(), pp1.z()),pp2.z()),pp3.z());
00313 
00314   }
00315   else
00316   {
00317     G4double xoffset = pTransform.NetTranslation().x() ;
00318     xMin    = xoffset + fXMin;
00319     xMax    = xoffset + fXMax;
00320     G4double yoffset = pTransform.NetTranslation().y() ;
00321     yMin    = yoffset + fYMin;
00322     yMax    = yoffset + fYMax;
00323     G4double zoffset = pTransform.NetTranslation().z() ;
00324     zMin    = zoffset + fZMin;
00325     zMax    = zoffset + fZMax;
00326   }
00327 
00328   if (pVoxelLimit.IsXLimited())
00329   {
00330     if ( (xMin > pVoxelLimit.GetMaxXExtent()+fTol) || 
00331          (xMax < pVoxelLimit.GetMinXExtent()-fTol)  )  { return false; }
00332     else
00333     {
00334       xMin = std::max(xMin, pVoxelLimit.GetMinXExtent());
00335       xMax = std::min(xMax, pVoxelLimit.GetMaxXExtent());
00336     }
00337   }
00338 
00339   if (pVoxelLimit.IsYLimited())
00340   {
00341     if ( (yMin > pVoxelLimit.GetMaxYExtent()+fTol) ||
00342          (yMax < pVoxelLimit.GetMinYExtent()-fTol)  )  { return false; }
00343     else
00344     {
00345       yMin = std::max(yMin, pVoxelLimit.GetMinYExtent());
00346       yMax = std::min(yMax, pVoxelLimit.GetMaxYExtent());
00347     }
00348     }
00349 
00350     if (pVoxelLimit.IsZLimited())
00351     {
00352       if ( (zMin > pVoxelLimit.GetMaxZExtent()+fTol) ||
00353            (zMax < pVoxelLimit.GetMinZExtent()-fTol)  )  { return false; }
00354     else
00355     {
00356       zMin = std::max(zMin, pVoxelLimit.GetMinZExtent());
00357       zMax = std::min(zMax, pVoxelLimit.GetMaxZExtent());
00358     }
00359   }
00360 
00361   switch (pAxis)
00362   {
00363     case kXAxis:
00364       pMin=xMin;
00365       pMax=xMax;
00366       break;
00367     case kYAxis:
00368       pMin=yMin;
00369       pMax=yMax;
00370       break;
00371     case kZAxis:
00372       pMin=zMin;
00373       pMax=zMax;
00374       break;
00375     default:
00376       break;
00377   }
00378 
00379   return true;
00380 } 
00381 
00383 //
00384 // Return whether point inside/outside/on surface, using tolerance
00385 
00386 EInside G4Tet::Inside(const G4ThreeVector& p) const
00387 {
00388   G4double r123, r134, r142, r234;
00389 
00390   // this is written to allow if-statement truncation so the outside test
00391   // (where most of the world is) can fail very quickly and efficiently
00392 
00393   if ( (r123=p.dot(fNormal123)-fCdotN123) > fTol ||
00394        (r134=p.dot(fNormal134)-fCdotN134) > fTol ||
00395        (r142=p.dot(fNormal142)-fCdotN142) > fTol ||
00396        (r234=p.dot(fNormal234)-fCdotN234) > fTol )
00397   {
00398     return kOutside; // at least one is out!
00399   }
00400   else if( (r123 < -fTol)&&(r134 < -fTol)&&(r142 < -fTol)&&(r234 < -fTol) )
00401   {
00402     return kInside; // all are definitively inside
00403   }
00404   else
00405   {
00406     return kSurface; // too close to tell
00407   }
00408 }
00409 
00411 //
00412 // Calculate side nearest to p, and return normal
00413 // If two sides are equidistant, normal of first side (x/y/z) 
00414 // encountered returned.
00415 // This assumes that we are looking from the inside!
00416 
00417 G4ThreeVector G4Tet::SurfaceNormal( const G4ThreeVector& p) const
00418 {
00419   G4double r123=std::fabs(p.dot(fNormal123)-fCdotN123);
00420   G4double r134=std::fabs(p.dot(fNormal134)-fCdotN134);
00421   G4double r142=std::fabs(p.dot(fNormal142)-fCdotN142);
00422   G4double r234=std::fabs(p.dot(fNormal234)-fCdotN234);
00423 
00424   static const G4double delta = 0.5*kCarTolerance;
00425   G4ThreeVector sumnorm(0., 0., 0.);
00426   G4int noSurfaces=0; 
00427 
00428   if (r123 <= delta)         
00429   {
00430      noSurfaces ++; 
00431      sumnorm= fNormal123; 
00432   }
00433 
00434   if (r134 <= delta)    
00435   {
00436      noSurfaces ++; 
00437      sumnorm += fNormal134; 
00438   }
00439  
00440   if (r142 <= delta)    
00441   {
00442      noSurfaces ++; 
00443      sumnorm += fNormal142;
00444   }
00445   if (r234 <= delta)    
00446   {
00447      noSurfaces ++; 
00448      sumnorm += fNormal234;
00449   }
00450   
00451   if( noSurfaces > 0 )
00452   { 
00453      if( noSurfaces == 1 )
00454      { 
00455        return sumnorm; 
00456      }
00457      else
00458      {
00459        return sumnorm.unit();
00460      }
00461   }
00462   else // Approximative Surface Normal
00463   {
00464 
00465     if( (r123<=r134) && (r123<=r142) && (r123<=r234) ) { return fNormal123; }
00466     else if ( (r134<=r142) && (r134<=r234) )           { return fNormal134; }
00467     else if (r142 <= r234)                             { return fNormal142; }
00468     return fNormal234;
00469   }
00470 }
00472 //
00473 // Calculate distance to box from an outside point
00474 // - return kInfinity if no intersection.
00475 // All this is very unrolled, for speed.
00476 
00477 G4double G4Tet::DistanceToIn(const G4ThreeVector& p,
00478                              const G4ThreeVector& v) const
00479 {
00480     G4ThreeVector vu(v.unit()), hp;
00481     G4double vdotn, t, tmin=kInfinity;
00482 
00483     G4double extraDistance=10.0*fTol; // a little ways into the solid
00484 
00485     vdotn=-vu.dot(fNormal123);
00486     if(vdotn > 1e-12)
00487     { // this is a candidate face, since it is pointing at us
00488       t=(p.dot(fNormal123)-fCdotN123)/vdotn; // #  distance to intersection
00489       if( (t>=-fTol) && (t<tmin) )
00490       { // if not true, we're going away from this face or it's not close
00491         hp=p+vu*(t+extraDistance); // a little beyond point of intersection
00492         if ( ( hp.dot(fNormal134)-fCdotN134 < 0.0 ) &&
00493              ( hp.dot(fNormal142)-fCdotN142 < 0.0 ) &&
00494              ( hp.dot(fNormal234)-fCdotN234 < 0.0 ) )
00495         {
00496           tmin=t;
00497         }
00498       }
00499     }
00500 
00501     vdotn=-vu.dot(fNormal134);
00502     if(vdotn > 1e-12)
00503     { // # this is a candidate face, since it is pointing at us
00504       t=(p.dot(fNormal134)-fCdotN134)/vdotn; // #  distance to intersection
00505       if( (t>=-fTol) && (t<tmin) )
00506       { // if not true, we're going away from this face
00507         hp=p+vu*(t+extraDistance); // a little beyond point of intersection
00508         if ( ( hp.dot(fNormal123)-fCdotN123 < 0.0 ) && 
00509              ( hp.dot(fNormal142)-fCdotN142 < 0.0 ) &&
00510              ( hp.dot(fNormal234)-fCdotN234 < 0.0 ) )
00511         {
00512           tmin=t;
00513         }
00514       }
00515     }
00516 
00517     vdotn=-vu.dot(fNormal142);
00518     if(vdotn > 1e-12)
00519     { // # this is a candidate face, since it is pointing at us
00520       t=(p.dot(fNormal142)-fCdotN142)/vdotn; // #  distance to intersection
00521       if( (t>=-fTol) && (t<tmin) )
00522       { // if not true, we're going away from this face
00523         hp=p+vu*(t+extraDistance); // a little beyond point of intersection
00524         if ( ( hp.dot(fNormal123)-fCdotN123 < 0.0 ) &&
00525              ( hp.dot(fNormal134)-fCdotN134 < 0.0 ) &&
00526              ( hp.dot(fNormal234)-fCdotN234 < 0.0 ) )
00527         {
00528           tmin=t;
00529         }
00530       }
00531     }
00532 
00533     vdotn=-vu.dot(fNormal234);
00534     if(vdotn > 1e-12)
00535     { // # this is a candidate face, since it is pointing at us
00536       t=(p.dot(fNormal234)-fCdotN234)/vdotn; // #  distance to intersection
00537       if( (t>=-fTol) && (t<tmin) )
00538       { // if not true, we're going away from this face
00539         hp=p+vu*(t+extraDistance); // a little beyond point of intersection
00540         if ( ( hp.dot(fNormal123)-fCdotN123 < 0.0 ) &&
00541              ( hp.dot(fNormal134)-fCdotN134 < 0.0 ) &&
00542              ( hp.dot(fNormal142)-fCdotN142 < 0.0 ) )
00543         {
00544           tmin=t;
00545         }
00546       }
00547     }
00548 
00549   return std::max(0.0,tmin);
00550 }
00551 
00553 // 
00554 // Approximate distance to tet.
00555 // returns distance to sphere centered on bounding box
00556 // - If inside return 0
00557 
00558 G4double G4Tet::DistanceToIn(const G4ThreeVector& p) const
00559 {
00560   G4double dd=(p-fMiddle).mag() - fMaxSize - fTol;
00561   return std::max(0.0, dd);
00562 }
00563 
00565 //
00566 // Calcluate distance to surface of box from inside
00567 // by calculating distances to box's x/y/z planes.
00568 // Smallest distance is exact distance to exiting.
00569 
00570 G4double G4Tet::DistanceToOut( const G4ThreeVector& p,const G4ThreeVector& v,
00571                                const G4bool calcNorm,
00572                                      G4bool *validNorm, G4ThreeVector *n) const
00573 {
00574     G4ThreeVector vu(v.unit());
00575     G4double t1=kInfinity,t2=kInfinity,t3=kInfinity,t4=kInfinity, vdotn, tt;
00576 
00577     vdotn=vu.dot(fNormal123);
00578     if(vdotn > 1e-12)  // #we're heading towards this face, so it is a candidate
00579     {
00580       t1=(fCdotN123-p.dot(fNormal123))/vdotn; // #  distance to intersection
00581     }
00582 
00583     vdotn=vu.dot(fNormal134);
00584     if(vdotn > 1e-12) // #we're heading towards this face, so it is a candidate
00585     {
00586       t2=(fCdotN134-p.dot(fNormal134))/vdotn; // #  distance to intersection
00587     }
00588 
00589     vdotn=vu.dot(fNormal142);
00590     if(vdotn > 1e-12) // #we're heading towards this face, so it is a candidate
00591     {
00592       t3=(fCdotN142-p.dot(fNormal142))/vdotn; // #  distance to intersection
00593     }
00594 
00595     vdotn=vu.dot(fNormal234);
00596     if(vdotn > 1e-12) // #we're heading towards this face, so it is a candidate
00597     {
00598       t4=(fCdotN234-p.dot(fNormal234))/vdotn; // #  distance to intersection
00599     }
00600 
00601     tt=std::min(std::min(std::min(t1,t2),t3),t4);
00602 
00603     if (warningFlag && (tt == kInfinity || tt < -fTol))
00604     {
00605       DumpInfo();
00606       std::ostringstream message;
00607       message << "No good intersection found or already outside!?" << G4endl
00608               << "p = " << p / mm << "mm" << G4endl
00609               << "v = " << v  << G4endl
00610               << "t1, t2, t3, t4 (mm) "
00611               << t1/mm << ", " << t2/mm << ", " << t3/mm << ", " << t4/mm;
00612       G4Exception("G4Tet::DistanceToOut(p,v,...)", "GeomSolids1002",
00613                   JustWarning, message);
00614       if(validNorm)
00615       {
00616         *validNorm=false; // flag normal as meaningless
00617       }
00618     }
00619     else if(calcNorm && n)
00620     {
00621       G4ThreeVector normal;
00622       if(tt==t1)        { normal=fNormal123; }
00623       else if (tt==t2)  { normal=fNormal134; }
00624       else if (tt==t3)  { normal=fNormal142; }
00625       else if (tt==t4)  { normal=fNormal234; }
00626       *n=normal;
00627       if(validNorm) { *validNorm=true; }
00628     }
00629 
00630     return std::max(tt,0.0); // avoid tt<0.0 by a tiny bit
00631                              // if we are right on a face
00632 }
00633 
00635 //
00636 // Calculate exact shortest distance to any boundary from inside
00637 // - If outside return 0
00638 
00639 G4double G4Tet::DistanceToOut(const G4ThreeVector& p) const
00640 {
00641   G4double t1,t2,t3,t4;
00642   t1=fCdotN123-p.dot(fNormal123); //  distance to plane, positive if inside
00643   t2=fCdotN134-p.dot(fNormal134); //  distance to plane
00644   t3=fCdotN142-p.dot(fNormal142); //  distance to plane
00645   t4=fCdotN234-p.dot(fNormal234); //  distance to plane
00646 
00647   // if any one of these is negative, we are outside,
00648   // so return zero in that case
00649 
00650   G4double tmin=std::min(std::min(std::min(t1,t2),t3),t4);
00651   return (tmin < fTol)? 0:tmin;
00652 }
00653 
00655 //
00656 // Create a List containing the transformed vertices
00657 // Note: Caller has deletion responsibility
00658 
00659 G4ThreeVectorList*
00660 G4Tet::CreateRotatedVertices(const G4AffineTransform& pTransform) const
00661 {
00662   G4ThreeVectorList* vertices = new G4ThreeVectorList();
00663 
00664   if (vertices)
00665   {
00666     vertices->reserve(4);
00667     G4ThreeVector vertex0(fAnchor);
00668     G4ThreeVector vertex1(fP2);
00669     G4ThreeVector vertex2(fP3);
00670     G4ThreeVector vertex3(fP4);
00671 
00672     vertices->push_back(pTransform.TransformPoint(vertex0));
00673     vertices->push_back(pTransform.TransformPoint(vertex1));
00674     vertices->push_back(pTransform.TransformPoint(vertex2));
00675     vertices->push_back(pTransform.TransformPoint(vertex3));
00676   }
00677   else
00678   {
00679     DumpInfo();
00680     G4Exception("G4Tet::CreateRotatedVertices()",
00681                 "GeomSolids0003", FatalException,
00682                 "Error in allocation of vertices. Out of memory !");
00683   }
00684   return vertices;
00685 }
00686 
00688 //
00689 // GetEntityType
00690 
00691 G4GeometryType G4Tet::GetEntityType() const
00692 {
00693   return G4String("G4Tet");
00694 }
00695 
00697 //
00698 // Make a clone of the object
00699 
00700 G4VSolid* G4Tet::Clone() const
00701 {
00702   return new G4Tet(*this);
00703 }
00704 
00706 //
00707 // Stream object contents to an output stream
00708 
00709 std::ostream& G4Tet::StreamInfo(std::ostream& os) const
00710 {
00711   G4int oldprc = os.precision(16);
00712   os << "-----------------------------------------------------------\n"
00713   << "    *** Dump for solid - " << GetName() << " ***\n"
00714   << "    ===================================================\n"
00715   << " Solid type: G4Tet\n"
00716   << " Parameters: \n"
00717   << "    anchor: " << fAnchor/mm << " mm \n"
00718   << "    p2: " << fP2/mm << " mm \n"
00719   << "    p3: " << fP3/mm << " mm \n"
00720   << "    p4: " << fP4/mm << " mm \n"
00721   << "    normal123: " << fNormal123 << " \n"
00722   << "    normal134: " << fNormal134 << " \n"
00723   << "    normal142: " << fNormal142 << " \n"
00724   << "    normal234: " << fNormal234 << " \n"
00725   << "-----------------------------------------------------------\n";
00726   os.precision(oldprc);
00727 
00728   return os;
00729 }
00730 
00731 
00733 //
00734 // GetPointOnFace
00735 //
00736 // Auxiliary method for get point on surface
00737 
00738 G4ThreeVector G4Tet::GetPointOnFace(G4ThreeVector p1, G4ThreeVector p2,
00739                                     G4ThreeVector p3, G4double& area) const
00740 {
00741   G4double lambda1,lambda2;
00742   G4ThreeVector v, w;
00743 
00744   v = p3 - p1;
00745   w = p1 - p2;
00746 
00747   lambda1 = RandFlat::shoot(0.,1.);
00748   lambda2 = RandFlat::shoot(0.,lambda1);
00749 
00750   area = 0.5*(v.cross(w)).mag();
00751 
00752   return (p2 + lambda1*w + lambda2*v);
00753 }
00754 
00756 //
00757 // GetPointOnSurface
00758 
00759 G4ThreeVector G4Tet::GetPointOnSurface() const
00760 {
00761   G4double chose,aOne,aTwo,aThree,aFour;
00762   G4ThreeVector p1, p2, p3, p4;
00763   
00764   p1 = GetPointOnFace(fAnchor,fP2,fP3,aOne);
00765   p2 = GetPointOnFace(fAnchor,fP4,fP3,aTwo);
00766   p3 = GetPointOnFace(fAnchor,fP4,fP2,aThree);
00767   p4 = GetPointOnFace(fP4,fP3,fP2,aFour);
00768   
00769   chose = RandFlat::shoot(0.,aOne+aTwo+aThree+aFour);
00770   if( (chose>=0.) && (chose <aOne) ) {return p1;}
00771   else if( (chose>=aOne) && (chose < aOne+aTwo) ) {return p2;}
00772   else if( (chose>=aOne+aTwo) && (chose<aOne+aTwo+aThree) ) {return p3;}
00773   return p4;
00774 }
00775 
00777 //
00778 // GetVertices
00779 
00780 std::vector<G4ThreeVector> G4Tet::GetVertices() const 
00781 {
00782   std::vector<G4ThreeVector> vertices(4);
00783   vertices[0] = fAnchor;
00784   vertices[1] = fP2;
00785   vertices[2] = fP3;
00786   vertices[3] = fP4;
00787 
00788   return vertices;
00789 }
00790 
00792 //
00793 // GetCubicVolume
00794 
00795 G4double G4Tet::GetCubicVolume()
00796 {
00797   return fCubicVolume;
00798 }
00799 
00801 //
00802 // GetSurfaceArea
00803 
00804 G4double G4Tet::GetSurfaceArea()
00805 {
00806   return fSurfaceArea;
00807 }
00808 
00809 // Methods for visualisation
00810 
00812 //
00813 // DescribeYourselfTo
00814 
00815 void G4Tet::DescribeYourselfTo (G4VGraphicsScene& scene) const 
00816 {
00817   scene.AddSolid (*this);
00818 }
00819 
00821 //
00822 // GetExtent
00823 
00824 G4VisExtent G4Tet::GetExtent() const 
00825 {
00826   return G4VisExtent (fXMin, fXMax, fYMin, fYMax, fZMin, fZMax);
00827 }
00828 
00830 //
00831 // CreatePolyhedron
00832 
00833 G4Polyhedron* G4Tet::CreatePolyhedron () const 
00834 {
00835   G4Polyhedron *ph=new G4Polyhedron;
00836   G4double xyz[4][3];
00837   const G4int faces[4][4]={{1,3,2,0},{1,4,3,0},{1,2,4,0},{2,3,4,0}};
00838   xyz[0][0]=fAnchor.x(); xyz[0][1]=fAnchor.y(); xyz[0][2]=fAnchor.z();
00839   xyz[1][0]=fP2.x(); xyz[1][1]=fP2.y(); xyz[1][2]=fP2.z();
00840   xyz[2][0]=fP3.x(); xyz[2][1]=fP3.y(); xyz[2][2]=fP3.z();
00841   xyz[3][0]=fP4.x(); xyz[3][1]=fP4.y(); xyz[3][2]=fP4.z();
00842 
00843   ph->createPolyhedron(4,4,xyz,faces);
00844 
00845   return ph;
00846 }
00847 
00849 //
00850 // CreateNURBS
00851 
00852 G4NURBS* G4Tet::CreateNURBS () const 
00853 {
00854   return new G4NURBSbox (fDx, fDy, fDz);
00855 }
00856 
00858 //
00859 // GetPolyhedron
00860 
00861 G4Polyhedron* G4Tet::GetPolyhedron () const
00862 {
00863   if (!fpPolyhedron ||
00864       fpPolyhedron->GetNumberOfRotationStepsAtTimeOfCreation() !=
00865       fpPolyhedron->GetNumberOfRotationSteps())
00866     {
00867       delete fpPolyhedron;
00868       fpPolyhedron = CreatePolyhedron();
00869     }
00870   return fpPolyhedron;
00871 }

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