Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4UnionSolid.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 G4UnionSolid
27//
28// 23.04.18 E.Tcherniaev: added extended BBox, yearly return in Inside()
29// 17.03.17 E.Tcherniaev: revision of SurfaceNormal()
30// 12.09.98 V.Grichine: first implementation
31// --------------------------------------------------------------------
32
33#include <sstream>
34
35#include "G4UnionSolid.hh"
36
37#include "G4SystemOfUnits.hh"
38#include "G4VoxelLimits.hh"
41
42#include "G4VGraphicsScene.hh"
43#include "G4Polyhedron.hh"
45
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 Init();
59}
60
61/////////////////////////////////////////////////////////////////////
62//
63// Constructor
64
66 G4VSolid* pSolidA ,
67 G4VSolid* pSolidB ,
68 G4RotationMatrix* rotMatrix,
69 const G4ThreeVector& transVector )
70 : G4BooleanSolid(pName,pSolidA,pSolidB,rotMatrix,transVector)
71
72{
73 Init();
74}
75
76///////////////////////////////////////////////////////////
77//
78// Constructor
79
81 G4VSolid* pSolidA ,
82 G4VSolid* pSolidB ,
83 const G4Transform3D& transform )
84 : G4BooleanSolid(pName,pSolidA,pSolidB,transform)
85{
86 Init();
87}
88
89//////////////////////////////////////////////////////////////////
90//
91// Fake default constructor - sets only member data and allocates memory
92// for usage restricted to object persistency.
93
96{
97}
98
99///////////////////////////////////////////////////////////
100//
101// Destructor
102
104{
105}
106
107///////////////////////////////////////////////////////////////
108//
109// Copy constructor
110
112 : G4BooleanSolid (rhs)
113{
114 fPMin = rhs.fPMin;
115 fPMax = rhs.fPMax;
116 halfCarTolerance=0.5*kCarTolerance;
117}
118
119///////////////////////////////////////////////////////////////
120//
121// Assignment operator
122
124{
125 // Check assignment to self
126 //
127 if (this == &rhs) { return *this; }
128
129 // Copy base class data
130 //
132
133 fPMin = rhs.fPMin;
134 fPMax = rhs.fPMax;
135 halfCarTolerance = rhs.halfCarTolerance;
136
137 return *this;
138}
139
140///////////////////////////////////////////////////////////////
141//
142// Initialisation
143
144void G4UnionSolid::Init()
145{
147 G4ThreeVector pmin, pmax;
148 BoundingLimits(pmin, pmax);
149 fPMin = pmin - pdelta;
150 fPMax = pmax + pdelta;
151 halfCarTolerance=0.5*kCarTolerance;
152}
153
154//////////////////////////////////////////////////////////////////////////
155//
156// Get bounding box
157
159 G4ThreeVector& pMax) const
160{
161 G4ThreeVector minA,maxA, minB,maxB;
162 fPtrSolidA->BoundingLimits(minA,maxA);
163 fPtrSolidB->BoundingLimits(minB,maxB);
164
165 pMin.set(std::min(minA.x(),minB.x()),
166 std::min(minA.y(),minB.y()),
167 std::min(minA.z(),minB.z()));
168
169 pMax.set(std::max(maxA.x(),maxB.x()),
170 std::max(maxA.y(),maxB.y()),
171 std::max(maxA.z(),maxB.z()));
172
173 // Check correctness of the bounding box
174 //
175 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
176 {
177 std::ostringstream message;
178 message << "Bad bounding box (min >= max) for solid: "
179 << GetName() << " !"
180 << "\npMin = " << pMin
181 << "\npMax = " << pMax;
182 G4Exception("G4UnionSolid::BoundingLimits()", "GeomMgt0001",
183 JustWarning, message);
184 DumpInfo();
185 }
186}
187
188//////////////////////////////////////////////////////////////////////////
189//
190// Calculate extent under transform and specified limit
191
192G4bool
194 const G4VoxelLimits& pVoxelLimit,
195 const G4AffineTransform& pTransform,
196 G4double& pMin,
197 G4double& pMax ) const
198{
199 G4bool touchesA, touchesB, out ;
200 G4double minA = kInfinity, minB = kInfinity,
201 maxA = -kInfinity, maxB = -kInfinity;
202
203 touchesA = fPtrSolidA->CalculateExtent( pAxis, pVoxelLimit,
204 pTransform, minA, maxA);
205 touchesB = fPtrSolidB->CalculateExtent( pAxis, pVoxelLimit,
206 pTransform, minB, maxB);
207 if( touchesA || touchesB )
208 {
209 pMin = std::min( minA, minB );
210 pMax = std::max( maxA, maxB );
211 out = true ;
212 }
213 else
214 {
215 out = false ;
216 }
217
218 return out ; // It exists in this slice if either one does.
219}
220
221/////////////////////////////////////////////////////
222//
223// Important comment: When solids A and B touch together along flat
224// surface the surface points will be considered as kSurface, while points
225// located around will correspond to kInside
226
228{
229 if (std::max(p.z()-fPMax.z(), fPMin.z()-p.z()) > 0) { return kOutside; }
230
231 EInside positionA = fPtrSolidA->Inside(p);
232 if (positionA == kInside) { return positionA; } // inside A
233 EInside positionB = fPtrSolidB->Inside(p);
234 if (positionA == kOutside) { return positionB; }
235
236 if (positionB == kInside) { return positionB; } // inside B
237 if (positionB == kOutside) { return positionA; } // surface A
238
239 // Both points are on surface
240 //
241 static const G4double rtol
243
244 return ((fPtrSolidA->SurfaceNormal(p) +
245 fPtrSolidB->SurfaceNormal(p)).mag2() < rtol) ? kInside : kSurface;
246}
247
248//////////////////////////////////////////////////////////////
249//
250// Get surface normal
251
254{
255 EInside positionA = fPtrSolidA->Inside(p);
256 EInside positionB = fPtrSolidB->Inside(p);
257
258 if (positionA == kSurface &&
259 positionB == kOutside) return fPtrSolidA->SurfaceNormal(p);
260
261 if (positionA == kOutside &&
262 positionB == kSurface) return fPtrSolidB->SurfaceNormal(p);
263
264 if (positionA == kSurface &&
265 positionB == kSurface)
266 {
267 if (Inside(p) == kSurface)
268 {
271 return (normalA + normalB).unit();
272 }
273 }
274#ifdef G4BOOLDEBUG
275 G4String surf[3] = { "OUTSIDE", "SURFACE", "INSIDE" };
276 std::ostringstream message;
277 G4int oldprc = message.precision(16);
278 message << "Invalid call of SurfaceNormal(p) for union solid: "
279 << GetName() << " !"
280 << "\nPoint p" << p << " is " << surf[Inside(p)] << " !!!";
281 message.precision(oldprc);
282 G4Exception("G4UnionSolid::SurfaceNormal()", "GeomMgt0001",
283 JustWarning, message);
284#endif
285 return fPtrSolidA->SurfaceNormal(p);
286}
287
288/////////////////////////////////////////////////////////////
289//
290// The same algorithm as in DistanceToIn(p)
291
294 const G4ThreeVector& v ) const
295{
296#ifdef G4BOOLDEBUG
297 if( Inside(p) == kInside )
298 {
299 G4cout << "WARNING - Invalid call in "
300 << "G4UnionSolid::DistanceToIn(p,v)" << G4endl
301 << " Point p is inside !" << G4endl;
302 G4cout << " p = " << p << G4endl;
303 G4cout << " v = " << v << G4endl;
304 G4cerr << "WARNING - Invalid call in "
305 << "G4UnionSolid::DistanceToIn(p,v)" << G4endl
306 << " Point p is inside !" << G4endl;
307 G4cerr << " p = " << p << G4endl;
308 G4cerr << " v = " << v << G4endl;
309 }
310#endif
311
312 return std::min(fPtrSolidA->DistanceToIn(p,v),
313 fPtrSolidB->DistanceToIn(p,v) ) ;
314}
315
316////////////////////////////////////////////////////////
317//
318// Approximate nearest distance from the point p to the union of
319// two solids
320
323{
324#ifdef G4BOOLDEBUG
325 if( Inside(p) == kInside )
326 {
327 G4cout << "WARNING - Invalid call in "
328 << "G4UnionSolid::DistanceToIn(p)" << G4endl
329 << " Point p is inside !" << G4endl;
330 G4cout << " p = " << p << G4endl;
331 G4cerr << "WARNING - Invalid call in "
332 << "G4UnionSolid::DistanceToIn(p)" << G4endl
333 << " Point p is inside !" << G4endl;
334 G4cerr << " p = " << p << G4endl;
335 }
336#endif
337 G4double distA = fPtrSolidA->DistanceToIn(p) ;
338 G4double distB = fPtrSolidB->DistanceToIn(p) ;
339 G4double safety = std::min(distA,distB) ;
340 if(safety < 0.0) safety = 0.0 ;
341 return safety ;
342}
343
344//////////////////////////////////////////////////////////
345//
346// The same algorithm as DistanceToOut(p)
347
350 const G4ThreeVector& v,
351 const G4bool calcNorm,
352 G4bool *validNorm,
353 G4ThreeVector *n ) const
354{
355 G4double dist = 0.0, disTmp = 0.0 ;
356 G4ThreeVector normTmp;
357 G4ThreeVector* nTmp = &normTmp;
358
359 if( Inside(p) == kOutside )
360 {
361#ifdef G4BOOLDEBUG
362 G4cout << "Position:" << G4endl << G4endl;
363 G4cout << "p.x() = " << p.x()/mm << " mm" << G4endl;
364 G4cout << "p.y() = " << p.y()/mm << " mm" << G4endl;
365 G4cout << "p.z() = " << p.z()/mm << " mm" << G4endl << G4endl;
366 G4cout << "Direction:" << G4endl << G4endl;
367 G4cout << "v.x() = " << v.x() << G4endl;
368 G4cout << "v.y() = " << v.y() << G4endl;
369 G4cout << "v.z() = " << v.z() << G4endl << G4endl;
370 G4cout << "WARNING - Invalid call in "
371 << "G4UnionSolid::DistanceToOut(p,v)" << G4endl
372 << " Point p is outside !" << G4endl;
373 G4cout << " p = " << p << G4endl;
374 G4cout << " v = " << v << G4endl;
375 G4cerr << "WARNING - Invalid call in "
376 << "G4UnionSolid::DistanceToOut(p,v)" << G4endl
377 << " Point p is outside !" << G4endl;
378 G4cerr << " p = " << p << G4endl;
379 G4cerr << " v = " << v << G4endl;
380#endif
381 }
382 else
383 {
384 EInside positionA = fPtrSolidA->Inside(p) ;
385
386 if( positionA != kOutside )
387 {
388 do // Loop checking, 13.08.2015, G.Cosmo
389 {
390 disTmp = fPtrSolidA->DistanceToOut(p+dist*v,v,calcNorm,
391 validNorm,nTmp);
392 dist += disTmp ;
393
394 if(fPtrSolidB->Inside(p+dist*v) != kOutside)
395 {
396 disTmp = fPtrSolidB->DistanceToOut(p+dist*v,v,calcNorm,
397 validNorm,nTmp);
398 dist += disTmp ;
399 }
400 }
401 while( (fPtrSolidA->Inside(p+dist*v) != kOutside)
402 && (disTmp > halfCarTolerance) );
403 }
404 else // if( positionB != kOutside )
405 {
406 do // Loop checking, 13.08.2015, G.Cosmo
407 {
408 disTmp = fPtrSolidB->DistanceToOut(p+dist*v,v,calcNorm,
409 validNorm,nTmp);
410 dist += disTmp ;
411
412 if(fPtrSolidA->Inside(p+dist*v) != kOutside)
413 {
414 disTmp = fPtrSolidA->DistanceToOut(p+dist*v,v,calcNorm,
415 validNorm,nTmp);
416 dist += disTmp ;
417 }
418 }
419 while( (fPtrSolidB->Inside(p+dist*v) != kOutside)
420 && (disTmp > halfCarTolerance) );
421 }
422 }
423 if( calcNorm )
424 {
425 *validNorm = false ;
426 *n = *nTmp ;
427 }
428 return dist ;
429}
430
431//////////////////////////////////////////////////////////////
432//
433// Inverted algorithm of DistanceToIn(p)
434
437{
438 G4double distout = 0.0;
439 if( Inside(p) == kOutside )
440 {
441#ifdef G4BOOLDEBUG
442 G4cout << "WARNING - Invalid call in "
443 << "G4UnionSolid::DistanceToOut(p)" << G4endl
444 << " Point p is outside !" << G4endl;
445 G4cout << " p = " << p << G4endl;
446 G4cerr << "WARNING - Invalid call in "
447 << "G4UnionSolid::DistanceToOut(p)" << G4endl
448 << " Point p is outside !" << G4endl;
449 G4cerr << " p = " << p << G4endl;
450#endif
451 }
452 else
453 {
454 EInside positionA = fPtrSolidA->Inside(p) ;
455 EInside positionB = fPtrSolidB->Inside(p) ;
456
457 // Is this equivalent ??
458 // if( ! ( (positionA == kOutside)) &&
459 // (positionB == kOutside)) )
460 if((positionA == kInside && positionB == kInside ) ||
461 (positionA == kInside && positionB == kSurface ) ||
462 (positionA == kSurface && positionB == kInside ) )
463 {
464 distout= std::max(fPtrSolidA->DistanceToOut(p),
466 }
467 else
468 {
469 if(positionA == kOutside)
470 {
471 distout= fPtrSolidB->DistanceToOut(p) ;
472 }
473 else
474 {
475 distout= fPtrSolidA->DistanceToOut(p) ;
476 }
477 }
478 }
479 return distout;
480}
481
482//////////////////////////////////////////////////////////////
483//
484// GetEntityType
485
487{
488 return G4String("G4UnionSolid");
489}
490
491//////////////////////////////////////////////////////////////////////////
492//
493// Make a clone of the object
494
496{
497 return new G4UnionSolid(*this);
498}
499
500//////////////////////////////////////////////////////////////
501//
502// ComputeDimensions
503
504void
506 const G4int,
507 const G4VPhysicalVolume* )
508{
509}
510
511/////////////////////////////////////////////////
512//
513// DescribeYourselfTo
514
515void
517{
518 scene.AddSolid (*this);
519}
520
521////////////////////////////////////////////////////
522//
523// CreatePolyhedron
524
527{
528 HepPolyhedronProcessor processor;
529 // Stack components and components of components recursively
530 // See G4BooleanSolid::StackPolyhedron
531 G4Polyhedron* top = StackPolyhedron(processor, this);
532 G4Polyhedron* result = new G4Polyhedron(*top);
533 if (processor.execute(*result)) { return result; }
534 else { return nullptr; }
535}
536
537
538////////////////////////////////////////////////////
539//
540// GetCubicVolume
541
543{
544 if( fCubicVolume != -1.0 ) {
545 return fCubicVolume;
546 }
547 G4double cubVolumeA = fPtrSolidA->GetCubicVolume();
548 G4double cubVolumeB = fPtrSolidB->GetCubicVolume();
549
550 G4ThreeVector bminA, bmaxA, bminB, bmaxB;
551 fPtrSolidA->BoundingLimits(bminA, bmaxA);
552 fPtrSolidB->BoundingLimits(bminB, bmaxB);
553
554 G4double intersection = 0.;
555 G4bool canIntersect =
556 bminA.x() < bmaxB.x() && bminA.y() < bmaxB.y() && bminA.z() < bmaxB.z() &&
557 bminB.x() < bmaxA.x() && bminB.y() < bmaxA.y() && bminB.z() < bmaxA.z();
558 if ( canIntersect )
559 {
560 G4IntersectionSolid intersectVol( "Temporary-Intersection-for-Union",
562 intersection = intersectVol.GetCubicVolume();
563 }
564
565 fCubicVolume = cubVolumeA + cubVolumeB - intersection;
566
567 return fCubicVolume;
568}
@ 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
Hep3Vector unit() const
double x() const
double y() const
void set(double x, double y, double z)
G4VSolid * fPtrSolidA
G4BooleanSolid & operator=(const G4BooleanSolid &rhs)
G4VSolid * fPtrSolidB
G4double fCubicVolume
G4Polyhedron * StackPolyhedron(HepPolyhedronProcessor &, const G4VSolid *) const
virtual G4double GetCubicVolume()
G4double GetRadialTolerance() const
static G4GeometryTolerance * GetInstance()
G4UnionSolid(const G4String &pName, G4VSolid *pSolidA, G4VSolid *pSolidB)
Definition: G4UnionSolid.cc:53
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const
virtual ~G4UnionSolid()
EInside Inside(const G4ThreeVector &p) const
G4UnionSolid & operator=(const G4UnionSolid &rhs)
G4GeometryType GetEntityType() const
G4Polyhedron * CreatePolyhedron() const
virtual G4double GetCubicVolume() final
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
void DescribeYourselfTo(G4VGraphicsScene &scene) const
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
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
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