61 const G4int    G4GenericTrap::fgkNofVertices = 8;
 
   62 const G4double G4GenericTrap::fgkTolerance = 1E-3;
 
   67                               const std::vector<G4TwoVector>&  vertices )
 
   85   G4String errorDescription = 
"InvalidSetup in \" ";
 
   86   errorDescription += 
name;
 
   87   errorDescription += 
"\""; 
 
   93   if ( 
G4int(vertices.size()) != fgkNofVertices )
 
   95     G4Exception(
"G4GenericTrap::G4GenericTrap()", 
"GeomSolids0002",
 
  103      G4Exception(
"G4GenericTrap::G4GenericTrap()", 
"GeomSolids0002",
 
  109   if(CheckOrder(vertices))
 
  111     for (
G4int i=0; i<fgkNofVertices; ++i) {fVertices.push_back(vertices[i]);}
 
  115     for (
G4int i=0; i <4; ++i) {fVertices.push_back(vertices[3-i]);}
 
  116     for (
G4int i=0; i <4; ++i) {fVertices.push_back(vertices[7-i]);}
 
  121   for (
G4int j=0; j < 2; j++)
 
  123     for (
G4int i=1; i<4; ++i)
 
  126       length = (fVertices[k]-fVertices[k-1]).mag();
 
  129         std::ostringstream message;
 
  130         message << 
"Length segment is too small." << 
G4endl 
  131                 << 
"Distance between " << fVertices[k-1] << 
" and " 
  132                 << fVertices[k] << 
" is only " << length << 
" mm !"; 
 
  133         G4Exception(
"G4GenericTrap::G4GenericTrap()", 
"GeomSolids1001",
 
  134                     JustWarning, message, 
"Vertices will be collapsed.");
 
  135         fVertices[k]=fVertices[k-1];
 
  142   for( 
G4int i=0; i<4; i++) { fTwist[i]=0.; }
 
  143   fIsTwisted = ComputeIsTwisted();
 
  153    if ( !fIsTwisted )  { fTessellatedSolid = CreateTessellatedSolid(); }
 
  162     halfCarTolerance(0.),
 
  166     fTessellatedSolid(0),
 
  176   for (
size_t i=0; i<4; ++i)  { fTwist[i]=0.; }
 
  184   delete fTessellatedSolid;
 
  191     fpPolyhedron(0), halfCarTolerance(rhs.halfCarTolerance),
 
  192     fDz(rhs.fDz), fVertices(rhs.fVertices),
 
  193     fIsTwisted(rhs.fIsTwisted), fTessellatedSolid(0),
 
  194     fMinBBoxVector(rhs.fMinBBoxVector), fMaxBBoxVector(rhs.fMaxBBoxVector),
 
  195     fVisSubdivisions(rhs.fVisSubdivisions),
 
  196     fSurfaceArea(rhs.fSurfaceArea), fCubicVolume(rhs.fCubicVolume) 
 
  198    for (
size_t i=0; i<4; ++i)  { fTwist[i] = rhs.fTwist[i]; }
 
  200    if (rhs.fTessellatedSolid && !fIsTwisted )
 
  201    { fTessellatedSolid = CreateTessellatedSolid(); } 
 
  211    if (
this == &rhs)  { 
return *
this; }
 
  219    fpPolyhedron = 0; halfCarTolerance = rhs.halfCarTolerance;
 
  220    fDz = rhs.fDz; fVertices = rhs.fVertices;
 
  221    fIsTwisted = rhs.fIsTwisted; fTessellatedSolid = 0;
 
  222    fMinBBoxVector = rhs.fMinBBoxVector; fMaxBBoxVector = rhs.fMaxBBoxVector;
 
  223    fVisSubdivisions = rhs.fVisSubdivisions;
 
  224    fSurfaceArea = rhs.fSurfaceArea; fCubicVolume = rhs.fCubicVolume;
 
  226    for (
size_t i=0; i<4; ++i)  { fTwist[i] = rhs.fTwist[i]; }
 
  228    if (rhs.fTessellatedSolid && !fIsTwisted )
 
  229    { 
delete fTessellatedSolid; fTessellatedSolid = CreateTessellatedSolid(); } 
 
  239                               const std::vector<G4TwoVector>& poly)
 const  
  245   for (
G4int i = 0; i < 4; i++)
 
  249     cross = (p.
x()-poly[i].x())*(poly[j].
y()-poly[i].y())-
 
  250             (p.
y()-poly[i].y())*(poly[j].
x()-poly[i].x());
 
  252     len2=(poly[i]-poly[j]).mag2();
 
  255       if(cross*cross<=len2*halfCarTolerance*halfCarTolerance)  
 
  259         if ( poly[i].
y() > poly[j].y() )  { k=i; l=j; }
 
  262         if ( poly[k].
x() != poly[l].x() )
 
  264           test = (p.
x()-poly[l].x())/(poly[k].
x()-poly[l].x())
 
  265                * (poly[k].
y()-poly[l].y())+poly[l].
y();
 
  274         if( (test>=(poly[l].
y()-halfCarTolerance))
 
  275          && (test<=(poly[k].
y()+halfCarTolerance)) )
 
  284       else if (cross<0.)  { 
return kOutside; }
 
  296     if ( (std::fabs(p.
x()-poly[0].x())+std::fabs(p.
y()-poly[0].y())) > halfCarTolerance )
 
  311    if ( fTessellatedSolid )
 
  313      return fTessellatedSolid->
Inside(p);
 
  318   std::vector<G4TwoVector> xy;
 
  320   if (std::fabs(p.
z()) <= fDz+halfCarTolerance)  
 
  325     for (
G4int i=0; i<4; i++)
 
  327       xy.push_back(fVertices[i+4]+cf*( fVertices[i]-fVertices[i+4]));
 
  330     innew=InsidePolygone(p,xy);
 
  334       if(std::fabs(p.
z()) > fDz-halfCarTolerance)  { innew=
kSurface; }
 
  348   if ( fTessellatedSolid )
 
  355                 p0, p1, p2, r1, r2, r3, r4;
 
  356   G4int noSurfaces = 0; 
 
  360   distz = fDz-std::fabs(p.
z());
 
  361   if (distz < halfCarTolerance)
 
  377   std:: vector<G4TwoVector> vertices;  
 
  379   for (
G4int i=0; i<4; i++)
 
  381     vertices.push_back(fVertices[i+4]+cf*(fVertices[i]-fVertices[i+4]));
 
  386   for (
G4int q=0; q<4; q++)
 
  409         p2=
G4ThreeVector(fVertices[(q+1)%4+4].
x(),fVertices[(q+1)%4+4].
y(),fDz);
 
  412     lnorm = (p1-p0).cross(p2-p0);
 
  413     lnorm = lnorm.
unit();
 
  414     if(zPlusSide)  { lnorm=-lnorm; }
 
  423         G4double proj=(p-p0).dot(p2-p0)/normP;
 
  424         if(proj<0)     { proj=0; }
 
  425         if(proj>normP) { proj=normP; }
 
  431         r1=r1+proj*(r2-r1)/normP;
 
  432         r3=r3+proj*(r4-r3)/normP;
 
  439     distxy=std::fabs((p0-p).dot(lnorm));
 
  440     if ( distxy<halfCarTolerance )
 
  446       sumnorm=sumnorm+lnorm;
 
  460   if ( noSurfaces == 0 )
 
  462     G4Exception(
"G4GenericTrap::SurfaceNormal(p)", 
"GeomSolids1002",
 
  467   else if ( noSurfaces == 1 )  { sumnorm = sumnorm; }
 
  468   else                         { sumnorm = sumnorm.unit(); }
 
  476                                             const G4int ipl )
 const 
  481   if ( fTessellatedSolid )
 
  497   u=fVertices[i+4]+cf*(fVertices[i]-fVertices[i+4]);
 
  498   v=fVertices[j+4]+cf*(fVertices[j]-fVertices[j+4]);
 
  504   if (std::fabs(distz)<halfCarTolerance)
 
  517     if ( std::fabs(p.
z()+fDz) > halfCarTolerance )
 
  526   lnorm=-(p1-p0).cross(p2-p0);
 
  527   if (distz>-halfCarTolerance)  { lnorm=-lnorm.
unit(); }
 
  528   else                          { lnorm=lnorm.
unit();  }
 
  537       G4double proj=(p-p0).dot(p2-p0)/normP;
 
  538       if (proj<0)     { proj=0; }
 
  539       if (proj>normP) { proj=normP; }
 
  545       r1=r1+proj*(r2-r1)/normP;
 
  546       r3=r3+proj*(r4-r3)/normP;
 
  560                                     const G4int ipl)
 const  
  572   xa=fVertices[ipl].x();
 
  573   ya=fVertices[ipl].y();
 
  574   xb=fVertices[ipl+4].x();
 
  575   yb=fVertices[ipl+4].y();
 
  578   xd=fVertices[4+j].x();
 
  579   yd=fVertices[4+j].y();
 
  596   G4double a = (dtx*v.
y()-dty*v.
x()+(tx1*ty2-tx2*ty1)*v.
z())*v.
z();
 
  597   G4double b = dxs*v.
y()-dys*v.
x()+(dtx*p.
y()-dty*p.
x()+ty2*xs1-ty1*xs2
 
  598              + tx1*ys2-tx2*ys1)*v.
z();
 
  610     if (q>-halfCarTolerance)
 
  612       if (q<halfCarTolerance)
 
  614         if (NormalToPlane(p,ipl).dot(v)<=0)
 
  617           { 
return kInfinity; }
 
  623       if (std::fabs(zi)<fDz)
 
  631         zi = (xp-
x1)*(xp-x2)+(yp-
y1)*(yp-y2);
 
  632         if (zi<=halfCarTolerance)  { 
return q; }
 
  640     if (a>0) { q=0.5*(-b-std::sqrt(d))/a; }
 
  641     else     { q=0.5*(-b+std::sqrt(d))/a; }
 
  645     if (q>-halfCarTolerance)
 
  647       if(q<halfCarTolerance)
 
  649         if (NormalToPlane(p,ipl).
dot(v)<=0)
 
  655           if (a>0) { q=0.5*(-b+std::sqrt(d))/a; }
 
  656           else     { q=0.5*(-b-std::sqrt(d))/a; }
 
  657           if (q<=halfCarTolerance) { 
return kInfinity; }
 
  663       if (std::fabs(zi)<fDz)
 
  671         zi = (xp-
x1)*(xp-x2)+(yp-
y1)*(yp-y2);
 
  672         if (zi<=halfCarTolerance)  { 
return q; }
 
  675     if (a>0)  { q=0.5*(-b+std::sqrt(d))/a; }
 
  676     else      { q=0.5*(-b-std::sqrt(d))/a; }
 
  680     if (q>-halfCarTolerance)
 
  682       if(q<halfCarTolerance)
 
  684         if (NormalToPlane(p,ipl).
dot(v)<=0)
 
  690           if (a>0) { q=0.5*(-b-std::sqrt(d))/a; }
 
  691           else     { q=0.5*(-b+std::sqrt(d))/a; }
 
  692           if (q<=halfCarTolerance)  { 
return kInfinity; }
 
  698       if (std::fabs(zi)<fDz)
 
  706         zi = (xp-
x1)*(xp-x2)+(yp-
y1)*(yp-y2);
 
  707         if (zi<=halfCarTolerance)  { 
return q; }
 
  720   if ( fTessellatedSolid )
 
  734     dist[i]=DistToPlane(p, v, i);  
 
  740   if (std::fabs(p.
z())>fDz-halfCarTolerance)
 
  747         dist[4] = (fDz-p.
z())/v.
z();
 
  751         dist[4] = (-fDz-p.
z())/v.
z();
 
  753       if (dist[4]<-halfCarTolerance)
 
  759         if(dist[4]<halfCarTolerance)
 
  763           if (n.
dot(v)<0) { dist[4]=0.; }
 
  764           else            { dist[4]=kInfinity; }
 
  774     if (dist[i] < distmin)  { distmin = dist[i]; }
 
  777   if (distmin<halfCarTolerance)  { distmin=0.; }
 
  789   if ( fTessellatedSolid )
 
  796   if(safz<0) { safz=0; }
 
  802   for (iseg=0; iseg<4; iseg++)
 
  804     safxy = SafetyToFace(p,iseg);
 
  805     if (safxy>safe)  { safe=safxy; }
 
  824   norm=NormalToPlane(p,iseg);
 
  825   safe = (p-p1).dot(norm); 
 
  846   if ( (std::fabs(xa-xc)+std::fabs(ya-yc)) < halfCarTolerance )
 
  848     xc=fVertices[j+4].x();
 
  849     yc=fVertices[j+4].y();
 
  855     if ( (std::fabs(xb-xc)+std::fabs(yb-yc)) < halfCarTolerance )
 
  862   G4double c=(xb-xa)*(yc-ya)-(xc-xa)*(yb-ya);
 
  868     t=-(a*p.
x()+b*p.
y()+c*p.
z()+
d)/t;
 
  870   if ( (t<halfCarTolerance) && (t>-halfCarTolerance) )
 
  895   if ( fTessellatedSolid )
 
  897     return fTessellatedSolid->
DistanceToOut(p, v, calcNorm, validNorm, n);
 
  902   G4bool lateral_cross = 
false;
 
  903   ESide side = kUndefined;
 
  905   if (calcNorm)  { *validNorm=
true; } 
 
  909     distmin=(-fDz-p.
z())/v.
z();
 
  916       distmin = (fDz-p.
z())/v.
z(); 
 
  919     else            { distmin = kInfinity; }
 
  926   for (
G4int ipl=0; ipl<4; ipl++)
 
  929     xa=fVertices[ipl].x();
 
  930     ya=fVertices[ipl].y();
 
  931     xb=fVertices[ipl+4].x();
 
  932     yb=fVertices[ipl+4].y();
 
  935     xd=fVertices[4+j].x();
 
  936     yd=fVertices[4+j].y();
 
  938     if ( ((std::fabs(xb-xd)+std::fabs(yb-yd))<halfCarTolerance)
 
  939       || ((std::fabs(xa-xc)+std::fabs(ya-yc))<halfCarTolerance) )
 
  941       G4double q=DistToTriangle(p,v,ipl) ;
 
  942       if ( (q>=0) && (q<distmin) )
 
  963     G4double a = (dtx*v.
y()-dty*v.
x()+(tx1*ty2-tx2*ty1)*v.
z())*v.
z();
 
  964     G4double b = dxs*v.
y()-dys*v.
x()+(dtx*p.
y()-dty*p.
x()+ty2*xs1-ty1*xs2
 
  965                + tx1*ys2-tx2*ys1)*v.
z();
 
  966     G4double c=dxs*p.
y()-dys*p.
x()+xs1*ys2-xs2*ys1;
 
  976       if ((q > -halfCarTolerance) && (q < distmin))
 
  978         if (q < halfCarTolerance)
 
  980           if (NormalToPlane(p,ipl).dot(v)<0.)  { 
continue; }
 
  991       if (a > 0)  { q=0.5*(-b-std::sqrt(d))/a; }
 
  992       else        { q=0.5*(-b+std::sqrt(d))/a; }
 
  996       if (q > -halfCarTolerance )
 
 1000           if(q < halfCarTolerance)
 
 1002             if (NormalToPlane(p,ipl).
dot(v)<0.)  
 
 1004               if (a > 0)  { q=0.5*(-b+std::sqrt(d))/a; }
 
 1005               else        { q=0.5*(-b-std::sqrt(d))/a; }
 
 1006               if (( q > halfCarTolerance) && (q < distmin))
 
 1009                 lateral_cross = 
true;
 
 1016           lateral_cross = 
true;
 
 1022         if (a > 0)  { q=0.5*(-b+std::sqrt(d))/a; }
 
 1023         else        { q=0.5*(-b-std::sqrt(d))/a; }
 
 1027         if ((q > -halfCarTolerance) && (q < distmin))
 
 1029           if (q < halfCarTolerance)
 
 1031             if (NormalToPlane(p,ipl).
dot(v)<0.)  
 
 1033               if (a > 0)  { q=0.5*(-b-std::sqrt(d))/a; }
 
 1034               else        { q=0.5*(-b+std::sqrt(d))/a; }
 
 1035               if ( ( q > halfCarTolerance) && (q < distmin) )
 
 1038                 lateral_cross = 
true;
 
 1045           lateral_cross = 
true;
 
 1059     if (v.z()>0.) { i=4; }
 
 1060     std::vector<G4TwoVector> xy;
 
 1061     for ( 
G4int j=0; j<4; j++)  { xy.push_back(fVertices[i+j]); }
 
 1065     if (InsidePolygone(pt,xy)==
kOutside)
 
 1076       if(v.z()>0) {side=kPZ;}
 
 1087         *n=NormalToPlane(pt,0);
 
 1090         *n=NormalToPlane(pt,1);
 
 1093         *n=NormalToPlane(pt,2);
 
 1096         *n=NormalToPlane(pt,3);
 
 1106         std::ostringstream message;
 
 1107         G4int oldprc = message.precision(16);
 
 1108         message << 
"Undefined side for valid surface normal to solid." << 
G4endl 
 1110                 << 
"  p.x() = "   << p.
x()/
mm << 
" mm" << 
G4endl 
 1111                 << 
"  p.y() = "   << p.
y()/
mm << 
" mm" << 
G4endl 
 1112                 << 
"  p.z() = "   << p.
z()/
mm << 
" mm" << 
G4endl 
 1113                 << 
"Direction:" << 
G4endl 
 1114                 << 
"  v.x() = "   << v.
x() << 
G4endl 
 1115                 << 
"  v.y() = "   << v.
y() << 
G4endl 
 1116                 << 
"  v.z() = "   << v.
z() << 
G4endl 
 1117                 << 
"Proposed distance :" << 
G4endl 
 1118                 << 
"  distmin = " << distmin/
mm << 
" mm";
 
 1119         message.precision(oldprc);
 
 1120         G4Exception(
"G4GenericTrap::DistanceToOut(p,v,..)",
 
 1126   if (distmin<halfCarTolerance)  { distmin=0.; }
 
 1137   if ( fTessellatedSolid )
 
 1144   if (safz<0) { safz = 0; }
 
 1149   for (
G4int iseg=0; iseg<4; iseg++)
 
 1151     safxy = std::fabs(SafetyToFace(p,iseg));
 
 1152     if (safxy < safe)  { safe = safxy; }
 
 1166   if ( fTessellatedSolid )
 
 1169                                               pTransform, pMin, pMax);
 
 1178   Dx = 0.5*(maxVec.
x()- minVec.
x());
 
 1179   Dy = 0.5*(maxVec.
y()- minVec.
y());
 
 1286     G4bool existsAfterClip=
false;
 
 1294     vertices=CreateRotatedVertices(pTransform);
 
 1299     if ( (pMin!=kInfinity) || (pMax!=-kInfinity) )
 
 1301       existsAfterClip=
true;
 
 1322         existsAfterClip=
true;
 
 1328     return existsAfterClip;
 
 1351     vertices->reserve(8);
 
 1372     G4Exception(
"G4GenericTrap::CreateRotatedVertices()", 
"FatalError",
 
 1396   G4int oldprc = os.precision(16);
 
 1397   os << 
"-----------------------------------------------------------\n" 
 1398      << 
"    *** Dump for solid - " << 
GetName() << 
" *** \n" 
 1399      << 
"    =================================================== \n" 
 1401      << 
"   half length Z: " << fDz/
mm << 
" mm \n" 
 1402      << 
"   list of vertices:\n";
 
 1404   for ( 
G4int i=0; i<fgkNofVertices; ++i )
 
 1406     os << std::setw(5) << 
"#" << i 
 
 1407        << 
"   vx = " << fVertices[i].x()/
mm << 
" mm"  
 1408        << 
"   vy = " << fVertices[i].y()/
mm << 
" mm" << 
G4endl;
 
 1410   os.precision(oldprc);
 
 1421   if ( fTessellatedSolid )
 
 1429   G4double rand,area,chose,cf,lambda0,lambda1,alfa,beta,zp;
 
 1432   std::vector<G4ThreeVector> vertices;
 
 1433   for (
G4int i=0; i<4;i++)
 
 1435     vertices.push_back(
G4ThreeVector(fVertices[i].
x(),fVertices[i].
y(),-fDz));
 
 1437   for (
G4int i=4; i<8;i++)
 
 1439     vertices.push_back(
G4ThreeVector(fVertices[i].
x(),fVertices[i].
y(),fDz));
 
 1444   G4double Surface0=GetFaceSurfaceArea(vertices[0],vertices[1],
 
 1445                                        vertices[2],vertices[3]);
 
 1446   G4double Surface1=GetFaceSurfaceArea(vertices[0],vertices[1],
 
 1447                                        vertices[5],vertices[4]);
 
 1448   G4double Surface2=GetFaceSurfaceArea(vertices[3],vertices[0],
 
 1449                                        vertices[4],vertices[7]);
 
 1450   G4double Surface3=GetFaceSurfaceArea(vertices[2],vertices[3],
 
 1451                                        vertices[7],vertices[6]);
 
 1452   G4double Surface4=GetFaceSurfaceArea(vertices[2],vertices[1],
 
 1453                                        vertices[5],vertices[6]);
 
 1454   G4double Surface5=GetFaceSurfaceArea(vertices[4],vertices[5],
 
 1455                                        vertices[6],vertices[7]);
 
 1457   area = Surface0+Surface1+Surface2+Surface3+Surface4+Surface5;
 
 1460   if ( ( chose < Surface0)
 
 1461     || ( chose > (Surface0+Surface1+Surface2+Surface3+Surface4)) )
 
 1465     if(chose < Surface0)
 
 1468       u = fVertices[ipl]; v = fVertices[j];
 
 1469       w = fVertices[(ipl+3)%4]; 
 
 1474       u = fVertices[ipl+4]; v = fVertices[j+4];
 
 1475       w = fVertices[(ipl+3)%4+4]; 
 
 1480     lambda0=alfa-lambda1;
 
 1483     v = u+lambda0*v+lambda1*w;
 
 1487     if (chose < Surface0+Surface1) { ipl=0; }
 
 1488     else if (chose < Surface0+Surface1+Surface2) { ipl=1; }
 
 1489     else if (chose < Surface0+Surface1+Surface2+Surface3) { ipl=2; }
 
 1493     cf = 0.5*(fDz-zp)/fDz;
 
 1494     u = fVertices[ipl+4]+cf*( fVertices[ipl]-fVertices[ipl+4]);
 
 1495     v = fVertices[j+4]+cf*(fVertices[j]-fVertices[j+4]); 
 
 1507   if(fCubicVolume != 0.) {;}
 
 1509   return fCubicVolume;
 
 1516   if(fSurfaceArea != 0.) {;}
 
 1519     std::vector<G4ThreeVector> vertices;
 
 1520     for (
G4int i=0; i<4;i++)
 
 1522       vertices.push_back(
G4ThreeVector(fVertices[i].
x(),fVertices[i].
y(),-fDz));
 
 1524     for (
G4int i=4; i<8;i++)
 
 1526       vertices.push_back(
G4ThreeVector(fVertices[i].
x(),fVertices[i].
y(),fDz));
 
 1531     G4double fSurface0=GetFaceSurfaceArea(vertices[0],vertices[1],
 
 1532                                           vertices[2],vertices[3]);
 
 1533     G4double fSurface1=GetFaceSurfaceArea(vertices[0],vertices[1],
 
 1534                                           vertices[5],vertices[4]);
 
 1535     G4double fSurface2=GetFaceSurfaceArea(vertices[3],vertices[0],
 
 1536                                           vertices[4],vertices[7]);
 
 1537     G4double fSurface3=GetFaceSurfaceArea(vertices[2],vertices[3],
 
 1538                                           vertices[7],vertices[6]);
 
 1539     G4double fSurface4=GetFaceSurfaceArea(vertices[2],vertices[1],
 
 1540                                           vertices[5],vertices[6]);
 
 1541     G4double fSurface5=GetFaceSurfaceArea(vertices[4],vertices[5],
 
 1542                                           vertices[6],vertices[7]);
 
 1548       fSurfaceArea = fSurface0+fSurface1+fSurface2
 
 1549                    + fSurface3+fSurface4+fSurface5;
 
 1556   return fSurfaceArea;
 
 1577   aOne = 0.5*Area.
mag();
 
 1580   aTwo = 0.5*Area.
mag();
 
 1587 G4bool G4GenericTrap::ComputeIsTwisted() 
 
 1594   G4int nv = fgkNofVertices/2;
 
 1596   for ( 
G4int i=0; i<4; i++ )
 
 1598     dx1 = fVertices[(i+1)%nv].
x()-fVertices[i].x();
 
 1599     dy1 = fVertices[(i+1)%nv].
y()-fVertices[i].y();
 
 1600     if ( (dx1 == 0) && (dy1 == 0) ) { 
continue; }
 
 1602     dx2 = fVertices[nv+(i+1)%nv].
x()-fVertices[nv+i].x();
 
 1603     dy2 = fVertices[nv+(i+1)%nv].
y()-fVertices[nv+i].y();
 
 1605     if ( dx2 == 0 && dy2 == 0 ) { 
continue; }
 
 1606     G4double twist_angle = std::fabs(dy1*dx2 - dx1*dy2);        
 
 1607     if ( twist_angle < fgkTolerance ) { 
continue; }
 
 1609     SetTwistAngle(i,twist_angle);
 
 1613     twist_angle = std::acos( (dx1*dx2 + dy1*dy2)
 
 1614                            / (std::sqrt(dx1*dx1+dy1*dy1)
 
 1615                              * std::sqrt(dx2*dx2+dy2*dy2)) );
 
 1619       std::ostringstream message;
 
 1620       message << 
"Twisted Angle is bigger than 90 degrees - " << 
GetName()
 
 1622               << 
"     Potential problem of malformed Solid !" << 
G4endl 
 1623               << 
"     TwistANGLE = " << twist_angle
 
 1624               << 
"*rad  for lateral plane N= " << i;
 
 1625       G4Exception(
"G4GenericTrap::ComputeIsTwisted()", 
"GeomSolids1002",
 
 1635 G4bool G4GenericTrap::CheckOrder(
const std::vector<G4TwoVector>& vertices)
 const 
 1640   G4bool clockwise_order=
true;
 
 1645   for (
G4int i=0; i<4; i++)
 
 1648     sum1 += vertices[i].x()*vertices[j].y() - vertices[j].x()*vertices[i].y();
 
 1649     sum2 += vertices[i+4].x()*vertices[j+4].y()
 
 1650           - vertices[j+4].x()*vertices[i+4].y();
 
 1652   if (sum1*sum2 < -fgkTolerance)
 
 1654      std::ostringstream message;
 
 1655      message << 
"Lower/upper faces defined with opposite clockwise - " 
 1657      G4Exception(
"G4GenericTrap::CheckOrder()", 
"GeomSolids0002",
 
 1661    if ((sum1 > 0.)||(sum2 > 0.))
 
 1663      std::ostringstream message;
 
 1664      message << 
"Vertices must be defined in clockwise XY planes - " 
 1666      G4Exception(
"G4GenericTrap::CheckOrder()", 
"GeomSolids1001",
 
 1668      clockwise_order = 
false;
 
 1673    G4bool illegal_cross = 
false;
 
 1674    illegal_cross = IsSegCrossingZ(vertices[0],vertices[4],
 
 1675                                   vertices[1],vertices[5]);
 
 1679      illegal_cross = IsSegCrossingZ(vertices[2],vertices[6],
 
 1680                                     vertices[3],vertices[7]);
 
 1685      illegal_cross = IsSegCrossing(vertices[0],vertices[1],
 
 1686                                    vertices[2],vertices[3]);
 
 1690      illegal_cross = IsSegCrossing(vertices[0],vertices[3],
 
 1691                                    vertices[1],vertices[2]);
 
 1695      illegal_cross = IsSegCrossing(vertices[4],vertices[5],
 
 1696                                    vertices[6],vertices[7]);
 
 1700      illegal_cross = IsSegCrossing(vertices[4],vertices[7],
 
 1701                                    vertices[5],vertices[6]);
 
 1706       std::ostringstream message;
 
 1707       message << 
"Malformed polygone with opposite sides - " << 
GetName();
 
 1708       G4Exception(
"G4GenericTrap::CheckOrderAndSetup()",
 
 1711    return clockwise_order;
 
 1716 void G4GenericTrap::ReorderVertices(std::vector<G4ThreeVector>& vertices)
 const 
 1720   std::vector<G4ThreeVector> oldVertices(vertices);
 
 1722   for ( 
G4int i=0; i < 
G4int(oldVertices.size()); ++i )
 
 1724     vertices[i] = oldVertices[oldVertices.size()-1-i];
 
 1738   G4double dx1,dx2,xm=0.,ym=0.,a1=0.,a2=0.,b1=0.,b2=0.;
 
 1742   if( std::fabs(dx1) < fgkTolerance )  { stand1 = 
true; }
 
 1743   if( std::fabs(dx2) < fgkTolerance )  { stand2 = 
true; }
 
 1746     a1 = (b.
x()*a.
y()-a.
x()*b.
y())/dx1;
 
 1751     a2 = (d.
x()*c.
y()-c.
x()*d.
y())/dx2;
 
 1754   if (stand1 && stand2)
 
 1758     if (std::fabs(a.
x()-c.
x())<fgkTolerance)
 
 1762        if ( ((c.
y()-a.
y())*(c.
y()-b.
y())<-fgkTolerance)
 
 1763          || ((d.
y()-a.
y())*(d.
y()-b.
y())<-fgkTolerance)
 
 1764          || ((a.
y()-c.
y())*(a.
y()-d.
y())<-fgkTolerance)
 
 1765          || ((b.
y()-c.
y())*(b.
y()-d.
y())<-fgkTolerance) )  { 
return true; }
 
 1788       if (std::fabs(b1-b2) < fgkTolerance)
 
 1792         if (std::fabs(c.
y()-(a1+b1*c.
x())) > fgkTolerance)  { 
return false; }
 
 1796         if ( ((c.
x()-a.
x())*(c.
x()-b.
x())<-fgkTolerance)
 
 1797           || ((d.
x()-a.
x())*(d.
x()-b.
x())<-fgkTolerance)
 
 1798           || ((a.
x()-c.
x())*(a.
x()-d.
x())<-fgkTolerance)
 
 1799           || ((b.
x()-c.
x())*(b.
x()-d.
x())<-fgkTolerance) )  { 
return true; }
 
 1803       xm = (a1-a2)/(b2-b1);
 
 1804       ym = (a1*b2-a2*b1)/(b2-b1);
 
 1810   G4double check = (xm-a.
x())*(xm-b.
x())+(ym-a.
y())*(ym-b.
y());
 
 1811   if (check > -fgkTolerance)  { 
return false; }
 
 1812   check = (xm-c.
x())*(xm-d.
x())+(ym-c.
y())*(ym-d.
y());
 
 1813   if (check > -fgkTolerance)  { 
return false; }
 
 1850   det = dv.
x()*v1.
y()*v2.
z()+dv.
y()*v1.
z()*v2.
x()
 
 1851       - dv.
x()*v1.
z()*v2.
y()-dv.
y()*v1.
x()*v2.
z();
 
 1855     temp1 = v1.
cross(v2);
 
 1856     temp2 = (p2-p1).cross(v2);
 
 1857     if (temp1.
dot(temp2) < 0)  { 
return false; } 
 
 1861     q = ((dv).cross(v2)).mag()/q;
 
 1871 G4GenericTrap::MakeDownFacet(
const std::vector<G4ThreeVector>& fromVertices, 
 
 1878   if ( (fromVertices[ind1] == fromVertices[ind2]) ||
 
 1879        (fromVertices[ind2] == fromVertices[ind3]) ||
 
 1880        (fromVertices[ind1] == fromVertices[ind3]) )  { 
return 0; }
 
 1882   std::vector<G4ThreeVector> vertices;
 
 1883   vertices.push_back(fromVertices[ind1]);
 
 1884   vertices.push_back(fromVertices[ind2]);
 
 1885   vertices.push_back(fromVertices[ind3]);
 
 1889   G4ThreeVector cross=(vertices[1]-vertices[0]).cross(vertices[2]-vertices[1]);
 
 1891   if ( cross.
z() > 0.0 )
 
 1895     std::ostringstream message;
 
 1896     message << 
"Vertices in wrong order - " << 
GetName();
 
 1897     G4Exception(
"G4GenericTrap::MakeDownFacet", 
"GeomSolids0002",
 
 1907 G4GenericTrap::MakeUpFacet(
const std::vector<G4ThreeVector>& fromVertices, 
 
 1915   if ( (fromVertices[ind1] == fromVertices[ind2]) ||
 
 1916        (fromVertices[ind2] == fromVertices[ind3]) ||
 
 1917        (fromVertices[ind1] == fromVertices[ind3]) )  { 
return 0; }
 
 1919   std::vector<G4ThreeVector> vertices;
 
 1920   vertices.push_back(fromVertices[ind1]);
 
 1921   vertices.push_back(fromVertices[ind2]);
 
 1922   vertices.push_back(fromVertices[ind3]);
 
 1926   G4ThreeVector cross=(vertices[1]-vertices[0]).cross(vertices[2]-vertices[1]);
 
 1928   if ( cross.
z() < 0.0 )
 
 1932     std::ostringstream message;
 
 1933     message << 
"Vertices in wrong order - " << 
GetName();
 
 1934     G4Exception(
"G4GenericTrap::MakeUpFacet", 
"GeomSolids0002",
 
 1944 G4GenericTrap::MakeSideFacet(
const G4ThreeVector& downVertex0,
 
 1952   if ( (downVertex0 == downVertex1) && (upVertex0 == upVertex1) )
 
 1957   if ( downVertex0 == downVertex1 )
 
 1962   if ( upVertex0 == upVertex1 )
 
 1977   G4int nv = fgkNofVertices/2;
 
 1978   std::vector<G4ThreeVector> downVertices;
 
 1979   for ( 
G4int i=0; i<nv; i++ )
 
 1982                                          fVertices[i].
y(), -fDz));
 
 1985   std::vector<G4ThreeVector> upVertices;
 
 1986   for ( 
G4int i=nv; i<2*nv; i++ )
 
 1989                                        fVertices[i].
y(), fDz));
 
 1995     = (downVertices[1]-downVertices[0]).cross(downVertices[2]-downVertices[1]);
 
 1997     = (upVertices[1]-upVertices[0]).cross(upVertices[2]-upVertices[1]);
 
 1998   if ( (cross.
z() > 0.0) || (cross1.
z() > 0.0) )
 
 2000     ReorderVertices(downVertices);
 
 2001     ReorderVertices(upVertices);
 
 2007   facet = MakeDownFacet(downVertices, 0, 1, 2);
 
 2008   if (facet)  { tessellatedSolid->
AddFacet( facet ); }
 
 2009   facet = MakeDownFacet(downVertices, 0, 2, 3);
 
 2010   if (facet)  { tessellatedSolid->
AddFacet( facet ); }
 
 2011   facet = MakeUpFacet(upVertices, 0, 2, 1);
 
 2012   if (facet)  { tessellatedSolid->
AddFacet( facet ); }
 
 2013   facet = MakeUpFacet(upVertices, 0, 3, 2);
 
 2014   if (facet)  { tessellatedSolid->
AddFacet( facet ); }
 
 2018   for ( 
G4int i = 0; i < nv; ++i )
 
 2020     G4int j = (i+1) % nv;
 
 2021     facet = MakeSideFacet(downVertices[j], downVertices[i], 
 
 2022                           upVertices[i], upVertices[j]);
 
 2024     if ( facet )  { tessellatedSolid->
AddFacet( facet ); }
 
 2029   return tessellatedSolid;
 
 2034 void G4GenericTrap::ComputeBBox() 
 
 2039   minX = maxX = fVertices[0].x();
 
 2040   minY = maxY = fVertices[0].y();
 
 2042   for (
G4int i=1; i< fgkNofVertices; i++)
 
 2044     if (minX>fVertices[i].
x()) { minX=fVertices[i].x(); }
 
 2045     if (maxX<fVertices[i].
x()) { maxX=fVertices[i].x(); }
 
 2046     if (minY>fVertices[i].
y()) { minY=fVertices[i].y(); }
 
 2047     if (maxY<fVertices[i].
y()) { maxY=fVertices[i].y(); }
 
 2059   if ( fTessellatedSolid )
 
 2081   if ( fTessellatedSolid )
 
 2097   if ( fTessellatedSolid )
 
 2106   Dx = 0.5*(maxVec.
x()- minVec.
x());
 
 2107   Dy = 0.5*(maxVec.
y()- minVec.
y());
 
 2118   if ( fTessellatedSolid )
 
 2128   size_t nVertices, nFacets;
 
 2130   G4int subdivisions=0;
 
 2153       Dx = 0.5*(maxVec.
x()- minVec.
y());
 
 2154       Dy = 0.5*(maxVec.
y()- minVec.
y());
 
 2155       if (Dy > Dx)  { Dx=Dy; }
 
 2157       subdivisions=8*
G4int(maxTwist/(Dx*Dx*Dx)*fDz);
 
 2158       if (subdivisions<4)  { subdivisions=4; }
 
 2159       if (subdivisions>30) { subdivisions=30; }
 
 2162   G4int sub4=4*subdivisions;
 
 2163   nVertices = 8+subdivisions*4;
 
 2164   nFacets = 6+subdivisions*4;
 
 2173                                         fVertices[i].
y(),-fDz));
 
 2175   for( i=0;i<subdivisions;i++)
 
 2177     for(
G4int j=0;j<4;j++)
 
 2179       G4TwoVector u=fVertices[j]+cf*(i+1)*( fVertices[j+4]-fVertices[j]);
 
 2186                                         fVertices[i].
y(),fDz));
 
 2192   for (i=0;i<subdivisions+1;i++)
 
 2195     polyhedron->
AddFacet(5+is,8+is,4+is,1+is);
 
 2196     polyhedron->
AddFacet(8+is,7+is,3+is,4+is);
 
 2197     polyhedron->
AddFacet(7+is,6+is,2+is,3+is);
 
 2198     polyhedron->
AddFacet(6+is,5+is,1+is,2+is); 
 
 2200   polyhedron->
AddFacet(5+sub4,6+sub4,7+sub4,8+sub4);  
 
void ClipCrossSection(G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const 
G4Polyhedron * GetPolyhedron() const 
G4int GetVisSubdivisions() const 
void SetSolidClosed(const G4bool t)
G4GenericTrap(const G4String &name, G4double halfZ, const std::vector< G4TwoVector > &vertices)
G4double GetMinYExtent() const 
CLHEP::Hep3Vector G4ThreeVector
G4VisExtent GetExtent() const 
virtual void DescribeYourselfTo(G4VGraphicsScene &scene) const 
double dot(const Hep3Vector &) const 
G4Polyhedron * fpPolyhedron
G4bool IsYLimited() const 
virtual G4VisExtent GetExtent() const 
virtual G4double GetCubicVolume()
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const 
G4GeometryType GetEntityType() const 
G4bool IsXLimited() const 
virtual void AddSolid(const G4Box &)=0
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const 
G4double GetSurfaceArea()
G4double GetMaxXExtent() const 
G4double GetMinZExtent() const 
G4ThreeVector GetPointOnSurface() const 
virtual G4double DistanceToOut(const G4ThreeVector &p) const 
G4GenericTrap & operator=(const G4GenericTrap &rhs)
virtual G4Polyhedron * CreatePolyhedron() const 
virtual EInside Inside(const G4ThreeVector &p) const 
std::ostream & StreamInfo(std::ostream &os) const 
G4bool AddFacet(G4VFacet *aFacet)
G4double GetTwistAngle(G4int index) const 
std::vector< G4ThreeVector > G4ThreeVectorList
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const 
G4double GetCubicVolume()
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const 
EInside Inside(const G4ThreeVector &p) const 
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
virtual G4Polyhedron * GetPolyhedron() const 
G4double GetMinXExtent() const 
G4Polyhedron * CreatePolyhedron() const 
G4double GetMaxZExtent() const 
static G4int GetNumberOfRotationSteps()
void AddFacet(const G4int iv1, const G4int iv2, const G4int iv3, const G4int iv4=0)
G4double GetMaxYExtent() const 
G4VSolid & operator=(const G4VSolid &rhs)
G4int GetNumberOfRotationStepsAtTimeOfCreation() const 
Hep3Vector cross(const Hep3Vector &) const 
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const 
G4double GetMaxExtent(const EAxis pAxis) const 
virtual G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const 
G4bool IsZLimited() const 
virtual G4double GetSurfaceArea()
virtual G4ThreeVector GetPointOnSurface() const 
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const 
void ClipBetweenSections(G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const 
G4double GetMinExtent(const EAxis pAxis) const 
void DescribeYourselfTo(G4VGraphicsScene &scene) const 
void AddVertex(const G4ThreeVector v)