50 #if !defined(G4GEOM_USE_UGENERICTRAP) 
   76 const G4int    G4GenericTrap::fgkNofVertices = 8;
 
   77 const G4double G4GenericTrap::fgkTolerance = 1E-3;
 
   82                               const std::vector<G4TwoVector>&  vertices )
 
   84     fRebuildPolyhedron(false),
 
  101   G4String errorDescription = 
"InvalidSetup in \" ";
 
  102   errorDescription += 
name;
 
  103   errorDescription += 
"\""; 
 
  109   if ( 
G4int(vertices.size()) != fgkNofVertices )
 
  111     G4Exception(
"G4GenericTrap::G4GenericTrap()", 
"GeomSolids0002",
 
  119      G4Exception(
"G4GenericTrap::G4GenericTrap()", 
"GeomSolids0002",
 
  125   if(CheckOrder(vertices))
 
  127     for (
G4int i=0; i<fgkNofVertices; ++i) {fVertices.push_back(vertices[i]);}
 
  131     for (
G4int i=0; i <4; ++i) {fVertices.push_back(vertices[3-i]);}
 
  132     for (
G4int i=0; i <4; ++i) {fVertices.push_back(vertices[7-i]);}
 
  137   for (
G4int j=0; j < 2; j++)
 
  139     for (
G4int i=1; i<4; ++i)
 
  142       length = (fVertices[k]-fVertices[k-1]).mag();
 
  145         std::ostringstream message;
 
  146         message << 
"Length segment is too small." << 
G4endl 
  147                 << 
"Distance between " << fVertices[k-1] << 
" and " 
  148                 << fVertices[k] << 
" is only " << length << 
" mm !"; 
 
  149         G4Exception(
"G4GenericTrap::G4GenericTrap()", 
"GeomSolids1001",
 
  150                     JustWarning, message, 
"Vertices will be collapsed.");
 
  151         fVertices[k]=fVertices[k-1];
 
  158   for( 
G4int i=0; i<4; i++) { fTwist[i]=0.; }
 
  159   fIsTwisted = ComputeIsTwisted();
 
  169    if ( !fIsTwisted )  { fTessellatedSolid = CreateTessellatedSolid(); }
 
  177     fRebuildPolyhedron(false),
 
  179     halfCarTolerance(0.),
 
  183     fTessellatedSolid(0),
 
  199   delete fTessellatedSolid;
 
  206     fRebuildPolyhedron(false), fpPolyhedron(0),
 
  207     halfCarTolerance(rhs.halfCarTolerance),
 
  208     fDz(rhs.fDz), fVertices(rhs.fVertices),
 
  209     fIsTwisted(rhs.fIsTwisted), fTessellatedSolid(0),
 
  210     fMinBBoxVector(rhs.fMinBBoxVector), fMaxBBoxVector(rhs.fMaxBBoxVector),
 
  211     fVisSubdivisions(rhs.fVisSubdivisions),
 
  212     fSurfaceArea(rhs.fSurfaceArea), fCubicVolume(rhs.fCubicVolume) 
 
  214    for (
size_t i=0; i<4; ++i)  { fTwist[i] = rhs.fTwist[i]; }
 
  216    if (rhs.fTessellatedSolid && !fIsTwisted )
 
  217    { fTessellatedSolid = CreateTessellatedSolid(); } 
 
  227    if (
this == &rhs)  { 
return *
this; }
 
  235    halfCarTolerance = rhs.halfCarTolerance;
 
  236    fDz = rhs.fDz; fVertices = rhs.fVertices;
 
  237    fIsTwisted = rhs.fIsTwisted; fTessellatedSolid = 0;
 
  238    fMinBBoxVector = rhs.fMinBBoxVector; fMaxBBoxVector = rhs.fMaxBBoxVector;
 
  239    fVisSubdivisions = rhs.fVisSubdivisions;
 
  240    fSurfaceArea = rhs.fSurfaceArea; fCubicVolume = rhs.fCubicVolume;
 
  242    for (
size_t i=0; i<4; ++i)  { fTwist[i] = rhs.fTwist[i]; }
 
  244    if (rhs.fTessellatedSolid && !fIsTwisted )
 
  245    { 
delete fTessellatedSolid; fTessellatedSolid = CreateTessellatedSolid(); } 
 
  257                               const std::vector<G4TwoVector>& poly)
 const  
  263   for (
G4int i = 0; i < 4; i++)
 
  267     cross = (p.
x()-poly[i].x())*(poly[j].y()-poly[i].y())-
 
  268             (p.
y()-poly[i].y())*(poly[j].x()-poly[i].x());
 
  270     len2=(poly[i]-poly[j]).mag2();
 
  273       if(cross*cross<=len2*halfCarTolerance*halfCarTolerance)  
 
  282         if (poly[j].x() > poly[i].x())
 
  291         if ( p.
x() > poly[iMax].x()+halfCarTolerance
 
  292           || p.
x() < poly[iMin].x()-halfCarTolerance )
 
  297         if (poly[j].y() > poly[i].y())
 
  307         if ( p.
y() > poly[iMax].y()+halfCarTolerance
 
  308           || p.
y() < poly[iMin].y()-halfCarTolerance )
 
  313         if ( poly[iMax].x() != poly[iMin].x() )
 
  315           test = (p.
x()-poly[iMin].x())/(poly[iMax].x()-poly[iMin].x())
 
  316                * (poly[iMax].y()-poly[iMin].y())+poly[iMin].y();
 
  325         if( (test>=(poly[iMin].y()-halfCarTolerance))
 
  326          && (test<=(poly[iMax].y()+halfCarTolerance)) )
 
  335       else if (cross<0.)  { 
return kOutside; }
 
  347     if ( (std::fabs(p.
x()-poly[0].x())+std::fabs(p.
y()-poly[0].y())) > halfCarTolerance )
 
  362    if ( fTessellatedSolid )
 
  364      return fTessellatedSolid->
Inside(p);
 
  369   std::vector<G4TwoVector> xy;
 
  371   if (std::fabs(p.
z()) <= fDz+halfCarTolerance)  
 
  376     for (
G4int i=0; i<4; i++)
 
  378       xy.push_back(fVertices[i+4]+cf*( fVertices[i]-fVertices[i+4]));
 
  381     innew=InsidePolygone(p,xy);
 
  385       if(std::fabs(p.
z()) > fDz-halfCarTolerance)  { innew=
kSurface; }
 
  399   if ( fTessellatedSolid )
 
  406                 p0, p1, p2, r1, r2, r3, r4;
 
  407   G4int noSurfaces = 0; 
 
  411   distz = fDz-std::fabs(p.
z());
 
  412   if (distz < halfCarTolerance)
 
  428   std:: vector<G4TwoVector> vertices;  
 
  430   for (
G4int i=0; i<4; i++)
 
  432     vertices.push_back(fVertices[i+4]+cf*(fVertices[i]-fVertices[i+4]));
 
  437   for (
G4int q=0; q<4; q++)
 
  448     p2=
G4ThreeVector(vertices[(q+1)%4].x(),vertices[(q+1)%4].y(),p.
z());
 
  456         p2=
G4ThreeVector(fVertices[(q+1)%4].x(),fVertices[(q+1)%4].y(),-fDz);
 
  460         p2=
G4ThreeVector(fVertices[(q+1)%4+4].x(),fVertices[(q+1)%4+4].y(),fDz);
 
  463     lnorm = (p1-p0).cross(p2-p0);
 
  464     lnorm = lnorm.
unit();
 
  465     if(zPlusSide)  { lnorm=-lnorm; }
 
  474         G4double proj=(p-p0).dot(p2-p0)/normP;
 
  475         if(proj<0)     { proj=0; }
 
  476         if(proj>normP) { proj=normP; }
 
  482         r1=r1+proj*(r2-r1)/normP;
 
  483         r3=r3+proj*(r4-r3)/normP;
 
  490     distxy=std::fabs((p0-p).dot(lnorm));
 
  491     if ( distxy<halfCarTolerance )
 
  497       sumnorm=sumnorm+lnorm;
 
  511   if ( noSurfaces == 0 )
 
  514     G4Exception(
"G4GenericTrap::SurfaceNormal(p)", 
"GeomSolids1002",
 
  520   else if ( noSurfaces == 1 )  { ; }
 
  521   else                         { sumnorm = sumnorm.unit(); }
 
  529                                             const G4int ipl )
 const 
  534   if ( fTessellatedSolid )
 
  550   u=fVertices[i+4]+cf*(fVertices[i]-fVertices[i+4]);
 
  551   v=fVertices[j+4]+cf*(fVertices[j]-fVertices[j+4]);
 
  557   if (std::fabs(distz)<halfCarTolerance)
 
  572     if ( std::fabs(p.
z()+fDz) > halfCarTolerance )
 
  581   lnorm=-(p1-p0).cross(p2-p0);
 
  582   if (distz>-halfCarTolerance)  { lnorm=-lnorm.
unit(); }
 
  583   else                          { lnorm=lnorm.
unit();  }
 
  592       G4double proj=(p-p0).dot(p2-p0)/normP;
 
  593       if (proj<0)     { proj=0; }
 
  594       if (proj>normP) { proj=normP; }
 
  600       r1=r1+proj*(r2-r1)/normP;
 
  601       r3=r3+proj*(r4-r3)/normP;
 
  615                                     const G4int ipl)
 const  
  627   xa=fVertices[ipl].x();
 
  628   ya=fVertices[ipl].y();
 
  629   xb=fVertices[ipl+4].x();
 
  630   yb=fVertices[ipl+4].y();
 
  633   xd=fVertices[4+j].x();
 
  634   yd=fVertices[4+j].y();
 
  651   G4double a = (dtx*v.
y()-dty*v.
x()+(tx1*ty2-tx2*ty1)*v.
z())*v.
z();
 
  652   G4double b = dxs*v.
y()-dys*v.
x()+(dtx*p.
y()-dty*p.
x()+ty2*xs1-ty1*xs2
 
  653              + tx1*ys2-tx2*ys1)*v.
z();
 
  654   G4double c=dxs*p.
y()-dys*p.
x()+xs1*ys2-xs2*ys1;
 
  665     if (q>-halfCarTolerance)
 
  667       if (q<halfCarTolerance)
 
  669         if (NormalToPlane(p,ipl).dot(v)<=0)
 
  678       if (std::fabs(zi)<fDz)
 
  686         zi = (xp-x1)*(xp-x2)+(yp-y1)*(yp-y2);
 
  687         if (zi<=halfCarTolerance)  { 
return q; }
 
  695     if (a>0) { q=0.5*(-b-std::sqrt(d))/a; }
 
  696     else     { q=0.5*(-b+std::sqrt(d))/a; }
 
  700     if (q>-halfCarTolerance)
 
  702       if(q<halfCarTolerance)
 
  704         if (NormalToPlane(p,ipl).
dot(v)<=0)
 
  710           if (a>0) { q=0.5*(-b+std::sqrt(d))/a; }
 
  711           else     { q=0.5*(-b-std::sqrt(d))/a; }
 
  712           if (q<=halfCarTolerance) { 
return kInfinity; }
 
  718       if (std::fabs(zi)<fDz)
 
  726         zi = (xp-x1)*(xp-x2)+(yp-y1)*(yp-y2);
 
  727         if (zi<=halfCarTolerance)  { 
return q; }
 
  730     if (a>0)  { q=0.5*(-b+std::sqrt(d))/a; }
 
  731     else      { q=0.5*(-b-std::sqrt(d))/a; }
 
  735     if (q>-halfCarTolerance)
 
  737       if(q<halfCarTolerance)
 
  739         if (NormalToPlane(p,ipl).
dot(v)<=0)
 
  745           if (a>0) { q=0.5*(-b-std::sqrt(d))/a; }
 
  746           else     { q=0.5*(-b+std::sqrt(d))/a; }
 
  747           if (q<=halfCarTolerance)  { 
return kInfinity; }
 
  753       if (std::fabs(zi)<fDz)
 
  761         zi = (xp-x1)*(xp-x2)+(yp-y1)*(yp-y2);
 
  762         if (zi<=halfCarTolerance)  { 
return q; }
 
  775   if ( fTessellatedSolid )
 
  789     dist[i]=DistToPlane(p, v, i);  
 
  795   if (std::fabs(p.
z())>fDz-halfCarTolerance)
 
  802         dist[4] = (fDz-p.
z())/v.
z();
 
  806         dist[4] = (-fDz-p.
z())/v.
z();
 
  808       if (dist[4]<-halfCarTolerance)
 
  814         if(dist[4]<halfCarTolerance)
 
  818           if (n.
dot(v)<0) { dist[4]=0.; }
 
  829     if (dist[i] < distmin)  { distmin = dist[i]; }
 
  832   if (distmin<halfCarTolerance)  { distmin=0.; }
 
  844   if ( fTessellatedSolid )
 
  851   if(safz<0) { safz=0; }
 
  857   for (iseg=0; iseg<4; iseg++)
 
  859     safxy = SafetyToFace(p,iseg);
 
  860     if (safxy>safe)  { safe=safxy; }
 
  877   p1=
G4ThreeVector(fVertices[iseg].x(),fVertices[iseg].y(),-fDz);
 
  879   norm=NormalToPlane(p,iseg);
 
  880   safe = (p-p1).dot(norm); 
 
  901   if ( (std::fabs(xa-xc)+std::fabs(ya-yc)) < halfCarTolerance )
 
  903     xc=fVertices[j+4].x();
 
  904     yc=fVertices[j+4].y();
 
  910     if ( (std::fabs(xb-xc)+std::fabs(yb-yc)) < halfCarTolerance )
 
  917   G4double c=(xb-xa)*(yc-ya)-(xc-xa)*(yb-ya);
 
  923     t=-(a*p.
x()+b*p.
y()+c*p.
z()+d)/t;
 
  925   if ( (t<halfCarTolerance) && (t>-halfCarTolerance) )
 
  950   if ( fTessellatedSolid )
 
  952     return fTessellatedSolid->
DistanceToOut(p, v, calcNorm, validNorm, n);
 
  957   G4bool lateral_cross = 
false;
 
  958   ESide side = kUndefined;
 
  960   if (calcNorm)  { *validNorm=
true; } 
 
  964     distmin=(-fDz-p.
z())/v.
z();
 
  971       distmin = (fDz-p.
z())/v.
z(); 
 
  981   for (
G4int ipl=0; ipl<4; ipl++)
 
  984     xa=fVertices[ipl].x();
 
  985     ya=fVertices[ipl].y();
 
  986     xb=fVertices[ipl+4].x();
 
  987     yb=fVertices[ipl+4].y();
 
  990     xd=fVertices[4+j].x();
 
  991     yd=fVertices[4+j].y();
 
  993     if ( ((std::fabs(xb-xd)+std::fabs(yb-yd))<halfCarTolerance)
 
  994       || ((std::fabs(xa-xc)+std::fabs(ya-yc))<halfCarTolerance) )
 
  996       G4double q=DistToTriangle(p,v,ipl) ;
 
  997       if ( (q>=0) && (q<distmin) )
 
 1018     G4double a = (dtx*v.
y()-dty*v.
x()+(tx1*ty2-tx2*ty1)*v.
z())*v.
z();
 
 1019     G4double b = dxs*v.
y()-dys*v.
x()+(dtx*p.
y()-dty*p.
x()+ty2*xs1-ty1*xs2
 
 1020                + tx1*ys2-tx2*ys1)*v.
z();
 
 1021     G4double c=dxs*p.
y()-dys*p.
x()+xs1*ys2-xs2*ys1;
 
 1031       if ((q > -halfCarTolerance) && (q < distmin))
 
 1033         if (q < halfCarTolerance)
 
 1035           if (NormalToPlane(p,ipl).dot(v)<0.)  { 
continue; }
 
 1046       if (a > 0)  { q=0.5*(-b-std::sqrt(d))/a; }
 
 1047       else        { q=0.5*(-b+std::sqrt(d))/a; }
 
 1051       if (q > -halfCarTolerance )
 
 1055           if(q < halfCarTolerance)
 
 1057             if (NormalToPlane(p,ipl).
dot(v)<0.)  
 
 1059               if (a > 0)  { q=0.5*(-b+std::sqrt(d))/a; }
 
 1060               else        { q=0.5*(-b-std::sqrt(d))/a; }
 
 1061               if (( q > halfCarTolerance) && (q < distmin))
 
 1064                 lateral_cross = 
true;
 
 1071           lateral_cross = 
true;
 
 1077         if (a > 0)  { q=0.5*(-b+std::sqrt(d))/a; }
 
 1078         else        { q=0.5*(-b-std::sqrt(d))/a; }
 
 1082         if ((q > -halfCarTolerance) && (q < distmin))
 
 1084           if (q < halfCarTolerance)
 
 1086             if (NormalToPlane(p,ipl).
dot(v)<0.)  
 
 1088               if (a > 0)  { q=0.5*(-b-std::sqrt(d))/a; }
 
 1089               else        { q=0.5*(-b+std::sqrt(d))/a; }
 
 1090               if ( ( q > halfCarTolerance) && (q < distmin) )
 
 1093                 lateral_cross = 
true;
 
 1100           lateral_cross = 
true;
 
 1114     if (v.z()>0.) { i=4; }
 
 1115     std::vector<G4TwoVector> xy;
 
 1116     for ( 
G4int j=0; j<4; j++)  { xy.push_back(fVertices[i+j]); }
 
 1120     if (InsidePolygone(pt,xy)==
kOutside)
 
 1131       if(v.z()>0) {side=kPZ;}
 
 1142         *n=NormalToPlane(pt,0);
 
 1145         *n=NormalToPlane(pt,1);
 
 1148         *n=NormalToPlane(pt,2);
 
 1151         *n=NormalToPlane(pt,3);
 
 1161         std::ostringstream message;
 
 1162         G4int oldprc = message.precision(16);
 
 1163         message << 
"Undefined side for valid surface normal to solid." << 
G4endl 
 1165                 << 
"  p.x() = "   << p.
x()/
mm << 
" mm" << 
G4endl 
 1166                 << 
"  p.y() = "   << p.
y()/
mm << 
" mm" << 
G4endl 
 1167                 << 
"  p.z() = "   << p.
z()/
mm << 
" mm" << 
G4endl 
 1168                 << 
"Direction:" << 
G4endl 
 1169                 << 
"  v.x() = "   << v.
x() << 
G4endl 
 1170                 << 
"  v.y() = "   << v.
y() << 
G4endl 
 1171                 << 
"  v.z() = "   << v.
z() << 
G4endl 
 1172                 << 
"Proposed distance :" << 
G4endl 
 1173                 << 
"  distmin = " << distmin/
mm << 
" mm";
 
 1174         message.precision(oldprc);
 
 1175         G4Exception(
"G4GenericTrap::DistanceToOut(p,v,..)",
 
 1181   if (distmin<halfCarTolerance)  { distmin=0.; }
 
 1192   if ( fTessellatedSolid )
 
 1199   if (safz<0) { safz = 0; }
 
 1204   for (
G4int iseg=0; iseg<4; iseg++)
 
 1206     safxy = std::fabs(SafetyToFace(p,iseg));
 
 1207     if (safxy < safe)  { safe = safxy; }
 
 1217   pMin = GetMinimumBBox();
 
 1218   pMax = GetMaximumBBox();
 
 1222   if (pMin.
x() >= pMax.
x() || pMin.
y() >= pMax.
y() || pMin.
z() >= pMax.
z())
 
 1224     std::ostringstream message;
 
 1225     message << 
"Bad bounding box (min >= max) for solid: " 
 1227             << 
"\npMin = " << pMin
 
 1228             << 
"\npMax = " << pMax;
 
 1249 #ifdef G4BBOX_EXTENT 
 1250   if (
true) 
return bbox.
CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
 
 1254     return exist = (pMin < pMax) ? 
true : 
false;
 
 1266   for (
G4int i=0; i<4; ++i)
 
 1270     baseA[2*i].set(va.
x(),va.
y(),-dz);
 
 1271     baseB[2*i].set(vb.
x(),vb.
y(), dz);
 
 1273   for (
G4int i=0; i<4; ++i)
 
 1275     G4int k1=2*i, k2=(2*i+2)%8;
 
 1276     G4double ax = (baseA[k2].x()-baseA[k1].x());
 
 1277     G4double ay = (baseA[k2].y()-baseA[k1].y());
 
 1278     G4double bx = (baseB[k2].x()-baseB[k1].x());
 
 1279     G4double by = (baseB[k2].y()-baseB[k1].y());
 
 1281     baseA[k1+1] = (znorm < 0.0) ? baseA[k2] : baseA[k1];
 
 1282     baseB[k1+1] = (znorm < 0.0) ? baseB[k1] : baseB[k2];
 
 1285   std::vector<const G4ThreeVectorList *> polygons(2);
 
 1286   polygons[0] = &baseA;
 
 1287   polygons[1] = &baseB;
 
 1312   G4int oldprc = os.precision(16);
 
 1313   os << 
"-----------------------------------------------------------\n" 
 1314      << 
"    *** Dump for solid - " << 
GetName() << 
" *** \n" 
 1315      << 
"    =================================================== \n" 
 1317      << 
"   half length Z: " << fDz/
mm << 
" mm \n" 
 1318      << 
"   list of vertices:\n";
 
 1320   for ( 
G4int i=0; i<fgkNofVertices; ++i )
 
 1322     os << std::setw(5) << 
"#" << i 
 
 1323        << 
"   vx = " << fVertices[i].x()/
mm << 
" mm"  
 1324        << 
"   vy = " << fVertices[i].y()/
mm << 
" mm" << 
G4endl;
 
 1326   os.precision(oldprc);
 
 1337   if ( fTessellatedSolid )
 
 1345   G4double rand,area,chose,cf,lambda0,lambda1,alfa,beta,zp;
 
 1348   std::vector<G4ThreeVector> vertices;
 
 1349   for (
G4int i=0; i<4;i++)
 
 1351     vertices.push_back(
G4ThreeVector(fVertices[i].x(),fVertices[i].y(),-fDz));
 
 1353   for (
G4int i=4; i<8;i++)
 
 1355     vertices.push_back(
G4ThreeVector(fVertices[i].x(),fVertices[i].y(),fDz));
 
 1360   G4double Surface0=GetFaceSurfaceArea(vertices[0],vertices[1],
 
 1361                                        vertices[2],vertices[3]);
 
 1362   G4double Surface1=GetFaceSurfaceArea(vertices[0],vertices[1],
 
 1363                                        vertices[5],vertices[4]);
 
 1364   G4double Surface2=GetFaceSurfaceArea(vertices[3],vertices[0],
 
 1365                                        vertices[4],vertices[7]);
 
 1366   G4double Surface3=GetFaceSurfaceArea(vertices[2],vertices[3],
 
 1367                                        vertices[7],vertices[6]);
 
 1368   G4double Surface4=GetFaceSurfaceArea(vertices[2],vertices[1],
 
 1369                                        vertices[5],vertices[6]);
 
 1370   G4double Surface5=GetFaceSurfaceArea(vertices[4],vertices[5],
 
 1371                                        vertices[6],vertices[7]);
 
 1373   area = Surface0+Surface1+Surface2+Surface3+Surface4+Surface5;
 
 1376   if ( ( chose < Surface0)
 
 1377     || ( chose > (Surface0+Surface1+Surface2+Surface3+Surface4)) )
 
 1381     if(chose < Surface0)
 
 1384       u = fVertices[ipl]; v = fVertices[j];
 
 1385       w = fVertices[(ipl+3)%4]; 
 
 1390       u = fVertices[ipl+4]; v = fVertices[j+4];
 
 1391       w = fVertices[(ipl+3)%4+4]; 
 
 1396     lambda0=alfa-lambda1;
 
 1399     v = u+lambda0*v+lambda1*w;
 
 1403     if (chose < Surface0+Surface1) { ipl=0; }
 
 1404     else if (chose < Surface0+Surface1+Surface2) { ipl=1; }
 
 1405     else if (chose < Surface0+Surface1+Surface2+Surface3) { ipl=2; }
 
 1409     cf = 0.5*(fDz-zp)/fDz;
 
 1410     u = fVertices[ipl+4]+cf*( fVertices[ipl]-fVertices[ipl+4]);
 
 1411     v = fVertices[j+4]+cf*(fVertices[j]-fVertices[j+4]); 
 
 1423   if (fSurfaceArea == 0.0) {
 
 1428       G4ThreeVector vertix0(fVertices[0].x(),fVertices[0].y(),-fDz);
 
 1429       G4ThreeVector vertix1(fVertices[1].x(),fVertices[1].y(),-fDz);
 
 1430       G4ThreeVector vertix2(fVertices[2].x(),fVertices[2].y(),-fDz);
 
 1431       G4ThreeVector vertix3(fVertices[3].x(),fVertices[3].y(),-fDz);
 
 1432       G4ThreeVector vertix4(fVertices[4].x(),fVertices[4].y(), fDz);
 
 1433       G4ThreeVector vertix5(fVertices[5].x(),fVertices[5].y(), fDz);
 
 1434       G4ThreeVector vertix6(fVertices[6].x(),fVertices[6].y(), fDz);
 
 1435       G4ThreeVector vertix7(fVertices[7].x(),fVertices[7].y(), fDz);
 
 1438       fSurfaceArea = GetFaceSurfaceArea(vertix0,vertix1,vertix2,vertix3)  
 
 1439                    + GetFaceSurfaceArea(vertix1,vertix0,vertix4,vertix5)  
 
 1440                    + GetFaceSurfaceArea(vertix2,vertix1,vertix5,vertix6)  
 
 1441                    + GetFaceSurfaceArea(vertix3,vertix2,vertix6,vertix7)  
 
 1442                    + GetFaceSurfaceArea(vertix0,vertix3,vertix7,vertix4)  
 
 1443                    + GetFaceSurfaceArea(vertix7,vertix6,vertix5,vertix4); 
 
 1446   return fSurfaceArea;
 
 1453   if (fCubicVolume == 0.0) {
 
 1458       G4ThreeVector vertix0(fVertices[0].x(),fVertices[0].y(),-fDz);
 
 1459       G4ThreeVector vertix1(fVertices[1].x(),fVertices[1].y(),-fDz);
 
 1460       G4ThreeVector vertix2(fVertices[2].x(),fVertices[2].y(),-fDz);
 
 1461       G4ThreeVector vertix3(fVertices[3].x(),fVertices[3].y(),-fDz);
 
 1462       G4ThreeVector vertix4(fVertices[4].x(),fVertices[4].y(), fDz);
 
 1463       G4ThreeVector vertix5(fVertices[5].x(),fVertices[5].y(), fDz);
 
 1464       G4ThreeVector vertix6(fVertices[6].x(),fVertices[6].y(), fDz);
 
 1465       G4ThreeVector vertix7(fVertices[7].x(),fVertices[7].y(), fDz);
 
 1468       fCubicVolume = GetFaceCubicVolume(vertix0,vertix1,vertix2,vertix3)  
 
 1469                    + GetFaceCubicVolume(vertix1,vertix0,vertix4,vertix5)  
 
 1470                    + GetFaceCubicVolume(vertix2,vertix1,vertix5,vertix6)  
 
 1471                    + GetFaceCubicVolume(vertix3,vertix2,vertix6,vertix7)  
 
 1472                    + GetFaceCubicVolume(vertix0,vertix3,vertix7,vertix4)  
 
 1473                    + GetFaceCubicVolume(vertix7,vertix6,vertix5,vertix4); 
 
 1476   return fCubicVolume;
 
 1487   return (((p2-p0).cross(p3-p1)).mag()) / 2.;
 
 1499   return (((p2-p0).cross(p3-p1)).dot(p0)) / 6.;
 
 1504 G4bool G4GenericTrap::ComputeIsTwisted() 
 
 1511   G4int nv = fgkNofVertices/2;
 
 1513   for ( 
G4int i=0; i<4; i++ )
 
 1515     dx1 = fVertices[(i+1)%nv].x()-fVertices[i].x();
 
 1516     dy1 = fVertices[(i+1)%nv].y()-fVertices[i].y();
 
 1517     if ( (dx1 == 0) && (dy1 == 0) ) { 
continue; }
 
 1519     dx2 = fVertices[nv+(i+1)%nv].x()-fVertices[nv+i].x();
 
 1520     dy2 = fVertices[nv+(i+1)%nv].y()-fVertices[nv+i].y();
 
 1522     if ( dx2 == 0 && dy2 == 0 ) { 
continue; }
 
 1523     G4double twist_angle = std::fabs(dy1*dx2 - dx1*dy2);        
 
 1524     if ( twist_angle < fgkTolerance ) { 
continue; }
 
 1526     SetTwistAngle(i,twist_angle);
 
 1530     twist_angle = std::acos( (dx1*dx2 + dy1*dy2)
 
 1531                            / (std::sqrt(dx1*dx1+dy1*dy1)
 
 1532                              * std::sqrt(dx2*dx2+dy2*dy2)) );
 
 1536       std::ostringstream message;
 
 1537       message << 
"Twisted Angle is bigger than 90 degrees - " << 
GetName()
 
 1539               << 
"     Potential problem of malformed Solid !" << 
G4endl 
 1540               << 
"     TwistANGLE = " << twist_angle
 
 1541               << 
"*rad  for lateral plane N= " << i;
 
 1542       G4Exception(
"G4GenericTrap::ComputeIsTwisted()", 
"GeomSolids1002",
 
 1552 G4bool G4GenericTrap::CheckOrder(
const std::vector<G4TwoVector>& vertices)
 const 
 1557   G4bool clockwise_order=
true;
 
 1562   for (
G4int i=0; i<4; i++)
 
 1565     sum1 += vertices[i].x()*vertices[j].y() - vertices[j].x()*vertices[i].y();
 
 1566     sum2 += vertices[i+4].x()*vertices[j+4].y()
 
 1567           - vertices[j+4].x()*vertices[i+4].y();
 
 1569   if (sum1*sum2 < -fgkTolerance)
 
 1571      std::ostringstream message;
 
 1572      message << 
"Lower/upper faces defined with opposite clockwise - " 
 1574      G4Exception(
"G4GenericTrap::CheckOrder()", 
"GeomSolids0002",
 
 1578    if ((sum1 > 0.)||(sum2 > 0.))
 
 1580      std::ostringstream message;
 
 1581      message << 
"Vertices must be defined in clockwise XY planes - " 
 1583      G4Exception(
"G4GenericTrap::CheckOrder()", 
"GeomSolids1001",
 
 1585      clockwise_order = 
false;
 
 1590    G4bool illegal_cross = 
false;
 
 1591    illegal_cross = IsSegCrossingZ(vertices[0],vertices[4],
 
 1592                                   vertices[1],vertices[5]);
 
 1596      illegal_cross = IsSegCrossingZ(vertices[2],vertices[6],
 
 1597                                     vertices[3],vertices[7]);
 
 1602      illegal_cross = IsSegCrossing(vertices[0],vertices[1],
 
 1603                                    vertices[2],vertices[3]);
 
 1607      illegal_cross = IsSegCrossing(vertices[0],vertices[3],
 
 1608                                    vertices[1],vertices[2]);
 
 1612      illegal_cross = IsSegCrossing(vertices[4],vertices[5],
 
 1613                                    vertices[6],vertices[7]);
 
 1617      illegal_cross = IsSegCrossing(vertices[4],vertices[7],
 
 1618                                    vertices[5],vertices[6]);
 
 1623       std::ostringstream message;
 
 1624       message << 
"Malformed polygone with opposite sides - " << 
GetName();
 
 1625       G4Exception(
"G4GenericTrap::CheckOrderAndSetup()",
 
 1628    return clockwise_order;
 
 1633 void G4GenericTrap::ReorderVertices(std::vector<G4ThreeVector>& vertices)
 const 
 1637   std::vector<G4ThreeVector> oldVertices(vertices);
 
 1639   for ( 
G4int i=0; i < 
G4int(oldVertices.size()); ++i )
 
 1641     vertices[i] = oldVertices[oldVertices.size()-1-i];
 
 1655   G4double dx1,dx2,xm=0.,ym=0.,a1=0.,a2=0.,b1=0.,b2=0.;
 
 1659   if( std::fabs(dx1) < fgkTolerance )  { stand1 = 
true; }
 
 1660   if( std::fabs(dx2) < fgkTolerance )  { stand2 = 
true; }
 
 1663     a1 = (b.
x()*a.
y()-a.
x()*b.
y())/dx1;
 
 1668     a2 = (d.
x()*c.
y()-c.
x()*d.
y())/dx2;
 
 1671   if (stand1 && stand2)
 
 1675     if (std::fabs(a.
x()-c.
x())<fgkTolerance)
 
 1679        if ( ((c.
y()-a.
y())*(c.
y()-b.
y())<-fgkTolerance)
 
 1680          || ((d.
y()-a.
y())*(d.
y()-b.
y())<-fgkTolerance)
 
 1681          || ((a.
y()-c.
y())*(a.
y()-d.
y())<-fgkTolerance)
 
 1682          || ((b.
y()-c.
y())*(b.
y()-d.
y())<-fgkTolerance) )  { 
return true; }
 
 1705       if (std::fabs(b1-b2) < fgkTolerance)
 
 1709         if (std::fabs(c.
y()-(a1+b1*c.
x())) > fgkTolerance)  { 
return false; }
 
 1713         if ( ((c.
x()-a.
x())*(c.
x()-b.
x())<-fgkTolerance)
 
 1714           || ((d.
x()-a.
x())*(d.
x()-b.
x())<-fgkTolerance)
 
 1715           || ((a.
x()-c.
x())*(a.
x()-d.
x())<-fgkTolerance)
 
 1716           || ((b.
x()-c.
x())*(b.
x()-d.
x())<-fgkTolerance) )  { 
return true; }
 
 1720       xm = (a1-a2)/(b2-b1);
 
 1721       ym = (a1*b2-a2*b1)/(b2-b1);
 
 1727   G4double check = (xm-a.
x())*(xm-b.
x())+(ym-a.
y())*(ym-b.
y());
 
 1728   if (check > -fgkTolerance)  { 
return false; }
 
 1729   check = (xm-c.
x())*(xm-d.
x())+(ym-c.
y())*(ym-d.
y());
 
 1730   if (check > -fgkTolerance)  { 
return false; }
 
 1763       (std::fabs((p4-p3).y()) < 
kCarTolerance ) )  { 
return false; }
 
 1767   det = dv.
x()*v1.
y()*v2.
z()+dv.
y()*v1.
z()*v2.
x()
 
 1768       - dv.
x()*v1.
z()*v2.
y()-dv.
y()*v1.
x()*v2.
z();
 
 1772     temp1 = v1.
cross(v2);
 
 1773     temp2 = (p2-p1).cross(v2);
 
 1774     if (temp1.
dot(temp2) < 0)  { 
return false; } 
 
 1778     q = ((dv).cross(v2)).mag()/q;
 
 1788 G4GenericTrap::MakeDownFacet(
const std::vector<G4ThreeVector>& fromVertices, 
 
 1795   if ( (fromVertices[ind1] == fromVertices[ind2]) ||
 
 1796        (fromVertices[ind2] == fromVertices[ind3]) ||
 
 1797        (fromVertices[ind1] == fromVertices[ind3]) )  { 
return 0; }
 
 1799   std::vector<G4ThreeVector> vertices;
 
 1800   vertices.push_back(fromVertices[ind1]);
 
 1801   vertices.push_back(fromVertices[ind2]);
 
 1802   vertices.push_back(fromVertices[ind3]);
 
 1806   G4ThreeVector cross=(vertices[1]-vertices[0]).cross(vertices[2]-vertices[1]);
 
 1808   if ( cross.
z() > 0.0 )
 
 1812     std::ostringstream message;
 
 1813     message << 
"Vertices in wrong order - " << 
GetName();
 
 1814     G4Exception(
"G4GenericTrap::MakeDownFacet", 
"GeomSolids0002",
 
 1824 G4GenericTrap::MakeUpFacet(
const std::vector<G4ThreeVector>& fromVertices, 
 
 1832   if ( (fromVertices[ind1] == fromVertices[ind2]) ||
 
 1833        (fromVertices[ind2] == fromVertices[ind3]) ||
 
 1834        (fromVertices[ind1] == fromVertices[ind3]) )  { 
return 0; }
 
 1836   std::vector<G4ThreeVector> vertices;
 
 1837   vertices.push_back(fromVertices[ind1]);
 
 1838   vertices.push_back(fromVertices[ind2]);
 
 1839   vertices.push_back(fromVertices[ind3]);
 
 1843   G4ThreeVector cross=(vertices[1]-vertices[0]).cross(vertices[2]-vertices[1]);
 
 1845   if ( cross.
z() < 0.0 )
 
 1849     std::ostringstream message;
 
 1850     message << 
"Vertices in wrong order - " << 
GetName();
 
 1851     G4Exception(
"G4GenericTrap::MakeUpFacet", 
"GeomSolids0002",
 
 1861 G4GenericTrap::MakeSideFacet(
const G4ThreeVector& downVertex0,
 
 1869   if ( (downVertex0 == downVertex1) && (upVertex0 == upVertex1) )
 
 1874   if ( downVertex0 == downVertex1 )
 
 1879   if ( upVertex0 == upVertex1 )
 
 1894   G4int nv = fgkNofVertices/2;
 
 1895   std::vector<G4ThreeVector> downVertices;
 
 1896   for ( 
G4int i=0; i<nv; i++ )
 
 1899                                          fVertices[i].y(), -fDz));
 
 1902   std::vector<G4ThreeVector> upVertices;
 
 1903   for ( 
G4int i=nv; i<2*nv; i++ )
 
 1906                                        fVertices[i].y(), fDz));
 
 1912     = (downVertices[1]-downVertices[0]).cross(downVertices[2]-downVertices[1]);
 
 1914     = (upVertices[1]-upVertices[0]).cross(upVertices[2]-upVertices[1]);
 
 1915   if ( (cross.
z() > 0.0) || (cross1.
z() > 0.0) )
 
 1917     ReorderVertices(downVertices);
 
 1918     ReorderVertices(upVertices);
 
 1924   facet = MakeDownFacet(downVertices, 0, 1, 2);
 
 1925   if (facet)  { tessellatedSolid->
AddFacet( facet ); }
 
 1926   facet = MakeDownFacet(downVertices, 0, 2, 3);
 
 1927   if (facet)  { tessellatedSolid->
AddFacet( facet ); }
 
 1928   facet = MakeUpFacet(upVertices, 0, 2, 1);
 
 1929   if (facet)  { tessellatedSolid->
AddFacet( facet ); }
 
 1930   facet = MakeUpFacet(upVertices, 0, 3, 2);
 
 1931   if (facet)  { tessellatedSolid->
AddFacet( facet ); }
 
 1935   for ( 
G4int i = 0; i < nv; ++i )
 
 1937     G4int j = (i+1) % nv;
 
 1938     facet = MakeSideFacet(downVertices[j], downVertices[i], 
 
 1939                           upVertices[i], upVertices[j]);
 
 1941     if ( facet )  { tessellatedSolid->
AddFacet( facet ); }
 
 1946   return tessellatedSolid;
 
 1951 void G4GenericTrap::ComputeBBox() 
 
 1956   minX = maxX = fVertices[0].x();
 
 1957   minY = maxY = fVertices[0].y();
 
 1959   for (
G4int i=1; i< fgkNofVertices; i++)
 
 1961     if (minX>fVertices[i].x()) { minX=fVertices[i].x(); }
 
 1962     if (maxX<fVertices[i].x()) { maxX=fVertices[i].x(); }
 
 1963     if (minY>fVertices[i].y()) { minY=fVertices[i].y(); }
 
 1964     if (maxY<fVertices[i].y()) { maxY=fVertices[i].y(); }
 
 1976   if ( fTessellatedSolid )
 
 2002   if ( fTessellatedSolid )
 
 2018   if ( fTessellatedSolid )
 
 2027                       minVec.
y(), maxVec.
y(),
 
 2028                       minVec.
z(), maxVec.
z()); 
 
 2037   if ( fTessellatedSolid )
 
 2047   size_t nVertices, nFacets;
 
 2049   G4int subdivisions=0;
 
 2072       Dx = 0.5*(maxVec.
x()- minVec.
y());
 
 2073       Dy = 0.5*(maxVec.
y()- minVec.
y());
 
 2074       if (Dy > Dx)  { Dx=Dy; }
 
 2076       subdivisions=8*
G4int(maxTwist/(Dx*Dx*Dx)*fDz);
 
 2077       if (subdivisions<4)  { subdivisions=4; }
 
 2078       if (subdivisions>30) { subdivisions=30; }
 
 2081   G4int sub4=4*subdivisions;
 
 2082   nVertices = 8+subdivisions*4;
 
 2083   nFacets = 6+subdivisions*4;
 
 2092                                         fVertices[i].y(),-fDz));
 
 2094   for( i=0;i<subdivisions;i++)
 
 2096     for(
G4int j=0;j<4;j++)
 
 2098       G4TwoVector u=fVertices[j]+cf*(i+1)*( fVertices[j+4]-fVertices[j]);
 
 2105                                         fVertices[i].y(),fDz));
 
 2111   for (i=0;i<subdivisions+1;i++)
 
 2114     polyhedron->
AddFacet(5+is,8+is,4+is,1+is);
 
 2115     polyhedron->
AddFacet(8+is,7+is,3+is,4+is);
 
 2116     polyhedron->
AddFacet(7+is,6+is,2+is,3+is);
 
 2117     polyhedron->
AddFacet(6+is,5+is,1+is,2+is); 
 
 2119   polyhedron->
AddFacet(5+sub4,6+sub4,7+sub4,8+sub4);  
 
G4Polyhedron * GetPolyhedron() const 
 
G4int GetVisSubdivisions() const 
 
void SetSolidClosed(const G4bool t)
 
static constexpr double mm
 
G4GenericTrap(const G4String &name, G4double halfZ, const std::vector< G4TwoVector > &vertices)
 
static const G4double kInfinity
 
CLHEP::Hep3Vector G4ThreeVector
 
G4VisExtent GetExtent() const 
 
virtual void DescribeYourselfTo(G4VGraphicsScene &scene) const 
 
double dot(const Hep3Vector &) const 
 
G4Polyhedron * fpPolyhedron
 
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 
 
G4TwoVector GetVertex(G4int index) const 
 
G4double GetZHalfLength() const 
 
virtual void AddSolid(const G4Box &)=0
 
#define G4MUTEX_INITIALIZER
 
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const 
 
G4double GetSurfaceArea()
 
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 
 
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimits, const G4Transform3D &pTransform3D, G4double &pMin, G4double &pMax) const 
 
G4bool BoundingBoxVsVoxelLimits(const EAxis pAxis, const G4VoxelLimits &pVoxelLimits, const G4Transform3D &pTransform3D, G4double &pMin, G4double &pMax) const 
 
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 
 
std::vector< G4ThreeVector > G4ThreeVectorList
 
void Extent(G4ThreeVector &pMin, G4ThreeVector &pMax) const 
 
void AddVertex(const G4ThreeVector &v)
 
EInside Inside(const G4ThreeVector &p) const 
 
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
 
virtual G4Polyhedron * GetPolyhedron() const 
 
G4Polyhedron * CreatePolyhedron() const 
 
static G4int GetNumberOfRotationSteps()
 
void AddFacet(const G4int iv1, const G4int iv2, const G4int iv3, const G4int iv4=0)
 
G4VSolid & operator=(const G4VSolid &rhs)
 
G4int GetNumberOfRotationStepsAtTimeOfCreation() const 
 
static constexpr double pi
 
Hep3Vector cross(const Hep3Vector &) const 
 
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const 
 
virtual G4double GetSurfaceArea()
 
virtual G4ThreeVector GetPointOnSurface() const 
 
G4bool fRebuildPolyhedron
 
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const 
 
void DescribeYourselfTo(G4VGraphicsScene &scene) const