Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4VCSGfaceted.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 // the GEANT4 collaboration.
27 //
28 // By copying, distributing or modifying the Program (or any work
29 // based on the Program) you indicate your acceptance of this statement,
30 // and all its terms.
31 //
32 // $Id: G4VCSGfaceted.cc 67011 2013-01-29 16:17:41Z gcosmo $
33 //
34 //
35 // --------------------------------------------------------------------
36 // GEANT 4 class source file
37 //
38 //
39 // G4VCSGfaceted.cc
40 //
41 // Implementation of the virtual class of a CSG type shape that is built
42 // entirely out of G4VCSGface faces.
43 //
44 // --------------------------------------------------------------------
45 
46 #include "G4VCSGfaceted.hh"
47 #include "G4VCSGface.hh"
48 #include "G4SolidExtentList.hh"
49 
50 #include "G4VoxelLimits.hh"
51 #include "G4AffineTransform.hh"
52 
53 #include "Randomize.hh"
54 
55 #include "G4Polyhedron.hh"
56 #include "G4VGraphicsScene.hh"
57 #include "G4NURBS.hh"
58 #include "G4NURBSbox.hh"
59 #include "G4VisExtent.hh"
60 
61 //
62 // Constructor
63 //
65  : G4VSolid(name),
66  numFace(0), faces(0), fCubicVolume(0.), fSurfaceArea(0.), fpPolyhedron(0),
67  fStatistics(1000000), fCubVolEpsilon(0.001), fAreaAccuracy(-1.)
68 {
69 }
70 
71 
72 //
73 // Fake default constructor - sets only member data and allocates memory
74 // for usage restricted to object persistency.
75 //
77  : G4VSolid(a),
78  numFace(0), faces(0), fCubicVolume(0.), fSurfaceArea(0.), fpPolyhedron(0),
79  fStatistics(1000000), fCubVolEpsilon(0.001), fAreaAccuracy(-1.)
80 {
81 }
82 
83 //
84 // Destructor
85 //
87 {
88  DeleteStuff();
89  delete fpPolyhedron;
90 }
91 
92 
93 //
94 // Copy constructor
95 //
97  : G4VSolid( source )
98 {
99  fStatistics = source.fStatistics;
100  fCubVolEpsilon = source.fCubVolEpsilon;
101  fAreaAccuracy = source.fAreaAccuracy;
102 
103  CopyStuff( source );
104 }
105 
106 
107 //
108 // Assignment operator
109 //
111 {
112  if (&source == this) { return *this; }
113 
114  // Copy base class data
115  //
116  G4VSolid::operator=(source);
117 
118  // Copy data
119  //
120  fStatistics = source.fStatistics;
121  fCubVolEpsilon = source.fCubVolEpsilon;
122  fAreaAccuracy = source.fAreaAccuracy;
123 
124  DeleteStuff();
125  CopyStuff( source );
126 
127  return *this;
128 }
129 
130 
131 //
132 // CopyStuff (protected)
133 //
134 // Copy the contents of source
135 //
137 {
138  numFace = source.numFace;
139  if (numFace == 0) { return; } // odd, but permissable?
140 
141  faces = new G4VCSGface*[numFace];
142 
143  G4VCSGface **face = faces,
144  **sourceFace = source.faces;
145  do
146  {
147  *face = (*sourceFace)->Clone();
148  } while( ++sourceFace, ++face < faces+numFace );
149  fCubicVolume = source.fCubicVolume;
150  fSurfaceArea = source.fSurfaceArea;
151  fpPolyhedron = 0;
152 }
153 
154 
155 //
156 // DeleteStuff (protected)
157 //
158 // Delete all allocated objects
159 //
161 {
162  if (numFace)
163  {
164  G4VCSGface **face = faces;
165  do
166  {
167  delete *face;
168  } while( ++face < faces + numFace );
169 
170  delete [] faces;
171  }
172 }
173 
174 
175 //
176 // CalculateExtent
177 //
179  const G4VoxelLimits &voxelLimit,
180  const G4AffineTransform &transform,
181  G4double &min,
182  G4double &max ) const
183 {
184  G4SolidExtentList extentList( axis, voxelLimit );
185 
186  //
187  // Loop over all faces, checking min/max extent as we go.
188  //
189  G4VCSGface **face = faces;
190  do
191  {
192  (*face)->CalculateExtent( axis, voxelLimit, transform, extentList );
193  } while( ++face < faces + numFace );
194 
195  //
196  // Return min/max value
197  //
198  return extentList.GetExtent( min, max );
199 }
200 
201 
202 //
203 // Inside
204 //
205 // It could be a good idea to override this virtual
206 // member to add first a simple test (such as spherical
207 // test or whatnot) and to call this version only if
208 // the simplier test fails.
209 //
211 {
212  EInside answer=kOutside;
213  G4VCSGface **face = faces;
214  G4double best = kInfinity;
215  do
216  {
217  G4double distance;
218  EInside result = (*face)->Inside( p, kCarTolerance/2, &distance );
219  if (result == kSurface) { return kSurface; }
220  if (distance < best)
221  {
222  best = distance;
223  answer = result;
224  }
225  } while( ++face < faces + numFace );
226 
227  return answer;
228 }
229 
230 
231 //
232 // SurfaceNormal
233 //
235 {
236  G4ThreeVector answer;
237  G4VCSGface **face = faces;
238  G4double best = kInfinity;
239  do
240  {
241  G4double distance;
242  G4ThreeVector normal = (*face)->Normal( p, &distance );
243  if (distance < best)
244  {
245  best = distance;
246  answer = normal;
247  }
248  } while( ++face < faces + numFace );
249 
250  return answer;
251 }
252 
253 
254 //
255 // DistanceToIn(p,v)
256 //
258  const G4ThreeVector &v ) const
259 {
260  G4double distance = kInfinity;
261  G4double distFromSurface = kInfinity;
262  G4VCSGface **face = faces;
263  G4VCSGface *bestFace = *face;
264  do
265  {
266  G4double faceDistance,
267  faceDistFromSurface;
268  G4ThreeVector faceNormal;
269  G4bool faceAllBehind;
270  if ((*face)->Intersect( p, v, false, kCarTolerance/2,
271  faceDistance, faceDistFromSurface,
272  faceNormal, faceAllBehind ) )
273  {
274  //
275  // Intersecting face
276  //
277  if (faceDistance < distance)
278  {
279  distance = faceDistance;
280  distFromSurface = faceDistFromSurface;
281  bestFace = *face;
282  if (distFromSurface <= 0) { return 0; }
283  }
284  }
285  } while( ++face < faces + numFace );
286 
287  if (distance < kInfinity && distFromSurface<kCarTolerance/2)
288  {
289  if (bestFace->Distance(p,false) < kCarTolerance/2) { distance = 0; }
290  }
291 
292  return distance;
293 }
294 
295 
296 //
297 // DistanceToIn(p)
298 //
300 {
301  return DistanceTo( p, false );
302 }
303 
304 
305 //
306 // DistanceToOut(p,v)
307 //
309  const G4ThreeVector &v,
310  const G4bool calcNorm,
311  G4bool *validNorm,
312  G4ThreeVector *n ) const
313 {
314  G4bool allBehind = true;
315  G4double distance = kInfinity;
316  G4double distFromSurface = kInfinity;
317  G4ThreeVector normal;
318 
319  G4VCSGface **face = faces;
320  G4VCSGface *bestFace = *face;
321  do
322  {
323  G4double faceDistance,
324  faceDistFromSurface;
325  G4ThreeVector faceNormal;
326  G4bool faceAllBehind;
327  if ((*face)->Intersect( p, v, true, kCarTolerance/2,
328  faceDistance, faceDistFromSurface,
329  faceNormal, faceAllBehind ) )
330  {
331  //
332  // Intersecting face
333  //
334  if ( (distance < kInfinity) || (!faceAllBehind) ) { allBehind = false; }
335  if (faceDistance < distance)
336  {
337  distance = faceDistance;
338  distFromSurface = faceDistFromSurface;
339  normal = faceNormal;
340  bestFace = *face;
341  if (distFromSurface <= 0) { break; }
342  }
343  }
344  } while( ++face < faces + numFace );
345 
346  if (distance < kInfinity)
347  {
348  if (distFromSurface <= 0)
349  {
350  distance = 0;
351  }
352  else if (distFromSurface<kCarTolerance/2)
353  {
354  if (bestFace->Distance(p,true) < kCarTolerance/2) { distance = 0; }
355  }
356 
357  if (calcNorm)
358  {
359  *validNorm = allBehind;
360  *n = normal;
361  }
362  }
363  else
364  {
365  if (Inside(p) == kSurface) { distance = 0; }
366  if (calcNorm) { *validNorm = false; }
367  }
368 
369  return distance;
370 }
371 
372 
373 //
374 // DistanceToOut(p)
375 //
377 {
378  return DistanceTo( p, true );
379 }
380 
381 
382 //
383 // DistanceTo
384 //
385 // Protected routine called by DistanceToIn and DistanceToOut
386 //
388  const G4bool outgoing ) const
389 {
390  G4VCSGface **face = faces;
391  G4double best = kInfinity;
392  do
393  {
394  G4double distance = (*face)->Distance( p, outgoing );
395  if (distance < best) { best = distance; }
396  } while( ++face < faces + numFace );
397 
398  return (best < 0.5*kCarTolerance) ? 0 : best;
399 }
400 
401 
402 //
403 // DescribeYourselfTo
404 //
406 {
407  scene.AddSolid( *this );
408 }
409 
410 
411 //
412 // GetExtent
413 //
414 // Define the sides of the box into which our solid instance would fit.
415 //
417 {
418  static const G4ThreeVector xMax(1,0,0), xMin(-1,0,0),
419  yMax(0,1,0), yMin(0,-1,0),
420  zMax(0,0,1), zMin(0,0,-1);
421  static const G4ThreeVector *axes[6] =
422  { &xMin, &xMax, &yMin, &yMax, &zMin, &zMax };
423 
424  G4double answers[6] =
425  {-kInfinity, -kInfinity, -kInfinity, -kInfinity, -kInfinity, -kInfinity};
426 
427  G4VCSGface **face = faces;
428  do
429  {
430  const G4ThreeVector **axis = axes+5 ;
431  G4double *answer = answers+5;
432  do
433  {
434  G4double testFace = (*face)->Extent( **axis );
435  if (testFace > *answer) { *answer = testFace; }
436  }
437  while( --axis, --answer >= answers );
438 
439  } while( ++face < faces + numFace );
440 
441  return G4VisExtent( -answers[0], answers[1],
442  -answers[2], answers[3],
443  -answers[4], answers[5] );
444 }
445 
446 
447 //
448 // GetEntityType
449 //
451 {
452  return G4String("G4CSGfaceted");
453 }
454 
455 
456 //
457 // Stream object contents to an output stream
458 //
459 std::ostream& G4VCSGfaceted::StreamInfo( std::ostream& os ) const
460 {
461  os << "-----------------------------------------------------------\n"
462  << " *** Dump for solid - " << GetName() << " ***\n"
463  << " ===================================================\n"
464  << " Solid type: G4VCSGfaceted\n"
465  << " Parameters: \n"
466  << " number of faces: " << numFace << "\n"
467  << "-----------------------------------------------------------\n";
468 
469  return os;
470 }
471 
472 
473 //
474 // GetCubVolStatistics
475 //
477 {
478  return fStatistics;
479 }
480 
481 
482 //
483 // GetCubVolEpsilon
484 //
486 {
487  return fCubVolEpsilon;
488 }
489 
490 
491 //
492 // SetCubVolStatistics
493 //
495 {
496  fCubicVolume=0.;
497  fStatistics=st;
498 }
499 
500 
501 //
502 // SetCubVolEpsilon
503 //
505 {
506  fCubicVolume=0.;
507  fCubVolEpsilon=ep;
508 }
509 
510 
511 //
512 // GetAreaStatistics
513 //
515 {
516  return fStatistics;
517 }
518 
519 
520 //
521 // GetAreaAccuracy
522 //
524 {
525  return fAreaAccuracy;
526 }
527 
528 
529 //
530 // SetAreaStatistics
531 //
533 {
534  fSurfaceArea=0.;
535  fStatistics=st;
536 }
537 
538 
539 //
540 // SetAreaAccuracy
541 //
543 {
544  fSurfaceArea=0.;
545  fAreaAccuracy=ep;
546 }
547 
548 
549 //
550 // GetCubicVolume
551 //
553 {
554  if(fCubicVolume != 0.) {;}
555  else { fCubicVolume = EstimateCubicVolume(fStatistics,fCubVolEpsilon); }
556  return fCubicVolume;
557 }
558 
559 
560 //
561 // GetSurfaceArea
562 //
564 {
565  if(fSurfaceArea != 0.) {;}
566  else { fSurfaceArea = EstimateSurfaceArea(fStatistics,fAreaAccuracy); }
567  return fSurfaceArea;
568 }
569 
570 
571 //
572 // GetPolyhedron
573 //
575 {
576  if (!fpPolyhedron ||
579  {
580  delete fpPolyhedron;
582  }
583  return fpPolyhedron;
584 }
585 
586 
587 //
588 // GetPointOnSurfaceGeneric proportional to Areas of faces
589 // in case of GenericPolycone or GenericPolyhedra
590 //
592 {
593  // Preparing variables
594  //
595  G4ThreeVector answer=G4ThreeVector(0.,0.,0.);
596  G4VCSGface **face = faces;
597  G4double area = 0;
598  G4int i;
599  std::vector<G4double> areas;
600 
601  // First step: calculate surface areas
602  //
603  do
604  {
605  G4double result = (*face)->SurfaceArea( );
606  areas.push_back(result);
607  area=area+result;
608  } while( ++face < faces + numFace );
609 
610  // Second Step: choose randomly one surface
611  //
612  G4VCSGface **face1 = faces;
613  G4double chose = area*G4UniformRand();
614  G4double Achose1, Achose2;
615  Achose1=0; Achose2=0.;
616  i=0;
617 
618  do
619  {
620  Achose2+=areas[i];
621  if(chose>=Achose1 && chose<Achose2)
622  {
623  G4ThreeVector point;
624  point= (*face1)->GetPointOnFace();
625  return point;
626  }
627  i++;
628  Achose1=Achose2;
629  } while( ++face1 < faces + numFace );
630 
631  return answer;
632 }