G4ExtrudedSolid.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 //
00027 // $Id: G4ExtrudedSolid.cc 69790 2013-05-15 12:39:10Z gcosmo $
00028 //
00029 //
00030 // --------------------------------------------------------------------
00031 // GEANT 4 class source file
00032 //
00033 // G4ExtrudedSolid.cc
00034 //
00035 // Author: Ivana Hrivnacova, IPN Orsay
00036 // --------------------------------------------------------------------
00037 
00038 #include <set>
00039 #include <algorithm>
00040 #include <cmath>
00041 #include <iomanip>
00042 
00043 #include "G4ExtrudedSolid.hh"
00044 #include "G4PhysicalConstants.hh"
00045 #include "G4SystemOfUnits.hh"
00046 #include "G4VFacet.hh"
00047 #include "G4TriangularFacet.hh"
00048 #include "G4QuadrangularFacet.hh"
00049 
00050 //_____________________________________________________________________________
00051 
00052 G4ExtrudedSolid::G4ExtrudedSolid( const G4String& pName,
00053                                         std::vector<G4TwoVector> polygon,
00054                                         std::vector<ZSection> zsections)
00055   : G4TessellatedSolid(pName),
00056     fNv(polygon.size()),
00057     fNz(zsections.size()),
00058     fPolygon(),
00059     fZSections(),
00060     fTriangles(),
00061     fIsConvex(false),
00062     fGeometryType("G4ExtrudedSolid")
00063     
00064 {
00065   // General constructor 
00066 
00067   // First check input parameters
00068 
00069   if ( fNv < 3 )
00070   {
00071     std::ostringstream message;
00072     message << "Number of polygon vertices < 3 - " << pName;
00073     G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0002",
00074                 FatalErrorInArgument, message);
00075   }
00076      
00077   if ( fNz < 2 )
00078   {
00079     std::ostringstream message;
00080     message << "Number of z-sides < 2 - " << pName;
00081     G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0002",
00082                 FatalErrorInArgument, message);
00083   }
00084      
00085   for ( G4int i=0; i<fNz-1; ++i ) 
00086   {
00087     if ( zsections[i].fZ > zsections[i+1].fZ ) 
00088     {
00089       std::ostringstream message;
00090       message << "Z-sections have to be ordered by z value (z0 < z1 < z2...) - "
00091               << pName;
00092       G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0002",
00093                   FatalErrorInArgument, message);
00094     }
00095     if ( std::fabs( zsections[i+1].fZ - zsections[i].fZ ) < kCarTolerance * 0.5 ) 
00096     {
00097       std::ostringstream message;
00098       message << "Z-sections with the same z position are not supported - "
00099               << pName;
00100       G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0001",
00101                   FatalException, message);
00102     }
00103   }  
00104   
00105   // Check if polygon vertices are defined clockwise
00106   // (the area is positive if polygon vertices are defined anti-clockwise)
00107   //
00108   G4double area = 0.;
00109   for ( G4int i=0; i<fNv; ++i ) {
00110     G4int j = i+1;
00111     if ( j == fNv ) j = 0;
00112     area += 0.5 * ( polygon[i].x()*polygon[j].y() - polygon[j].x()*polygon[i].y());
00113   }
00114  
00115   // Copy polygon
00116   //
00117   if  ( area < 0. ) {   
00118     // Polygon vertices are defined clockwise, we just copy the polygon       
00119     for ( G4int i=0; i<fNv; ++i ) { fPolygon.push_back(polygon[i]); }
00120   }
00121   else {
00122     // Polygon vertices are defined anti-clockwise, we revert them
00123     //G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids1001",
00124     //            JustWarning, 
00125     //            "Polygon vertices defined anti-clockwise, reverting polygon");      
00126     for ( G4int i=0; i<fNv; ++i ) { fPolygon.push_back(polygon[fNv-i-1]); }
00127   }
00128     
00129   
00130   // Copy z-sections
00131   //
00132   for ( G4int i=0; i<fNz; ++i ) { fZSections.push_back(zsections[i]); }
00133     
00134 
00135   G4bool result = MakeFacets();
00136   if (!result)
00137   {   
00138     std::ostringstream message;
00139     message << "Making facets failed - " << pName;
00140     G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0003",
00141                 FatalException, message);
00142   }
00143   fIsConvex = IsConvex();
00144 
00145   
00146   ComputeProjectionParameters();
00147 }
00148 
00149 //_____________________________________________________________________________
00150 
00151 G4ExtrudedSolid::G4ExtrudedSolid( const G4String& pName,
00152                                         std::vector<G4TwoVector> polygon,       
00153                                         G4double dz,
00154                                         G4TwoVector off1, G4double scale1,
00155                                         G4TwoVector off2, G4double scale2 )
00156   : G4TessellatedSolid(pName),
00157     fNv(polygon.size()),
00158     fNz(2),
00159     fPolygon(),
00160     fZSections(),
00161     fTriangles(),
00162     fIsConvex(false),
00163     fGeometryType("G4ExtrudedSolid")
00164     
00165 {
00166   // Special constructor for solid with 2 z-sections
00167 
00168   // First check input parameters
00169   //
00170   if ( fNv < 3 )
00171   {
00172     std::ostringstream message;
00173     message << "Number of polygon vertices < 3 - " << pName;
00174     G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0002",
00175                 FatalErrorInArgument, message);
00176   }
00177      
00178   // Check if polygon vertices are defined clockwise
00179   // (the area is positive if polygon vertices are defined anti-clockwise)
00180   
00181   G4double area = 0.;
00182   for ( G4int i=0; i<fNv; ++i )
00183   {
00184     G4int j = i+1;
00185     if ( j == fNv )  { j = 0; }
00186     area += 0.5 * ( polygon[i].x()*polygon[j].y()
00187                   - polygon[j].x()*polygon[i].y());
00188   }
00189  
00190   // Copy polygon
00191   //
00192   if  ( area < 0. )
00193   {   
00194     // Polygon vertices are defined clockwise, we just copy the polygon       
00195     for ( G4int i=0; i<fNv; ++i ) { fPolygon.push_back(polygon[i]); }
00196   }
00197   else
00198   {
00199     // Polygon vertices are defined anti-clockwise, we revert them
00200     //G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids1001",
00201     //            JustWarning, 
00202     //            "Polygon vertices defined anti-clockwise, reverting polygon");      
00203     for ( G4int i=0; i<fNv; ++i ) { fPolygon.push_back(polygon[fNv-i-1]); }
00204   }
00205   
00206   // Copy z-sections
00207   //
00208   fZSections.push_back(ZSection(-dz, off1, scale1));
00209   fZSections.push_back(ZSection( dz, off2, scale2));
00210     
00211   G4bool result = MakeFacets();
00212   if (!result)
00213   {   
00214     std::ostringstream message;
00215     message << "Making facets failed - " << pName;
00216     G4Exception("G4ExtrudedSolid::G4ExtrudedSolid()", "GeomSolids0003",
00217                 FatalException, message);
00218   }
00219   fIsConvex = IsConvex();
00220 
00221   ComputeProjectionParameters();
00222 }
00223 
00224 //_____________________________________________________________________________
00225 
00226 G4ExtrudedSolid::G4ExtrudedSolid( __void__& a )
00227   : G4TessellatedSolid(a), fNv(0), fNz(0), fPolygon(), fZSections(),
00228     fTriangles(), fIsConvex(false), fGeometryType("G4ExtrudedSolid")
00229 {
00230   // Fake default constructor - sets only member data and allocates memory
00231   //                            for usage restricted to object persistency.
00232 }
00233 
00234 //_____________________________________________________________________________
00235 
00236 G4ExtrudedSolid::G4ExtrudedSolid(const G4ExtrudedSolid& rhs)
00237   : G4TessellatedSolid(rhs), fNv(rhs.fNv), fNz(rhs.fNz),
00238     fPolygon(rhs.fPolygon), fZSections(rhs.fZSections),
00239     fTriangles(rhs.fTriangles), fIsConvex(rhs.fIsConvex),
00240     fGeometryType(rhs.fGeometryType), fKScales(rhs.fKScales),
00241     fScale0s(rhs.fScale0s), fKOffsets(rhs.fKOffsets), fOffset0s(rhs.fOffset0s)
00242 {
00243 }
00244 
00245 
00246 //_____________________________________________________________________________
00247 
00248 G4ExtrudedSolid& G4ExtrudedSolid::operator = (const G4ExtrudedSolid& rhs) 
00249 {
00250    // Check assignment to self
00251    //
00252    if (this == &rhs)  { return *this; }
00253 
00254    // Copy base class data
00255    //
00256    G4TessellatedSolid::operator=(rhs);
00257 
00258    // Copy data
00259    //
00260    fNv = rhs.fNv; fNz = rhs.fNz;
00261    fPolygon = rhs.fPolygon; fZSections = rhs.fZSections;
00262    fTriangles = rhs.fTriangles; fIsConvex = rhs.fIsConvex;
00263    fGeometryType = rhs.fGeometryType; fKScales = rhs.fKScales;
00264    fScale0s = rhs.fScale0s; fKOffsets = rhs.fKOffsets;
00265    fOffset0s = rhs.fOffset0s;
00266 
00267    return *this;
00268 }
00269 
00270 //_____________________________________________________________________________
00271 
00272 G4ExtrudedSolid::~G4ExtrudedSolid()
00273 {
00274   // Destructor
00275 }
00276 
00277 //_____________________________________________________________________________
00278 
00279 void G4ExtrudedSolid::ComputeProjectionParameters()
00280 {
00281   // Compute parameters for point projections p(z) 
00282   // to the polygon scale & offset:
00283   // scale(z) = k*z + scale0
00284   // offset(z) = l*z + offset0
00285   // p(z) = scale(z)*p0 + offset(z)  
00286   // p0 = (p(z) - offset(z))/scale(z);
00287   //  
00288 
00289   for ( G4int iz=0; iz<fNz-1; ++iz) 
00290   {
00291     G4double z1      = fZSections[iz].fZ;
00292     G4double z2      = fZSections[iz+1].fZ;
00293     G4double scale1  = fZSections[iz].fScale;
00294     G4double scale2  = fZSections[iz+1].fScale;
00295     G4TwoVector off1 = fZSections[iz].fOffset;
00296     G4TwoVector off2 = fZSections[iz+1].fOffset;
00297     
00298     G4double kscale = (scale2 - scale1)/(z2 - z1);
00299     G4double scale0 =  scale2 - kscale*(z2 - z1)/2.0; 
00300     G4TwoVector koff = (off2 - off1)/(z2 - z1);
00301     G4TwoVector off0 =  off2 - koff*(z2 - z1)/2.0; 
00302 
00303     fKScales.push_back(kscale);
00304     fScale0s.push_back(scale0);
00305     fKOffsets.push_back(koff);
00306     fOffset0s.push_back(off0);
00307   }  
00308 }
00309 
00310 
00311 //_____________________________________________________________________________
00312 
00313 G4ThreeVector G4ExtrudedSolid::GetVertex(G4int iz, G4int ind) const
00314 {
00315   // Shift and scale vertices
00316 
00317   return G4ThreeVector( fPolygon[ind].x() * fZSections[iz].fScale
00318                       + fZSections[iz].fOffset.x(), 
00319                         fPolygon[ind].y() * fZSections[iz].fScale
00320                       + fZSections[iz].fOffset.y(), fZSections[iz].fZ);
00321 }       
00322 
00323 //_____________________________________________________________________________
00324 
00325 
00326 G4TwoVector G4ExtrudedSolid::ProjectPoint(const G4ThreeVector& point) const
00327 {
00328   // Project point in the polygon scale
00329   // scale(z) = k*z + scale0
00330   // offset(z) = l*z + offset0
00331   // p(z) = scale(z)*p0 + offset(z)  
00332   // p0 = (p(z) - offset(z))/scale(z);
00333   
00334   // Select projection (z-segment of the solid) according to p.z()
00335   //
00336   G4int iz = 0;
00337   while ( point.z() > fZSections[iz+1].fZ && iz < fNz-2 ) { ++iz; }
00338   
00339   G4double z0 = ( fZSections[iz+1].fZ + fZSections[iz].fZ )/2.0;
00340   G4TwoVector p2(point.x(), point.y());
00341   G4double pscale  = fKScales[iz]*(point.z()-z0) + fScale0s[iz];
00342   G4TwoVector poffset = fKOffsets[iz]*(point.z()-z0) + fOffset0s[iz];
00343   
00344   // G4cout << point << " projected to " 
00345   //        << iz << "-th z-segment polygon as "
00346   //        << (p2 - poffset)/pscale << G4endl;
00347 
00348   // pscale is always >0 as it is an interpolation between two
00349   // positive scale values
00350   //
00351   return (p2 - poffset)/pscale;
00352 }  
00353 
00354 //_____________________________________________________________________________
00355 
00356 G4bool G4ExtrudedSolid::IsSameLine(G4TwoVector p,  
00357                                    G4TwoVector l1, G4TwoVector l2) const
00358 {
00359   // Return true if p is on the line through l1, l2 
00360 
00361   if ( l1.x() == l2.x() )
00362   {
00363     return std::fabs(p.x() - l1.x()) < kCarTolerance * 0.5; 
00364   }
00365    G4double  slope= ((l2.y() - l1.y())/(l2.x() - l1.x())); 
00366    G4double predy= l1.y() +  slope *(p.x() - l1.x()); 
00367    G4double dy= p.y() - predy; 
00368 
00369    // Calculate perpendicular distance
00370    //
00371    // G4double perpD= std::fabs(dy) / std::sqrt( 1 + slope * slope ); 
00372    // G4bool   simpleComp= (perpD<0.5*kCarTolerance); 
00373 
00374    // Check perpendicular distance vs tolerance 'directly'
00375    //
00376    const G4double tol= 0.5 * kCarTolerance ; 
00377    G4bool    squareComp=  (dy*dy < (1+slope*slope) * tol * tol); 
00378 
00379    // return  simpleComp; 
00380    return squareComp;
00381 }                    
00382 
00383 //_____________________________________________________________________________
00384 
00385 G4bool G4ExtrudedSolid::IsSameLineSegment(G4TwoVector p,  
00386                                    G4TwoVector l1, G4TwoVector l2) const
00387 {
00388   // Return true if p is on the line through l1, l2 and lies between
00389   // l1 and l2 
00390 
00391   if ( p.x() < std::min(l1.x(), l2.x()) - kCarTolerance * 0.5 || 
00392        p.x() > std::max(l1.x(), l2.x()) + kCarTolerance * 0.5 ||
00393        p.y() < std::min(l1.y(), l2.y()) - kCarTolerance * 0.5 || 
00394        p.y() > std::max(l1.y(), l2.y()) + kCarTolerance * 0.5 )
00395   {
00396     return false;
00397   }
00398 
00399   return IsSameLine(p, l1, l2);
00400 }
00401 
00402 //_____________________________________________________________________________
00403 
00404 G4bool G4ExtrudedSolid::IsSameSide(G4TwoVector p1, G4TwoVector p2, 
00405                                    G4TwoVector l1, G4TwoVector l2) const
00406 {
00407   // Return true if p1 and p2 are on the same side of the line through l1, l2 
00408 
00409   return   ( (p1.x() - l1.x()) * (l2.y() - l1.y())
00410          - (l2.x() - l1.x()) * (p1.y() - l1.y()) )
00411          * ( (p2.x() - l1.x()) * (l2.y() - l1.y())
00412          - (l2.x() - l1.x()) * (p2.y() - l1.y()) ) > 0;
00413 }       
00414 
00415 //_____________________________________________________________________________
00416 
00417 G4bool G4ExtrudedSolid::IsPointInside(G4TwoVector a, G4TwoVector b,
00418                                       G4TwoVector c, G4TwoVector p) const
00419 {
00420   // Return true if p is inside of triangle abc or on its edges, 
00421   // else returns false 
00422 
00423   // Check extent first
00424   //
00425   if ( ( p.x() < a.x() && p.x() < b.x() && p.x() < c.x() ) || 
00426        ( p.x() > a.x() && p.x() > b.x() && p.x() > c.x() ) || 
00427        ( p.y() < a.y() && p.y() < b.y() && p.y() < c.y() ) || 
00428        ( p.y() > a.y() && p.y() > b.y() && p.y() > c.y() ) ) return false;
00429   
00430   G4bool inside 
00431     = IsSameSide(p, a, b, c)
00432       && IsSameSide(p, b, a, c)
00433       && IsSameSide(p, c, a, b);
00434 
00435   G4bool onEdge
00436     = IsSameLineSegment(p, a, b)
00437       || IsSameLineSegment(p, b, c)
00438       || IsSameLineSegment(p, c, a);
00439       
00440   return inside || onEdge;    
00441 }     
00442 
00443 //_____________________________________________________________________________
00444 
00445 G4double 
00446 G4ExtrudedSolid::GetAngle(G4TwoVector po, G4TwoVector pa, G4TwoVector pb) const
00447 {
00448   // Return the angle of the vertex in po
00449 
00450   G4TwoVector t1 = pa - po;
00451   G4TwoVector t2 = pb - po;
00452   
00453   G4double result = (std::atan2(t1.y(), t1.x()) - std::atan2(t2.y(), t2.x()));
00454 
00455   if ( result < 0 ) result += 2*pi;
00456 
00457   return result;
00458 }
00459 
00460 //_____________________________________________________________________________
00461 
00462 G4VFacet*
00463 G4ExtrudedSolid::MakeDownFacet(G4int ind1, G4int ind2, G4int ind3) const
00464 {
00465   // Create a triangular facet from the polygon points given by indices
00466   // forming the down side ( the normal goes in -z)
00467 
00468   std::vector<G4ThreeVector> vertices;
00469   vertices.push_back(GetVertex(0, ind1));
00470   vertices.push_back(GetVertex(0, ind2));
00471   vertices.push_back(GetVertex(0, ind3));
00472   
00473   // first vertex most left
00474   //
00475   G4ThreeVector cross 
00476     = (vertices[1]-vertices[0]).cross(vertices[2]-vertices[1]);
00477   
00478   if ( cross.z() > 0.0 )
00479   {
00480     // vertices ardered clock wise has to be reordered
00481 
00482     // G4cout << "G4ExtrudedSolid::MakeDownFacet: reordering vertices " 
00483     //        << ind1 << ", " << ind2 << ", " << ind3 << G4endl; 
00484 
00485     G4ThreeVector tmp = vertices[1];
00486     vertices[1] = vertices[2];
00487     vertices[2] = tmp;
00488   }
00489   
00490   return new G4TriangularFacet(vertices[0], vertices[1],
00491                                vertices[2], ABSOLUTE);
00492 }      
00493 
00494 //_____________________________________________________________________________
00495 
00496 G4VFacet*
00497 G4ExtrudedSolid::MakeUpFacet(G4int ind1, G4int ind2, G4int ind3) const      
00498 {
00499   // Creates a triangular facet from the polygon points given by indices
00500   // forming the upper side ( z>0 )
00501 
00502   std::vector<G4ThreeVector> vertices;
00503   vertices.push_back(GetVertex(fNz-1, ind1));
00504   vertices.push_back(GetVertex(fNz-1, ind2));
00505   vertices.push_back(GetVertex(fNz-1, ind3));
00506   
00507   // first vertex most left
00508   //
00509   G4ThreeVector cross 
00510     = (vertices[1]-vertices[0]).cross(vertices[2]-vertices[1]);
00511   
00512   if ( cross.z() < 0.0 )
00513   {
00514     // vertices ordered clock wise has to be reordered
00515 
00516     // G4cout << "G4ExtrudedSolid::MakeUpFacet: reordering vertices " 
00517     //        << ind1 << ", " << ind2 << ", " << ind3 << G4endl; 
00518 
00519     G4ThreeVector tmp = vertices[1];
00520     vertices[1] = vertices[2];
00521     vertices[2] = tmp;
00522   }
00523   
00524   return new G4TriangularFacet(vertices[0], vertices[1],
00525                                vertices[2], ABSOLUTE);
00526 }      
00527 
00528 //_____________________________________________________________________________
00529 
00530 G4bool G4ExtrudedSolid::AddGeneralPolygonFacets()
00531 {
00532   // Decompose polygonal sides in triangular facets
00533 
00534   typedef std::pair < G4TwoVector, G4int > Vertex;
00535 
00536   // Fill one more vector
00537   //
00538   std::vector< Vertex > verticesToBeDone;
00539   for ( G4int i=0; i<fNv; ++i )
00540   {
00541     verticesToBeDone.push_back(Vertex(fPolygon[i], i));
00542   }
00543   std::vector< Vertex > ears;
00544   
00545   std::vector< Vertex >::iterator c1 = verticesToBeDone.begin();
00546   std::vector< Vertex >::iterator c2 = c1+1;  
00547   std::vector< Vertex >::iterator c3 = c1+2;  
00548   while ( verticesToBeDone.size()>2 )
00549   {
00550 
00551     // G4cout << "Looking at triangle : "
00552     //        << c1->second << "  " << c2->second
00553     //        << "  " << c3->second << G4endl;  
00554 
00555     // skip concave vertices
00556     //
00557     G4double angle = GetAngle(c2->first, c3->first, c1->first);
00558     //G4cout << "angle " << angle  << G4endl;
00559 
00560     G4int counter = 0;
00561     while ( angle > pi )
00562     {
00563       // G4cout << "Skipping concave vertex " << c2->second << G4endl;
00564 
00565       // try next three consecutive vertices
00566       //
00567       c1 = c2;
00568       c2 = c3;
00569       ++c3; 
00570       if ( c3 == verticesToBeDone.end() ) { c3 = verticesToBeDone.begin(); }
00571 
00572       // G4cout << "Looking at triangle : "
00573       //        << c1->second << "  " << c2->second
00574       //        << "  " << c3->second << G4endl; 
00575       
00576       angle = GetAngle(c2->first, c3->first, c1->first); 
00577       //G4cout << "angle " << angle  << G4endl;
00578       
00579       counter++;
00580       
00581       if ( counter > fNv) {
00582         G4Exception("G4ExtrudedSolid::AddGeneralPolygonFacets",
00583                     "GeomSolids0003", FatalException,
00584                     "Triangularisation has failed.");
00585         break;
00586       }  
00587     }
00588 
00589     G4bool good = true;
00590     std::vector< Vertex >::iterator it;
00591     for ( it=verticesToBeDone.begin(); it != verticesToBeDone.end(); ++it )
00592     {
00593       // skip vertices of tested triangle
00594       //
00595       if ( it == c1 || it == c2 || it == c3 ) { continue; }
00596 
00597       if ( IsPointInside(c1->first, c2->first, c3->first, it->first) )
00598       {
00599         // G4cout << "Point " << it->second << " is inside" << G4endl;
00600         good = false;
00601 
00602         // try next three consecutive vertices
00603         //
00604         c1 = c2;
00605         c2 = c3;
00606         ++c3; 
00607         if ( c3 == verticesToBeDone.end() ) { c3 = verticesToBeDone.begin(); }
00608         break;
00609       }
00610       // else 
00611       //   { G4cout << "Point " << it->second << " is outside" << G4endl; }
00612     }
00613     if ( good )
00614     {
00615       // all points are outside triangle, we can make a facet
00616 
00617       // G4cout << "Found triangle : "
00618       //        << c1->second << "  " << c2->second
00619       //        << "  " << c3->second << G4endl;  
00620 
00621       G4bool result;
00622       result = AddFacet( MakeDownFacet(c1->second, c2->second, c3->second) );
00623       if ( ! result ) { return false; }
00624 
00625       result = AddFacet( MakeUpFacet(c1->second, c2->second, c3->second) );
00626       if ( ! result ) { return false; }
00627 
00628       std::vector<G4int> triangle(3);
00629       triangle[0] = c1->second;
00630       triangle[1] = c2->second;
00631       triangle[2] = c3->second;
00632       fTriangles.push_back(triangle);
00633 
00634       // remove the ear point from verticesToBeDone
00635       //
00636       verticesToBeDone.erase(c2);
00637       c1 = verticesToBeDone.begin();
00638       c2 = c1+1;  
00639       c3 = c1+2;  
00640     } 
00641   }
00642   return true;
00643 }
00644 
00645 //_____________________________________________________________________________
00646 
00647 G4bool G4ExtrudedSolid::MakeFacets()
00648 {
00649   // Define facets
00650 
00651   G4bool good;
00652   
00653   // Decomposition of polygonal sides in the facets
00654   //
00655   if ( fNv == 3 )
00656   {
00657     good = AddFacet( new G4TriangularFacet( GetVertex(0, 0), GetVertex(0, 1),
00658                                             GetVertex(0, 2), ABSOLUTE) );
00659     if ( ! good ) { return false; }
00660 
00661     good = AddFacet( new G4TriangularFacet( GetVertex(fNz-1, 2), GetVertex(fNz-1, 1),
00662                                             GetVertex(fNz-1, 0), ABSOLUTE) );
00663     if ( ! good ) { return false; }
00664     
00665     std::vector<G4int> triangle(3);
00666     triangle[0] = 0;
00667     triangle[1] = 1;
00668     triangle[2] = 2;
00669     fTriangles.push_back(triangle);
00670   }
00671   
00672   else if ( fNv == 4 )
00673   {
00674     good = AddFacet( new G4QuadrangularFacet( GetVertex(0, 0),GetVertex(0, 1),
00675                                               GetVertex(0, 2),GetVertex(0, 3),
00676                                               ABSOLUTE) );
00677     if ( ! good ) { return false; }
00678 
00679     good = AddFacet( new G4QuadrangularFacet( GetVertex(fNz-1, 3), GetVertex(fNz-1, 2), 
00680                                               GetVertex(fNz-1, 1), GetVertex(fNz-1, 0),
00681                                               ABSOLUTE) );
00682     if ( ! good ) { return false; }
00683 
00684     std::vector<G4int> triangle1(3);
00685     triangle1[0] = 0;
00686     triangle1[1] = 1;
00687     triangle1[2] = 2;
00688     fTriangles.push_back(triangle1);
00689 
00690     std::vector<G4int> triangle2(3);
00691     triangle2[0] = 0;
00692     triangle2[1] = 2;
00693     triangle2[2] = 3;
00694     fTriangles.push_back(triangle2);
00695   }  
00696   else
00697   {
00698     good = AddGeneralPolygonFacets();
00699     if ( ! good ) { return false; }
00700   }
00701     
00702   // The quadrangular sides
00703   //
00704   for ( G4int iz = 0; iz < fNz-1; ++iz ) 
00705   {
00706     for ( G4int i = 0; i < fNv; ++i )
00707     {
00708       G4int j = (i+1) % fNv;
00709       good = AddFacet( new G4QuadrangularFacet
00710                         ( GetVertex(iz, j), GetVertex(iz, i), 
00711                           GetVertex(iz+1, i), GetVertex(iz+1, j), ABSOLUTE) );
00712       if ( ! good ) { return false; }
00713     }
00714   }  
00715 
00716   SetSolidClosed(true);
00717 
00718   return good;
00719 }
00720 
00721 //_____________________________________________________________________________
00722 
00723 G4bool G4ExtrudedSolid::IsConvex() const
00724 {
00725   // Get polygon convexity (polygon is convex if all vertex angles are < pi )
00726 
00727   for ( G4int i=0; i< fNv; ++i )
00728   {
00729     G4int j = ( i + 1 ) % fNv;
00730     G4int k = ( i + 2 ) % fNv;
00731     G4TwoVector v1 = fPolygon[i]-fPolygon[j];
00732     G4TwoVector v2 = fPolygon[k]-fPolygon[j];
00733     G4double dphi = v2.phi() - v1.phi();
00734     if ( dphi < 0. )  { dphi += 2.*pi; }
00735     
00736     if ( dphi >= pi ) { return false; }
00737   }
00738   
00739   return true;
00740 }     
00741 
00742 //_____________________________________________________________________________
00743 
00744 G4GeometryType G4ExtrudedSolid::GetEntityType () const
00745 {
00746   // Return entity type
00747 
00748   return fGeometryType;
00749 }
00750 
00751 //_____________________________________________________________________________
00752 
00753 G4VSolid* G4ExtrudedSolid::Clone() const
00754 {
00755   return new G4ExtrudedSolid(*this);
00756 }
00757 
00758 //_____________________________________________________________________________
00759 
00760 EInside G4ExtrudedSolid::Inside (const G4ThreeVector &p) const
00761 {
00762   // Override the base class function  as it fails in case of concave polygon.
00763   // Project the point in the original polygon scale and check if it is inside
00764   // for each triangle.
00765 
00766   // Check first if outside extent
00767   //
00768   if ( p.x() < GetMinXExtent() - kCarTolerance * 0.5 ||
00769        p.x() > GetMaxXExtent() + kCarTolerance * 0.5 ||
00770        p.y() < GetMinYExtent() - kCarTolerance * 0.5 ||
00771        p.y() > GetMaxYExtent() + kCarTolerance * 0.5 ||
00772        p.z() < GetMinZExtent() - kCarTolerance * 0.5 ||
00773        p.z() > GetMaxZExtent() + kCarTolerance * 0.5 )
00774   {
00775     // G4cout << "G4ExtrudedSolid::Outside extent: " << p << G4endl;
00776     return kOutside;
00777   }  
00778 
00779   // Project point p(z) to the polygon scale p0
00780   //
00781   G4TwoVector pscaled = ProjectPoint(p);
00782   
00783   // Check if on surface of polygon
00784   //
00785   for ( G4int i=0; i<fNv; ++i )
00786   {
00787     G4int j = (i+1) % fNv;
00788     if ( IsSameLineSegment(pscaled, fPolygon[i], fPolygon[j]) )
00789     {
00790       // G4cout << "G4ExtrudedSolid::Inside return Surface (on polygon) "
00791       //        << G4endl;
00792 
00793       return kSurface;
00794     }  
00795   }   
00796 
00797   // Now check if inside triangles
00798   //
00799   std::vector< std::vector<G4int> >::const_iterator it = fTriangles.begin();
00800   G4bool inside = false;
00801   do
00802   {
00803     if ( IsPointInside(fPolygon[(*it)[0]], fPolygon[(*it)[1]],
00804                        fPolygon[(*it)[2]], pscaled) )  { inside = true; }
00805     ++it;
00806   } while ( (inside == false) && (it != fTriangles.end()) );
00807   
00808   if ( inside )
00809   {
00810     // Check if on surface of z sides
00811     //
00812     if ( std::fabs( p.z() - fZSections[0].fZ ) < kCarTolerance * 0.5 ||
00813          std::fabs( p.z() - fZSections[fNz-1].fZ ) < kCarTolerance * 0.5 )
00814     {
00815       // G4cout << "G4ExtrudedSolid::Inside return Surface (on z side)"
00816       //        << G4endl;
00817 
00818       return kSurface;
00819     }  
00820   
00821     // G4cout << "G4ExtrudedSolid::Inside return Inside" << G4endl;
00822 
00823     return kInside;
00824   }  
00825                             
00826   // G4cout << "G4ExtrudedSolid::Inside return Outside " << G4endl;
00827 
00828   return kOutside; 
00829 }  
00830 
00831 //_____________________________________________________________________________
00832 
00833 G4double G4ExtrudedSolid::DistanceToOut (const G4ThreeVector &p,
00834                                          const G4ThreeVector &v,
00835                                          const G4bool calcNorm,
00836                                                G4bool *validNorm,
00837                                                G4ThreeVector *n) const
00838 {
00839   // Override the base class function to redefine validNorm
00840   // (the solid can be concave) 
00841 
00842   G4double distOut =
00843     G4TessellatedSolid::DistanceToOut(p, v, calcNorm, validNorm, n);
00844   if (validNorm) { *validNorm = fIsConvex; }
00845 
00846   return distOut;
00847 }
00848 
00849 
00850 //_____________________________________________________________________________
00851 
00852 G4double G4ExtrudedSolid::DistanceToOut (const G4ThreeVector &p) const
00853 {
00854   // Override the overloaded base class function
00855 
00856   return G4TessellatedSolid::DistanceToOut(p);
00857 }
00858 
00859 //_____________________________________________________________________________
00860 
00861 std::ostream& G4ExtrudedSolid::StreamInfo(std::ostream &os) const
00862 {
00863   G4int oldprc = os.precision(16);
00864   os << "-----------------------------------------------------------\n"
00865      << "    *** Dump for solid - " << GetName() << " ***\n"
00866      << "    ===================================================\n"
00867      << " Solid geometry type: " << fGeometryType  << G4endl;
00868 
00869   if ( fIsConvex) 
00870     { os << " Convex polygon; list of vertices:" << G4endl; }
00871   else  
00872     { os << " Concave polygon; list of vertices:" << G4endl; }
00873   
00874   for ( G4int i=0; i<fNv; ++i )
00875   {
00876     os << std::setw(5) << "#" << i 
00877        << "   vx = " << fPolygon[i].x()/mm << " mm" 
00878        << "   vy = " << fPolygon[i].y()/mm << " mm" << G4endl;
00879   }
00880   
00881   os << " Sections:" << G4endl;
00882   for ( G4int iz=0; iz<fNz; ++iz ) 
00883   {
00884     os << "   z = "   << fZSections[iz].fZ/mm          << " mm  "
00885        << "  x0= "    << fZSections[iz].fOffset.x()/mm << " mm  "
00886        << "  y0= "    << fZSections[iz].fOffset.y()/mm << " mm  " 
00887        << "  scale= " << fZSections[iz].fScale << G4endl;
00888   }     
00889 
00890 /*
00891   // Triangles (for debugging)
00892   os << G4endl; 
00893   os << " Triangles:" << G4endl;
00894   os << " Triangle #   vertex1   vertex2   vertex3" << G4endl;
00895 
00896   G4int counter = 0;
00897   std::vector< std::vector<G4int> >::const_iterator it;
00898   for ( it = fTriangles.begin(); it != fTriangles.end(); it++ ) {
00899      std::vector<G4int> triangle = *it;
00900      os << std::setw(10) << counter++ 
00901         << std::setw(10) << triangle[0] << std::setw(10)  << triangle[1]  << std::setw(10)  << triangle[2] 
00902         << G4endl;
00903   }          
00904 */
00905   os.precision(oldprc);
00906 
00907   return os;
00908 }  

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