Geant4.10
 All Data Structures 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 66356 2012-12-18 09:02:32Z 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 "G4VisExtent.hh"
58 
59 //
60 // Constructor
61 //
63  : G4VSolid(name),
64  numFace(0), faces(0), fCubicVolume(0.), fSurfaceArea(0.), fpPolyhedron(0),
65  fStatistics(1000000), fCubVolEpsilon(0.001), fAreaAccuracy(-1.)
66 {
67 }
68 
69 
70 //
71 // Fake default constructor - sets only member data and allocates memory
72 // for usage restricted to object persistency.
73 //
75  : G4VSolid(a),
76  numFace(0), faces(0), fCubicVolume(0.), fSurfaceArea(0.), fpPolyhedron(0),
77  fStatistics(1000000), fCubVolEpsilon(0.001), fAreaAccuracy(-1.)
78 {
79 }
80 
81 //
82 // Destructor
83 //
85 {
86  DeleteStuff();
87  delete fpPolyhedron;
88 }
89 
90 
91 //
92 // Copy constructor
93 //
95  : G4VSolid( source )
96 {
97  fStatistics = source.fStatistics;
98  fCubVolEpsilon = source.fCubVolEpsilon;
99  fAreaAccuracy = source.fAreaAccuracy;
100 
101  CopyStuff( source );
102 }
103 
104 
105 //
106 // Assignment operator
107 //
109 {
110  if (&source == this) { return *this; }
111 
112  // Copy base class data
113  //
114  G4VSolid::operator=(source);
115 
116  // Copy data
117  //
118  fStatistics = source.fStatistics;
119  fCubVolEpsilon = source.fCubVolEpsilon;
120  fAreaAccuracy = source.fAreaAccuracy;
121 
122  DeleteStuff();
123  CopyStuff( source );
124 
125  return *this;
126 }
127 
128 
129 //
130 // CopyStuff (protected)
131 //
132 // Copy the contents of source
133 //
135 {
136  numFace = source.numFace;
137  if (numFace == 0) { return; } // odd, but permissable?
138 
139  faces = new G4VCSGface*[numFace];
140 
141  G4VCSGface **face = faces,
142  **sourceFace = source.faces;
143  do
144  {
145  *face = (*sourceFace)->Clone();
146  } while( ++sourceFace, ++face < faces+numFace );
147  fCubicVolume = source.fCubicVolume;
148  fSurfaceArea = source.fSurfaceArea;
149  fpPolyhedron = 0;
150 }
151 
152 
153 //
154 // DeleteStuff (protected)
155 //
156 // Delete all allocated objects
157 //
159 {
160  if (numFace)
161  {
162  G4VCSGface **face = faces;
163  do
164  {
165  delete *face;
166  } while( ++face < faces + numFace );
167 
168  delete [] faces;
169  }
170 }
171 
172 
173 //
174 // CalculateExtent
175 //
177  const G4VoxelLimits &voxelLimit,
178  const G4AffineTransform &transform,
179  G4double &min,
180  G4double &max ) const
181 {
182  G4SolidExtentList extentList( axis, voxelLimit );
183 
184  //
185  // Loop over all faces, checking min/max extent as we go.
186  //
187  G4VCSGface **face = faces;
188  do
189  {
190  (*face)->CalculateExtent( axis, voxelLimit, transform, extentList );
191  } while( ++face < faces + numFace );
192 
193  //
194  // Return min/max value
195  //
196  return extentList.GetExtent( min, max );
197 }
198 
199 
200 //
201 // Inside
202 //
203 // It could be a good idea to override this virtual
204 // member to add first a simple test (such as spherical
205 // test or whatnot) and to call this version only if
206 // the simplier test fails.
207 //
209 {
210  EInside answer=kOutside;
211  G4VCSGface **face = faces;
212  G4double best = kInfinity;
213  do
214  {
215  G4double distance;
216  EInside result = (*face)->Inside( p, kCarTolerance/2, &distance );
217  if (result == kSurface) { return kSurface; }
218  if (distance < best)
219  {
220  best = distance;
221  answer = result;
222  }
223  } while( ++face < faces + numFace );
224 
225  return answer;
226 }
227 
228 
229 //
230 // SurfaceNormal
231 //
233 {
234  G4ThreeVector answer;
235  G4VCSGface **face = faces;
236  G4double best = kInfinity;
237  do
238  {
239  G4double distance;
240  G4ThreeVector normal = (*face)->Normal( p, &distance );
241  if (distance < best)
242  {
243  best = distance;
244  answer = normal;
245  }
246  } while( ++face < faces + numFace );
247 
248  return answer;
249 }
250 
251 
252 //
253 // DistanceToIn(p,v)
254 //
256  const G4ThreeVector &v ) const
257 {
258  G4double distance = kInfinity;
259  G4double distFromSurface = kInfinity;
260  G4VCSGface **face = faces;
261  G4VCSGface *bestFace = *face;
262  do
263  {
264  G4double faceDistance,
265  faceDistFromSurface;
266  G4ThreeVector faceNormal;
267  G4bool faceAllBehind;
268  if ((*face)->Intersect( p, v, false, kCarTolerance/2,
269  faceDistance, faceDistFromSurface,
270  faceNormal, faceAllBehind ) )
271  {
272  //
273  // Intersecting face
274  //
275  if (faceDistance < distance)
276  {
277  distance = faceDistance;
278  distFromSurface = faceDistFromSurface;
279  bestFace = *face;
280  if (distFromSurface <= 0) { return 0; }
281  }
282  }
283  } while( ++face < faces + numFace );
284 
285  if (distance < kInfinity && distFromSurface<kCarTolerance/2)
286  {
287  if (bestFace->Distance(p,false) < kCarTolerance/2) { distance = 0; }
288  }
289 
290  return distance;
291 }
292 
293 
294 //
295 // DistanceToIn(p)
296 //
298 {
299  return DistanceTo( p, false );
300 }
301 
302 
303 //
304 // DistanceToOut(p,v)
305 //
307  const G4ThreeVector &v,
308  const G4bool calcNorm,
309  G4bool *validNorm,
310  G4ThreeVector *n ) const
311 {
312  G4bool allBehind = true;
313  G4double distance = kInfinity;
314  G4double distFromSurface = kInfinity;
315  G4ThreeVector normal;
316 
317  G4VCSGface **face = faces;
318  G4VCSGface *bestFace = *face;
319  do
320  {
321  G4double faceDistance,
322  faceDistFromSurface;
323  G4ThreeVector faceNormal;
324  G4bool faceAllBehind;
325  if ((*face)->Intersect( p, v, true, kCarTolerance/2,
326  faceDistance, faceDistFromSurface,
327  faceNormal, faceAllBehind ) )
328  {
329  //
330  // Intersecting face
331  //
332  if ( (distance < kInfinity) || (!faceAllBehind) ) { allBehind = false; }
333  if (faceDistance < distance)
334  {
335  distance = faceDistance;
336  distFromSurface = faceDistFromSurface;
337  normal = faceNormal;
338  bestFace = *face;
339  if (distFromSurface <= 0) { break; }
340  }
341  }
342  } while( ++face < faces + numFace );
343 
344  if (distance < kInfinity)
345  {
346  if (distFromSurface <= 0)
347  {
348  distance = 0;
349  }
350  else if (distFromSurface<kCarTolerance/2)
351  {
352  if (bestFace->Distance(p,true) < kCarTolerance/2) { distance = 0; }
353  }
354 
355  if (calcNorm)
356  {
357  *validNorm = allBehind;
358  *n = normal;
359  }
360  }
361  else
362  {
363  if (Inside(p) == kSurface) { distance = 0; }
364  if (calcNorm) { *validNorm = false; }
365  }
366 
367  return distance;
368 }
369 
370 
371 //
372 // DistanceToOut(p)
373 //
375 {
376  return DistanceTo( p, true );
377 }
378 
379 
380 //
381 // DistanceTo
382 //
383 // Protected routine called by DistanceToIn and DistanceToOut
384 //
386  const G4bool outgoing ) const
387 {
388  G4VCSGface **face = faces;
389  G4double best = kInfinity;
390  do
391  {
392  G4double distance = (*face)->Distance( p, outgoing );
393  if (distance < best) { best = distance; }
394  } while( ++face < faces + numFace );
395 
396  return (best < 0.5*kCarTolerance) ? 0 : best;
397 }
398 
399 
400 //
401 // DescribeYourselfTo
402 //
404 {
405  scene.AddSolid( *this );
406 }
407 
408 
409 //
410 // GetExtent
411 //
412 // Define the sides of the box into which our solid instance would fit.
413 //
415 {
416  static const G4ThreeVector xMax(1,0,0), xMin(-1,0,0),
417  yMax(0,1,0), yMin(0,-1,0),
418  zMax(0,0,1), zMin(0,0,-1);
419  static const G4ThreeVector *axes[6] =
420  { &xMin, &xMax, &yMin, &yMax, &zMin, &zMax };
421 
422  G4double answers[6] =
423  {-kInfinity, -kInfinity, -kInfinity, -kInfinity, -kInfinity, -kInfinity};
424 
425  G4VCSGface **face = faces;
426  do
427  {
428  const G4ThreeVector **axis = axes+5 ;
429  G4double *answer = answers+5;
430  do
431  {
432  G4double testFace = (*face)->Extent( **axis );
433  if (testFace > *answer) { *answer = testFace; }
434  }
435  while( --axis, --answer >= answers );
436 
437  } while( ++face < faces + numFace );
438 
439  return G4VisExtent( -answers[0], answers[1],
440  -answers[2], answers[3],
441  -answers[4], answers[5] );
442 }
443 
444 
445 //
446 // GetEntityType
447 //
449 {
450  return G4String("G4CSGfaceted");
451 }
452 
453 
454 //
455 // Stream object contents to an output stream
456 //
457 std::ostream& G4VCSGfaceted::StreamInfo( std::ostream& os ) const
458 {
459  os << "-----------------------------------------------------------\n"
460  << " *** Dump for solid - " << GetName() << " ***\n"
461  << " ===================================================\n"
462  << " Solid type: G4VCSGfaceted\n"
463  << " Parameters: \n"
464  << " number of faces: " << numFace << "\n"
465  << "-----------------------------------------------------------\n";
466 
467  return os;
468 }
469 
470 
471 //
472 // GetCubVolStatistics
473 //
475 {
476  return fStatistics;
477 }
478 
479 
480 //
481 // GetCubVolEpsilon
482 //
484 {
485  return fCubVolEpsilon;
486 }
487 
488 
489 //
490 // SetCubVolStatistics
491 //
493 {
494  fCubicVolume=0.;
495  fStatistics=st;
496 }
497 
498 
499 //
500 // SetCubVolEpsilon
501 //
503 {
504  fCubicVolume=0.;
505  fCubVolEpsilon=ep;
506 }
507 
508 
509 //
510 // GetAreaStatistics
511 //
513 {
514  return fStatistics;
515 }
516 
517 
518 //
519 // GetAreaAccuracy
520 //
522 {
523  return fAreaAccuracy;
524 }
525 
526 
527 //
528 // SetAreaStatistics
529 //
531 {
532  fSurfaceArea=0.;
533  fStatistics=st;
534 }
535 
536 
537 //
538 // SetAreaAccuracy
539 //
541 {
542  fSurfaceArea=0.;
543  fAreaAccuracy=ep;
544 }
545 
546 
547 //
548 // GetCubicVolume
549 //
551 {
552  if(fCubicVolume != 0.) {;}
553  else { fCubicVolume = EstimateCubicVolume(fStatistics,fCubVolEpsilon); }
554  return fCubicVolume;
555 }
556 
557 
558 //
559 // GetSurfaceArea
560 //
562 {
563  if(fSurfaceArea != 0.) {;}
564  else { fSurfaceArea = EstimateSurfaceArea(fStatistics,fAreaAccuracy); }
565  return fSurfaceArea;
566 }
567 
568 
569 //
570 // GetPolyhedron
571 //
573 {
574  if (!fpPolyhedron ||
577  {
578  delete fpPolyhedron;
580  }
581  return fpPolyhedron;
582 }
583 
584 
585 //
586 // GetPointOnSurfaceGeneric proportional to Areas of faces
587 // in case of GenericPolycone or GenericPolyhedra
588 //
590 {
591  // Preparing variables
592  //
593  G4ThreeVector answer=G4ThreeVector(0.,0.,0.);
594  G4VCSGface **face = faces;
595  G4double area = 0;
596  G4int i;
597  std::vector<G4double> areas;
598 
599  // First step: calculate surface areas
600  //
601  do
602  {
603  G4double result = (*face)->SurfaceArea( );
604  areas.push_back(result);
605  area=area+result;
606  } while( ++face < faces + numFace );
607 
608  // Second Step: choose randomly one surface
609  //
610  G4VCSGface **face1 = faces;
611  G4double chose = area*G4UniformRand();
612  G4double Achose1, Achose2;
613  Achose1=0; Achose2=0.;
614  i=0;
615 
616  do
617  {
618  Achose2+=areas[i];
619  if(chose>=Achose1 && chose<Achose2)
620  {
621  G4ThreeVector point;
622  point= (*face1)->GetPointOnFace();
623  return point;
624  }
625  i++;
626  Achose1=Achose2;
627  } while( ++face1 < faces + numFace );
628 
629  return answer;
630 }
G4String GetName() const
virtual void DescribeYourselfTo(G4VGraphicsScene &scene) const
G4double fSurfaceArea
void SetAreaStatistics(G4int st)
G4Polyhedron * fpPolyhedron
CLHEP::Hep3Vector G4ThreeVector
virtual G4VCSGface * Clone()=0
const char * p
Definition: xmltok.h:285
G4int GetCubVolStatistics() const
G4double EstimateCubicVolume(G4int nStat, G4double epsilon) const
Definition: G4VSolid.cc:203
G4bool GetExtent(G4double &min, G4double &max) const
void SetCubVolStatistics(G4int st)
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
virtual ~G4VCSGfaceted()
const XML_Char * name
G4VCSGfaceted(const G4String &name)
void CopyStuff(const G4VCSGfaceted &source)
virtual void AddSolid(const G4Box &)=0
int G4int
Definition: G4Types.hh:78
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
virtual G4double SurfaceArea()=0
G4ThreeVector GetPointOnSurfaceGeneric() const
virtual G4double GetSurfaceArea()
G4double fCubicVolume
#define G4UniformRand()
Definition: Randomize.hh:87
virtual G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
virtual G4Polyhedron * CreatePolyhedron() const =0
G4double GetAreaAccuracy() const
bool G4bool
Definition: G4Types.hh:79
virtual G4VisExtent GetExtent() const
const G4VCSGfaceted & operator=(const G4VCSGfaceted &source)
const G4int n
G4double EstimateSurfaceArea(G4int nStat, G4double ell) const
Definition: G4VSolid.cc:263
virtual G4GeometryType GetEntityType() const
G4VCSGface ** faces
virtual std::ostream & StreamInfo(std::ostream &os) const
static G4int GetNumberOfRotationSteps()
T max(const T t1, const T t2)
brief Return the largest of the two arguments
EInside
Definition: geomdefs.hh:58
EAxis
Definition: geomdefs.hh:54
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
virtual G4double DistanceTo(const G4ThreeVector &p, const G4bool outgoing) const
virtual EInside Inside(const G4ThreeVector &p, G4double tolerance, G4double *bestDistance)=0
void SetAreaAccuracy(G4double ep)
G4VSolid & operator=(const G4VSolid &rhs)
Definition: G4VSolid.cc:110
G4double kCarTolerance
Definition: G4VSolid.hh:305
G4int GetNumberOfRotationStepsAtTimeOfCreation() const
G4int GetAreaStatistics() const
double G4double
Definition: G4Types.hh:76
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const
virtual EInside Inside(const G4ThreeVector &p) const
virtual void CalculateExtent(const EAxis axis, const G4VoxelLimits &voxelLimit, const G4AffineTransform &tranform, G4SolidExtentList &extentList)=0
virtual G4double GetCubicVolume()
G4double GetCubVolEpsilon() const
void SetCubVolEpsilon(G4double ep)
virtual G4double Distance(const G4ThreeVector &p, G4bool outgoing)=0
virtual G4Polyhedron * GetPolyhedron() const