G4FCylindricalSurface Class Reference

#include <G4FCylindricalSurface.hh>

Inheritance diagram for G4FCylindricalSurface:

G4Surface G4STEPEntity

Public Member Functions

 G4FCylindricalSurface ()
 G4FCylindricalSurface (const G4Point3D &o, const G4Vector3D &a, G4double r, G4double l)
virtual ~G4FCylindricalSurface ()
G4int operator== (const G4FCylindricalSurface &c) const
virtual G4Vector3D SurfaceNormal (const G4Point3D &p) const
virtual G4int Inside (const G4Vector3D &x) const
G4String GetEntityType () const
G4int Intersect (const G4Ray &)
virtual G4double HowNear (const G4Vector3D &x) const
void CalcBBox ()
virtual const char * NameOf () const
virtual void PrintOn (std::ostream &os=G4cout) const
virtual G4int WithinBoundary (const G4Vector3D &x) const
virtual G4double Scale () const
virtual G4double Area () const
virtual void resize (G4double r, G4double l)
G4double GetLength () const
G4Vector3D GetAxis () const
G4double GetRadius () const
void SetRadius (G4double r)
void InitValues ()

Protected Attributes

G4Axis2Placement3D Position
G4double radius
G4double length

Detailed Description

Definition at line 64 of file G4FCylindricalSurface.hh.


Constructor & Destructor Documentation

G4FCylindricalSurface::G4FCylindricalSurface (  ) 

Definition at line 41 of file G4FCylindricalSurface.cc.

00042   : radius(0.), length(1.)
00043 {
00044 }

G4FCylindricalSurface::G4FCylindricalSurface ( const G4Point3D o,
const G4Vector3D a,
G4double  r,
G4double  l 
)

Definition at line 52 of file G4FCylindricalSurface.cc.

References G4endl, G4Exception(), G4Axis2Placement3D::Init(), JustWarning, length, G4Surface::origin, Position, and radius.

00057 { 
00058   //  make a G4FCylindricalSurface with origin o, axis a, 
00059   //  radius r, and length l
00060   G4Vector3D dir(1,1,1);
00061   Position.Init(dir, a, o);
00062 
00063   origin = o;
00064   radius = r;
00065   
00066   //  Require length to be positive or zero
00067   if ( l >= 0.0 )
00068     length = l;
00069   else 
00070   {
00071     std::ostringstream message;
00072     message << "Negative length." << G4endl
00073             << "Default length of 0.0 is used.";    
00074     G4Exception("G4FCylindricalSurface::G4FCylindricalSurface()",
00075                 "GeomSolids1001", JustWarning, message);
00076 
00077     length = 0.0;
00078   }
00079 
00080   //  Require radius to be non-negative (i.e., allow zero)
00081   if ( r >= 0.0 )
00082     radius = r;
00083   else 
00084   {
00085     std::ostringstream message;
00086     message << "Negative radius." << G4endl
00087             << "Default value of 0.0 is used.";    
00088     G4Exception("G4FCylindricalSurface::G4FCylindricalSurface()",
00089                 "GeomSolids1001", JustWarning, message);
00090 
00091     radius = 0.0;
00092   }
00093 }

G4FCylindricalSurface::~G4FCylindricalSurface (  )  [virtual]

Definition at line 47 of file G4FCylindricalSurface.cc.

00048 {
00049 }


Member Function Documentation

G4double G4FCylindricalSurface::Area (  )  const [virtual]

Definition at line 111 of file G4FCylindricalSurface.cc.

References length, G4INCL::Math::pi, and radius.

00112 {
00113   return ( 2.0 * pi * radius * length );
00114 }

void G4FCylindricalSurface::CalcBBox (  )  [virtual]

Reimplemented from G4Surface.

Definition at line 119 of file G4FCylindricalSurface.cc.

References G4Surface::bbox, G4BoundingBox3D::Extend(), G4Axis2Placement3D::GetAxis(), G4Axis2Placement3D::GetLocation(), G4BoundingBox3D::Init(), G4Surface::kCarTolerance, length, PINFINITY(), Position, and radius.

00120 {
00121   // Finds the bounds of the surface iow
00122   // calculates the bounds for a bounding box
00123   // to the surface. The bounding box is used
00124   // for a preliminary check of intersection.
00125   G4Point3D Max = G4Point3D(-PINFINITY);
00126   G4Point3D Min = G4Point3D( PINFINITY);
00127 
00128   G4Point3D Tmp; 
00129   G4Point3D Origin    = Position.GetLocation();
00130   G4Point3D EndOrigin = G4Point3D( Origin + (length*Position.GetAxis()) );
00131   G4Point3D Radius(radius, radius, 0);
00132  
00133   // Default BBox
00134   G4Point3D Tolerance(kCarTolerance, kCarTolerance, kCarTolerance);
00135   G4Point3D BoxMin(Origin-Tolerance);
00136   G4Point3D BoxMax(Origin+Tolerance);
00137 
00138   bbox = new G4BoundingBox3D();
00139   bbox->Init(BoxMin, BoxMax);
00140 
00141 
00142   Tmp = (Origin - Radius);
00143   bbox->Extend(Tmp);
00144   
00145   Tmp = Origin + Radius;
00146   bbox->Extend(Tmp);
00147 
00148   Tmp = EndOrigin - Radius;
00149   bbox->Extend(Tmp);
00150   
00151   Tmp = EndOrigin + Radius;
00152   bbox->Extend(Tmp);
00153 }

G4Vector3D G4FCylindricalSurface::GetAxis (  )  const [inline]

Definition at line 59 of file G4FCylindricalSurface.icc.

References G4Axis2Placement3D::GetAxis(), and Position.

00060 {
00061   return Position.GetAxis();
00062 }

G4String G4FCylindricalSurface::GetEntityType (  )  const [inline, virtual]

Reimplemented from G4Surface.

Definition at line 47 of file G4FCylindricalSurface.icc.

00048 {
00049   return G4String("Cylindrical_Surface");
00050 }

G4double G4FCylindricalSurface::GetLength (  )  const [inline]

Definition at line 53 of file G4FCylindricalSurface.icc.

References length.

00054 { 
00055   return length; 
00056 }

G4double G4FCylindricalSurface::GetRadius (  )  const [inline]

Definition at line 65 of file G4FCylindricalSurface.icc.

References radius.

00066 {
00067   return radius;
00068 }

G4double G4FCylindricalSurface::HowNear ( const G4Vector3D x  )  const [virtual]

Reimplemented from G4Surface.

Definition at line 241 of file G4FCylindricalSurface.cc.

References length, G4Surface::origin, and radius.

Referenced by Inside().

00242 {
00243   // Shortest distance from the point x to the G4FCylindricalSurface.
00244   // The distance will be always positive 
00245 
00246   G4double   hownear;
00247 
00248   G4Vector3D upcorner = G4Vector3D ( radius, 0 , origin.z()+length);
00249   G4Vector3D downcorner = G4Vector3D ( radius, 0 , origin.z());
00250   G4Vector3D xd;  
00251   
00252   xd = G4Vector3D ( std::sqrt ( x.x()*x.x() + x.y()*x.y() ) , 0 , x.z() );
00253     
00254   
00255   G4double Zinter = (xd.z()) ;
00256   
00257   if ( ((Zinter >= downcorner.z()) && (Zinter <=upcorner.z())) ) {
00258     hownear = std::fabs( radius - xd.x() );
00259   } else {
00260     hownear = std::min ( (xd-upcorner).mag() , (xd-downcorner).mag() );
00261   }
00262 
00263   return hownear;
00264 }

void G4FCylindricalSurface::InitValues (  ) 

G4int G4FCylindricalSurface::Inside ( const G4Vector3D x  )  const [virtual]

Definition at line 309 of file G4FCylindricalSurface.cc.

References HowNear(), and G4Surface::kCarTolerance.

Referenced by Intersect().

00310 { 
00311   //  Return 0 if point x is outside G4CylindricalSurface, 1 if Inside.
00312   //  Outside means that the distance to the G4CylindricalSurface would 
00313   //  be negative.
00314   //  Use the HowNear function to calculate this distance.
00315   if ( HowNear( x ) >= -0.5*kCarTolerance )
00316     return 1;
00317   else
00318     return 0; 
00319 }

G4int G4FCylindricalSurface::Intersect ( const G4Ray  )  [virtual]

Reimplemented from G4Surface.

Definition at line 156 of file G4FCylindricalSurface.cc.

References G4Surface::closest_hit, G4Surface::distance, G4Axis2Placement3D::GetAxis(), G4Surface::GetBBox(), G4Ray::GetDir(), G4Axis2Placement3D::GetLocation(), G4Ray::GetStart(), Inside(), G4Surface::kCarTolerance, PINFINITY(), Position, and radius.

00157 {
00158   // This function count the number of intersections of a 
00159   // bounded cylindrical surface by a ray.
00160   // At first, calculates the intersections with the infinite 
00161   // cylindrical surfsace. After, count the intersections within the
00162   // finite cylindrical surface boundaries, and set "distance" to the 
00163   // closest distance from the start point to the nearest intersection
00164   // If the point is on the surface it returns or the intersection with
00165   // the opposite surface or kInfinity
00166 
00167   // If no intersection is founded, set distance = kInfinity and
00168   // return 0
00169 
00170   distance    = kInfinity;
00171   closest_hit = PINFINITY;
00172 
00173   // origin and direction of the ray
00174   G4Point3D  x    = ry.GetStart();
00175   G4Vector3D dhat = ry.GetDir();
00176 
00177   // cylinder axis 
00178   G4Vector3D ahat  = Position.GetAxis();
00179  
00180   //  array of solutions in distance along the ray
00181   G4double sol[2];
00182   sol[0]=-1.0;
00183   sol[1]=-1.0;
00184 
00185   // calculate the two intersections (quadratic equation)   
00186   G4Vector3D gamma =  G4Vector3D( x - Position.GetLocation() );
00187   
00188   G4double ga = gamma * ahat;
00189   G4double da = dhat * ahat;
00190 
00191   G4double A = da * da - dhat * dhat;
00192   G4double B = 2 * ( -gamma * dhat + ga * da );
00193   G4double C =  -gamma * gamma + ga * ga  + radius * radius ;
00194 
00195   G4double radical = B * B  -  4.0 * A * C; 
00196 
00197   if ( radical < 0.0 ) 
00198     // no intersection
00199     return 0;
00200   else 
00201   {
00202     G4double root = std::sqrt( radical );
00203     sol[0] = ( - B + root ) / ( 2. * A );
00204     sol[1] = ( - B - root ) / ( 2. * A );
00205   }
00206   
00207   // validity of the solutions
00208   // the hit point must be into the bounding box of the cylindrical surface
00209   G4Point3D p0 = G4Point3D( x + sol[0]*dhat );
00210   G4Point3D p1 = G4Point3D( x + sol[1]*dhat );
00211 
00212   if( !GetBBox()->Inside(p0) )
00213     sol[0] = kInfinity;
00214 
00215   if( !GetBBox()->Inside(p1) )
00216     sol[1] = kInfinity;
00217   
00218   //  now loop over each positive solution, keeping the first one (smallest
00219   //  distance along the Ray) which is within the boundary of the sub-shape
00220   G4int nbinter = 0;
00221   distance = kInfinity;
00222 
00223   for ( G4int i = 0; i < 2; i++ ) 
00224   {  
00225     if(sol[i] < kInfinity) {
00226       if ( sol[i] >= kCarTolerance*0.5 ) {
00227         nbinter ++;
00228         // real intersection
00229         // set the distance if it is the smallest
00230         if( distance > sol[i]*sol[i]) {
00231           distance = sol[i]*sol[i];
00232         }
00233       }
00234     }    
00235   }
00236 
00237   return nbinter;
00238 }

const char * G4FCylindricalSurface::NameOf (  )  const [virtual]

Definition at line 96 of file G4FCylindricalSurface.cc.

00097 {
00098   return "G4FCylindricalSurface"; 
00099 }

G4int G4FCylindricalSurface::operator== ( const G4FCylindricalSurface c  )  const [inline]

Definition at line 38 of file G4FCylindricalSurface.icc.

References length, G4Surface::origin, Position, and radius.

00039 {
00040   return ( origin == c.origin && 
00041            Position == c.Position   && 
00042            radius == c.radius && 
00043            length == c.length    );
00044 }

void G4FCylindricalSurface::PrintOn ( std::ostream &  os = G4cout  )  const [virtual]

Definition at line 102 of file G4FCylindricalSurface.cc.

References G4Axis2Placement3D::GetAxis(), length, G4Surface::origin, Position, and radius.

00103 { 
00104   os << "G4FCylindricalSurface with origin: " << origin << "\t"
00105      << "and axis: " << Position.GetAxis() << "\n"
00106      << "\t radius: " << radius << "\t and length: "
00107      << length << "\n";
00108 }

void G4FCylindricalSurface::resize ( G4double  r,
G4double  l 
) [virtual]

Definition at line 322 of file G4FCylindricalSurface.cc.

References G4endl, G4Exception(), JustWarning, length, and radius.

00323 {
00324   //  Resize a G4FCylindricalSurface to a new radius r and new length l
00325   //  Require radius to be non-negative
00326   if ( r >= 0.0 )
00327     radius = r;
00328   else 
00329   {
00330     std::ostringstream message;
00331     message << "Negative radius." << G4endl
00332             << "Original value of " << radius << " is retained.";    
00333     G4Exception("G4FCylindricalSurface::resize()",
00334                 "GeomSolids1001", JustWarning, message);
00335   }
00336 
00337   //  Require length to be positive
00338   if ( l > 0.0 )
00339     length = l;
00340   else 
00341   {
00342     std::ostringstream message;
00343     message << "Negative or zero length." << G4endl
00344             << "Original value of " << length << " is retained.";    
00345     G4Exception("G4FCylindricalSurface::resize()",
00346                 "GeomSolids1001", JustWarning, message);
00347   }
00348 }

G4double G4FCylindricalSurface::Scale (  )  const [virtual]

Definition at line 278 of file G4FCylindricalSurface.cc.

References length, and radius.

00279 {
00280   //  Returns the radius of a G4FCylindricalSurface unless it is zero,
00281   //  in which case returns the length.
00282   //  Used for Scale-invariant tests of surface thickness.
00283   if ( radius == 0.0 )
00284     return length;
00285   else
00286     return radius;
00287 }

void G4FCylindricalSurface::SetRadius ( G4double  r  ) 

G4Vector3D G4FCylindricalSurface::SurfaceNormal ( const G4Point3D p  )  const [virtual]

Implements G4Surface.

Definition at line 290 of file G4FCylindricalSurface.cc.

References G4Axis2Placement3D::GetAxis(), G4Axis2Placement3D::GetLocation(), CLHEP::detail::n, Position, and G4Surface::sameSense.

00291 {
00292   //  return the Normal unit vector to the G4CylindricalSurface at a point 
00293   //  p on (or nearly on) the G4CylindricalSurface
00294   
00295   G4Vector3D n = G4Vector3D( ( p - Position.GetLocation() ) - 
00296                            ( ( p - Position.GetLocation()) *
00297                                Position.GetAxis() ) * Position.GetAxis() );
00298   G4double nmag = n.mag();
00299   
00300   if ( nmag != 0.0 )
00301     n = n * (1/nmag);
00302 
00303   if( !sameSense )
00304     n = -n;
00305   
00306   return n;
00307 }

G4int G4FCylindricalSurface::WithinBoundary ( const G4Vector3D x  )  const [virtual]

Definition at line 266 of file G4FCylindricalSurface.cc.

References G4Axis2Placement3D::GetAxis(), G4Axis2Placement3D::GetLocation(), length, and Position.

00267 {
00268   //  return 1 if point x is within the boundaries of the G4FCylindricalSurface
00269   //  return 0 otherwise (assume it is on the cylinder)
00270   if ( std::fabs( ( x - Position.GetLocation()) * Position.GetAxis() )
00271        <= 0.5 * length )
00272     return 1;
00273   else
00274     return 0;
00275 }


Field Documentation

G4double G4FCylindricalSurface::length [protected]

Definition at line 154 of file G4FCylindricalSurface.hh.

Referenced by Area(), CalcBBox(), G4FCylindricalSurface(), GetLength(), HowNear(), operator==(), PrintOn(), resize(), Scale(), and WithinBoundary().

G4Axis2Placement3D G4FCylindricalSurface::Position [protected]

Definition at line 152 of file G4FCylindricalSurface.hh.

Referenced by CalcBBox(), G4FCylindricalSurface(), GetAxis(), Intersect(), operator==(), PrintOn(), SurfaceNormal(), and WithinBoundary().

G4double G4FCylindricalSurface::radius [protected]

Definition at line 153 of file G4FCylindricalSurface.hh.

Referenced by Area(), CalcBBox(), G4FCylindricalSurface(), GetRadius(), HowNear(), Intersect(), operator==(), PrintOn(), resize(), and Scale().


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:51:59 2013 for Geant4 by  doxygen 1.4.7