Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4BREPSolid.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$
28 //
29 // ----------------------------------------------------------------------
30 // GEANT 4 class source file
31 //
32 // G4BREPSolid.cc
33 //
34 // ----------------------------------------------------------------------
35 
36 #include "G4BREPSolid.hh"
37 #include "G4VoxelLimits.hh"
38 #include "G4AffineTransform.hh"
39 #include "G4VGraphicsScene.hh"
40 #include "G4Polyhedron.hh"
41 #include "G4NURBSbox.hh"
42 #include "G4BoundingBox3D.hh"
43 #include "G4FPlane.hh"
44 #include "G4BSplineSurface.hh"
45 #include "G4ToroidalSurface.hh"
46 #include "G4SphericalSurface.hh"
47 
51 
53  : G4VSolid(name),
54  Box(0), Convex(0), AxisBox(0), PlaneSolid(0), place(0), bbox(0),
55  intersectionDistance(kInfinity), active(1), startInside(0),
56  nb_of_surfaces(0), SurfaceVec(0), RealDist(0.), solidname(name), Id(0),
57  fStatistics(1000000), fCubVolEpsilon(0.001), fAreaAccuracy(-1.),
58  fCubicVolume(0.), fSurfaceArea(0.), fpPolyhedron(0)
59 {
60  static G4bool warn=true;
61  if (warn)
62  {
63  G4cout
64  << "--------------------------------------------------------" << G4endl
65  << "WARNING: BREPS classes are being dismissed. They will |" << G4endl
66  << " be removed starting from next Geant4 major |" << G4endl
67  << " release. Please, consider switching to adopt |" << G4endl
68  << " correspondent CSG or specific primitives. |" << G4endl
69  << "--------------------------------------------------------"
70  << G4endl;
71  warn = false;
72  }
73 }
74 
76  G4Surface** srfVec ,
77  G4int numberOfSrfs )
78  : G4VSolid(name),
79  Box(0), Convex(0), AxisBox(0), PlaneSolid(0), place(0), bbox(0),
80  intersectionDistance(kInfinity), active(1), startInside(0),
81  nb_of_surfaces(numberOfSrfs), SurfaceVec(srfVec), RealDist(0.),
82  solidname(name), Id(0),
83  fStatistics(1000000), fCubVolEpsilon(0.001), fAreaAccuracy(-1.),
84  fCubicVolume(0.), fSurfaceArea(0.), fpPolyhedron(0)
85 {
86  static G4bool warn=true;
87  if (warn)
88  {
89  G4cout
90  << "--------------------------------------------------------" << G4endl
91  << "WARNING: BREPS classes are being dismissed. They will |" << G4endl
92  << " be removed starting from next Geant4 major |" << G4endl
93  << " release. Please, consider switching to adopt |" << G4endl
94  << " correspondent CSG or specific primitives. |" << G4endl
95  << "--------------------------------------------------------"
96  << G4endl;
97  warn = false;
98  }
99  Initialize();
100 }
101 
103  : G4VSolid(a),
104  Box(0), Convex(0), AxisBox(0), PlaneSolid(0), place(0), bbox(0),
105  intersectionDistance(kInfinity), active(1), startInside(0),
106  nb_of_surfaces(0), SurfaceVec(0), RealDist(0.), solidname(""), Id(0),
107  fStatistics(1000000), fCubVolEpsilon(0.001), fAreaAccuracy(-1.),
108  fCubicVolume(0.), fSurfaceArea(0.), fpPolyhedron(0)
109 {
110 }
111 
113 {
114  delete place;
115  delete bbox;
116 
117  for(G4int a=0;a<nb_of_surfaces;++a)
118  { delete SurfaceVec[a]; }
119 
120  delete [] SurfaceVec;
121  delete fpPolyhedron;
122 }
123 
125  : G4VSolid(rhs), Box(rhs.Box), Convex(rhs.Convex), AxisBox(rhs.AxisBox),
126  PlaneSolid(rhs.PlaneSolid), place(0), bbox(0),
127  intersectionDistance(rhs.intersectionDistance), active(rhs.active),
128  startInside(rhs.startInside), nb_of_surfaces(rhs.nb_of_surfaces),
129  intersection_point(rhs.intersection_point), SurfaceVec(0),
130  RealDist(rhs.RealDist), solidname(rhs.solidname), Id(rhs.Id),
131  fStatistics(rhs.fStatistics), fCubVolEpsilon(rhs.fCubVolEpsilon),
132  fAreaAccuracy(rhs.fAreaAccuracy), fCubicVolume(rhs.fCubicVolume),
133  fSurfaceArea(rhs.fSurfaceArea), fpPolyhedron(0)
134 {
135  Initialize();
136 }
137 
139 {
140  // Check assignment to self
141  //
142  if (this == &rhs) { return *this; }
143 
144  // Copy base class data
145  //
146  G4VSolid::operator=(rhs);
147 
148  // Copy data
149  //
150  Box= rhs.Box; Convex= rhs.Convex; AxisBox= rhs.AxisBox;
151  PlaneSolid= rhs.PlaneSolid;
155  RealDist= rhs.RealDist; solidname= rhs.solidname; Id= rhs.Id;
156  fStatistics= rhs.fStatistics; fCubVolEpsilon= rhs.fCubVolEpsilon;
157  fAreaAccuracy= rhs.fAreaAccuracy; fCubicVolume= rhs.fCubicVolume;
158  fSurfaceArea= rhs.fSurfaceArea;
159 
160  delete place; delete bbox; place= 0; bbox= 0;
161  for(G4int a=0;a<nb_of_surfaces;++a)
162  { delete SurfaceVec[a]; }
163  delete [] SurfaceVec; SurfaceVec= 0;
164  delete fpPolyhedron; fpPolyhedron= 0;
165 
166  Initialize();
167 
168  return *this;
169 }
170 
172 {
173  if(active)
174  {
175  // Compute bounding box for solids and surfaces
176  // Convert concave planes to convex
177  //
178  ShortestDistance= kInfinity;
179  if (!SurfaceVec) { return; }
180 
181  IsBox();
183 
184  if(!Box || !AxisBox)
185  IsConvex();
186 
187  CalcBBoxes();
188  }
189 }
190 
192 {
193  return "Closed_Shell";
194 }
195 
197 {
198  return new G4BREPSolid(*this);
199 }
200 
201 void G4BREPSolid::Reset() const
202 {
203  ((G4BREPSolid*)this)->active=1;
204  ((G4BREPSolid*)this)->intersectionDistance=kInfinity;
205  ((G4BREPSolid*)this)->startInside=0;
206 
207  for(register G4int a=0;a<nb_of_surfaces;a++)
208  SurfaceVec[a]->Reset();
209 
210  ShortestDistance = kInfinity;
211 }
212 
214 {
215  if(!PlaneSolid)
216  return; // All faces must be planar
217 
218  Convex=1;
219 
220  // Checks that the normals of the surfaces point outwards.
221  // If not, turns the Normal to point out.
222 
223  // Loop through each face and check the G4Vector3D of the Normal
224  //
225  G4Surface* srf;
226  G4Point3D V;
227 
228  G4int SrfNum = 0;
229  G4double YValue=0;
230  G4Point3D Pt;
231 
232  G4int a, b;
233  for(a=0; a<nb_of_surfaces; a++)
234  {
235  // Find vertex point containing extreme y value
236  //
237  srf = SurfaceVec[a];
238  G4int Points = srf->GetNumberOfPoints();
239 
240  for(b =0; b<Points; b++)
241  {
242  Pt = (G4Point3D)srf->GetPoint(b);
243  if(YValue < Pt.y())
244  {
245  YValue = Pt.y();
246  SrfNum = a; // Save srf number
247  }
248  }
249  }
250 
251  // Move the selected face to the first in the List
252  //
253  srf = SurfaceVec[SrfNum];
254 
255  // Start handling the surfaces in order and compare
256  // the neighbouring ones and turn their normals if they
257  // point inwards
258  //
259  G4Point3D Pt1;
260  G4Point3D Pt2;
261  G4Point3D Pt3;
262  G4Point3D Pt4;
263 
264  G4Vector3D N1;
265  G4Vector3D N2;
266  G4Vector3D N3;
267  G4Vector3D N4;
268 
269  G4int* ConnectedList = new G4int[nb_of_surfaces];
270 
271  for(a=0; a<nb_of_surfaces; a++)
272  ConnectedList[a]=0;
273 
274  G4Surface* ConnectedSrf;
275 
276  for(a=0; a<nb_of_surfaces-1; a++)
277  {
278  if(ConnectedList[a] == 0)
279  break;
280  else
281  ConnectedList[a]=1;
282 
283  srf = SurfaceVec[a];
284  G4int SrfPoints = srf->GetNumberOfPoints();
285  N1 = (srf->Norm())->GetDir();
286 
287  for(b=a+1; b<nb_of_surfaces; b++)
288  {
289  if(ConnectedList[b] == 1)
290  break;
291  else
292  ConnectedList[b]=1;
293 
294  // Get next in List
295  //
296  ConnectedSrf = SurfaceVec[b];
297 
298  // Check if it is connected to srf by looping through the points.
299  //
300  G4int ConnSrfPoints = ConnectedSrf->GetNumberOfPoints();
301 
302  for(G4int c=0;c<SrfPoints;c++)
303  {
304  Pt1 = srf->GetPoint(c);
305 
306  for(G4int d=0;d<ConnSrfPoints;d++)
307  {
308  // Find common points
309  //
310  Pt2 = (ConnectedSrf)->GetPoint(d);
311 
312  if( Pt1 == Pt2 )
313  {
314  // Common point found. Compare normals.
315  //
316  N2 = ((ConnectedSrf)->Norm())->GetDir();
317 
318  // Check cross product.
319  //
320  G4Vector3D CP1 = G4Vector3D( N1.cross(N2) );
321  G4double CrossProd1 = CP1.x()+CP1.y()+CP1.z();
322 
323  // Create the other normals
324  //
325  if(c==0)
326  Pt3 = srf->GetPoint(c+1);
327  else
328  Pt3 = srf->GetPoint(0);
329 
330  N3 = (Pt1-Pt3);
331 
332  if(d==0)
333  Pt4 = (ConnectedSrf)->GetPoint(d+1);
334  else
335  Pt4 = (ConnectedSrf)->GetPoint(0);
336 
337  N4 = (Pt1-Pt4);
338 
339  G4Vector3D CP2 = G4Vector3D( N3.cross(N4) );
340  G4double CrossProd2 = CP2.x()+CP2.y()+CP2.z();
341 
342  G4cout << "\nCroosProd2: " << CrossProd2;
343 
344  if( (CrossProd1 < 0 && CrossProd2 < 0) ||
345  (CrossProd1 > 0 && CrossProd2 > 0) )
346  {
347  // Turn Normal
348  //
349  (ConnectedSrf)->Norm()
350  ->SetDir(-1 * (ConnectedSrf)->Norm()->GetDir());
351 
352  // Take the CrossProd1 again as the other Normal was turned.
353  //
354  CP1 = N1.cross(N2);
355  CrossProd1 = CP1.x()+CP1.y()+CP1.z();
356  }
357  if(CrossProd1 > 0)
358  Convex=0;
359  }
360  }
361  }
362  }
363  }
364 
365  delete []ConnectedList;
366 }
367 
368 G4int G4BREPSolid::IsBox()
369 {
370  // This is done by checking that the solid consists of 6 planes.
371  // Then the type is checked to be planar face by face.
372  // For each G4Plane the Normal is computed. The dot product
373  // of one face Normal and each other face Normal is computed.
374  // One result should be 1 and the rest 0 in order to the solid
375  // to be a box.
376 
377  Box=0;
378  G4Surface* srf1, *srf2;
379  register G4int a;
380 
381  // Compute the Normal for the planes
382  //
383  for(a=0; a < nb_of_surfaces;a++)
384  {
385  srf1 = SurfaceVec[a];
386 
387  if(srf1->MyType()==1)
388  (srf1)->Project(); // Compute the projection
389  else
390  {
391  PlaneSolid=0;
392  return 0;
393  }
394  }
395 
396  // Check that all faces are planar
397  //
398  for(a=0; a < nb_of_surfaces;a++)
399  {
400  srf1 = SurfaceVec[a];
401 
402  if (srf1->MyType()!=1)
403  return 0;
404  }
405 
406  PlaneSolid = 1;
407 
408  // Check that the amount of faces is correct
409  //
410  if(nb_of_surfaces!=6) return 0;
411 
412  G4Point3D Pt;
413  G4int Points;
414  G4int Sides=0;
415  G4int Opposite=0;
416 
417  srf1 = SurfaceVec[0];
418  Points = (srf1)->GetNumberOfPoints();
419 
420  if(Points!=4)
421  return 0;
422 
423  G4Vector3D Normal1 = (srf1->Norm())->GetDir();
424  G4double Result;
425 
426  for(G4int b=1; b < nb_of_surfaces;b++)
427  {
428  srf2 = SurfaceVec[b];
429  G4Vector3D Normal2 = ((srf2)->Norm())->GetDir();
430  Result = std::fabs(Normal1 * Normal2);
431 
432  if((Result != 0) && (Result != 1))
433  return 0;
434  else
435  {
436  if(!(G4int)Result)
437  Sides++;
438  else
439  if(((G4int)Result) == 1)
440  Opposite++;
441  }
442  }
443 
444  if((Opposite != 1) && (Sides != nb_of_surfaces-2))
445  return 0;
446 
447  G4Vector3D x_axis(1,0,0);
448  G4Vector3D y_axis(0,1,0);
449 
450  if(((std::fabs(x_axis * Normal1) == 1) && (std::fabs(y_axis * Normal1) == 0)) ||
451  ((std::fabs(x_axis * Normal1) == 0) && (std::fabs(y_axis * Normal1) == 1)) ||
452  ((std::fabs(x_axis * Normal1) == 0) && (std::fabs(y_axis * Normal1) == 0)))
453  AxisBox=1;
454  else
455  Box=1;
456 
457  return 1;
458 }
459 
461 {
462  if (!PlaneSolid) { return 0; } // All faces must be planar
463 
464  // This is not robust. There can be concave solids
465  // where the concavity comes for example from three triangles.
466 
467  // Additional checking 20.8. For each face the connecting faces are
468  // found and the cross product computed between the face and each
469  // connecting face. If the result changes value at any point the
470  // solid is concave.
471 
472  G4Surface* Srf=0;
473  G4Surface* ConnectedSrf=0;
474  G4int Result;
475  Convex = 1;
476 
477  G4int a, b, c, d;
478  for(a=0;a<nb_of_surfaces;a++)
479  {
480  Srf = SurfaceVec[a];
481 
482  // Primary test. Test wether any one of the faces
483  // is concave -> solid is concave. This is not enough to
484  // distinguish all the cases of concavity.
485  //
486  Result = Srf->IsConvex();
487 
488  if (Result != -1)
489  {
490  Convex = 0;
491  return 0;
492  }
493  }
494 
495  Srf = SurfaceVec[0];
496 
497  G4int ConnectingPoints=0;
498  G4Point3D Pt1, Pt2;
499  G4Vector3D N1, N2;
500 
501  // L. Broglia
502  // The number of connecting points can be
503  // (nb_of_surfaces-1) * nb_of_surfaces (loop a & loop b)
504 
505  // G4int* ConnectedList = new G4int[nb_of_surfaces];
506  const G4int maxCNum = (nb_of_surfaces-1)*nb_of_surfaces;
507  G4int* ConnectedList = new G4int[maxCNum];
508 
509  for(a=0; a<maxCNum; a++)
510  {
511  ConnectedList[a]=0;
512  }
513 
514  G4int Connections=0;
515 
516  for(a=0; a<nb_of_surfaces-1; a++)
517  {
518  Srf = SurfaceVec[a];
519  G4int SrfPoints = Srf->GetNumberOfPoints();
520  Result=0;
521 
522  for(b=0; b<nb_of_surfaces; b++)
523  {
524  if (b==a) { continue; }
525 
526  // Get next in List
527  //
528  ConnectedSrf = SurfaceVec[b];
529 
530  // Check if it is connected to Srf by looping through the points.
531  //
532  G4int ConnSrfPoints = ConnectedSrf->GetNumberOfPoints();
533 
534  for(c=0; c<SrfPoints; c++)
535  {
536  const G4Point3D& Pts1 =Srf->GetPoint(c);
537 
538  for(d=0; d<ConnSrfPoints; d++)
539  {
540  // Find common points
541  //
542  const G4Point3D& Pts2 = ConnectedSrf->GetPoint(d);
543  if (Pts1 == Pts2) { ConnectingPoints++; }
544  }
545  if (ConnectingPoints > 0) { break; }
546  }
547 
548  if (ConnectingPoints > 0)
549  {
550  Connections++;
551  ConnectedList[Connections]=b;
552  }
553  ConnectingPoints=0;
554  }
555  }
556 
557  // If connected, check for concavity.
558  // Get surfaces from ConnectedList and compare their normals
559  //
560  for(c=0; c<Connections; c++)
561  {
562  G4int Left=0;
563  G4int Right =0;
564  G4int tmp = ConnectedList[c];
565 
566  Srf = SurfaceVec[tmp];
567  ConnectedSrf = SurfaceVec[tmp+1];
568 
569  // Get normals
570  //
571  N1 = Srf->Norm()->GetDir();
572  N2 = ConnectedSrf->Norm()->GetDir();
573 
574  // Check cross product
575  //
576  G4Vector3D CP = G4Vector3D( N1.cross(N2) );
577  G4double CrossProd = CP.x()+CP.y()+CP.z();
578  if (CrossProd > 0) { Left++; }
579  if (CrossProd < 0) { Right++; }
580  if (Left&&Right)
581  {
582  Convex = 0;
583  delete [] ConnectedList;
584  return 0;
585  }
586  Connections=0;
587  }
588 
589  Convex=1;
590 
591  delete [] ConnectedList;
592 
593  return 1;
594 }
595 
597  const G4VoxelLimits& pVoxelLimit,
598  const G4AffineTransform& pTransform,
599  G4double& pMin, G4double& pMax) const
600 {
601  G4Point3D Min = bbox->GetBoxMin();
602  G4Point3D Max = bbox->GetBoxMax();
603 
604  if (!pTransform.IsRotated())
605  {
606  // Special case handling for unrotated boxes
607  // Compute x/y/z mins and maxs respecting limits, with early returns
608  // if outside limits. Then switch() on pAxis
609  //
610  G4double xoffset,xMin,xMax;
611  G4double yoffset,yMin,yMax;
612  G4double zoffset,zMin,zMax;
613 
614  xoffset=pTransform.NetTranslation().x();
615  xMin=xoffset+Min.x();
616  xMax=xoffset+Max.x();
617  if (pVoxelLimit.IsXLimited())
618  {
619  if (xMin>pVoxelLimit.GetMaxXExtent()
620  ||xMax<pVoxelLimit.GetMinXExtent())
621  {
622  return false;
623  }
624  else
625  {
626  if (xMin<pVoxelLimit.GetMinXExtent())
627  {
628  xMin=pVoxelLimit.GetMinXExtent();
629  }
630  if (xMax>pVoxelLimit.GetMaxXExtent())
631  {
632  xMax=pVoxelLimit.GetMaxXExtent();
633  }
634  }
635  }
636 
637  yoffset=pTransform.NetTranslation().y();
638  yMin=yoffset+Min.y();
639  yMax=yoffset+Max.y();
640  if (pVoxelLimit.IsYLimited())
641  {
642  if (yMin>pVoxelLimit.GetMaxYExtent()
643  ||yMax<pVoxelLimit.GetMinYExtent())
644  {
645  return false;
646  }
647  else
648  {
649  if (yMin<pVoxelLimit.GetMinYExtent())
650  {
651  yMin=pVoxelLimit.GetMinYExtent();
652  }
653  if (yMax>pVoxelLimit.GetMaxYExtent())
654  {
655  yMax=pVoxelLimit.GetMaxYExtent();
656  }
657  }
658  }
659 
660  zoffset=pTransform.NetTranslation().z();
661  zMin=zoffset+Min.z();
662  zMax=zoffset+Max.z();
663  if (pVoxelLimit.IsZLimited())
664  {
665  if (zMin>pVoxelLimit.GetMaxZExtent()
666  ||zMax<pVoxelLimit.GetMinZExtent())
667  {
668  return false;
669  }
670  else
671  {
672  if (zMin<pVoxelLimit.GetMinZExtent())
673  {
674  zMin=pVoxelLimit.GetMinZExtent();
675  }
676  if (zMax>pVoxelLimit.GetMaxZExtent())
677  {
678  zMax=pVoxelLimit.GetMaxZExtent();
679  }
680  }
681  }
682 
683  switch (pAxis)
684  {
685  case kXAxis:
686  pMin=xMin;
687  pMax=xMax;
688  break;
689  case kYAxis:
690  pMin=yMin;
691  pMax=yMax;
692  break;
693  case kZAxis:
694  pMin=zMin;
695  pMax=zMax;
696  break;
697  default:
698  break;
699  }
700  pMin-=kCarTolerance;
701  pMax+=kCarTolerance;
702 
703  return true;
704  }
705  else
706  {
707  // General rotated case - create and clip mesh to boundaries
708 
709  G4bool existsAfterClip=false;
710  G4ThreeVectorList *vertices;
711 
712  pMin=+kInfinity;
713  pMax=-kInfinity;
714 
715  // Calculate rotated vertex coordinates
716  //
717  vertices=CreateRotatedVertices(pTransform);
718  ClipCrossSection(vertices,0,pVoxelLimit,pAxis,pMin,pMax);
719  ClipCrossSection(vertices,4,pVoxelLimit,pAxis,pMin,pMax);
720  ClipBetweenSections(vertices,0,pVoxelLimit,pAxis,pMin,pMax);
721 
722  if ( (pMin!=kInfinity) || (pMax!=-kInfinity) )
723  {
724  existsAfterClip=true;
725 
726  // Add 2*tolerance to avoid precision troubles
727  //
728  pMin-=kCarTolerance;
729  pMax+=kCarTolerance;
730  }
731  else
732  {
733  // Check for case where completely enveloping clipping volume.
734  // If point inside then we are confident that the solid completely
735  // envelopes the clipping volume. Hence set min/max extents according
736  // to clipping volume extents along the specified axis.
737  //
738  G4ThreeVector clipCentre(
739  (pVoxelLimit.GetMinXExtent()+pVoxelLimit.GetMaxXExtent())*0.5,
740  (pVoxelLimit.GetMinYExtent()+pVoxelLimit.GetMaxYExtent())*0.5,
741  (pVoxelLimit.GetMinZExtent()+pVoxelLimit.GetMaxZExtent())*0.5);
742 
743  if (Inside(pTransform.Inverse().TransformPoint(clipCentre))!=kOutside)
744  {
745  existsAfterClip=true;
746  pMin=pVoxelLimit.GetMinExtent(pAxis);
747  pMax=pVoxelLimit.GetMaxExtent(pAxis);
748  }
749  }
750  delete vertices;
751  return existsAfterClip;
752  }
753 }
754 
757 {
758  G4Point3D Min = bbox->GetBoxMin();
759  G4Point3D Max = bbox->GetBoxMax();
760 
761  G4ThreeVectorList *vertices;
762  vertices=new G4ThreeVectorList();
763 
764  if (vertices)
765  {
766  vertices->reserve(8);
767  G4ThreeVector vertex0(Min.x(),Min.y(),Min.z());
768  G4ThreeVector vertex1(Max.x(),Min.y(),Min.z());
769  G4ThreeVector vertex2(Max.x(),Max.y(),Min.z());
770  G4ThreeVector vertex3(Min.x(),Max.y(),Min.z());
771  G4ThreeVector vertex4(Min.x(),Min.y(),Max.z());
772  G4ThreeVector vertex5(Max.x(),Min.y(),Max.z());
773  G4ThreeVector vertex6(Max.x(),Max.y(),Max.z());
774  G4ThreeVector vertex7(Min.x(),Max.y(),Max.z());
775 
776  vertices->push_back(pTransform.TransformPoint(vertex0));
777  vertices->push_back(pTransform.TransformPoint(vertex1));
778  vertices->push_back(pTransform.TransformPoint(vertex2));
779  vertices->push_back(pTransform.TransformPoint(vertex3));
780  vertices->push_back(pTransform.TransformPoint(vertex4));
781  vertices->push_back(pTransform.TransformPoint(vertex5));
782  vertices->push_back(pTransform.TransformPoint(vertex6));
783  vertices->push_back(pTransform.TransformPoint(vertex7));
784  }
785  else
786  {
787  G4Exception("G4BREPSolid::CreateRotatedVertices()", "GeomSolids0003",
788  FatalException, "Out of memory - Cannot allocate vertices!");
789  }
790  return vertices;
791 }
792 
793 EInside G4BREPSolid::Inside(register const G4ThreeVector& Pt) const
794 {
795  // This function finds if the point Pt is inside,
796  // outside or on the surface of the solid
797 
798  const G4double sqrHalfTolerance = kCarTolerance*kCarTolerance*0.25;
799 
800  G4Vector3D v(1, 0, 0.01);
801  G4Vector3D Pttmp(Pt);
802  G4Vector3D Vtmp(v);
803  G4Ray r(Pttmp, Vtmp);
804 
805  // Check if point is inside the PCone bounding box
806  //
807  if( !GetBBox()->Inside(Pttmp) )
808  return kOutside;
809 
810  // Set the surfaces to active again
811  //
812  Reset();
813 
814  // Test if the bounding box of each surface is intersected
815  // by the ray. If not, the surface become deactive.
816  //
818 
819  G4int hits=0, samehit=0;
820 
821  for(G4int a=0; a < nb_of_surfaces; a++)
822  {
823  if(SurfaceVec[a]->IsActive())
824  {
825  // Count the number of intersections. If this number is odd,
826  // the start of the ray is inside the volume bounded by the surfaces,
827  // so increment the number of intersection by 1 if the point is not
828  // on the surface and if this intersection was not found before.
829  //
830  if( (SurfaceVec[a]->Intersect(r)) & 1 )
831  {
832  // Test if the point is on the surface
833  //
834  if(SurfaceVec[a]->GetDistance() < sqrHalfTolerance)
835  return kSurface;
836 
837  // Test if this intersection was found before
838  //
839  for(G4int i=0; i<a; i++)
840  if(SurfaceVec[a]->GetDistance() == SurfaceVec[i]->GetDistance())
841  {
842  samehit++;
843  break;
844  }
845 
846  // Count the number of surfaces intersected by the ray
847  //
848  if(!samehit)
849  hits++;
850  }
851  }
852  }
853 
854  // If the number of surfaces intersected is odd,
855  // the point is inside the solid
856  //
857  if(hits&1)
858  return kInside;
859  else
860  return kOutside;
861 }
862 
864 {
865  // This function calculates the normal of the surface at a point on the
866  // surface. If the point is not on the surface the result is undefined.
867  // Note : the sense of the normal depends on the sense of the surface.
868 
869  const G4double sqrHalfTolerance = kCarTolerance*kCarTolerance*0.25;
870  G4int iplane;
871 
872  // Find on which surface the point is
873  //
874  for(iplane = 0; iplane < nb_of_surfaces; iplane++)
875  {
876  if(SurfaceVec[iplane]->HowNear(Pt) < sqrHalfTolerance)
877  // the point is on this surface
878  break;
879  }
880 
881  // Calculate the normal at this point
882  //
884 
885  return norm.unit();
886 }
887 
889 {
890  // Calculates the shortest distance ("safety") from a point
891  // outside the solid to any boundary of this solid.
892  // Return 0 if the point is already inside.
893 
894  G4double *dists = new G4double[nb_of_surfaces];
895  G4int a;
896 
897  // Set the surfaces to active again
898  //
899  Reset();
900 
901  // Compute the shortest distance of the point to each surface.
902  // Be careful : it's a signed value
903  //
904  for(a=0; a< nb_of_surfaces; a++)
905  dists[a] = SurfaceVec[a]->HowNear(Pt);
906 
907  G4double Dist = kInfinity;
908 
909  // If dists[] is positive, the point is outside, so take the shortest of
910  // the shortest positive distances dists[] can be equal to 0 : point on
911  // a surface.
912  // ( Problem with the G4FPlane : there is no inside and no outside...
913  // So, to test if the point is inside to return 0, utilize the Inside()
914  // function. But I don't know if it is really needed because dToIn is
915  // called only if the point is outside )
916  //
917  for(a = 0; a < nb_of_surfaces; a++)
918  if( std::fabs(Dist) > std::fabs(dists[a]) )
919  //if( dists[a] >= 0)
920  Dist = dists[a];
921 
922  delete[] dists;
923 
924  if(Dist == kInfinity)
925  return 0; // the point is inside the solid or on a surface
926  else
927  return std::fabs(Dist);
928 }
929 
931  register const G4ThreeVector& V ) const
932 {
933  // Calculates the distance from a point outside the solid
934  // to the solid's boundary along a specified direction vector.
935  //
936  // Note : Intersections with boundaries less than the tolerance must be
937  // ignored if the direction is away from the boundary.
938 
939  G4int a;
940 
941  // Set the surfaces to active again
942  //
943  Reset();
944 
945  const G4double sqrHalfTolerance = kCarTolerance*kCarTolerance*0.25;
946  G4Vector3D Pttmp(Pt);
947  G4Vector3D Vtmp(V);
948  G4Ray r(Pttmp, Vtmp);
949 
950  // Test if the bounding box of each surface is intersected
951  // by the ray. If not, the surface become deactive.
952  //
954 
955  ShortestDistance = kInfinity;
956 
957  for(a=0; a< nb_of_surfaces; a++)
958  {
959  if( SurfaceVec[a]->IsActive() )
960  {
961  // Test if the ray intersects the surface
962  //
963  if( SurfaceVec[a]->Intersect(r) )
964  {
965  G4double surfDistance = SurfaceVec[a]->GetDistance();
966 
967  // If more than 1 surface is intersected, take the nearest one
968  //
969  if( surfDistance < ShortestDistance )
970  {
971  if( surfDistance > sqrHalfTolerance )
972  {
973  ShortestDistance = surfDistance;
974  }
975  else
976  {
977  // The point is within the boundary. It is ignored it if
978  // the direction is away from the boundary
979  //
980  G4Vector3D Norm = SurfaceVec[a]->SurfaceNormal(Pttmp);
981 
982  if( (Norm * Vtmp) < 0 )
983  {
984  ShortestDistance = surfDistance;
985  }
986  }
987  }
988  }
989  }
990  }
991 
992  // Be careful !
993  // SurfaceVec->Distance is in fact the squared distance
994  //
995  if(ShortestDistance != kInfinity)
996  return std::sqrt(ShortestDistance);
997  else
998  return kInfinity; // No intersection
999 }
1000 
1002  register const G4ThreeVector& D,
1003  const G4bool,
1004  G4bool *validNorm,
1005  G4ThreeVector* ) const
1006 {
1007  // Calculates the distance from a point inside the solid to the solid's
1008  // boundary along a specified direction vector.
1009  // Returns 0 if the point is already outside.
1010  //
1011  // Note : If the shortest distance to a boundary is less than the tolerance,
1012  // it is ignored. This allows for a point within a tolerant boundary
1013  // to leave immediately.
1014 
1015  // Set the surfaces to active again
1016  //
1017  Reset();
1018 
1019  const G4double sqrHalfTolerance = kCarTolerance*kCarTolerance*0.25;
1020  G4Vector3D Ptv = P;
1021  G4int a;
1022 
1023  if(validNorm)
1024  *validNorm=false;
1025 
1026  G4Vector3D Pttmp(Ptv);
1027  G4Vector3D Vtmp(D);
1028 
1029  G4Ray r(Pttmp, Vtmp);
1030 
1031  // Test if the bounding box of each surface is intersected
1032  // by the ray. If not, the surface become deactive.
1033  //
1034  TestSurfaceBBoxes(r);
1035 
1036  ShortestDistance = kInfinity;
1037 
1038  for(a=0; a< nb_of_surfaces; a++)
1039  {
1040  if(SurfaceVec[a]->IsActive())
1041  {
1042  // Test if the ray intersect the surface
1043  //
1044  if( (SurfaceVec[a]->Intersect(r)) )
1045  {
1046  // If more than 1 surface is intersected, take the nearest one
1047  //
1048  G4double surfDistance = SurfaceVec[a]->GetDistance();
1049  if( surfDistance < ShortestDistance )
1050  {
1051  if( surfDistance > sqrHalfTolerance )
1052  {
1053  ShortestDistance = surfDistance;
1054  }
1055  else
1056  {
1057  // The point is within the boundary: ignore it
1058  }
1059  }
1060  }
1061  }
1062  }
1063 
1064  // Be careful !
1065  // SurfaceVec->Distance is in fact the squared distance
1066  //
1067  if(ShortestDistance != kInfinity)
1068  return std::sqrt(ShortestDistance);
1069  else
1070  return 0.0; // No intersection is found, the point is outside
1071 }
1072 
1074 {
1075  // Calculates the shortest distance ("safety") from a point
1076  // inside the solid to any boundary of this solid.
1077  // Returns 0 if the point is already outside.
1078 
1079  G4double *dists = new G4double[nb_of_surfaces];
1080  G4int a;
1081 
1082  // Set the surfaces to active again
1083  //
1084  Reset();
1085 
1086  // Compute the shortest distance of the point to each surfaces
1087  // Be careful : it's a signed value
1088  //
1089  for(a=0; a< nb_of_surfaces; a++)
1090  dists[a] = SurfaceVec[a]->HowNear(Pt);
1091 
1092  G4double Dist = kInfinity;
1093 
1094  // If dists[] is negative, the point is inside so take the shortest of the
1095  // shortest negative distances dists[] can be equal to 0 : point on a
1096  // surface
1097  // ( Problem with the G4FPlane : there is no inside and no outside...
1098  // So, to test if the point is outside to return 0, utilize the Inside()
1099  // function. But I don`t know if it is really needed because dToOut is
1100  // called only if the point is inside )
1101  //
1102  for(a = 0; a < nb_of_surfaces; a++)
1103  if( std::fabs(Dist) > std::fabs(dists[a]) )
1104  //if( dists[a] <= 0)
1105  Dist = dists[a];
1106 
1107  delete[] dists;
1108 
1109  if(Dist == kInfinity)
1110  return 0; // The point is ouside the solid or on a surface
1111  else
1112  return std::fabs(Dist);
1113 }
1114 
1116 {
1117  scene.AddSolid (*this);
1118 }
1119 
1121 {
1122  // Approximate implementation, just a box ...
1123 
1124  G4Point3D Min = bbox->GetBoxMin();
1125  G4Point3D Max = bbox->GetBoxMax();
1126 
1127  return new G4PolyhedronBox (Max.x(), Max.y(), Max.z());
1128 }
1129 
1131 {
1132  // Approximate implementation, just a box ...
1133 
1134  G4Point3D Min = bbox->GetBoxMin();
1135  G4Point3D Max = bbox->GetBoxMax();
1136 
1137  return new G4NURBSbox (Max.x(), Max.y(), Max.z());
1138 }
1139 
1141 {
1142  // First initialization. Calculates the bounding boxes
1143  // for the surfaces and for the solid.
1144 
1145  G4Surface* srf;
1146  G4Point3D min, max;
1147 
1148  if(active)
1149  {
1150  min = PINFINITY;
1151  max = -PINFINITY;
1152 
1153  for(G4int a = 0;a < nb_of_surfaces;a++)
1154  {
1155  // Get first in List
1156  //
1157  srf = SurfaceVec[a];
1158 /*
1159  G4int convex=1;
1160  G4int concavepoint=-1;
1161 
1162  if (srf->MyType() == 1)
1163  {
1164  concavepoint = srf->IsConvex();
1165  convex = srf->GetConvex();
1166  }
1167 
1168  // Make bbox for face
1169  //
1170  if(convex && Concavepoint==-1)
1171 */
1172  {
1173  srf->CalcBBox();
1174  G4Point3D box_min = srf->GetBBox()->GetBoxMin();
1175  G4Point3D box_max = srf->GetBBox()->GetBoxMax();
1176 
1177  // Find max and min of face bboxes to make solids bbox.
1178 
1179  // replace by Extend
1180  // max < box_max
1181  //
1182  if(max.x() < box_max.x()) max.setX(box_max.x());
1183  if(max.y() < box_max.y()) max.setY(box_max.y());
1184  if(max.z() < box_max.z()) max.setZ(box_max.z());
1185 
1186  // min > box_min
1187  //
1188  if(min.x() > box_min.x()) min.setX(box_min.x());
1189  if(min.y() > box_min.y()) min.setY(box_min.y());
1190  if(min.z() > box_min.z()) min.setZ(box_min.z());
1191  }
1192  }
1193  bbox = new G4BoundingBox3D(min, max);
1194  return;
1195  }
1196  G4Exception("G4BREPSolid::CalcBBoxes()", "GeomSolids1002",
1197  JustWarning, "No bbox calculated for solid.");
1198 }
1199 
1200 void G4BREPSolid::RemoveHiddenFaces(register const G4Ray& rayref,
1201  G4int In ) const
1202 {
1203  // Deactivates the planar faces that are on the "back" side of a solid.
1204  // B-splines are not handled by this function. Also cases where the ray
1205  // starting point is Inside the bbox of the solid are ignored as we don't
1206  // know if the starting point is Inside the actual solid except for
1207  // axis-oriented box-like solids.
1208 
1209  register G4Surface* srf;
1210  register const G4Vector3D& RayDir = rayref.GetDir();
1211  register G4double Result;
1212  G4int a;
1213 
1214  // In all other cases the ray starting point is outside the solid
1215  //
1216  if(!In)
1217  for(a=0; a<nb_of_surfaces; a++)
1218  {
1219  // Deactivates the solids faces that are hidden
1220  //
1221  srf = SurfaceVec[a];
1222  if(srf->MyType()==1)
1223  {
1224  const G4Vector3D& Normal = (srf->Norm())->GetDir();
1225  Result = (RayDir * Normal);
1226 
1227  if( Result >= 0 )
1228  srf->Deactivate();
1229  }
1230  }
1231  else
1232  for(a=0; a<nb_of_surfaces; a++)
1233  {
1234  // Deactivates the AxisBox type solids faces whose normals
1235  // point in the G4Vector3D opposite to the rays G4Vector3D
1236  // i.e. are behind the ray starting point as in this case the
1237  // ray starts from Inside the solid.
1238  //
1239  srf = SurfaceVec[a];
1240  if(srf->MyType()==1)
1241  {
1242  const G4Vector3D& Normal = (srf->Norm())->GetDir();
1243  Result = (RayDir * Normal);
1244 
1245  if( Result < 0 )
1246  srf->Deactivate();
1247  }
1248  }
1249 }
1250 
1251 void G4BREPSolid::TestSurfaceBBoxes(register const G4Ray& rayref) const
1252 {
1253  register G4Surface* srf;
1254  G4int active_srfs = nb_of_surfaces;
1255 
1256  // Do the bbox tests to all surfaces in List
1257  // for planar faces the intersection is instead evaluated.
1258  //
1259  G4int intersection=0;
1260 
1261  for(G4int a=0;a<nb_of_surfaces;a++)
1262  {
1263  // Get first in List
1264  //
1265  srf = SurfaceVec[a];
1266 
1267  if(srf->IsActive())
1268  {
1269  // Get type
1270  //
1271  if(srf->MyType() != 1) // 1 == planar face
1272  {
1273  if(srf->GetBBox()->Test(rayref))
1274  srf->SetDistance(bbox->GetDistance());
1275  else
1276  {
1277  // Test failed. Flag as inactive.
1278  //
1279  srf->Deactivate();
1280  active_srfs--;
1281  }
1282  }
1283  else
1284  {
1285  // Type was convex planar face
1286  intersection = srf->Intersect(rayref);
1287 
1288  if(!intersection)
1289  active_srfs--;
1290  }
1291  }
1292  else
1293  active_srfs--;
1294  }
1295 
1296  if(!active_srfs) Active(0);
1297 }
1298 
1299 
1300 G4int G4BREPSolid::Intersect(register const G4Ray& rayref) const
1301 {
1302  // Gets the roughly calculated closest intersection point for
1303  // a b_spline & accurate point for others.
1304 
1305  const G4double sqrHalfTolerance = kCarTolerance*kCarTolerance*0.25;
1306 
1307  register G4Surface* srf;
1308  G4double HitDistance = -1;
1309  const G4Point3D& RayStart = rayref.GetStart();
1310  const G4Point3D& RayDir = rayref.GetDir();
1311 
1312  G4int result=1;
1313 
1314  // Sort List of active surfaces according to
1315  // bbox distances to ray starting point
1316  //
1317  QuickSort(SurfaceVec, 0, nb_of_surfaces-1);
1318  G4int Number=0;
1319 
1320  // Start handling active surfaces in order
1321  //
1322  for(register G4int a=0;a<nb_of_surfaces;a++)
1323  {
1324  srf = SurfaceVec[a];
1325 
1326  if(srf->IsActive())
1327  {
1328  result = srf->Intersect(rayref);
1329  if(result)
1330  {
1331  // Get the evaluated point on the surface
1332  //
1333  const G4Point3D& closest_point = srf->GetClosestHit();
1334 
1335  // Test for DistanceToIn(pt, vec)
1336  // if d = 0 and vec.norm > 0, do not see the surface
1337  //
1338  if( !( (srf->GetDistance() < sqrHalfTolerance) ||
1339  (RayDir.dot(srf->SurfaceNormal(closest_point)) > 0) ) )
1340  {
1341 
1342  if(srf->MyType()==1)
1343  HitDistance = srf->GetDistance();
1344  else
1345  {
1346  // Check if the evaluated point is in front of the
1347  // bbox of the next surface.
1348  //
1349  HitDistance = RayStart.distance2(closest_point);
1350  }
1351  }
1352  }
1353  else // No hit
1354  {
1355  srf->Deactivate();
1356  }
1357  }
1358  Number++;
1359  }
1360 
1361  if(HitDistance < 0)
1362  return 0;
1363 
1364  QuickSort(SurfaceVec, 0, nb_of_surfaces-1);
1365 
1366  if(!(SurfaceVec[0]->IsActive()))
1367  return 0;
1368 
1369  ((G4BREPSolid*)this)->intersection_point = SurfaceVec[0]->GetClosestHit();
1370  bbox->SetDistance(HitDistance);
1371 
1372  return 1;
1373 }
1374 
1375 G4int G4BREPSolid::FinalEvaluation(register const G4Ray& rayref,
1376  G4int ToIn ) const
1377 {
1378  const G4double sqrHalfTolerance = kCarTolerance*kCarTolerance*0.25;
1379  register G4Surface* srf;
1380  G4double Dist=0;
1381 
1382  ((G4BREPSolid*)this)->intersectionDistance = kInfinity;
1383 
1384  for(register G4int a=0;a<nb_of_surfaces;a++)
1385  {
1386  srf = SurfaceVec[a];
1387 
1388  if(srf->IsActive())
1389  {
1390  const G4Point3D& srf_intersection = srf->Evaluation(rayref);
1391 
1392  // Compute hit point distance from ray starting point
1393  //
1394  if(srf->MyType() != 1)
1395  {
1396  G4Point3D start = rayref.GetStart();
1397  Dist = srf_intersection.distance2(start);
1398  }
1399  else
1400  Dist = srf->GetDistance();
1401 
1402  // Skip point wich are on the surface i.e. within tolerance of the
1403  // surface. Special handling for DistanceToIn & reflections
1404  //
1405  if(Dist < sqrHalfTolerance)
1406  {
1407  if(ToIn)
1408  {
1409  const G4Vector3D& Dir = rayref.GetDir();
1410  const G4Point3D& Hit = srf->GetClosestHit();
1411  const G4Vector3D& Norm = srf->SurfaceNormal(Hit);
1412 
1413  if(( Dir * Norm ) >= 0)
1414  {
1415  Dist = kInfinity;
1416  srf->Deactivate();
1417  }
1418 
1419  // else continue with the distance, even though < tolerance
1420  }
1421  else
1422  {
1423  Dist = kInfinity;
1424  srf->Deactivate();
1425  }
1426  }
1427 
1428  // If more than one surfaces are evaluated till the
1429  // final stage, only the closest point is taken
1430  //
1431  if(Dist < intersectionDistance)
1432  {
1433  // Check that Hit is in the direction of the ray
1434  // from the starting point
1435  //
1436  const G4Point3D& Pt = rayref.GetStart();
1437  const G4Vector3D& Dir = rayref.GetDir();
1438 
1439  G4Point3D TestPoint = (0.00001*Dir) + Pt;
1440  G4double TestDistance = srf_intersection.distance2(TestPoint);
1441 
1442  if(TestDistance > Dist)
1443  {
1444  // Hit behind ray starting point, no intersection
1445  //
1446  Dist = kInfinity;
1447  srf->Deactivate();
1448  }
1449  else
1450  {
1451  ((G4BREPSolid*)this)->intersectionDistance = Dist;
1452  ((G4BREPSolid*)this)->intersection_point = srf_intersection;
1453  }
1454 
1455  // Check that the intersection is closer than the
1456  // next surfaces approximated point
1457  //
1458  if(srf->IsActive())
1459  {
1460  if(a+1<nb_of_surfaces)
1461  {
1462  const G4Point3D& Hit = srf->GetClosestHit();
1463  const G4Vector3D& Norm = srf->SurfaceNormal(Hit);
1464 
1465  // L. Broglia
1466  //if(( Dir * Norm ) >= 0)
1467  if(( Dir * Norm ) < 0)
1468  {
1469  Dist = kInfinity;
1470  srf->Deactivate();
1471  }
1472 
1473  // else continue with the distance, even though < tolerance
1474 
1475  ShortestDistance = Dist;
1476  }
1477  else
1478  {
1479  ShortestDistance = Dist;
1480  return 1;
1481  }
1482  }
1483  }
1484  }
1485  else // if srf NOT active
1486  {
1487  /* if(intersectionDistance < kInfinity)
1488  return 1;
1489  return 0;*/
1490  }
1491  }
1492  if(intersectionDistance < kInfinity)
1493  return 1;
1494 
1495  return 0;
1496 }
1497 
1499 {
1500  G4Point3D scope;
1501  G4Point3D Max = bbox->GetBoxMax();
1502  G4Point3D Min = bbox->GetBoxMin();
1503 
1504  scope.setX(std::fabs(Max.x()) - std::fabs(Min.x()));
1505  scope.setY(std::fabs(Max.y()) - std::fabs(Min.y()));
1506  scope.setZ(std::fabs(Max.z()) - std::fabs(Min.z()));
1507 
1508  return scope;
1509 }
1510 
1511 std::ostream& G4BREPSolid::StreamInfo(std::ostream& os) const
1512 {
1513  os << "-----------------------------------------------------------\n"
1514  << " *** Dump for solid - " << GetName() << " ***\n"
1515  << " ===================================================\n"
1516  << " Solid type: " << GetEntityType() << "\n"
1517  << " Parameters: \n"
1518  << " Number of solids: " << NumberOfSolids << "\n"
1519  << "-----------------------------------------------------------\n";
1520 
1521  return os;
1522 }
1523 
1525 {
1526  if (!fpPolyhedron ||
1527  fpPolyhedron->GetNumberOfRotationStepsAtTimeOfCreation() !=
1528  fpPolyhedron->GetNumberOfRotationSteps())
1529  {
1530  delete fpPolyhedron;
1531  fpPolyhedron = CreatePolyhedron();
1532  }
1533  return fpPolyhedron;
1534 }