Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4SubtractionSolid.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// Implementation of methods for the class G4IntersectionSolid
27//
28// 22.07.11 T.Nikitina: added detection of infinite loop in DistanceToIn(p,v)
29// 19.10.98 V.Grichine: new algorithm of DistanceToIn(p,v)
30// 14.10.98 V.Grichine: implementation of the first version
31// --------------------------------------------------------------------
32
33#include "G4SubtractionSolid.hh"
34
35#include "G4SystemOfUnits.hh"
36#include "G4VoxelLimits.hh"
39
40#include "G4VGraphicsScene.hh"
41#include "G4Polyhedron.hh"
43
45
46#include <sstream>
47
48///////////////////////////////////////////////////////////////////
49//
50// Transfer all data members to G4BooleanSolid which is responsible
51// for them. pName will be in turn sent to G4VSolid
52
54 G4VSolid* pSolidA ,
55 G4VSolid* pSolidB )
56 : G4BooleanSolid(pName,pSolidA,pSolidB)
57{
58}
59
60///////////////////////////////////////////////////////////////
61//
62// Constructor
63
65 G4VSolid* pSolidA ,
66 G4VSolid* pSolidB ,
67 G4RotationMatrix* rotMatrix,
68 const G4ThreeVector& transVector )
69 : G4BooleanSolid(pName,pSolidA,pSolidB,rotMatrix,transVector)
70{
71}
72
73///////////////////////////////////////////////////////////////
74//
75// Constructor
76
78 G4VSolid* pSolidA ,
79 G4VSolid* pSolidB ,
80 const G4Transform3D& transform )
81 : G4BooleanSolid(pName,pSolidA,pSolidB,transform)
82{
83}
84
85//////////////////////////////////////////////////////////////////
86//
87// Fake default constructor - sets only member data and allocates memory
88// for usage restricted to object persistency.
89
92{
93}
94
95///////////////////////////////////////////////////////////////
96//
97// Destructor
98
100{
101}
102
103///////////////////////////////////////////////////////////////
104//
105// Copy constructor
106
108 : G4BooleanSolid (rhs)
109{
110}
111
112///////////////////////////////////////////////////////////////
113//
114// Assignment operator
115
118{
119 // Check assignment to self
120 //
121 if (this == &rhs) { return *this; }
122
123 // Copy base class data
124 //
126
127 return *this;
128}
129
130//////////////////////////////////////////////////////////////////////////
131//
132// Get bounding box
133
134void
136 G4ThreeVector& pMax) const
137{
138 // Since it is unclear how the shape of the first solid will be changed
139 // after subtraction, just return its original bounding box.
140 //
141 fPtrSolidA->BoundingLimits(pMin,pMax);
142
143 // Check correctness of the bounding box
144 //
145 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
146 {
147 std::ostringstream message;
148 message << "Bad bounding box (min >= max) for solid: "
149 << GetName() << " !"
150 << "\npMin = " << pMin
151 << "\npMax = " << pMax;
152 G4Exception("G4SubtractionSolid::BoundingLimits()", "GeomMgt0001",
153 JustWarning, message);
154 DumpInfo();
155 }
156}
157
158//////////////////////////////////////////////////////////////////////////
159//
160// Calculate extent under transform and specified limit
161
162G4bool
164 const G4VoxelLimits& pVoxelLimit,
165 const G4AffineTransform& pTransform,
166 G4double& pMin,
167 G4double& pMax ) const
168{
169 // Since we cannot be sure how much the second solid subtracts
170 // from the first, we must use the first solid's extent!
171
172 return fPtrSolidA->CalculateExtent( pAxis, pVoxelLimit,
173 pTransform, pMin, pMax );
174}
175
176/////////////////////////////////////////////////////
177//
178// Touching ? Empty subtraction ?
179
181{
182 EInside positionA = fPtrSolidA->Inside(p);
183 if (positionA == kOutside) return positionA; // outside A
184
185 EInside positionB = fPtrSolidB->Inside(p);
186 if (positionB == kOutside) return positionA;
187
188 if (positionB == kInside) return kOutside;
189 if (positionA == kInside) return kSurface; // surface B
190
191 // Point is on both surfaces
192 //
193 static const G4double rtol = 1000*kCarTolerance;
194
195 return ((fPtrSolidA->SurfaceNormal(p) -
196 fPtrSolidB->SurfaceNormal(p)).mag2() > rtol) ? kSurface : kOutside;
197}
198
199//////////////////////////////////////////////////////////////
200//
201// SurfaceNormal
202
205{
206 G4ThreeVector normal;
207
208 EInside InsideA = fPtrSolidA->Inside(p);
209 EInside InsideB = fPtrSolidB->Inside(p);
210
211 if( InsideA == kOutside )
212 {
213#ifdef G4BOOLDEBUG
214 G4cout << "WARNING - Invalid call [1] in "
215 << "G4SubtractionSolid::SurfaceNormal(p)" << G4endl
216 << " Point p is outside !" << G4endl;
217 G4cout << " p = " << p << G4endl;
218 G4cerr << "WARNING - Invalid call [1] in "
219 << "G4SubtractionSolid::SurfaceNormal(p)" << G4endl
220 << " Point p is outside !" << G4endl;
221 G4cerr << " p = " << p << G4endl;
222#endif
223 normal = fPtrSolidA->SurfaceNormal(p) ;
224 }
225 else if( InsideA == kSurface &&
226 InsideB != kInside )
227 {
228 normal = fPtrSolidA->SurfaceNormal(p) ;
229 }
230 else if( InsideA == kInside &&
231 InsideB != kOutside )
232 {
233 normal = -fPtrSolidB->SurfaceNormal(p) ;
234 }
235 else
236 {
238 {
239 normal = fPtrSolidA->SurfaceNormal(p) ;
240 }
241 else
242 {
243 normal = -fPtrSolidB->SurfaceNormal(p) ;
244 }
245#ifdef G4BOOLDEBUG
246 if(Inside(p) == kInside)
247 {
248 G4cout << "WARNING - Invalid call [2] in "
249 << "G4SubtractionSolid::SurfaceNormal(p)" << G4endl
250 << " Point p is inside !" << G4endl;
251 G4cout << " p = " << p << G4endl;
252 G4cerr << "WARNING - Invalid call [2] in "
253 << "G4SubtractionSolid::SurfaceNormal(p)" << G4endl
254 << " Point p is inside !" << G4endl;
255 G4cerr << " p = " << p << G4endl;
256 }
257#endif
258 }
259 return normal;
260}
261
262/////////////////////////////////////////////////////////////
263//
264// The same algorithm as in DistanceToIn(p)
265
268 const G4ThreeVector& v ) const
269{
270 G4double dist = 0.0, dist2 = 0.0, disTmp = 0.0;
271
272#ifdef G4BOOLDEBUG
273 if( Inside(p) == kInside )
274 {
275 G4cout << "WARNING - Invalid call in "
276 << "G4SubtractionSolid::DistanceToIn(p,v)" << G4endl
277 << " Point p is inside !" << G4endl;
278 G4cout << " p = " << p << G4endl;
279 G4cout << " v = " << v << G4endl;
280 G4cerr << "WARNING - Invalid call in "
281 << "G4SubtractionSolid::DistanceToIn(p,v)" << G4endl
282 << " Point p is inside !" << G4endl;
283 G4cerr << " p = " << p << G4endl;
284 G4cerr << " v = " << v << G4endl;
285 }
286#endif
287
288 // if( // ( fPtrSolidA->Inside(p) != kOutside) && // case1:p in both A&B
289 if ( fPtrSolidB->Inside(p) != kOutside ) // start: out of B
290 {
291 dist = fPtrSolidB->DistanceToOut(p,v) ; // ,calcNorm,validNorm,n) ;
292
293 if( fPtrSolidA->Inside(p+dist*v) != kInside )
294 {
295 G4int count1=0;
296 do // Loop checking, 13.08.2015, G.Cosmo
297 {
298 disTmp = fPtrSolidA->DistanceToIn(p+dist*v,v) ;
299
300 if(disTmp == kInfinity)
301 {
302 return kInfinity ;
303 }
304 dist += disTmp ;
305
306 if( Inside(p+dist*v) == kOutside )
307 {
308 disTmp = fPtrSolidB->DistanceToOut(p+dist*v,v) ;
309 dist2 = dist+disTmp;
310 if (dist == dist2) { return dist; } // no progress
311 dist = dist2 ;
312 ++count1;
313 if( count1 > 1000 ) // Infinite loop detected
314 {
315 G4String nameB = fPtrSolidB->GetName();
316 if(fPtrSolidB->GetEntityType()=="G4DisplacedSolid")
317 {
318 nameB = (dynamic_cast<G4DisplacedSolid*>(fPtrSolidB))
319 ->GetConstituentMovedSolid()->GetName();
320 }
321 std::ostringstream message;
322 message << "Illegal condition caused by solids: "
323 << fPtrSolidA->GetName() << " and " << nameB << G4endl;
324 message.precision(16);
325 message << "Looping detected in point " << p+dist*v
326 << ", from original point " << p
327 << " and direction " << v << G4endl
328 << "Computed candidate distance: " << dist << "*mm. ";
329 message.precision(6);
330 DumpInfo();
331 G4Exception("G4SubtractionSolid::DistanceToIn(p,v)",
332 "GeomSolids1001", JustWarning, message,
333 "Returning candidate distance.");
334 return dist;
335 }
336 }
337 }
338 while( Inside(p+dist*v) == kOutside ) ;
339 }
340 }
341 else // p outside A, start in A
342 {
343 dist = fPtrSolidA->DistanceToIn(p,v) ;
344
345 if( dist == kInfinity ) // past A, hence past A\B
346 {
347 return kInfinity ;
348 }
349 else
350 {
351 G4int count2=0;
352 while( Inside(p+dist*v) == kOutside ) // pushing loop
353 {
354 disTmp = fPtrSolidB->DistanceToOut(p+dist*v,v) ;
355 dist += disTmp ;
356
357 if( Inside(p+dist*v) == kOutside )
358 {
359 disTmp = fPtrSolidA->DistanceToIn(p+dist*v,v) ;
360
361 if(disTmp == kInfinity) // past A, hence past A\B
362 {
363 return kInfinity ;
364 }
365 dist2 = dist+disTmp;
366 if (dist == dist2) { return dist; } // no progress
367 dist = dist2 ;
368 ++count2;
369 if( count2 > 1000 ) // Infinite loop detected
370 {
371 G4String nameB = fPtrSolidB->GetName();
372 if(fPtrSolidB->GetEntityType()=="G4DisplacedSolid")
373 {
374 nameB = (dynamic_cast<G4DisplacedSolid*>(fPtrSolidB))
375 ->GetConstituentMovedSolid()->GetName();
376 }
377 std::ostringstream message;
378 message << "Illegal condition caused by solids: "
379 << fPtrSolidA->GetName() << " and " << nameB << G4endl;
380 message.precision(16);
381 message << "Looping detected in point " << p+dist*v
382 << ", from original point " << p
383 << " and direction " << v << G4endl
384 << "Computed candidate distance: " << dist << "*mm. ";
385 message.precision(6);
386 DumpInfo();
387 G4Exception("G4SubtractionSolid::DistanceToIn(p,v)",
388 "GeomSolids1001", JustWarning, message,
389 "Returning candidate distance.");
390 return dist;
391 }
392 }
393 } // Loop checking, 13.08.2015, G.Cosmo
394 }
395 }
396
397 return dist ;
398}
399
400////////////////////////////////////////////////////////
401//
402// Approximate nearest distance from the point p to the intersection of
403// two solids. It is usually underestimated from the point of view of
404// isotropic safety
405
408{
409 G4double dist = 0.0;
410
411#ifdef G4BOOLDEBUG
412 if( Inside(p) == kInside )
413 {
414 G4cout << "WARNING - Invalid call in "
415 << "G4SubtractionSolid::DistanceToIn(p)" << G4endl
416 << " Point p is inside !" << G4endl;
417 G4cout << " p = " << p << G4endl;
418 G4cerr << "WARNING - Invalid call in "
419 << "G4SubtractionSolid::DistanceToIn(p)" << G4endl
420 << " Point p is inside !" << G4endl;
421 G4cerr << " p = " << p << G4endl;
422 }
423#endif
424
425 if( ( fPtrSolidA->Inside(p) != kOutside) && // case 1
426 ( fPtrSolidB->Inside(p) != kOutside) )
427 {
428 dist = fPtrSolidB->DistanceToOut(p);
429 }
430 else
431 {
432 dist = fPtrSolidA->DistanceToIn(p);
433 }
434
435 return dist;
436}
437
438//////////////////////////////////////////////////////////
439//
440// The same algorithm as DistanceToOut(p)
441
444 const G4ThreeVector& v,
445 const G4bool calcNorm,
446 G4bool *validNorm,
447 G4ThreeVector *n ) const
448{
449#ifdef G4BOOLDEBUG
450 if( Inside(p) == kOutside )
451 {
452 G4cout << "Position:" << G4endl << G4endl;
453 G4cout << "p.x() = " << p.x()/mm << " mm" << G4endl;
454 G4cout << "p.y() = " << p.y()/mm << " mm" << G4endl;
455 G4cout << "p.z() = " << p.z()/mm << " mm" << G4endl << G4endl;
456 G4cout << "Direction:" << G4endl << G4endl;
457 G4cout << "v.x() = " << v.x() << G4endl;
458 G4cout << "v.y() = " << v.y() << G4endl;
459 G4cout << "v.z() = " << v.z() << G4endl << G4endl;
460 G4cout << "WARNING - Invalid call in "
461 << "G4SubtractionSolid::DistanceToOut(p,v)" << G4endl
462 << " Point p is outside !" << G4endl;
463 G4cout << " p = " << p << G4endl;
464 G4cout << " v = " << v << G4endl;
465 G4cerr << "WARNING - Invalid call in "
466 << "G4SubtractionSolid::DistanceToOut(p,v)" << G4endl
467 << " Point p is outside !" << G4endl;
468 G4cerr << " p = " << p << G4endl;
469 G4cerr << " v = " << v << G4endl;
470 }
471#endif
472
473 G4double distout;
474 G4double distA = fPtrSolidA->DistanceToOut(p,v,calcNorm,validNorm,n) ;
475 G4double distB = fPtrSolidB->DistanceToIn(p,v) ;
476 if(distB < distA)
477 {
478 if(calcNorm)
479 {
480 *n = -(fPtrSolidB->SurfaceNormal(p+distB*v)) ;
481 *validNorm = false ;
482 }
483 distout= distB ;
484 }
485 else
486 {
487 distout= distA ;
488 }
489 return distout;
490}
491
492//////////////////////////////////////////////////////////////
493//
494// Inverted algorithm of DistanceToIn(p)
495
498{
499 G4double dist=0.0;
500
501 if( Inside(p) == kOutside )
502 {
503#ifdef G4BOOLDEBUG
504 G4cout << "WARNING - Invalid call in "
505 << "G4SubtractionSolid::DistanceToOut(p)" << G4endl
506 << " Point p is outside" << G4endl;
507 G4cout << " p = " << p << G4endl;
508 G4cerr << "WARNING - Invalid call in "
509 << "G4SubtractionSolid::DistanceToOut(p)" << G4endl
510 << " Point p is outside" << G4endl;
511 G4cerr << " p = " << p << G4endl;
512#endif
513 }
514 else
515 {
516 dist= std::min(fPtrSolidA->DistanceToOut(p),
517 fPtrSolidB->DistanceToIn(p) ) ;
518 }
519 return dist;
520}
521
522//////////////////////////////////////////////////////////////
523//
524//
525
527{
528 return G4String("G4SubtractionSolid");
529}
530
531//////////////////////////////////////////////////////////////////////////
532//
533// Make a clone of the object
534
536{
537 return new G4SubtractionSolid(*this);
538}
539
540//////////////////////////////////////////////////////////////
541//
542// ComputeDimensions
543
544void
546 const G4int,
547 const G4VPhysicalVolume* )
548{
549}
550
551/////////////////////////////////////////////////
552//
553// DescribeYourselfTo
554
555void
557{
558 scene.AddSolid (*this);
559}
560
561////////////////////////////////////////////////////
562//
563// CreatePolyhedron
564
567{
568 HepPolyhedronProcessor processor;
569 // Stack components and components of components recursively
570 // See G4BooleanSolid::StackPolyhedron
571 G4Polyhedron* top = StackPolyhedron(processor, this);
572 G4Polyhedron* result = new G4Polyhedron(*top);
573 if (processor.execute(*result)) { return result; }
574 else { return nullptr; }
575}
576
577////////////////////////////////////////////////////
578//
579// GetCubicVolume
580//
581
583{
584 if( fCubicVolume != -1.0 ) {
585 return fCubicVolume;
586 }
587
588 G4double cubVolumeA = fPtrSolidA->GetCubicVolume();
589
590 G4ThreeVector bminA, bmaxA, bminB, bmaxB;
591 fPtrSolidA->BoundingLimits(bminA, bmaxA);
592 fPtrSolidB->BoundingLimits(bminB, bmaxB);
593 G4double intersection = 0.;
594 G4bool canIntersect =
595 bminA.x() < bmaxB.x() && bminA.y() < bmaxB.y() && bminA.z() < bmaxB.z() &&
596 bminB.x() < bmaxA.x() && bminB.y() < bmaxA.y() && bminB.z() < bmaxA.z();
597 if ( canIntersect )
598 {
599 G4IntersectionSolid intersectVol( "Temporary-Intersection-for-Union",
601 intersection = intersectVol.GetCubicVolume();
602 }
603
604 fCubicVolume = cubVolumeA - intersection;
605 if (fCubicVolume < 0.01*cubVolumeA) fCubicVolume = G4VSolid::GetCubicVolume();
606
607 return fCubicVolume;
608}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double z() const
double x() const
double y() const
G4VSolid * fPtrSolidA
G4BooleanSolid & operator=(const G4BooleanSolid &rhs)
G4VSolid * fPtrSolidB
G4double fCubicVolume
G4Polyhedron * StackPolyhedron(HepPolyhedronProcessor &, const G4VSolid *) const
virtual G4double GetCubicVolume()
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
virtual G4double GetCubicVolume() final
G4SubtractionSolid & operator=(const G4SubtractionSolid &rhs)
void DescribeYourselfTo(G4VGraphicsScene &scene) const
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
G4Polyhedron * CreatePolyhedron() const
G4GeometryType GetEntityType() const
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const
G4SubtractionSolid(const G4String &pName, G4VSolid *pSolidA, G4VSolid *pSolidB)
EInside Inside(const G4ThreeVector &p) const
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
G4VSolid * Clone() const
virtual void AddSolid(const G4Box &)=0
G4String GetName() const
virtual G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const =0
virtual EInside Inside(const G4ThreeVector &p) const =0
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const =0
void DumpInfo() const
G4double kCarTolerance
Definition: G4VSolid.hh:299
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const =0
virtual void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition: G4VSolid.cc:665
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0
virtual G4double GetCubicVolume()
Definition: G4VSolid.cc:188
virtual G4GeometryType GetEntityType() const =0
bool execute(HepPolyhedron &)
EAxis
Definition: geomdefs.hh:54
EInside
Definition: geomdefs.hh:67
@ kInside
Definition: geomdefs.hh:70
@ kOutside
Definition: geomdefs.hh:68
@ kSurface
Definition: geomdefs.hh:69