Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4GenericTrap.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 //
27 // $Id: G4GenericTrap.cc 67011 2013-01-29 16:17:41Z gcosmo $
28 //
29 //
30 // --------------------------------------------------------------------
31 // GEANT 4 class source file
32 //
33 // G4GenericTrap.cc
34 //
35 // Authors:
36 // Tatiana Nikitina, CERN; Ivana Hrivnacova, IPN Orsay
37 // Adapted from Root Arb8 implementation by Andrei Gheata, CERN
38 //
39 // History :
40 // 04 August 2011 T.Nikitina Add SetReferences() and InvertFacets()
41 // to CreatePolyhedron() for Visualisation of Boolean
42 // --------------------------------------------------------------------
43 
44 #include <iomanip>
45 
46 #include "G4GenericTrap.hh"
47 #include "G4PhysicalConstants.hh"
48 #include "G4SystemOfUnits.hh"
49 #include "G4TessellatedSolid.hh"
50 #include "G4TriangularFacet.hh"
51 #include "G4QuadrangularFacet.hh"
52 #include "G4VoxelLimits.hh"
53 #include "G4AffineTransform.hh"
54 #include "Randomize.hh"
55 
56 #include "G4VGraphicsScene.hh"
57 #include "G4Polyhedron.hh"
58 #include "G4PolyhedronArbitrary.hh"
59 #include "G4NURBS.hh"
60 #include "G4NURBSbox.hh"
61 #include "G4VisExtent.hh"
62 
63 const G4int G4GenericTrap::fgkNofVertices = 8;
64 const G4double G4GenericTrap::fgkTolerance = 1E-3;
65 
66 // --------------------------------------------------------------------
67 
69  const std::vector<G4TwoVector>& vertices )
70  : G4VSolid(name),
71  fpPolyhedron(0),
72  fDz(halfZ),
73  fVertices(),
74  fIsTwisted(false),
75  fTessellatedSolid(0),
76  fMinBBoxVector(G4ThreeVector(0,0,0)),
77  fMaxBBoxVector(G4ThreeVector(0,0,0)),
78  fVisSubdivisions(0),
79  fSurfaceArea(0.),
80  fCubicVolume(0.)
81 
82 {
83  // General constructor
84  const G4double min_length=5*1.e-6;
85  G4double length = 0.;
86  G4int k=0;
87  G4String errorDescription = "InvalidSetup in \" ";
88  errorDescription += name;
89  errorDescription += "\"";
90 
91  // Check vertices size
92 
93  if ( G4int(vertices.size()) != fgkNofVertices )
94  {
95  G4Exception("G4GenericTrap::G4GenericTrap()", "GeomSolids0002",
96  FatalErrorInArgument, "Number of vertices != 8");
97  }
98 
99  // Check dZ
100  //
101  if (halfZ < kCarTolerance)
102  {
103  G4Exception("G4GenericTrap::G4GenericTrap()", "GeomSolids0002",
104  FatalErrorInArgument, "dZ is too small or negative");
105  }
106 
107  // Check Ordering and Copy vertices
108  //
109  if(CheckOrder(vertices))
110  {
111  for (G4int i=0; i<fgkNofVertices; ++i) {fVertices.push_back(vertices[i]);}
112  }
113  else
114  {
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]);}
117  }
118 
119  // Check length of segments and Adjust
120  //
121  for (G4int j=0; j < 2; j++)
122  {
123  for (G4int i=1; i<4; ++i)
124  {
125  k = j*4+i;
126  length = (fVertices[k]-fVertices[k-1]).mag();
127  if ( ( length < min_length) && ( length > kCarTolerance ) )
128  {
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];
136  }
137  }
138  }
139 
140  // Compute Twist
141  //
142  for( G4int i=0; i<4; i++) { fTwist[i]=0.; }
143  fIsTwisted = ComputeIsTwisted();
144 
145  // Compute Bounding Box
146  //
147  ComputeBBox();
148 
149  // If not twisted - create tessellated solid
150  // (an alternative implementation for testing)
151  //
152 #ifdef G4TESS_TEST
153  if ( !fIsTwisted ) { fTessellatedSolid = CreateTessellatedSolid(); }
154 #endif
155 }
156 
157 // --------------------------------------------------------------------
158 
160  : G4VSolid(a),
161  fpPolyhedron(0),
162  fDz(0.),
163  fVertices(),
164  fIsTwisted(false),
165  fTessellatedSolid(0),
166  fMinBBoxVector(G4ThreeVector(0,0,0)),
167  fMaxBBoxVector(G4ThreeVector(0,0,0)),
168  fVisSubdivisions(0),
169  fSurfaceArea(0.),
170  fCubicVolume(0.)
171 {
172  // Fake default constructor - sets only member data and allocates memory
173  // for usage restricted to object persistency.
174 
175  for (size_t i=0; i<4; ++i) { fTwist[i]=0.; }
176 }
177 
178 // --------------------------------------------------------------------
179 
181 {
182  // Destructor
183  delete fTessellatedSolid;
184 }
185 
186 // --------------------------------------------------------------------
187 
189  : G4VSolid(rhs),
190  fpPolyhedron(0), fDz(rhs.fDz), fVertices(rhs.fVertices),
191  fIsTwisted(rhs.fIsTwisted), fTessellatedSolid(0),
192  fMinBBoxVector(rhs.fMinBBoxVector), fMaxBBoxVector(rhs.fMaxBBoxVector),
193  fVisSubdivisions(rhs.fVisSubdivisions),
194  fSurfaceArea(rhs.fSurfaceArea), fCubicVolume(rhs.fCubicVolume)
195 {
196  for (size_t i=0; i<4; ++i) { fTwist[i] = rhs.fTwist[i]; }
197 #ifdef G4TESS_TEST
198  if (rhs.fTessellatedSolid && !fIsTwisted )
199  { fTessellatedSolid = CreateTessellatedSolid(); }
200 #endif
201 }
202 
203 // --------------------------------------------------------------------
204 
206 {
207  // Check assignment to self
208  //
209  if (this == &rhs) { return *this; }
210 
211  // Copy base class data
212  //
213  G4VSolid::operator=(rhs);
214 
215  // Copy data
216  //
217  fpPolyhedron = 0; fDz = rhs.fDz; fVertices = rhs.fVertices;
218  fIsTwisted = rhs.fIsTwisted; fTessellatedSolid = 0;
219  fMinBBoxVector = rhs.fMinBBoxVector; fMaxBBoxVector = rhs.fMaxBBoxVector;
220  fVisSubdivisions = rhs.fVisSubdivisions;
221  fSurfaceArea = rhs.fSurfaceArea; fCubicVolume = rhs.fCubicVolume;
222 
223  for (size_t i=0; i<4; ++i) { fTwist[i] = rhs.fTwist[i]; }
224 #ifdef G4TESS_TEST
225  if (rhs.fTessellatedSolid && !fIsTwisted )
226  { delete fTessellatedSolid; fTessellatedSolid = CreateTessellatedSolid(); }
227 #endif
228 
229  return *this;
230 }
231 
232 // --------------------------------------------------------------------
233 
234 EInside
235 G4GenericTrap::InsidePolygone(const G4ThreeVector& p,
236  const std::vector<G4TwoVector>& poly) const
237 {
238  static const G4double halfCarTolerance=kCarTolerance*0.5;
239  EInside in = kInside;
240  G4double cross, len2;
241  G4int count=0;
242 
243  for (G4int i = 0; i < 4; i++)
244  {
245  G4int j = (i+1) % 4;
246 
247  cross = (p.x()-poly[i].x())*(poly[j].y()-poly[i].y())-
248  (p.y()-poly[i].y())*(poly[j].x()-poly[i].x());
249 
250  len2=(poly[i]-poly[j]).mag2();
251  if (len2 > kCarTolerance)
252  {
253  if(cross*cross<=len2*halfCarTolerance*halfCarTolerance) // Surface check
254  {
255  G4double test;
256  G4int k,l;
257  if ( poly[i].y() > poly[j].y() ) { k=i; l=j; }
258  else { k=j; l=i; }
259 
260  if ( poly[k].x() != poly[l].x() )
261  {
262  test = (p.x()-poly[l].x())/(poly[k].x()-poly[l].x())
263  * (poly[k].y()-poly[l].y())+poly[l].y();
264  }
265  else
266  {
267  test = p.y();
268  }
269 
270  // Check if point is Inside Segment
271  //
272  if( (test>=(poly[l].y()-halfCarTolerance))
273  && (test<=(poly[k].y()+halfCarTolerance)) )
274  {
275  return kSurface;
276  }
277  else
278  {
279  return kOutside;
280  }
281  }
282  else if (cross<0.) { return kOutside; }
283  }
284  else
285  {
286  count++;
287  }
288  }
289 
290  // All collapsed vertices, Tet like
291  //
292  if(count==4)
293  {
294  if ( (std::fabs(p.x()-poly[0].x())+std::fabs(p.y()-poly[0].y())) > halfCarTolerance )
295  {
296  in=kOutside;
297  }
298  }
299  return in;
300 }
301 
302 // --------------------------------------------------------------------
303 
305 {
306  // Test if point is inside this shape
307 
308 #ifdef G4TESS_TEST
309  if ( fTessellatedSolid )
310  {
311  return fTessellatedSolid->Inside(p);
312  }
313 #endif
314 
315  static const G4double halfCarTolerance=kCarTolerance*0.5;
316  EInside innew=kOutside;
317  std::vector<G4TwoVector> xy;
318 
319  if (std::fabs(p.z()) <= fDz+halfCarTolerance) // First check Z range
320  {
321  // Compute intersection between Z plane containing point and the shape
322  //
323  G4double cf = 0.5*(fDz-p.z())/fDz;
324  for (G4int i=0; i<4; i++)
325  {
326  xy.push_back(fVertices[i+4]+cf*( fVertices[i]-fVertices[i+4]));
327  }
328 
329  innew=InsidePolygone(p,xy);
330 
331  if( (innew==kInside) || (innew==kSurface) )
332  {
333  if(std::fabs(p.z()) > fDz-halfCarTolerance) { innew=kSurface; }
334  }
335  }
336  return innew;
337 }
338 
339 // --------------------------------------------------------------------
340 
342 {
343  // Calculate side nearest to p, and return normal
344  // If two sides are equidistant, sum of the Normal is returned
345 
346 #ifdef G4TESS_TEST
347  if ( fTessellatedSolid )
348  {
349  return fTessellatedSolid->SurfaceNormal(p);
350  }
351 #endif
352 
353  static const G4double halfCarTolerance=kCarTolerance*0.5;
354 
355  G4ThreeVector lnorm, sumnorm(0.,0.,0.), apprnorm(0.,0.,1.),
356  p0, p1, p2, r1, r2, r3, r4;
357  G4int noSurfaces = 0;
358  G4double distxy,distz;
359  G4bool zPlusSide=false;
360 
361  distz = fDz-std::fabs(p.z());
362  if (distz < halfCarTolerance)
363  {
364  if(p.z()>0)
365  {
366  zPlusSide=true;
367  sumnorm=G4ThreeVector(0,0,1);
368  }
369  else
370  {
371  sumnorm=G4ThreeVector(0,0,-1);
372  }
373  noSurfaces ++;
374  }
375 
376  // Check lateral planes
377  //
378  std:: vector<G4TwoVector> vertices;
379  G4double cf = 0.5*(fDz-p.z())/fDz;
380  for (G4int i=0; i<4; i++)
381  {
382  vertices.push_back(fVertices[i+4]+cf*(fVertices[i]-fVertices[i+4]));
383  }
384 
385  // Compute distance for lateral planes
386  //
387  for (G4int q=0; q<4; q++)
388  {
389  p0=G4ThreeVector(vertices[q].x(),vertices[q].y(),p.z());
390  if(zPlusSide)
391  {
392  p1=G4ThreeVector(fVertices[q].x(),fVertices[q].y(),-fDz);
393  }
394  else
395  {
396  p1=G4ThreeVector(fVertices[q+4].x(),fVertices[q+4].y(),fDz);
397  }
398  p2=G4ThreeVector(vertices[(q+1)%4].x(),vertices[(q+1)%4].y(),p.z());
399 
400  // Collapsed vertices
401  //
402  if ( (p2-p0).mag2() < kCarTolerance )
403  {
404  if ( std::fabs(p.z()+fDz) > kCarTolerance )
405  {
406  p2=G4ThreeVector(fVertices[(q+1)%4].x(),fVertices[(q+1)%4].y(),-fDz);
407  }
408  else
409  {
410  p2=G4ThreeVector(fVertices[(q+1)%4+4].x(),fVertices[(q+1)%4+4].y(),fDz);
411  }
412  }
413  lnorm = (p1-p0).cross(p2-p0);
414  lnorm = lnorm.unit();
415  if(zPlusSide) { lnorm=-lnorm; }
416 
417  // Adjust Normal for Twisted Surface
418  //
419  if ( (fIsTwisted) && (GetTwistAngle(q)!=0) )
420  {
421  G4double normP=(p2-p0).mag();
422  if(normP)
423  {
424  G4double proj=(p-p0).dot(p2-p0)/normP;
425  if(proj<0) { proj=0; }
426  if(proj>normP) { proj=normP; }
427  G4int j=(q+1)%4;
428  r1=G4ThreeVector(fVertices[q+4].x(),fVertices[q+4].y(),fDz);
429  r2=G4ThreeVector(fVertices[j+4].x(),fVertices[j+4].y(),fDz);
430  r3=G4ThreeVector(fVertices[q].x(),fVertices[q].y(),-fDz);
431  r4=G4ThreeVector(fVertices[j].x(),fVertices[j].y(),-fDz);
432  r1=r1+proj*(r2-r1)/normP;
433  r3=r3+proj*(r4-r3)/normP;
434  r2=r1-r3;
435  r4=r2.cross(p2-p0); r4=r4.unit();
436  lnorm=r4;
437  }
438  } // End if fIsTwisted
439 
440  distxy=std::fabs((p0-p).dot(lnorm));
441  if ( distxy<halfCarTolerance )
442  {
443  noSurfaces ++;
444 
445  // Negative sign for Normal is taken for Outside Normal
446  //
447  sumnorm=sumnorm+lnorm;
448  }
449 
450  // For ApproxSurfaceNormal
451  //
452  if (distxy<distz)
453  {
454  distz=distxy;
455  apprnorm=lnorm;
456  }
457  } // End for loop
458 
459  // Calculate final Normal, add Normal in the Corners and Touching Sides
460  //
461  if ( noSurfaces == 0 )
462  {
463  G4Exception("G4GenericTrap::SurfaceNormal(p)", "GeomSolids1002",
464  JustWarning, "Point p is not on surface !?" );
465  sumnorm=apprnorm;
466  // Add Approximative Surface Normal Calculation?
467  }
468  else if ( noSurfaces == 1 ) { sumnorm = sumnorm; }
469  else { sumnorm = sumnorm.unit(); }
470 
471  return sumnorm ;
472 }
473 
474 // --------------------------------------------------------------------
475 
476 G4ThreeVector G4GenericTrap::NormalToPlane( const G4ThreeVector& p,
477  const G4int ipl ) const
478 {
479  // Return normal to given lateral plane ipl
480 
481 #ifdef G4TESS_TEST
482  if ( fTessellatedSolid )
483  {
484  return fTessellatedSolid->SurfaceNormal(p);
485  }
486 #endif
487 
488  static const G4double halfCarTolerance=kCarTolerance*0.5;
489  G4ThreeVector lnorm, norm(0.,0.,0.), p0,p1,p2;
490 
491  G4double distz = fDz-p.z();
492  G4int i=ipl; // current plane index
493 
494  G4TwoVector u,v;
495  G4ThreeVector r1,r2,r3,r4;
496  G4double cf = 0.5*(fDz-p.z())/fDz;
497  G4int j=(i+1)%4;
498 
499  u=fVertices[i+4]+cf*(fVertices[i]-fVertices[i+4]);
500  v=fVertices[j+4]+cf*(fVertices[j]-fVertices[j+4]);
501 
502  // Compute cross product
503  //
504  p0=G4ThreeVector(u.x(),u.y(),p.z());
505 
506  if (std::fabs(distz)<halfCarTolerance)
507  {
508  p1=G4ThreeVector(fVertices[i].x(),fVertices[i].y(),-fDz);distz=-1;}
509  else
510  {
511  p1=G4ThreeVector(fVertices[i+4].x(),fVertices[i+4].y(),fDz);
512  }
513  p2=G4ThreeVector(v.x(),v.y(),p.z());
514 
515  // Collapsed vertices
516  //
517  if ( (p2-p0).mag2() < kCarTolerance )
518  {
519  if ( std::fabs(p.z()+fDz) > halfCarTolerance )
520  {
521  p2=G4ThreeVector(fVertices[j].x(),fVertices[j].y(),-fDz);
522  }
523  else
524  {
525  p2=G4ThreeVector(fVertices[j+4].x(),fVertices[j+4].y(),fDz);
526  }
527  }
528  lnorm=-(p1-p0).cross(p2-p0);
529  if (distz>-halfCarTolerance) { lnorm=-lnorm.unit(); }
530  else { lnorm=lnorm.unit(); }
531 
532  // Adjust Normal for Twisted Surface
533  //
534  if( (fIsTwisted) && (GetTwistAngle(ipl)!=0) )
535  {
536  G4double normP=(p2-p0).mag();
537  if(normP)
538  {
539  G4double proj=(p-p0).dot(p2-p0)/normP;
540  if (proj<0) { proj=0; }
541  if (proj>normP) { proj=normP; }
542 
543  r1=G4ThreeVector(fVertices[i+4].x(),fVertices[i+4].y(),fDz);
544  r2=G4ThreeVector(fVertices[j+4].x(),fVertices[j+4].y(),fDz);
545  r3=G4ThreeVector(fVertices[i].x(),fVertices[i].y(),-fDz);
546  r4=G4ThreeVector(fVertices[j].x(),fVertices[j].y(),-fDz);
547  r1=r1+proj*(r2-r1)/normP;
548  r3=r3+proj*(r4-r3)/normP;
549  r2=r1-r3;
550  r4=r2.cross(p2-p0);r4=r4.unit();
551  lnorm=r4;
552  }
553  } // End if fIsTwisted
554 
555  return lnorm;
556 }
557 
558 // --------------------------------------------------------------------
559 
560 G4double G4GenericTrap::DistToPlane(const G4ThreeVector& p,
561  const G4ThreeVector& v,
562  const G4int ipl) const
563 {
564  // Computes distance to plane ipl :
565  // ipl=0 : points 0,4,1,5
566  // ipl=1 : points 1,5,2,6
567  // ipl=2 : points 2,6,3,7
568  // ipl=3 : points 3,7,0,4
569 
570  static const G4double halfCarTolerance=0.5*kCarTolerance;
571  G4double xa,xb,xc,xd,ya,yb,yc,yd;
572 
573  G4int j = (ipl+1)%4;
574 
575  xa=fVertices[ipl].x();
576  ya=fVertices[ipl].y();
577  xb=fVertices[ipl+4].x();
578  yb=fVertices[ipl+4].y();
579  xc=fVertices[j].x();
580  yc=fVertices[j].y();
581  xd=fVertices[4+j].x();
582  yd=fVertices[4+j].y();
583 
584  G4double dz2 =0.5/fDz;
585  G4double tx1 =dz2*(xb-xa);
586  G4double ty1 =dz2*(yb-ya);
587  G4double tx2 =dz2*(xd-xc);
588  G4double ty2 =dz2*(yd-yc);
589  G4double dzp =fDz+p.z();
590  G4double xs1 =xa+tx1*dzp;
591  G4double ys1 =ya+ty1*dzp;
592  G4double xs2 =xc+tx2*dzp;
593  G4double ys2 =yc+ty2*dzp;
594  G4double dxs =xs2-xs1;
595  G4double dys =ys2-ys1;
596  G4double dtx =tx2-tx1;
597  G4double dty =ty2-ty1;
598 
599  G4double a = (dtx*v.y()-dty*v.x()+(tx1*ty2-tx2*ty1)*v.z())*v.z();
600  G4double b = dxs*v.y()-dys*v.x()+(dtx*p.y()-dty*p.x()+ty2*xs1-ty1*xs2
601  + tx1*ys2-tx2*ys1)*v.z();
602  G4double c=dxs*p.y()-dys*p.x()+xs1*ys2-xs2*ys1;
603  G4double q=kInfinity;
604  G4double x1,x2,y1,y2,xp,yp,zi;
605 
606  if (std::fabs(a)<kCarTolerance)
607  {
608  if (std::fabs(b)<kCarTolerance) { return kInfinity; }
609  q=-c/b;
610 
611  // Check if Point is on the Surface
612 
613  if (q>-halfCarTolerance)
614  {
615  if (q<halfCarTolerance)
616  {
617  if (NormalToPlane(p,ipl).dot(v)<=0)
618  { if(Inside(p) != kOutside) { return 0.; } }
619  else
620  { return kInfinity; }
621  }
622 
623  // Check the Intersection
624  //
625  zi=p.z()+q*v.z();
626  if (std::fabs(zi)<fDz)
627  {
628  x1=xs1+tx1*v.z()*q;
629  x2=xs2+tx2*v.z()*q;
630  xp=p.x()+q*v.x();
631  y1=ys1+ty1*v.z()*q;
632  y2=ys2+ty2*v.z()*q;
633  yp=p.y()+q*v.y();
634  zi = (xp-x1)*(xp-x2)+(yp-y1)*(yp-y2);
635  if (zi<=halfCarTolerance) { return q; }
636  }
637  }
638  return kInfinity;
639  }
640  G4double d=b*b-4*a*c;
641  if (d>=0)
642  {
643  if (a>0) { q=0.5*(-b-std::sqrt(d))/a; }
644  else { q=0.5*(-b+std::sqrt(d))/a; }
645 
646  // Check if Point is on the Surface
647  //
648  if (q>-halfCarTolerance)
649  {
650  if(q<halfCarTolerance)
651  {
652  if (NormalToPlane(p,ipl).dot(v)<=0)
653  {
654  if(Inside(p)!= kOutside) { return 0.; }
655  }
656  else // Check second root; return kInfinity
657  {
658  if (a>0) { q=0.5*(-b+std::sqrt(d))/a; }
659  else { q=0.5*(-b-std::sqrt(d))/a; }
660  if (q<=halfCarTolerance) { return kInfinity; }
661  }
662  }
663  // Check the Intersection
664  //
665  zi=p.z()+q*v.z();
666  if (std::fabs(zi)<fDz)
667  {
668  x1=xs1+tx1*v.z()*q;
669  x2=xs2+tx2*v.z()*q;
670  xp=p.x()+q*v.x();
671  y1=ys1+ty1*v.z()*q;
672  y2=ys2+ty2*v.z()*q;
673  yp=p.y()+q*v.y();
674  zi = (xp-x1)*(xp-x2)+(yp-y1)*(yp-y2);
675  if (zi<=halfCarTolerance) { return q; }
676  }
677  }
678  if (a>0) { q=0.5*(-b+std::sqrt(d))/a; }
679  else { q=0.5*(-b-std::sqrt(d))/a; }
680 
681  // Check if Point is on the Surface
682  //
683  if (q>-halfCarTolerance)
684  {
685  if(q<halfCarTolerance)
686  {
687  if (NormalToPlane(p,ipl).dot(v)<=0)
688  {
689  if(Inside(p) != kOutside) { return 0.; }
690  }
691  else // Check second root; return kInfinity.
692  {
693  if (a>0) { q=0.5*(-b-std::sqrt(d))/a; }
694  else { q=0.5*(-b+std::sqrt(d))/a; }
695  if (q<=halfCarTolerance) { return kInfinity; }
696  }
697  }
698  // Check the Intersection
699  //
700  zi=p.z()+q*v.z();
701  if (std::fabs(zi)<fDz)
702  {
703  x1=xs1+tx1*v.z()*q;
704  x2=xs2+tx2*v.z()*q;
705  xp=p.x()+q*v.x();
706  y1=ys1+ty1*v.z()*q;
707  y2=ys2+ty2*v.z()*q;
708  yp=p.y()+q*v.y();
709  zi = (xp-x1)*(xp-x2)+(yp-y1)*(yp-y2);
710  if (zi<=halfCarTolerance) { return q; }
711  }
712  }
713  }
714  return kInfinity;
715 }
716 
717 // --------------------------------------------------------------------
718 
720  const G4ThreeVector& v) const
721 {
722 #ifdef G4TESS_TEST
723  if ( fTessellatedSolid )
724  {
725  return fTessellatedSolid->DistanceToIn(p, v);
726  }
727 #endif
728 
729  static const G4double halfCarTolerance=kCarTolerance*0.5;
730 
731  G4double dist[5];
733 
734  // Check lateral faces
735  //
736  G4int i;
737  for (i=0; i<4; i++)
738  {
739  dist[i]=DistToPlane(p, v, i);
740  }
741 
742  // Check Z planes
743  //
744  dist[4]=kInfinity;
745  if (std::fabs(p.z())>fDz-halfCarTolerance)
746  {
747  if (v.z())
748  {
750  if (p.z()>0)
751  {
752  dist[4] = (fDz-p.z())/v.z();
753  }
754  else
755  {
756  dist[4] = (-fDz-p.z())/v.z();
757  }
758  if (dist[4]<-halfCarTolerance)
759  {
760  dist[4]=kInfinity;
761  }
762  else
763  {
764  if(dist[4]<halfCarTolerance)
765  {
766  if(p.z()>0) { n=G4ThreeVector(0,0,1); }
767  else { n=G4ThreeVector(0,0,-1); }
768  if (n.dot(v)<0) { dist[4]=0.; }
769  else { dist[4]=kInfinity; }
770  }
771  pt=p+dist[4]*v;
772  if (Inside(pt)==kOutside) { dist[4]=kInfinity; }
773  }
774  }
775  }
776  G4double distmin = dist[0];
777  for (i=1;i<5;i++)
778  {
779  if (dist[i] < distmin) { distmin = dist[i]; }
780  }
781 
782  if (distmin<halfCarTolerance) { distmin=0.; }
783 
784  return distmin;
785 }
786 
787 // --------------------------------------------------------------------
788 
790 {
791  // Computes the closest distance from given point to this shape
792 
793 #ifdef G4TESS_TEST
794  if ( fTessellatedSolid )
795  {
796  return fTessellatedSolid->DistanceToIn(p);
797  }
798 #endif
799 
800  G4double safz = std::fabs(p.z())-fDz;
801  if(safz<0) { safz=0; }
802 
803  G4int iseg;
804  G4double safe = safz;
805  G4double safxy = safz;
806 
807  for (iseg=0; iseg<4; iseg++)
808  {
809  safxy = SafetyToFace(p,iseg);
810  if (safxy>safe) { safe=safxy; }
811  }
812 
813  return safe;
814 }
815 
816 // --------------------------------------------------------------------
817 
818 G4double
819 G4GenericTrap::SafetyToFace(const G4ThreeVector& p, const G4int iseg) const
820 {
821  // Estimate distance to lateral plane defined by segment iseg in range [0,3]
822  // Might be negative: plane seen only from inside
823 
824  G4ThreeVector p1,norm;
825  G4double safe;
826 
827  p1=G4ThreeVector(fVertices[iseg].x(),fVertices[iseg].y(),-fDz);
828 
829  norm=NormalToPlane(p,iseg);
830  safe = (p-p1).dot(norm); // Can be negative
831 
832  return safe;
833 }
834 
835 // --------------------------------------------------------------------
836 
837 G4double
838 G4GenericTrap::DistToTriangle(const G4ThreeVector& p,
839  const G4ThreeVector& v, const G4int ipl) const
840 {
841  static const G4double halfCarTolerance=kCarTolerance*0.5;
842 
843  G4double xa=fVertices[ipl].x();
844  G4double ya=fVertices[ipl].y();
845  G4double xb=fVertices[ipl+4].x();
846  G4double yb=fVertices[ipl+4].y();
847  G4int j=(ipl+1)%4;
848  G4double xc=fVertices[j].x();
849  G4double yc=fVertices[j].y();
850  G4double zab=2*fDz;
851  G4double zac=0;
852 
853  if ( (std::fabs(xa-xc)+std::fabs(ya-yc)) < halfCarTolerance )
854  {
855  xc=fVertices[j+4].x();
856  yc=fVertices[j+4].y();
857  zac=2*fDz;
858  zab=2*fDz;
859 
860  //Line case
861  //
862  if ( (std::fabs(xb-xc)+std::fabs(yb-yc)) < halfCarTolerance )
863  {
864  return kInfinity;
865  }
866  }
867  G4double a=(yb-ya)*zac-(yc-ya)*zab;
868  G4double b=(xc-xa)*zab-(xb-xa)*zac;
869  G4double c=(xb-xa)*(yc-ya)-(xc-xa)*(yb-ya);
870  G4double d=-xa*a-ya*b+fDz*c;
871  G4double t=a*v.x()+b*v.y()+c*v.z();
872 
873  if (t!=0)
874  {
875  t=-(a*p.x()+b*p.y()+c*p.z()+d)/t;
876  }
877  if ( (t<halfCarTolerance) && (t>-halfCarTolerance) )
878  {
879  if (NormalToPlane(p,ipl).dot(v)<kCarTolerance)
880  {
881  t=kInfinity;
882  }
883  else
884  {
885  t=0;
886  }
887  }
888  if (Inside(p+v*t) != kSurface) { t=kInfinity; }
889 
890  return t;
891 }
892 
893 // --------------------------------------------------------------------
894 
896  const G4ThreeVector& v,
897  const G4bool calcNorm,
898  G4bool* validNorm,
899  G4ThreeVector* n) const
900 {
901 #ifdef G4TESS_TEST
902  if ( fTessellatedSolid )
903  {
904  return fTessellatedSolid->DistanceToOut(p, v, calcNorm, validNorm, n);
905  }
906 #endif
907 
908  static const G4double halfCarTolerance=kCarTolerance*0.5;
909 
910  G4double distmin;
911  G4bool lateral_cross = false;
912  ESide side = kUndefined;
913 
914  if (calcNorm) { *validNorm=true; } // All normals are valid
915 
916  if (v.z() < 0)
917  {
918  distmin=(-fDz-p.z())/v.z();
919  if (calcNorm) { side=kMZ; *n=G4ThreeVector(0,0,-1); }
920  }
921  else
922  {
923  if (v.z() > 0)
924  {
925  distmin = (fDz-p.z())/v.z();
926  if (calcNorm) { side=kPZ; *n=G4ThreeVector(0,0,1); }
927  }
928  else { distmin = kInfinity; }
929  }
930 
931  G4double dz2 =0.5/fDz;
932  G4double xa,xb,xc,xd;
933  G4double ya,yb,yc,yd;
934 
935  for (G4int ipl=0; ipl<4; ipl++)
936  {
937  G4int j = (ipl+1)%4;
938  xa=fVertices[ipl].x();
939  ya=fVertices[ipl].y();
940  xb=fVertices[ipl+4].x();
941  yb=fVertices[ipl+4].y();
942  xc=fVertices[j].x();
943  yc=fVertices[j].y();
944  xd=fVertices[4+j].x();
945  yd=fVertices[4+j].y();
946 
947  if ( ((std::fabs(xb-xd)+std::fabs(yb-yd))<halfCarTolerance)
948  || ((std::fabs(xa-xc)+std::fabs(ya-yc))<halfCarTolerance) )
949  {
950  G4double q=DistToTriangle(p,v,ipl) ;
951  if ( (q>=0) && (q<distmin) )
952  {
953  distmin=q;
954  lateral_cross=true;
955  side=ESide(ipl+1);
956  }
957  continue;
958  }
959  G4double tx1 =dz2*(xb-xa);
960  G4double ty1 =dz2*(yb-ya);
961  G4double tx2 =dz2*(xd-xc);
962  G4double ty2 =dz2*(yd-yc);
963  G4double dzp =fDz+p.z();
964  G4double xs1 =xa+tx1*dzp;
965  G4double ys1 =ya+ty1*dzp;
966  G4double xs2 =xc+tx2*dzp;
967  G4double ys2 =yc+ty2*dzp;
968  G4double dxs =xs2-xs1;
969  G4double dys =ys2-ys1;
970  G4double dtx =tx2-tx1;
971  G4double dty =ty2-ty1;
972  G4double a = (dtx*v.y()-dty*v.x()+(tx1*ty2-tx2*ty1)*v.z())*v.z();
973  G4double b = dxs*v.y()-dys*v.x()+(dtx*p.y()-dty*p.x()+ty2*xs1-ty1*xs2
974  + tx1*ys2-tx2*ys1)*v.z();
975  G4double c=dxs*p.y()-dys*p.x()+xs1*ys2-xs2*ys1;
976  G4double q=kInfinity;
977 
978  if (std::fabs(a) < kCarTolerance)
979  {
980  if (std::fabs(b) < kCarTolerance) { continue; }
981  q=-c/b;
982 
983  // Check for Point on the Surface
984  //
985  if ((q > -halfCarTolerance) && (q < distmin))
986  {
987  if (q < halfCarTolerance)
988  {
989  if (NormalToPlane(p,ipl).dot(v)<0.) { continue; }
990  }
991  distmin =q;
992  lateral_cross=true;
993  side=ESide(ipl+1);
994  }
995  continue;
996  }
997  G4double d=b*b-4*a*c;
998  if (d >= 0.)
999  {
1000  if (a > 0) { q=0.5*(-b-std::sqrt(d))/a; }
1001  else { q=0.5*(-b+std::sqrt(d))/a; }
1002 
1003  // Check for Point on the Surface
1004  //
1005  if (q > -halfCarTolerance )
1006  {
1007  if (q < distmin)
1008  {
1009  if(q < halfCarTolerance)
1010  {
1011  if (NormalToPlane(p,ipl).dot(v)<0.) // Check second root
1012  {
1013  if (a > 0) { q=0.5*(-b+std::sqrt(d))/a; }
1014  else { q=0.5*(-b-std::sqrt(d))/a; }
1015  if (( q > halfCarTolerance) && (q < distmin))
1016  {
1017  distmin=q;
1018  lateral_cross = true;
1019  side=ESide(ipl+1);
1020  }
1021  continue;
1022  }
1023  }
1024  distmin = q;
1025  lateral_cross = true;
1026  side=ESide(ipl+1);
1027  }
1028  }
1029  else
1030  {
1031  if (a > 0) { q=0.5*(-b+std::sqrt(d))/a; }
1032  else { q=0.5*(-b-std::sqrt(d))/a; }
1033 
1034  // Check for Point on the Surface
1035  //
1036  if ((q > -halfCarTolerance) && (q < distmin))
1037  {
1038  if (q < halfCarTolerance)
1039  {
1040  if (NormalToPlane(p,ipl).dot(v)<0.) // Check second root
1041  {
1042  if (a > 0) { q=0.5*(-b-std::sqrt(d))/a; }
1043  else { q=0.5*(-b+std::sqrt(d))/a; }
1044  if ( ( q > halfCarTolerance) && (q < distmin) )
1045  {
1046  distmin=q;
1047  lateral_cross = true;
1048  side=ESide(ipl+1);
1049  }
1050  continue;
1051  }
1052  }
1053  distmin =q;
1054  lateral_cross = true;
1055  side=ESide(ipl+1);
1056  }
1057  }
1058  }
1059  }
1060  if (!lateral_cross) // Make sure that track crosses the top or bottom
1061  {
1062  if (distmin >= kInfinity) { distmin=kCarTolerance; }
1063  G4ThreeVector pt=p+distmin*v;
1064 
1065  // Check if propagated point is in the polygon
1066  //
1067  G4int i=0;
1068  if (v.z()>0.) { i=4; }
1069  std::vector<G4TwoVector> xy;
1070  for ( G4int j=0; j<4; j++) { xy.push_back(fVertices[i+j]); }
1071 
1072  // Check Inside
1073  //
1074  if (InsidePolygone(pt,xy)==kOutside)
1075  {
1076  if(calcNorm)
1077  {
1078  if (v.z()>0) {side= kPZ; *n = G4ThreeVector(0,0,1);}
1079  else { side=kMZ; *n = G4ThreeVector(0,0,-1);}
1080  }
1081  return 0.;
1082  }
1083  else
1084  {
1085  if(v.z()>0) {side=kPZ;}
1086  else {side=kMZ;}
1087  }
1088  }
1089 
1090  if (calcNorm)
1091  {
1092  G4ThreeVector pt=p+v*distmin;
1093  switch (side)
1094  {
1095  case kXY0:
1096  *n=NormalToPlane(pt,0);
1097  break;
1098  case kXY1:
1099  *n=NormalToPlane(pt,1);
1100  break;
1101  case kXY2:
1102  *n=NormalToPlane(pt,2);
1103  break;
1104  case kXY3:
1105  *n=NormalToPlane(pt,3);
1106  break;
1107  case kPZ:
1108  *n=G4ThreeVector(0,0,1);
1109  break;
1110  case kMZ:
1111  *n=G4ThreeVector(0,0,-1);
1112  break;
1113  default:
1114  DumpInfo();
1115  std::ostringstream message;
1116  G4int oldprc = message.precision(16);
1117  message << "Undefined side for valid surface normal to solid." << G4endl
1118  << "Position:" << G4endl
1119  << " p.x() = " << p.x()/mm << " mm" << G4endl
1120  << " p.y() = " << p.y()/mm << " mm" << G4endl
1121  << " p.z() = " << p.z()/mm << " mm" << G4endl
1122  << "Direction:" << G4endl
1123  << " v.x() = " << v.x() << G4endl
1124  << " v.y() = " << v.y() << G4endl
1125  << " v.z() = " << v.z() << G4endl
1126  << "Proposed distance :" << G4endl
1127  << " distmin = " << distmin/mm << " mm";
1128  message.precision(oldprc);
1129  G4Exception("G4GenericTrap::DistanceToOut(p,v,..)",
1130  "GeomSolids1002", JustWarning, message);
1131  break;
1132  }
1133  }
1134 
1135  if (distmin<halfCarTolerance) { distmin=0.; }
1136 
1137  return distmin;
1138 }
1139 
1140 // --------------------------------------------------------------------
1141 
1143 {
1144 
1145 #ifdef G4TESS_TEST
1146  if ( fTessellatedSolid )
1147  {
1148  return fTessellatedSolid->DistanceToOut(p);
1149  }
1150 #endif
1151 
1152  G4double safz = fDz-std::fabs(p.z());
1153  if (safz<0) { safz = 0; }
1154 
1155  G4double safe = safz;
1156  G4double safxy = safz;
1157 
1158  for (G4int iseg=0; iseg<4; iseg++)
1159  {
1160  safxy = std::fabs(SafetyToFace(p,iseg));
1161  if (safxy < safe) { safe = safxy; }
1162  }
1163 
1164  return safe;
1165 }
1166 
1167 // --------------------------------------------------------------------
1168 
1170  const G4VoxelLimits& pVoxelLimit,
1171  const G4AffineTransform& pTransform,
1172  G4double& pMin, G4double& pMax) const
1173 {
1174 #ifdef G4TESS_TEST
1175  if ( fTessellatedSolid )
1176  {
1177  return fTessellatedSolid->CalculateExtent(pAxis, pVoxelLimit,
1178  pTransform, pMin, pMax);
1179  }
1180 #endif
1181 
1182  // Computes bounding vectors for a shape
1183  //
1184  G4double Dx,Dy;
1185  G4ThreeVector minVec = GetMinimumBBox();
1186  G4ThreeVector maxVec = GetMaximumBBox();
1187  Dx = 0.5*(maxVec.x()- minVec.x());
1188  Dy = 0.5*(maxVec.y()- minVec.y());
1189 
1190  if (!pTransform.IsRotated())
1191  {
1192  // Special case handling for unrotated shapes
1193  // Compute x/y/z mins and maxs respecting limits, with early returns
1194  // if outside limits. Then switch() on pAxis
1195  //
1196  G4double xoffset,xMin,xMax;
1197  G4double yoffset,yMin,yMax;
1198  G4double zoffset,zMin,zMax;
1199 
1200  xoffset=pTransform.NetTranslation().x();
1201  xMin=xoffset-Dx;
1202  xMax=xoffset+Dx;
1203  if (pVoxelLimit.IsXLimited())
1204  {
1205  if ( (xMin>pVoxelLimit.GetMaxXExtent()+kCarTolerance)
1206  || (xMax<pVoxelLimit.GetMinXExtent()-kCarTolerance) )
1207  {
1208  return false;
1209  }
1210  else
1211  {
1212  if (xMin<pVoxelLimit.GetMinXExtent())
1213  {
1214  xMin=pVoxelLimit.GetMinXExtent();
1215  }
1216  if (xMax>pVoxelLimit.GetMaxXExtent())
1217  {
1218  xMax=pVoxelLimit.GetMaxXExtent();
1219  }
1220  }
1221  }
1222 
1223  yoffset=pTransform.NetTranslation().y();
1224  yMin=yoffset-Dy;
1225  yMax=yoffset+Dy;
1226  if (pVoxelLimit.IsYLimited())
1227  {
1228  if ( (yMin>pVoxelLimit.GetMaxYExtent()+kCarTolerance)
1229  || (yMax<pVoxelLimit.GetMinYExtent()-kCarTolerance) )
1230  {
1231  return false;
1232  }
1233  else
1234  {
1235  if (yMin<pVoxelLimit.GetMinYExtent())
1236  {
1237  yMin=pVoxelLimit.GetMinYExtent();
1238  }
1239  if (yMax>pVoxelLimit.GetMaxYExtent())
1240  {
1241  yMax=pVoxelLimit.GetMaxYExtent();
1242  }
1243  }
1244  }
1245 
1246  zoffset=pTransform.NetTranslation().z();
1247  zMin=zoffset-fDz;
1248  zMax=zoffset+fDz;
1249  if (pVoxelLimit.IsZLimited())
1250  {
1251  if ( (zMin>pVoxelLimit.GetMaxZExtent()+kCarTolerance)
1252  || (zMax<pVoxelLimit.GetMinZExtent()-kCarTolerance) )
1253  {
1254  return false;
1255  }
1256  else
1257  {
1258  if (zMin<pVoxelLimit.GetMinZExtent())
1259  {
1260  zMin=pVoxelLimit.GetMinZExtent();
1261  }
1262  if (zMax>pVoxelLimit.GetMaxZExtent())
1263  {
1264  zMax=pVoxelLimit.GetMaxZExtent();
1265  }
1266  }
1267  }
1268 
1269  switch (pAxis)
1270  {
1271  case kXAxis:
1272  pMin = xMin;
1273  pMax = xMax;
1274  break;
1275  case kYAxis:
1276  pMin = yMin;
1277  pMax = yMax;
1278  break;
1279  case kZAxis:
1280  pMin = zMin;
1281  pMax = zMax;
1282  break;
1283  default:
1284  break;
1285  }
1286  pMin-=kCarTolerance;
1287  pMax+=kCarTolerance;
1288 
1289  return true;
1290  }
1291  else
1292  {
1293  // General rotated case - create and clip mesh to boundaries
1294 
1295  G4bool existsAfterClip=false;
1296  G4ThreeVectorList *vertices;
1297 
1298  pMin=+kInfinity;
1299  pMax=-kInfinity;
1300 
1301  // Calculate rotated vertex coordinates
1302  //
1303  vertices=CreateRotatedVertices(pTransform);
1304  ClipCrossSection(vertices,0,pVoxelLimit,pAxis,pMin,pMax);
1305  ClipCrossSection(vertices,4,pVoxelLimit,pAxis,pMin,pMax);
1306  ClipBetweenSections(vertices,0,pVoxelLimit,pAxis,pMin,pMax);
1307 
1308  if ( (pMin!=kInfinity) || (pMax!=-kInfinity) )
1309  {
1310  existsAfterClip=true;
1311 
1312  // Add 2*tolerance to avoid precision troubles
1313  //
1314  pMin-=kCarTolerance;
1315  pMax+=kCarTolerance;
1316  }
1317  else
1318  {
1319  // Check for case where completely enveloping clipping volume.
1320  // If point inside then we are confident that the solid completely
1321  // envelopes the clipping volume. Hence set min/max extents according
1322  // to clipping volume extents along the specified axis.
1323  //
1324  G4ThreeVector clipCentre(
1325  (pVoxelLimit.GetMinXExtent()+pVoxelLimit.GetMaxXExtent())*0.5,
1326  (pVoxelLimit.GetMinYExtent()+pVoxelLimit.GetMaxYExtent())*0.5,
1327  (pVoxelLimit.GetMinZExtent()+pVoxelLimit.GetMaxZExtent())*0.5);
1328 
1329  if (Inside(pTransform.Inverse().TransformPoint(clipCentre))!=kOutside)
1330  {
1331  existsAfterClip=true;
1332  pMin=pVoxelLimit.GetMinExtent(pAxis);
1333  pMax=pVoxelLimit.GetMaxExtent(pAxis);
1334  }
1335  }
1336  delete vertices;
1337  return existsAfterClip;
1338  }
1339 }
1340 
1341 // --------------------------------------------------------------------
1342 
1344 G4GenericTrap::CreateRotatedVertices(const G4AffineTransform& pTransform) const
1345 {
1346  // Create a List containing the transformed vertices
1347  // Ordering [0-3] -fDz cross section
1348  // [4-7] +fDz cross section such that [0] is below [4],
1349  // [1] below [5] etc.
1350  // Note: caller has deletion responsibility
1351 
1352  G4ThreeVector Min = GetMinimumBBox();
1353  G4ThreeVector Max = GetMaximumBBox();
1354 
1355  G4ThreeVectorList *vertices;
1356  vertices=new G4ThreeVectorList();
1357 
1358  if (vertices)
1359  {
1360  vertices->reserve(8);
1361  G4ThreeVector vertex0(Min.x(),Min.y(),Min.z());
1362  G4ThreeVector vertex1(Max.x(),Min.y(),Min.z());
1363  G4ThreeVector vertex2(Max.x(),Max.y(),Min.z());
1364  G4ThreeVector vertex3(Min.x(),Max.y(),Min.z());
1365  G4ThreeVector vertex4(Min.x(),Min.y(),Max.z());
1366  G4ThreeVector vertex5(Max.x(),Min.y(),Max.z());
1367  G4ThreeVector vertex6(Max.x(),Max.y(),Max.z());
1368  G4ThreeVector vertex7(Min.x(),Max.y(),Max.z());
1369 
1370  vertices->push_back(pTransform.TransformPoint(vertex0));
1371  vertices->push_back(pTransform.TransformPoint(vertex1));
1372  vertices->push_back(pTransform.TransformPoint(vertex2));
1373  vertices->push_back(pTransform.TransformPoint(vertex3));
1374  vertices->push_back(pTransform.TransformPoint(vertex4));
1375  vertices->push_back(pTransform.TransformPoint(vertex5));
1376  vertices->push_back(pTransform.TransformPoint(vertex6));
1377  vertices->push_back(pTransform.TransformPoint(vertex7));
1378  }
1379  else
1380  {
1381  G4Exception("G4GenericTrap::CreateRotatedVertices()", "FatalError",
1382  FatalException, "Out of memory - Cannot allocate vertices!");
1383  }
1384  return vertices;
1385 }
1386 
1387 // --------------------------------------------------------------------
1388 
1390 {
1391  return G4String("G4GenericTrap");
1392 }
1393 
1394 // --------------------------------------------------------------------
1395 
1397 {
1398  return new G4GenericTrap(*this);
1399 }
1400 
1401 // --------------------------------------------------------------------
1402 
1403 std::ostream& G4GenericTrap::StreamInfo(std::ostream& os) const
1404 {
1405  G4int oldprc = os.precision(16);
1406  os << "-----------------------------------------------------------\n"
1407  << " *** Dump for solid - " << GetName() << " *** \n"
1408  << " =================================================== \n"
1409  << " Solid geometry type: " << GetEntityType() << G4endl
1410  << " half length Z: " << fDz/mm << " mm \n"
1411  << " list of vertices:\n";
1412 
1413  for ( G4int i=0; i<fgkNofVertices; ++i )
1414  {
1415  os << std::setw(5) << "#" << i
1416  << " vx = " << fVertices[i].x()/mm << " mm"
1417  << " vy = " << fVertices[i].y()/mm << " mm" << G4endl;
1418  }
1419  os.precision(oldprc);
1420 
1421  return os;
1422 }
1423 
1424 // --------------------------------------------------------------------
1425 
1427 {
1428 
1429 #ifdef G4TESS_TEST
1430  if ( fTessellatedSolid )
1431  {
1432  return fTessellatedSolid->GetPointOnSurface();
1433  }
1434 #endif
1435 
1436  G4ThreeVector point;
1437  G4TwoVector u,v,w;
1438  G4double rand,area,chose,cf,lambda0,lambda1,alfa,beta,zp;
1439  G4int ipl,j;
1440 
1441  std::vector<G4ThreeVector> vertices;
1442  for (G4int i=0; i<4;i++)
1443  {
1444  vertices.push_back(G4ThreeVector(fVertices[i].x(),fVertices[i].y(),-fDz));
1445  }
1446  for (G4int i=4; i<8;i++)
1447  {
1448  vertices.push_back(G4ThreeVector(fVertices[i].x(),fVertices[i].y(),fDz));
1449  }
1450 
1451  // Surface Area of Planes(only estimation for twisted)
1452  //
1453  G4double Surface0=GetFaceSurfaceArea(vertices[0],vertices[1],
1454  vertices[2],vertices[3]);//-fDz plane
1455  G4double Surface1=GetFaceSurfaceArea(vertices[0],vertices[1],
1456  vertices[5],vertices[4]);// Lat plane
1457  G4double Surface2=GetFaceSurfaceArea(vertices[3],vertices[0],
1458  vertices[4],vertices[7]);// Lat plane
1459  G4double Surface3=GetFaceSurfaceArea(vertices[2],vertices[3],
1460  vertices[7],vertices[6]);// Lat plane
1461  G4double Surface4=GetFaceSurfaceArea(vertices[2],vertices[1],
1462  vertices[5],vertices[6]);// Lat plane
1463  G4double Surface5=GetFaceSurfaceArea(vertices[4],vertices[5],
1464  vertices[6],vertices[7]);// fDz plane
1465  rand = G4UniformRand();
1466  area = Surface0+Surface1+Surface2+Surface3+Surface4+Surface5;
1467  chose = rand*area;
1468 
1469  if ( ( chose < Surface0)
1470  || ( chose > (Surface0+Surface1+Surface2+Surface3+Surface4)) )
1471  { // fDz or -fDz Plane
1472  ipl = G4int(G4UniformRand()*4);
1473  j = (ipl+1)%4;
1474  if(chose < Surface0)
1475  {
1476  zp = -fDz;
1477  u = fVertices[ipl]; v = fVertices[j];
1478  w = fVertices[(ipl+3)%4];
1479  }
1480  else
1481  {
1482  zp = fDz;
1483  u = fVertices[ipl+4]; v = fVertices[j+4];
1484  w = fVertices[(ipl+3)%4+4];
1485  }
1486  alfa = G4UniformRand();
1487  beta = G4UniformRand();
1488  lambda1=alfa*beta;
1489  lambda0=alfa-lambda1;
1490  v = v-u;
1491  w = w-u;
1492  v = u+lambda0*v+lambda1*w;
1493  }
1494  else // Lateral Plane Twisted or Not
1495  {
1496  if (chose < Surface0+Surface1) { ipl=0; }
1497  else if (chose < Surface0+Surface1+Surface2) { ipl=1; }
1498  else if (chose < Surface0+Surface1+Surface2+Surface3) { ipl=2; }
1499  else { ipl=3; }
1500  j = (ipl+1)%4;
1501  zp = -fDz+G4UniformRand()*2*fDz;
1502  cf = 0.5*(fDz-zp)/fDz;
1503  u = fVertices[ipl+4]+cf*( fVertices[ipl]-fVertices[ipl+4]);
1504  v = fVertices[j+4]+cf*(fVertices[j]-fVertices[j+4]);
1505  v = u+(v-u)*G4UniformRand();
1506  }
1507  point=G4ThreeVector(v.x(),v.y(),zp);
1508 
1509  return point;
1510 }
1511 
1512 // --------------------------------------------------------------------
1513 
1515 {
1516  if(fCubicVolume != 0.) {;}
1517  else { fCubicVolume = G4VSolid::GetCubicVolume(); }
1518  return fCubicVolume;
1519 }
1520 
1521 // --------------------------------------------------------------------
1522 
1524 {
1525  if(fSurfaceArea != 0.) {;}
1526  else
1527  {
1528  std::vector<G4ThreeVector> vertices;
1529  for (G4int i=0; i<4;i++)
1530  {
1531  vertices.push_back(G4ThreeVector(fVertices[i].x(),fVertices[i].y(),-fDz));
1532  }
1533  for (G4int i=4; i<8;i++)
1534  {
1535  vertices.push_back(G4ThreeVector(fVertices[i].x(),fVertices[i].y(),fDz));
1536  }
1537 
1538  // Surface Area of Planes(only estimation for twisted)
1539  //
1540  G4double fSurface0=GetFaceSurfaceArea(vertices[0],vertices[1],
1541  vertices[2],vertices[3]);//-fDz plane
1542  G4double fSurface1=GetFaceSurfaceArea(vertices[0],vertices[1],
1543  vertices[5],vertices[4]);// Lat plane
1544  G4double fSurface2=GetFaceSurfaceArea(vertices[3],vertices[0],
1545  vertices[4],vertices[7]);// Lat plane
1546  G4double fSurface3=GetFaceSurfaceArea(vertices[2],vertices[3],
1547  vertices[7],vertices[6]);// Lat plane
1548  G4double fSurface4=GetFaceSurfaceArea(vertices[2],vertices[1],
1549  vertices[5],vertices[6]);// Lat plane
1550  G4double fSurface5=GetFaceSurfaceArea(vertices[4],vertices[5],
1551  vertices[6],vertices[7]);// fDz plane
1552 
1553  // Total Surface Area
1554  //
1555  if(!fIsTwisted)
1556  {
1557  fSurfaceArea = fSurface0+fSurface1+fSurface2
1558  + fSurface3+fSurface4+fSurface5;
1559  }
1560  else
1561  {
1562  fSurfaceArea = G4VSolid::GetSurfaceArea();
1563  }
1564  }
1565  return fSurfaceArea;
1566 }
1567 
1568 // --------------------------------------------------------------------
1569 
1570 G4double G4GenericTrap::GetFaceSurfaceArea(const G4ThreeVector& p0,
1571  const G4ThreeVector& p1,
1572  const G4ThreeVector& p2,
1573  const G4ThreeVector& p3) const
1574 {
1575  // Auxiliary method for Get Surface Area of Face
1576 
1577  G4double aOne, aTwo;
1578  G4ThreeVector t, u, v, w, Area, normal;
1579 
1580  t = p2 - p1;
1581  u = p0 - p1;
1582  v = p2 - p3;
1583  w = p0 - p3;
1584 
1585  Area = w.cross(v);
1586  aOne = 0.5*Area.mag();
1587 
1588  Area = t.cross(u);
1589  aTwo = 0.5*Area.mag();
1590 
1591  return aOne + aTwo;
1592 }
1593 
1594 // --------------------------------------------------------------------
1595 
1596 G4bool G4GenericTrap::ComputeIsTwisted()
1597 {
1598  // Computes tangents of twist angles (angles between projections on XY plane
1599  // of corresponding -dz +dz edges).
1600 
1601  G4bool twisted = false;
1602  G4double dx1, dy1, dx2, dy2;
1603  G4int nv = fgkNofVertices/2;
1604 
1605  for ( G4int i=0; i<4; i++ )
1606  {
1607  dx1 = fVertices[(i+1)%nv].x()-fVertices[i].x();
1608  dy1 = fVertices[(i+1)%nv].y()-fVertices[i].y();
1609  if ( (dx1 == 0) && (dy1 == 0) ) { continue; }
1610 
1611  dx2 = fVertices[nv+(i+1)%nv].x()-fVertices[nv+i].x();
1612  dy2 = fVertices[nv+(i+1)%nv].y()-fVertices[nv+i].y();
1613 
1614  if ( dx2 == 0 && dy2 == 0 ) { continue; }
1615  G4double twist_angle = std::fabs(dy1*dx2 - dx1*dy2);
1616  if ( twist_angle < fgkTolerance ) { continue; }
1617  twisted = true;
1618  SetTwistAngle(i,twist_angle);
1619 
1620  // Check on big angles, potentially navigation problem
1621 
1622  twist_angle = std::acos( (dx1*dx2 + dy1*dy2)
1623  / (std::sqrt(dx1*dx1+dy1*dy1)
1624  * std::sqrt(dx2*dx2+dy2*dy2)) );
1625 
1626  if ( std::fabs(twist_angle) > 0.5*pi+kCarTolerance )
1627  {
1628  std::ostringstream message;
1629  message << "Twisted Angle is bigger than 90 degrees - " << GetName()
1630  << G4endl
1631  << " Potential problem of malformed Solid !" << G4endl
1632  << " TwistANGLE = " << twist_angle
1633  << "*rad for lateral plane N= " << i;
1634  G4Exception("G4GenericTrap::ComputeIsTwisted()", "GeomSolids1002",
1635  JustWarning, message);
1636  }
1637  }
1638 
1639  return twisted;
1640 }
1641 
1642 // --------------------------------------------------------------------
1643 
1644 G4bool G4GenericTrap::CheckOrder(const std::vector<G4TwoVector>& vertices) const
1645 {
1646  // Test if the vertices are in a clockwise order, if not reorder them.
1647  // Also test if they're well defined without crossing opposite segments
1648 
1649  G4bool clockwise_order=true;
1650  G4double sum1 = 0.;
1651  G4double sum2 = 0.;
1652  G4int j;
1653 
1654  for (G4int i=0; i<4; i++)
1655  {
1656  j = (i+1)%4;
1657  sum1 += vertices[i].x()*vertices[j].y() - vertices[j].x()*vertices[i].y();
1658  sum2 += vertices[i+4].x()*vertices[j+4].y()
1659  - vertices[j+4].x()*vertices[i+4].y();
1660  }
1661  if (sum1*sum2 < -fgkTolerance)
1662  {
1663  std::ostringstream message;
1664  message << "Lower/upper faces defined with opposite clockwise - "
1665  << GetName();
1666  G4Exception("G4GenericTrap::CheckOrder()", "GeomSolids0002",
1667  FatalException, message);
1668  }
1669 
1670  if ((sum1 > 0.)||(sum2 > 0.))
1671  {
1672  std::ostringstream message;
1673  message << "Vertices must be defined in clockwise XY planes - "
1674  << GetName();
1675  G4Exception("G4GenericTrap::CheckOrder()", "GeomSolids1001",
1676  JustWarning,message, "Re-ordering...");
1677  clockwise_order = false;
1678  }
1679 
1680  // Check for illegal crossings
1681  //
1682  G4bool illegal_cross = false;
1683  illegal_cross = IsSegCrossingZ(vertices[0],vertices[4],
1684  vertices[1],vertices[5]);
1685 
1686  if (!illegal_cross)
1687  {
1688  illegal_cross = IsSegCrossingZ(vertices[2],vertices[6],
1689  vertices[3],vertices[7]);
1690  }
1691  // +/- dZ planes
1692  if (!illegal_cross)
1693  {
1694  illegal_cross = IsSegCrossing(vertices[0],vertices[1],
1695  vertices[2],vertices[3]);
1696  }
1697  if (!illegal_cross)
1698  {
1699  illegal_cross = IsSegCrossing(vertices[0],vertices[3],
1700  vertices[1],vertices[2]);
1701  }
1702  if (!illegal_cross)
1703  {
1704  illegal_cross = IsSegCrossing(vertices[4],vertices[5],
1705  vertices[6],vertices[7]);
1706  }
1707  if (!illegal_cross)
1708  {
1709  illegal_cross = IsSegCrossing(vertices[4],vertices[7],
1710  vertices[5],vertices[6]);
1711  }
1712 
1713  if (illegal_cross)
1714  {
1715  std::ostringstream message;
1716  message << "Malformed polygone with opposite sides - " << GetName();
1717  G4Exception("G4GenericTrap::CheckOrderAndSetup()",
1718  "GeomSolids0002", FatalException, message);
1719  }
1720  return clockwise_order;
1721 }
1722 
1723 // --------------------------------------------------------------------
1724 
1725 void G4GenericTrap::ReorderVertices(std::vector<G4ThreeVector>& vertices) const
1726 {
1727  // Reorder the vector of vertices
1728 
1729  std::vector<G4ThreeVector> oldVertices(vertices);
1730 
1731  for ( G4int i=0; i < G4int(oldVertices.size()); ++i )
1732  {
1733  vertices[i] = oldVertices[oldVertices.size()-1-i];
1734  }
1735 }
1736 
1737 // --------------------------------------------------------------------
1738 
1739 G4bool
1740 G4GenericTrap::IsSegCrossing(const G4TwoVector& a, const G4TwoVector& b,
1741  const G4TwoVector& c, const G4TwoVector& d) const
1742 {
1743  // Check if segments [A,B] and [C,D] are crossing
1744 
1745  G4bool stand1 = false;
1746  G4bool stand2 = false;
1747  G4double dx1,dx2,xm=0.,ym=0.,a1=0.,a2=0.,b1=0.,b2=0.;
1748  dx1=(b-a).x();
1749  dx2=(d-c).x();
1750 
1751  if( std::fabs(dx1) < fgkTolerance ) { stand1 = true; }
1752  if( std::fabs(dx2) < fgkTolerance ) { stand2 = true; }
1753  if (!stand1)
1754  {
1755  a1 = (b.x()*a.y()-a.x()*b.y())/dx1;
1756  b1 = (b-a).y()/dx1;
1757  }
1758  if (!stand2)
1759  {
1760  a2 = (d.x()*c.y()-c.x()*d.y())/dx2;
1761  b2 = (d-c).y()/dx2;
1762  }
1763  if (stand1 && stand2)
1764  {
1765  // Segments parallel and vertical
1766  //
1767  if (std::fabs(a.x()-c.x())<fgkTolerance)
1768  {
1769  // Check if segments are overlapping
1770  //
1771  if ( ((c.y()-a.y())*(c.y()-b.y())<-fgkTolerance)
1772  || ((d.y()-a.y())*(d.y()-b.y())<-fgkTolerance)
1773  || ((a.y()-c.y())*(a.y()-d.y())<-fgkTolerance)
1774  || ((b.y()-c.y())*(b.y()-d.y())<-fgkTolerance) ) { return true; }
1775 
1776  return false;
1777  }
1778  // Different x values
1779  //
1780  return false;
1781  }
1782 
1783  if (stand1) // First segment vertical
1784  {
1785  xm = a.x();
1786  ym = a2+b2*xm;
1787  }
1788  else
1789  {
1790  if (stand2) // Second segment vertical
1791  {
1792  xm = c.x();
1793  ym = a1+b1*xm;
1794  }
1795  else // Normal crossing
1796  {
1797  if (std::fabs(b1-b2) < fgkTolerance)
1798  {
1799  // Parallel segments, are they aligned
1800  //
1801  if (std::fabs(c.y()-(a1+b1*c.x())) > fgkTolerance) { return false; }
1802 
1803  // Aligned segments, are they overlapping
1804  //
1805  if ( ((c.x()-a.x())*(c.x()-b.x())<-fgkTolerance)
1806  || ((d.x()-a.x())*(d.x()-b.x())<-fgkTolerance)
1807  || ((a.x()-c.x())*(a.x()-d.x())<-fgkTolerance)
1808  || ((b.x()-c.x())*(b.x()-d.x())<-fgkTolerance) ) { return true; }
1809 
1810  return false;
1811  }
1812  xm = (a1-a2)/(b2-b1);
1813  ym = (a1*b2-a2*b1)/(b2-b1);
1814  }
1815  }
1816 
1817  // Check if crossing point is both between A,B and C,D
1818  //
1819  G4double check = (xm-a.x())*(xm-b.x())+(ym-a.y())*(ym-b.y());
1820  if (check > -fgkTolerance) { return false; }
1821  check = (xm-c.x())*(xm-d.x())+(ym-c.y())*(ym-d.y());
1822  if (check > -fgkTolerance) { return false; }
1823 
1824  return true;
1825 }
1826 
1827 // --------------------------------------------------------------------
1828 
1829 G4bool
1830 G4GenericTrap::IsSegCrossingZ(const G4TwoVector& a, const G4TwoVector& b,
1831  const G4TwoVector& c, const G4TwoVector& d) const
1832 {
1833  // Check if segments [A,B] and [C,D] are crossing when
1834  // A and C are on -dZ and B and D are on +dZ
1835 
1836  // Calculate the Intersection point between two lines in 3D
1837  //
1838  G4ThreeVector temp1,temp2;
1839  G4ThreeVector v1,v2,p1,p2,p3,p4,dv;
1840  G4double q,det;
1841  p1=G4ThreeVector(a.x(),a.y(),-fDz);
1842  p2=G4ThreeVector(c.x(),c.y(),-fDz);
1843  p3=G4ThreeVector(b.x(),b.y(),fDz);
1844  p4=G4ThreeVector(d.x(),d.y(),fDz);
1845  v1=p3-p1;
1846  v2=p4-p2;
1847  dv=p2-p1;
1848 
1849  // In case of Collapsed Vertices No crossing
1850  //
1851  if( (std::fabs(dv.x()) < kCarTolerance )&&
1852  (std::fabs(dv.y()) < kCarTolerance ) ) { return false; }
1853 
1854  if( (std::fabs((p4-p3).x()) < kCarTolerance )&&
1855  (std::fabs((p4-p3).y()) < kCarTolerance ) ) { return false; }
1856 
1857  // First estimate if Intersection is possible( if det is 0)
1858  //
1859  det = dv.x()*v1.y()*v2.z()+dv.y()*v1.z()*v2.x()
1860  - dv.x()*v1.z()*v2.y()-dv.y()*v1.x()*v2.z();
1861 
1862  if (std::fabs(det)<kCarTolerance) //Intersection
1863  {
1864  temp1 = v1.cross(v2);
1865  temp2 = (p2-p1).cross(v2);
1866  if (temp1.dot(temp2) < 0) { return false; } // intersection negative
1867  q = temp1.mag();
1868 
1869  if ( q < kCarTolerance ) { return false; } // parallel lines
1870  q = ((dv).cross(v2)).mag()/q;
1871 
1872  if(q < 1.-kCarTolerance) { return true; }
1873  }
1874  return false;
1875 }
1876 
1877 // --------------------------------------------------------------------
1878 
1879 G4VFacet*
1880 G4GenericTrap::MakeDownFacet(const std::vector<G4ThreeVector>& fromVertices,
1881  G4int ind1, G4int ind2, G4int ind3) const
1882 {
1883  // Create a triangular facet from the polygon points given by indices
1884  // forming the down side ( the normal goes in -z)
1885  // Do not create facet if 2 vertices are the same
1886 
1887  if ( (fromVertices[ind1] == fromVertices[ind2]) ||
1888  (fromVertices[ind2] == fromVertices[ind3]) ||
1889  (fromVertices[ind1] == fromVertices[ind3]) ) { return 0; }
1890 
1891  std::vector<G4ThreeVector> vertices;
1892  vertices.push_back(fromVertices[ind1]);
1893  vertices.push_back(fromVertices[ind2]);
1894  vertices.push_back(fromVertices[ind3]);
1895 
1896  // first vertex most left
1897  //
1898  G4ThreeVector cross=(vertices[1]-vertices[0]).cross(vertices[2]-vertices[1]);
1899 
1900  if ( cross.z() > 0.0 )
1901  {
1902  // Should not happen, as vertices should have been reordered at this stage
1903 
1904  std::ostringstream message;
1905  message << "Vertices in wrong order - " << GetName();
1906  G4Exception("G4GenericTrap::MakeDownFacet", "GeomSolids0002",
1907  FatalException, message);
1908  }
1909 
1910  return new G4TriangularFacet(vertices[0], vertices[1], vertices[2], ABSOLUTE);
1911 }
1912 
1913 // --------------------------------------------------------------------
1914 
1915 G4VFacet*
1916 G4GenericTrap::MakeUpFacet(const std::vector<G4ThreeVector>& fromVertices,
1917  G4int ind1, G4int ind2, G4int ind3) const
1918 {
1919  // Create a triangular facet from the polygon points given by indices
1920  // forming the upper side ( z>0 )
1921 
1922  // Do not create facet if 2 vertices are the same
1923  //
1924  if ( (fromVertices[ind1] == fromVertices[ind2]) ||
1925  (fromVertices[ind2] == fromVertices[ind3]) ||
1926  (fromVertices[ind1] == fromVertices[ind3]) ) { return 0; }
1927 
1928  std::vector<G4ThreeVector> vertices;
1929  vertices.push_back(fromVertices[ind1]);
1930  vertices.push_back(fromVertices[ind2]);
1931  vertices.push_back(fromVertices[ind3]);
1932 
1933  // First vertex most left
1934  //
1935  G4ThreeVector cross=(vertices[1]-vertices[0]).cross(vertices[2]-vertices[1]);
1936 
1937  if ( cross.z() < 0.0 )
1938  {
1939  // Should not happen, as vertices should have been reordered at this stage
1940 
1941  std::ostringstream message;
1942  message << "Vertices in wrong order - " << GetName();
1943  G4Exception("G4GenericTrap::MakeUpFacet", "GeomSolids0002",
1944  FatalException, message);
1945  }
1946 
1947  return new G4TriangularFacet(vertices[0], vertices[1], vertices[2], ABSOLUTE);
1948 }
1949 
1950 // --------------------------------------------------------------------
1951 
1952 G4VFacet*
1953 G4GenericTrap::MakeSideFacet(const G4ThreeVector& downVertex0,
1954  const G4ThreeVector& downVertex1,
1955  const G4ThreeVector& upVertex1,
1956  const G4ThreeVector& upVertex0) const
1957 {
1958  // Creates a triangular facet from the polygon points given by indices
1959  // forming the upper side ( z>0 )
1960 
1961  if ( (downVertex0 == downVertex1) && (upVertex0 == upVertex1) )
1962  {
1963  return 0;
1964  }
1965 
1966  if ( downVertex0 == downVertex1 )
1967  {
1968  return new G4TriangularFacet(downVertex0, upVertex1, upVertex0, ABSOLUTE);
1969  }
1970 
1971  if ( upVertex0 == upVertex1 )
1972  {
1973  return new G4TriangularFacet(downVertex0, downVertex1, upVertex0, ABSOLUTE);
1974  }
1975 
1976  return new G4QuadrangularFacet(downVertex0, downVertex1,
1977  upVertex1, upVertex0, ABSOLUTE);
1978 }
1979 
1980 // --------------------------------------------------------------------
1981 
1982 G4TessellatedSolid* G4GenericTrap::CreateTessellatedSolid() const
1983 {
1984  // 3D vertices
1985  //
1986  G4int nv = fgkNofVertices/2;
1987  std::vector<G4ThreeVector> downVertices;
1988  for ( G4int i=0; i<nv; i++ )
1989  {
1990  downVertices.push_back(G4ThreeVector(fVertices[i].x(),
1991  fVertices[i].y(), -fDz));
1992  }
1993 
1994  std::vector<G4ThreeVector> upVertices;
1995  for ( G4int i=nv; i<2*nv; i++ )
1996  {
1997  upVertices.push_back(G4ThreeVector(fVertices[i].x(),
1998  fVertices[i].y(), fDz));
1999  }
2000 
2001  // Reorder vertices if they are not ordered anti-clock wise
2002  //
2003  G4ThreeVector cross
2004  = (downVertices[1]-downVertices[0]).cross(downVertices[2]-downVertices[1]);
2005  G4ThreeVector cross1
2006  = (upVertices[1]-upVertices[0]).cross(upVertices[2]-upVertices[1]);
2007  if ( (cross.z() > 0.0) || (cross1.z() > 0.0) )
2008  {
2009  ReorderVertices(downVertices);
2010  ReorderVertices(upVertices);
2011  }
2012 
2013  G4TessellatedSolid* tessellatedSolid = new G4TessellatedSolid(GetName());
2014 
2015  G4VFacet* facet = 0;
2016  facet = MakeDownFacet(downVertices, 0, 1, 2);
2017  if (facet) { tessellatedSolid->AddFacet( facet ); }
2018  facet = MakeDownFacet(downVertices, 0, 2, 3);
2019  if (facet) { tessellatedSolid->AddFacet( facet ); }
2020  facet = MakeUpFacet(upVertices, 0, 2, 1);
2021  if (facet) { tessellatedSolid->AddFacet( facet ); }
2022  facet = MakeUpFacet(upVertices, 0, 3, 2);
2023  if (facet) { tessellatedSolid->AddFacet( facet ); }
2024 
2025  // The quadrangular sides
2026  //
2027  for ( G4int i = 0; i < nv; ++i )
2028  {
2029  G4int j = (i+1) % nv;
2030  facet = MakeSideFacet(downVertices[j], downVertices[i],
2031  upVertices[i], upVertices[j]);
2032 
2033  if ( facet ) { tessellatedSolid->AddFacet( facet ); }
2034  }
2035 
2036  tessellatedSolid->SetSolidClosed(true);
2037 
2038  return tessellatedSolid;
2039 }
2040 
2041 // --------------------------------------------------------------------
2042 
2043 void G4GenericTrap::ComputeBBox()
2044 {
2045  // Computes bounding box for a shape.
2046 
2047  G4double minX, maxX, minY, maxY;
2048  minX = maxX = fVertices[0].x();
2049  minY = maxY = fVertices[0].y();
2050 
2051  for (G4int i=1; i< fgkNofVertices; i++)
2052  {
2053  if (minX>fVertices[i].x()) { minX=fVertices[i].x(); }
2054  if (maxX<fVertices[i].x()) { maxX=fVertices[i].x(); }
2055  if (minY>fVertices[i].y()) { minY=fVertices[i].y(); }
2056  if (maxY<fVertices[i].y()) { maxY=fVertices[i].y(); }
2057  }
2058  fMinBBoxVector = G4ThreeVector(minX,minY,-fDz);
2059  fMaxBBoxVector = G4ThreeVector(maxX,maxY, fDz);
2060 }
2061 
2062 // --------------------------------------------------------------------
2063 
2065 {
2066 
2067 #ifdef G4TESS_TEST
2068  if ( fTessellatedSolid )
2069  {
2070  return fTessellatedSolid->GetPolyhedron();
2071  }
2072 #endif
2073 
2074  if ( (!fpPolyhedron)
2077  {
2078  delete fpPolyhedron;
2080  }
2081  return fpPolyhedron;
2082 }
2083 
2084 // --------------------------------------------------------------------
2085 
2087 {
2088 
2089 #ifdef G4TESS_TEST
2090  if ( fTessellatedSolid )
2091  {
2092  return fTessellatedSolid->DescribeYourselfTo(scene);
2093  }
2094 #endif
2095 
2096  scene.AddSolid(*this);
2097 }
2098 
2099 // --------------------------------------------------------------------
2100 
2102 {
2103  // Computes bounding vectors for the shape
2104 
2105 #ifdef G4TESS_TEST
2106  if ( fTessellatedSolid )
2107  {
2108  return fTessellatedSolid->GetExtent();
2109  }
2110 #endif
2111 
2112  G4double Dx,Dy;
2113  G4ThreeVector minVec = GetMinimumBBox();
2114  G4ThreeVector maxVec = GetMaximumBBox();
2115  Dx = 0.5*(maxVec.x()- minVec.x());
2116  Dy = 0.5*(maxVec.y()- minVec.y());
2117 
2118  return G4VisExtent (-Dx, Dx, -Dy, Dy, -fDz, fDz);
2119 }
2120 
2121 // --------------------------------------------------------------------
2122 
2124 {
2125 
2126 #ifdef G4TESS_TEST
2127  if ( fTessellatedSolid )
2128  {
2129  return fTessellatedSolid->CreatePolyhedron();
2130  }
2131 #endif
2132 
2133  // Approximation of Twisted Side
2134  // Construct extra Points, if Twisted Side
2135  //
2136  G4PolyhedronArbitrary* polyhedron;
2137  size_t nVertices, nFacets;
2138 
2139  G4int subdivisions=0;
2140  G4int i;
2141  if(fIsTwisted)
2142  {
2143  if ( GetVisSubdivisions()!= 0 )
2144  {
2145  subdivisions=GetVisSubdivisions();
2146  }
2147  else
2148  {
2149  // Estimation of Number of Subdivisions for smooth visualisation
2150  //
2151  G4double maxTwist=0.;
2152  for(i=0; i<4; i++)
2153  {
2154  if(GetTwistAngle(i)>maxTwist) { maxTwist=GetTwistAngle(i); }
2155  }
2156 
2157  // Computes bounding vectors for the shape
2158  //
2159  G4double Dx,Dy;
2160  G4ThreeVector minVec = GetMinimumBBox();
2161  G4ThreeVector maxVec = GetMaximumBBox();
2162  Dx = 0.5*(maxVec.x()- minVec.y());
2163  Dy = 0.5*(maxVec.y()- minVec.y());
2164  if (Dy > Dx) { Dx=Dy; }
2165 
2166  subdivisions=8*G4int(maxTwist/(Dx*Dx*Dx)*fDz);
2167  if (subdivisions<4) { subdivisions=4; }
2168  if (subdivisions>30) { subdivisions=30; }
2169  }
2170  }
2171  G4int sub4=4*subdivisions;
2172  nVertices = 8+subdivisions*4;
2173  nFacets = 6+subdivisions*4;
2174  G4double cf=1./(subdivisions+1);
2175  polyhedron = new G4PolyhedronArbitrary (nVertices, nFacets);
2176 
2177  // Add Vertex
2178  //
2179  for (i=0;i<4;i++)
2180  {
2181  polyhedron->AddVertex(G4ThreeVector(fVertices[i].x(),
2182  fVertices[i].y(),-fDz));
2183  }
2184  for( i=0;i<subdivisions;i++)
2185  {
2186  for(G4int j=0;j<4;j++)
2187  {
2188  G4TwoVector u=fVertices[j]+cf*(i+1)*( fVertices[j+4]-fVertices[j]);
2189  polyhedron->AddVertex(G4ThreeVector(u.x(),u.y(),-fDz+cf*2*fDz*(i+1)));
2190  }
2191  }
2192  for (i=4;i<8;i++)
2193  {
2194  polyhedron->AddVertex(G4ThreeVector(fVertices[i].x(),
2195  fVertices[i].y(),fDz));
2196  }
2197 
2198  // Add Facets
2199  //
2200  polyhedron->AddFacet(1,4,3,2); //Z-plane
2201  for (i=0;i<subdivisions+1;i++)
2202  {
2203  G4int is=i*4;
2204  polyhedron->AddFacet(5+is,8+is,4+is,1+is);
2205  polyhedron->AddFacet(8+is,7+is,3+is,4+is);
2206  polyhedron->AddFacet(7+is,6+is,2+is,3+is);
2207  polyhedron->AddFacet(6+is,5+is,1+is,2+is);
2208  }
2209  polyhedron->AddFacet(5+sub4,6+sub4,7+sub4,8+sub4); //Z-plane
2210 
2211  polyhedron->SetReferences();
2212  polyhedron->InvertFacets();
2213 
2214  return (G4Polyhedron*) polyhedron;
2215 }
2216 
2217 // --------------------------------------------------------------------
2218 
2220 {
2221 #ifdef G4TESS_TEST
2222  if ( fTessellatedSolid )
2223  {
2224  return fTessellatedSolid->CreateNURBS();
2225  }
2226 #endif
2227 
2228  // Computes bounding vectors for the shape
2229  //
2230  G4double Dx,Dy;
2231  G4ThreeVector minVec = GetMinimumBBox();
2232  G4ThreeVector maxVec = GetMaximumBBox();
2233  Dx = 0.5*(maxVec.x()- minVec.y());
2234  Dy = 0.5*(maxVec.y()- minVec.y());
2235 
2236  return new G4NURBSbox (Dx, Dy, fDz);
2237 }
2238 
2239 // --------------------------------------------------------------------