Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4VSolid.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 // class G4VSolid
30 //
31 // Implementation for solid base class
32 //
33 // History:
34 //
35 // 06.12.02 V.Grichine, restored original conditions in ClipPolygon()
36 // 10.05.02 V.Grichine, ClipPolygon(): clip only other axis and limited voxels
37 // 15.04.02 V.Grichine, bug fixed in ClipPolygon(): clip only one axis
38 // 13.03.02 V.Grichine, cosmetics of voxel limit functions
39 // 15.11.00 D.Williams, V.Grichine, fix in CalculateClippedPolygonExtent()
40 // 10.07.95 P.Kent, Added == operator, solid Store entry
41 // 30.06.95 P.Kent, Created.
42 // --------------------------------------------------------------------
43 
44 #include "G4VSolid.hh"
45 #include "G4SolidStore.hh"
46 #include "globals.hh"
47 #include "Randomize.hh"
48 #include "G4GeometryTolerance.hh"
49 
50 #include "G4VoxelLimits.hh"
51 #include "G4AffineTransform.hh"
52 #include "G4VisExtent.hh"
53 
55 //
56 // Constructor
57 // - Copies name
58 // - Add ourselves to solid Store
59 
61  : fshapeName(name)
62 {
64 
65  // Register to store
66  //
68 }
69 
71 //
72 // Copy constructor
73 //
74 
76  : kCarTolerance(rhs.kCarTolerance), fshapeName(rhs.fshapeName)
77 {
78  // Register to store
79  //
81 }
82 
84 //
85 // Fake default constructor - sets only member data and allocates memory
86 // for usage restricted to object persistency.
87 //
88 G4VSolid::G4VSolid( __void__& )
89  : fshapeName("")
90 {
91  // Register to store
92  //
94 }
95 
97 //
98 // Destructor (virtual)
99 // - Remove ourselves from solid Store
100 
102 {
104 }
105 
107 //
108 // Assignment operator
109 
111 {
112  // Check assignment to self
113  //
114  if (this == &rhs) { return *this; }
115 
116  // Copy data
117  //
119  fshapeName = rhs.fshapeName;
120 
121  return *this;
122 }
123 
125 //
126 // Streaming operator dumping solid contents
127 
128 std::ostream& operator<< ( std::ostream& os, const G4VSolid& e )
129 {
130  return e.StreamInfo(os);
131 }
132 
134 //
135 // Throw exception if ComputeDimensions called for illegal derived class
136 
138  const G4int,
139  const G4VPhysicalVolume*)
140 {
141  std::ostringstream message;
142  message << "Illegal call to G4VSolid::ComputeDimensions()" << G4endl
143  << "Method not overloaded by derived class !";
144  G4Exception("G4VSolid::ComputeDimensions()", "GeomMgt0003",
145  FatalException, message);
146 }
147 
149 //
150 // Throw exception (warning) for solids not implementing the method
151 
153 {
154  std::ostringstream message;
155  message << "Not implemented for solid: "
156  << this->GetEntityType() << " !" << G4endl
157  << "Returning origin.";
158  G4Exception("G4VSolid::GetPointOnSurface()", "GeomMgt1001",
159  JustWarning, message);
160  return G4ThreeVector(0,0,0);
161 }
162 
164 //
165 // Dummy implementations ...
166 
168 { return 0; }
169 
171 { return 0; }
172 
174 { return 0; }
175 
177 { return 0; }
178 
180 //
181 // Returns an estimation of the solid volume in internal units.
182 // The number of statistics and error accuracy is fixed.
183 // This method may be overloaded by derived classes to compute the
184 // exact geometrical quantity for solids where this is possible.
185 // or anyway to cache the computed value.
186 // This implementation does NOT cache the computed value.
187 
189 {
190  G4int cubVolStatistics = 1000000;
191  G4double cubVolEpsilon = 0.001;
192  return EstimateCubicVolume(cubVolStatistics, cubVolEpsilon);
193 }
194 
196 //
197 // Calculate cubic volume based on Inside() method.
198 // Accuracy is limited by the second argument or the statistics
199 // expressed by the first argument.
200 // Implementation is courtesy of Vasiliki Despoina Mitsou,
201 // University of Athens.
202 
204 {
205  G4int iInside=0;
206  G4double px,py,pz,minX,maxX,minY,maxY,minZ,maxZ,volume;
208  EInside in;
209 
210  // values needed for CalculateExtent signature
211 
212  G4VoxelLimits limit; // Unlimited
213  G4AffineTransform origin;
214 
215  // min max extents of pSolid along X,Y,Z
216 
217  this->CalculateExtent(kXAxis,limit,origin,minX,maxX);
218  this->CalculateExtent(kYAxis,limit,origin,minY,maxY);
219  this->CalculateExtent(kZAxis,limit,origin,minZ,maxZ);
220 
221  // limits
222 
223  if(nStat < 100) nStat = 100;
224  if(epsilon > 0.01) epsilon = 0.01;
225 
226  for(G4int i = 0; i < nStat; i++ )
227  {
228  px = minX+(maxX-minX)*G4UniformRand();
229  py = minY+(maxY-minY)*G4UniformRand();
230  pz = minZ+(maxZ-minZ)*G4UniformRand();
231  p = G4ThreeVector(px,py,pz);
232  in = this->Inside(p);
233  if(in != kOutside) iInside++;
234  }
235  volume = (maxX-minX)*(maxY-minY)*(maxZ-minZ)*iInside/nStat;
236  return volume;
237 }
238 
240 //
241 // Returns an estimation of the solid surface area in internal units.
242 // The number of statistics and error accuracy is fixed.
243 // This method may be overloaded by derived classes to compute the
244 // exact geometrical quantity for solids where this is possible.
245 // or anyway to cache the computed value.
246 // This implementation does NOT cache the computed value.
247 
249 {
250  G4int stat = 1000000;
251  G4double ell = -1.;
252  return EstimateSurfaceArea(stat,ell);
253 }
254 
256 //
257 // Estimate surface area based on Inside(), DistanceToIn(), and
258 // DistanceToOut() methods. Accuracy is limited by the statistics
259 // defined by the first argument. Implemented by Mikhail Kosov.
260 
262 {
263  G4int inside=0;
264  G4double px,py,pz,minX,maxX,minY,maxY,minZ,maxZ,surf;
266  EInside in;
267 
268  // values needed for CalculateExtent signature
269 
270  G4VoxelLimits limit; // Unlimited
271  G4AffineTransform origin;
272 
273  // min max extents of pSolid along X,Y,Z
274 
275  this->CalculateExtent(kXAxis,limit,origin,minX,maxX);
276  this->CalculateExtent(kYAxis,limit,origin,minY,maxY);
277  this->CalculateExtent(kZAxis,limit,origin,minZ,maxZ);
278 
279  // limits
280 
281  if(nStat < 100) { nStat = 100; }
282 
283  G4double dX=maxX-minX;
284  G4double dY=maxY-minY;
285  G4double dZ=maxZ-minZ;
286  if(ell<=0.) // Automatic definition of skin thickness
287  {
288  G4double minval=dX;
289  if(dY<dX) { minval=dY; }
290  if(dZ<minval) { minval=dZ; }
291  ell=.01*minval;
292  }
293 
294  G4double dd=2*ell;
295  minX-=ell; minY-=ell; minZ-=ell; dX+=dd; dY+=dd; dZ+=dd;
296 
297  for(G4int i = 0; i < nStat; i++ )
298  {
299  px = minX+dX*G4UniformRand();
300  py = minY+dY*G4UniformRand();
301  pz = minZ+dZ*G4UniformRand();
302  p = G4ThreeVector(px,py,pz);
303  in = this->Inside(p);
304  if(in != kOutside)
305  {
306  if (DistanceToOut(p)<ell) { inside++; }
307  }
308  else if(DistanceToIn(p)<ell) { inside++; }
309  }
310  // @@ The conformal correction can be upgraded
311  surf = dX*dY*dZ*inside/dd/nStat;
312  return surf;
313 }
314 
316 //
317 // Returns a pointer of a dynamically allocated copy of the solid.
318 // Returns NULL pointer with warning in case the concrete solid does not
319 // implement this method. The caller has responsibility for ownership.
320 //
321 
323 {
324  std::ostringstream message;
325  message << "Clone() method not implemented for type: "
326  << GetEntityType() << "!" << G4endl
327  << "Returning NULL pointer!";
328  G4Exception("G4VSolid::Clone()", "GeomMgt1001", JustWarning, message);
329  return 0;
330 }
331 
333 //
334 // Calculate the maximum and minimum extents of the polygon described
335 // by the vertices: pSectionIndex->pSectionIndex+1->
336 // pSectionIndex+2->pSectionIndex+3->pSectionIndex
337 // in the List pVertices
338 //
339 // If the minimum is <pMin pMin is set to the new minimum
340 // If the maximum is >pMax pMax is set to the new maximum
341 //
342 // No modifications are made to pVertices
343 //
344 
346  const G4int pSectionIndex,
347  const G4VoxelLimits& pVoxelLimit,
348  const EAxis pAxis,
349  G4double& pMin, G4double& pMax) const
350 {
351 
352  G4ThreeVectorList polygon;
353  polygon.reserve(4);
354  polygon.push_back((*pVertices)[pSectionIndex]);
355  polygon.push_back((*pVertices)[pSectionIndex+1]);
356  polygon.push_back((*pVertices)[pSectionIndex+2]);
357  polygon.push_back((*pVertices)[pSectionIndex+3]);
358  // G4cout<<"ClipCrossSection: 0-1-2-3"<<G4endl;
359  CalculateClippedPolygonExtent(polygon,pVoxelLimit,pAxis,pMin,pMax);
360  return;
361 }
362 
364 //
365 // Calculate the maximum and minimum extents of the polygons
366 // joining the CrossSections at pSectionIndex->pSectionIndex+3 and
367 // pSectionIndex+4->pSectionIndex7
368 //
369 // in the List pVertices, within the boundaries of the voxel limits pVoxelLimit
370 //
371 // If the minimum is <pMin pMin is set to the new minimum
372 // If the maximum is >pMax pMax is set to the new maximum
373 //
374 // No modifications are made to pVertices
375 
377  const G4int pSectionIndex,
378  const G4VoxelLimits& pVoxelLimit,
379  const EAxis pAxis,
380  G4double& pMin, G4double& pMax) const
381 {
382  G4ThreeVectorList polygon;
383  polygon.reserve(4);
384  polygon.push_back((*pVertices)[pSectionIndex]);
385  polygon.push_back((*pVertices)[pSectionIndex+4]);
386  polygon.push_back((*pVertices)[pSectionIndex+5]);
387  polygon.push_back((*pVertices)[pSectionIndex+1]);
388  // G4cout<<"ClipBetweenSections: 0-4-5-1"<<G4endl;
389  CalculateClippedPolygonExtent(polygon,pVoxelLimit,pAxis,pMin,pMax);
390  polygon.clear();
391 
392  polygon.push_back((*pVertices)[pSectionIndex+1]);
393  polygon.push_back((*pVertices)[pSectionIndex+5]);
394  polygon.push_back((*pVertices)[pSectionIndex+6]);
395  polygon.push_back((*pVertices)[pSectionIndex+2]);
396  // G4cout<<"ClipBetweenSections: 1-5-6-2"<<G4endl;
397  CalculateClippedPolygonExtent(polygon,pVoxelLimit,pAxis,pMin,pMax);
398  polygon.clear();
399 
400  polygon.push_back((*pVertices)[pSectionIndex+2]);
401  polygon.push_back((*pVertices)[pSectionIndex+6]);
402  polygon.push_back((*pVertices)[pSectionIndex+7]);
403  polygon.push_back((*pVertices)[pSectionIndex+3]);
404  // G4cout<<"ClipBetweenSections: 2-6-7-3"<<G4endl;
405  CalculateClippedPolygonExtent(polygon,pVoxelLimit,pAxis,pMin,pMax);
406  polygon.clear();
407 
408  polygon.push_back((*pVertices)[pSectionIndex+3]);
409  polygon.push_back((*pVertices)[pSectionIndex+7]);
410  polygon.push_back((*pVertices)[pSectionIndex+4]);
411  polygon.push_back((*pVertices)[pSectionIndex]);
412  // G4cout<<"ClipBetweenSections: 3-7-4-0"<<G4endl;
413  CalculateClippedPolygonExtent(polygon,pVoxelLimit,pAxis,pMin,pMax);
414  return;
415 }
416 
417 
419 //
420 // Calculate the maximum and minimum extents of the convex polygon pPolygon
421 // along the axis pAxis, within the limits pVoxelLimit
422 //
423 
424 void
426  const G4VoxelLimits& pVoxelLimit,
427  const EAxis pAxis,
428  G4double& pMin,
429  G4double& pMax) const
430 {
431  G4int noLeft,i;
432  G4double component;
433  /*
434  G4cout<<G4endl;
435  for(i = 0 ; i < pPolygon.size() ; i++ )
436  {
437  G4cout << i << "\t"
438  << "p.x = " << pPolygon[i].operator()(pAxis) << "\t"
439  // << "p.y = " << pPolygon[i].y() << "\t"
440  // << "p.z = " << pPolygon[i].z() << "\t"
441  << G4endl;
442  }
443  G4cout<<G4endl;
444  */
445  ClipPolygon(pPolygon,pVoxelLimit,pAxis);
446  noLeft = pPolygon.size();
447 
448  if ( noLeft )
449  {
450  // G4cout<<G4endl;
451  for (i=0;i<noLeft;i++)
452  {
453  component = pPolygon[i].operator()(pAxis);
454  // G4cout <<i<<"\t"<<component<<G4endl;
455 
456  if (component < pMin)
457  {
458  // G4cout <<i<<"\t"<<"Pmin = "<<component<<G4endl;
459  pMin = component;
460  }
461  if (component > pMax)
462  {
463  // G4cout <<i<<"\t"<<"PMax = "<<component<<G4endl;
464  pMax = component;
465  }
466  }
467  // G4cout<<G4endl;
468  }
469  // G4cout<<"pMin = "<<pMin<<"\t"<<"pMax = "<<pMax<<G4endl;
470 }
471 
473 //
474 // Clip the convex polygon described by the vertices at
475 // pSectionIndex ->pSectionIndex+3 within pVertices to the limits pVoxelLimit
476 //
477 // Set pMin to the smallest
478 //
479 // Calculate the extent of the polygon along pAxis, when clipped to the
480 // limits pVoxelLimit. If the polygon exists after clippin, set pMin to
481 // the polygon's minimum extent along the axis if <pMin, and set pMax to
482 // the polygon's maximum extent along the axis if >pMax.
483 //
484 // The polygon is described by a set of vectors, where each vector represents
485 // a vertex, so that the polygon is described by the vertex sequence:
486 // 0th->1st 1st->2nd 2nd->... nth->0th
487 //
488 // Modifications to the polygon are made
489 //
490 // NOTE: Execessive copying during clipping
491 
493  const G4VoxelLimits& pVoxelLimit,
494  const EAxis ) const
495 {
496  G4ThreeVectorList outputPolygon;
497 
498  if ( pVoxelLimit.IsLimited() )
499  {
500  if (pVoxelLimit.IsXLimited() ) // && pAxis != kXAxis)
501  {
502  G4VoxelLimits simpleLimit1;
503  simpleLimit1.AddLimit(kXAxis,pVoxelLimit.GetMinXExtent(),kInfinity);
504  // G4cout<<"MinXExtent()"<<G4endl;
505  ClipPolygonToSimpleLimits(pPolygon,outputPolygon,simpleLimit1);
506 
507  pPolygon.clear();
508 
509  if ( !outputPolygon.size() ) return;
510 
511  G4VoxelLimits simpleLimit2;
512  // G4cout<<"MaxXExtent()"<<G4endl;
513  simpleLimit2.AddLimit(kXAxis,-kInfinity,pVoxelLimit.GetMaxXExtent());
514  ClipPolygonToSimpleLimits(outputPolygon,pPolygon,simpleLimit2);
515 
516  if ( !pPolygon.size() ) return;
517  else outputPolygon.clear();
518  }
519  if ( pVoxelLimit.IsYLimited() ) // && pAxis != kYAxis)
520  {
521  G4VoxelLimits simpleLimit1;
522  simpleLimit1.AddLimit(kYAxis,pVoxelLimit.GetMinYExtent(),kInfinity);
523  ClipPolygonToSimpleLimits(pPolygon,outputPolygon,simpleLimit1);
524 
525  // Must always clear pPolygon - for clip to simpleLimit2 and in case of
526  // early exit
527 
528  pPolygon.clear();
529 
530  if ( !outputPolygon.size() ) return;
531 
532  G4VoxelLimits simpleLimit2;
533  simpleLimit2.AddLimit(kYAxis,-kInfinity,pVoxelLimit.GetMaxYExtent());
534  ClipPolygonToSimpleLimits(outputPolygon,pPolygon,simpleLimit2);
535 
536  if ( !pPolygon.size() ) return;
537  else outputPolygon.clear();
538  }
539  if ( pVoxelLimit.IsZLimited() ) // && pAxis != kZAxis)
540  {
541  G4VoxelLimits simpleLimit1;
542  simpleLimit1.AddLimit(kZAxis,pVoxelLimit.GetMinZExtent(),kInfinity);
543  ClipPolygonToSimpleLimits(pPolygon,outputPolygon,simpleLimit1);
544 
545  // Must always clear pPolygon - for clip to simpleLimit2 and in case of
546  // early exit
547 
548  pPolygon.clear();
549 
550  if ( !outputPolygon.size() ) return;
551 
552  G4VoxelLimits simpleLimit2;
553  simpleLimit2.AddLimit(kZAxis,-kInfinity,pVoxelLimit.GetMaxZExtent());
554  ClipPolygonToSimpleLimits(outputPolygon,pPolygon,simpleLimit2);
555 
556  // Return after final clip - no cleanup
557  }
558  }
559 }
560 
562 //
563 // pVoxelLimits must be only limited along one axis, and either the maximum
564 // along the axis must be +kInfinity, or the minimum -kInfinity
565 
566 void
567 G4VSolid::ClipPolygonToSimpleLimits( G4ThreeVectorList& pPolygon,
568  G4ThreeVectorList& outputPolygon,
569  const G4VoxelLimits& pVoxelLimit ) const
570 {
571  G4int i;
572  G4int noVertices=pPolygon.size();
573  G4ThreeVector vEnd,vStart;
574 
575  for (i = 0 ; i < noVertices ; i++ )
576  {
577  vStart = pPolygon[i];
578  // G4cout << "i = " << i << G4endl;
579  if ( i == noVertices-1 ) vEnd = pPolygon[0];
580  else vEnd = pPolygon[i+1];
581 
582  if ( pVoxelLimit.Inside(vStart) )
583  {
584  if (pVoxelLimit.Inside(vEnd))
585  {
586  // vStart and vEnd inside -> output end point
587  //
588  outputPolygon.push_back(vEnd);
589  }
590  else
591  {
592  // vStart inside, vEnd outside -> output crossing point
593  //
594  // G4cout << "vStart inside, vEnd outside" << G4endl;
595  pVoxelLimit.ClipToLimits(vStart,vEnd);
596  outputPolygon.push_back(vEnd);
597  }
598  }
599  else
600  {
601  if (pVoxelLimit.Inside(vEnd))
602  {
603  // vStart outside, vEnd inside -> output inside section
604  //
605  // G4cout << "vStart outside, vEnd inside" << G4endl;
606  pVoxelLimit.ClipToLimits(vStart,vEnd);
607  outputPolygon.push_back(vStart);
608  outputPolygon.push_back(vEnd);
609  }
610  else // Both point outside -> no output
611  {
612  // outputPolygon.push_back(vStart);
613  // outputPolygon.push_back(vEnd);
614  }
615  }
616  }
617 }
618 
620 {
621  G4VisExtent extent;
622  G4VoxelLimits voxelLimits; // Defaults to "infinite" limits.
623  G4AffineTransform affineTransform;
624  G4double vmin, vmax;
625  CalculateExtent(kXAxis,voxelLimits,affineTransform,vmin,vmax);
626  extent.SetXmin (vmin);
627  extent.SetXmax (vmax);
628  CalculateExtent(kYAxis,voxelLimits,affineTransform,vmin,vmax);
629  extent.SetYmin (vmin);
630  extent.SetYmax (vmax);
631  CalculateExtent(kZAxis,voxelLimits,affineTransform,vmin,vmax);
632  extent.SetZmin (vmin);
633  extent.SetZmax (vmax);
634  return extent;
635 }
636 
638 {
639  return 0;
640 }
641 
643 {
644  return 0;
645 }
646 
648 {
649  return 0;
650 }