Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4IntersectionSolid.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// 12.09.98 V.Grichine: first implementation
29// --------------------------------------------------------------------
30
31#include <sstream>
32
34
35#include "G4SystemOfUnits.hh"
36#include "G4VoxelLimits.hh"
38
39#include "G4VGraphicsScene.hh"
40#include "G4Polyhedron.hh"
43
44//////////////////////////////////////////////////////////////////////////
45//
46// Transfer all data members to G4BooleanSolid which is responsible
47// for them. pName will be in turn sent to G4VSolid
48//
49
51 G4VSolid* pSolidA ,
52 G4VSolid* pSolidB )
53 : G4BooleanSolid(pName,pSolidA,pSolidB)
54{
55}
56
57//////////////////////////////////////////////////////////////////////////
58//
59
61 G4VSolid* pSolidA,
62 G4VSolid* pSolidB,
63 G4RotationMatrix* rotMatrix,
64 const G4ThreeVector& transVector )
65 : G4BooleanSolid(pName,pSolidA,pSolidB,rotMatrix,transVector)
66{
67}
68
69//////////////////////////////////////////////////////////////////////////
70//
71//
72
74 G4VSolid* pSolidA,
75 G4VSolid* pSolidB,
76 const G4Transform3D& transform )
77 : G4BooleanSolid(pName,pSolidA,pSolidB,transform)
78{
79}
80
81//////////////////////////////////////////////////////////////////////////
82//
83// Fake default constructor - sets only member data and allocates memory
84// for usage restricted to object persistency.
85
90
91//////////////////////////////////////////////////////////////////////////
92//
93// Copy constructor
94
96
97//////////////////////////////////////////////////////////////////////////
98//
99// Assignment operator
100
103{
104 // Check assignment to self
105 //
106 if (this == &rhs) { return *this; }
107
108 // Copy base class data
109 //
111
112 return *this;
113}
114
115//////////////////////////////////////////////////////////////////////////
116//
117// Get bounding box
118
119void
121 G4ThreeVector& pMax) const
122{
123 G4ThreeVector minA,maxA, minB,maxB;
124 fPtrSolidA->BoundingLimits(minA,maxA);
125 fPtrSolidB->BoundingLimits(minB,maxB);
126
127 pMin.set(std::max(minA.x(),minB.x()),
128 std::max(minA.y(),minB.y()),
129 std::max(minA.z(),minB.z()));
130
131 pMax.set(std::min(maxA.x(),maxB.x()),
132 std::min(maxA.y(),maxB.y()),
133 std::min(maxA.z(),maxB.z()));
134
135 // Check correctness of the bounding box
136 //
137 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
138 {
139 std::ostringstream message;
140 message << "Bad bounding box (min >= max) for solid: "
141 << GetName() << " !"
142 << "\npMin = " << pMin
143 << "\npMax = " << pMax;
144 G4Exception("G4IntersectionSolid::BoundingLimits()", "GeomMgt0001",
145 JustWarning, message);
146 DumpInfo();
147 }
148}
149
150//////////////////////////////////////////////////////////////////////////
151//
152// Calculate extent under transform and specified limit
153
154G4bool
156 const G4VoxelLimits& pVoxelLimit,
157 const G4AffineTransform& pTransform,
158 G4double& pMin,
159 G4double& pMax) const
160{
161 G4bool retA, retB, out;
162 G4double minA, minB, maxA, maxB;
163
164 retA = fPtrSolidA
165 ->CalculateExtent( pAxis, pVoxelLimit, pTransform, minA, maxA);
166 retB = fPtrSolidB
167 ->CalculateExtent( pAxis, pVoxelLimit, pTransform, minB, maxB);
168
169 if( retA && retB )
170 {
171 pMin = std::max( minA, minB );
172 pMax = std::min( maxA, maxB );
173 out = (pMax > pMin); // true;
174 }
175 else
176 {
177 out = false;
178 }
179
180 return out; // It exists in this slice only if both exist in it.
181}
182
183//////////////////////////////////////////////////////////////////////////
184//
185// Touching ? Empty intersection ?
186
188{
189 EInside positionA = fPtrSolidA->Inside(p);
190 if(positionA == kOutside) { return positionA; } // outside A
191
192 EInside positionB = fPtrSolidB->Inside(p);
193 if(positionA == kInside) { return positionB; }
194
195 if(positionB == kOutside) { return positionB; } // outside B
196
197 return kSurface; // surface A & B
198}
199
200//////////////////////////////////////////////////////////////////////////
201//
202
205{
206 G4ThreeVector normal;
207 EInside insideA, insideB;
208
209 insideA = fPtrSolidA->Inside(p);
210 insideB = fPtrSolidB->Inside(p);
211
212#ifdef G4BOOLDEBUG
213 if( (insideA == kOutside) || (insideB == kOutside) )
214 {
215 G4cout << "WARNING - Invalid call in "
216 << "G4IntersectionSolid::SurfaceNormal(p)" << G4endl
217 << " Point p is outside !" << G4endl;
218 G4cout << " p = " << p << G4endl;
219 G4cerr << "WARNING - Invalid call in "
220 << "G4IntersectionSolid::SurfaceNormal(p)" << G4endl
221 << " Point p is outside !" << G4endl;
222 G4cerr << " p = " << p << G4endl;
223 }
224#endif
225
226 // On the surface of both is difficult ... treat it like on A now!
227 //
228 if( insideA == kSurface )
229 {
230 normal = fPtrSolidA->SurfaceNormal(p) ;
231 }
232 else if( insideB == kSurface )
233 {
234 normal = fPtrSolidB->SurfaceNormal(p) ;
235 }
236 else // We are on neither surface, so we should generate an exception
237 {
238 if(fPtrSolidA->DistanceToOut(p) <= fPtrSolidB->DistanceToOut(p) )
239 {
240 normal= fPtrSolidA->SurfaceNormal(p) ;
241 }
242 else
243 {
244 normal= fPtrSolidB->SurfaceNormal(p) ;
245 }
246#ifdef G4BOOLDEBUG
247 G4cout << "WARNING - Invalid call in "
248 << "G4IntersectionSolid::SurfaceNormal(p)" << G4endl
249 << " Point p is out of surface !" << G4endl;
250 G4cout << " p = " << p << G4endl;
251 G4cerr << "WARNING - Invalid call in "
252 << "G4IntersectionSolid::SurfaceNormal(p)" << G4endl
253 << " Point p is out of surface !" << G4endl;
254 G4cerr << " p = " << p << G4endl;
255#endif
256 }
257
258 return normal;
259}
260
261//////////////////////////////////////////////////////////////////////////
262//
263// The same algorithm as in DistanceToIn(p)
264
267 const G4ThreeVector& v ) const
268{
269 G4double dist = 0.0;
270 if( Inside(p) == kInside )
271 {
272#ifdef G4BOOLDEBUG
273 G4cout << "WARNING - Invalid call in "
274 << "G4IntersectionSolid::DistanceToIn(p,v)" << G4endl
275 << " Point p is inside !" << G4endl;
276 G4cout << " p = " << p << G4endl;
277 G4cout << " v = " << v << G4endl;
278 G4cerr << "WARNING - Invalid call in "
279 << "G4IntersectionSolid::DistanceToIn(p,v)" << G4endl
280 << " Point p is inside !" << G4endl;
281 G4cerr << " p = " << p << G4endl;
282 G4cerr << " v = " << v << G4endl;
283#endif
284 }
285 else // if( Inside(p) == kSurface )
286 {
287 EInside wA = fPtrSolidA->Inside(p);
288 EInside wB = fPtrSolidB->Inside(p);
289
290 G4ThreeVector pA = p, pB = p;
291 G4double dA = 0., dA1=0., dA2=0.;
292 G4double dB = 0., dB1=0., dB2=0.;
293 G4bool doA = true, doB = true;
294
295 static const std::size_t max_trials=10000;
296 for (std::size_t trial=0; trial<max_trials; ++trial)
297 {
298 if(doA)
299 {
300 // find next valid range for A
301
302 dA1 = 0.;
303
304 if( wA != kInside )
305 {
306 dA1 = fPtrSolidA->DistanceToIn(pA, v);
307
308 if( dA1 == kInfinity ) { return kInfinity; }
309
310 pA += dA1*v;
311 }
312 dA2 = dA1 + fPtrSolidA->DistanceToOut(pA, v);
313 }
314 dA1 += dA;
315 dA2 += dA;
316
317 if(doB)
318 {
319 // find next valid range for B
320
321 dB1 = 0.;
322 if(wB != kInside)
323 {
324 dB1 = fPtrSolidB->DistanceToIn(pB, v);
325
326 if(dB1 == kInfinity) { return kInfinity; }
327
328 pB += dB1*v;
329 }
330 dB2 = dB1 + fPtrSolidB->DistanceToOut(pB, v);
331 }
332 dB1 += dB;
333 dB2 += dB;
334
335 // check if they overlap
336
337 if( dA1 < dB1 )
338 {
339 if( dB1 < dA2 ) { return dB1; }
340
341 dA = dA2;
342 pA = p + dA*v; // continue from here
343 wA = kSurface;
344 doA = true;
345 doB = false;
346 }
347 else
348 {
349 if( dA1 < dB2 ) { return dA1; }
350
351 dB = dB2;
352 pB = p + dB*v; // continue from here
353 wB = kSurface;
354 doB = true;
355 doA = false;
356 }
357 }
358 }
359#ifdef G4BOOLDEBUG
360 G4Exception("G4IntersectionSolid::DistanceToIn(p,v)",
361 "GeomSolids0001", JustWarning,
362 "Reached maximum number of iterations! Returning zero.");
363#endif
364 return dist ;
365}
366
367//////////////////////////////////////////////////////////////////////////
368//
369// Approximate nearest distance from the point p to the intersection of
370// two solids
371
374{
375#ifdef G4BOOLDEBUG
376 if( Inside(p) == kInside )
377 {
378 G4cout << "WARNING - Invalid call in "
379 << "G4IntersectionSolid::DistanceToIn(p)" << G4endl
380 << " Point p is inside !" << G4endl;
381 G4cout << " p = " << p << G4endl;
382 G4cerr << "WARNING - Invalid call in "
383 << "G4IntersectionSolid::DistanceToIn(p)" << G4endl
384 << " Point p is inside !" << G4endl;
385 G4cerr << " p = " << p << G4endl;
386 }
387#endif
388 EInside sideA = fPtrSolidA->Inside(p) ;
389 EInside sideB = fPtrSolidB->Inside(p) ;
390 G4double dist=0.0 ;
391
392 if( sideA != kInside && sideB != kOutside )
393 {
394 dist = fPtrSolidA->DistanceToIn(p) ;
395 }
396 else
397 {
398 if( sideB != kInside && sideA != kOutside )
399 {
400 dist = fPtrSolidB->DistanceToIn(p) ;
401 }
402 else
403 {
404 dist = std::min(fPtrSolidA->DistanceToIn(p),
405 fPtrSolidB->DistanceToIn(p) ) ;
406 }
407 }
408 return dist ;
409}
410
411//////////////////////////////////////////////////////////////////////////
412//
413// The same algorithm as DistanceToOut(p)
414
417 const G4ThreeVector& v,
418 const G4bool calcNorm,
419 G4bool* validNorm,
420 G4ThreeVector* n ) const
421{
422 G4bool validNormA, validNormB;
423 G4ThreeVector nA, nB;
424
425#ifdef G4BOOLDEBUG
426 if( Inside(p) == kOutside )
427 {
428 G4cout << "Position:" << G4endl << G4endl;
429 G4cout << "p.x() = " << p.x()/mm << " mm" << G4endl;
430 G4cout << "p.y() = " << p.y()/mm << " mm" << G4endl;
431 G4cout << "p.z() = " << p.z()/mm << " mm" << G4endl << G4endl;
432 G4cout << "Direction:" << G4endl << G4endl;
433 G4cout << "v.x() = " << v.x() << G4endl;
434 G4cout << "v.y() = " << v.y() << G4endl;
435 G4cout << "v.z() = " << v.z() << G4endl << G4endl;
436 G4cout << "WARNING - Invalid call in "
437 << "G4IntersectionSolid::DistanceToOut(p,v)" << G4endl
438 << " Point p is outside !" << G4endl;
439 G4cout << " p = " << p << G4endl;
440 G4cout << " v = " << v << G4endl;
441 G4cerr << "WARNING - Invalid call in "
442 << "G4IntersectionSolid::DistanceToOut(p,v)" << G4endl
443 << " Point p is outside !" << G4endl;
444 G4cerr << " p = " << p << G4endl;
445 G4cerr << " v = " << v << G4endl;
446 }
447#endif
448 G4double distA = fPtrSolidA->DistanceToOut(p,v,calcNorm,&validNormA,&nA) ;
449 G4double distB = fPtrSolidB->DistanceToOut(p,v,calcNorm,&validNormB,&nB) ;
450
451 G4double dist = std::min(distA,distB) ;
452
453 if( calcNorm )
454 {
455 if ( distA < distB )
456 {
457 *validNorm = validNormA;
458 *n = nA;
459 }
460 else
461 {
462 *validNorm = validNormB;
463 *n = nB;
464 }
465 }
466
467 return dist ;
468}
469
470//////////////////////////////////////////////////////////////////////////
471//
472// Inverted algorithm of DistanceToIn(p)
473
476{
477#ifdef G4BOOLDEBUG
478 if( Inside(p) == kOutside )
479 {
480 G4cout << "WARNING - Invalid call in "
481 << "G4IntersectionSolid::DistanceToOut(p)" << G4endl
482 << " Point p is outside !" << G4endl;
483 G4cout << " p = " << p << G4endl;
484 G4cerr << "WARNING - Invalid call in "
485 << "G4IntersectionSolid::DistanceToOut(p)" << G4endl
486 << " Point p is outside !" << G4endl;
487 G4cerr << " p = " << p << G4endl;
488 }
489#endif
490
491 return std::min(fPtrSolidA->DistanceToOut(p),
492 fPtrSolidB->DistanceToOut(p) ) ;
493
494}
495
496//////////////////////////////////////////////////////////////////////////
497//
498// ComputeDimensions
499
500void
502 const G4int,
503 const G4VPhysicalVolume* )
504{
505 DumpInfo();
506 G4Exception("G4IntersectionSolid::ComputeDimensions()",
507 "GeomSolids0001", FatalException,
508 "Method not applicable in this context!");
509}
510
511//////////////////////////////////////////////////////////////////////////
512//
513// GetEntityType
514
516{
517 return {"G4IntersectionSolid"};
518}
519
520//////////////////////////////////////////////////////////////////////////
521//
522// Make a clone of the object
523
525{
526 return new G4IntersectionSolid(*this);
527}
528
529//////////////////////////////////////////////////////////////////////////
530//
531// DescribeYourselfTo
532
533void
535{
536 scene.AddSolid (*this);
537}
538
539//////////////////////////////////////////////////////////////////////////
540//
541// CreatePolyhedron
542
545{
546 if (fExternalBoolProcessor == nullptr)
547 {
548 HepPolyhedronProcessor processor;
549 // Stack components and components of components recursively
550 // See G4BooleanSolid::StackPolyhedron
551 G4Polyhedron* top = StackPolyhedron(processor, this);
552 auto result = new G4Polyhedron(*top);
553 if (processor.execute(*result))
554 {
555 return result;
556 }
557 return nullptr;
558 }
559 else
560 {
561 return fExternalBoolProcessor->Process(this);
562 }
563}
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
CLHEP::HepRotation G4RotationMatrix
CLHEP::Hep3Vector G4ThreeVector
HepGeom::Transform3D G4Transform3D
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4String G4GeometryType
Definition G4VSolid.hh:70
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
double z() const
double x() const
double y() const
void set(double x, double y, double z)
G4AffineTransform is a class for geometric affine transformations. It supports efficient arbitrary ro...
G4VSolid * fPtrSolidA
G4BooleanSolid & operator=(const G4BooleanSolid &rhs)
G4VSolid * fPtrSolidB
G4BooleanSolid(const G4String &pName, G4VSolid *pSolidA, G4VSolid *pSolidB)
static G4VBooleanProcessor * fExternalBoolProcessor
G4Polyhedron * StackPolyhedron(HepPolyhedronProcessor &, const G4VSolid *) const
G4IntersectionSolid is a solid describing the Boolean intersection of two solids.
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const override
G4GeometryType GetEntityType() const override
G4Polyhedron * CreatePolyhedron() const override
G4IntersectionSolid & operator=(const G4IntersectionSolid &rhs)
G4IntersectionSolid(const G4String &pName, G4VSolid *pSolidA, G4VSolid *pSolidB)
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const override
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const override
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const override
void DescribeYourselfTo(G4VGraphicsScene &scene) const override
EInside Inside(const G4ThreeVector &p) const override
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const override
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep) override
G4VSolid * Clone() const override
virtual void AddSolid(const G4Box &)=0
G4VPVParameterisation ia an abstract base class for Parameterisation, able to compute the transformat...
G4VPhysicalVolume is an abstract base class for the representation of a positioned volume....
G4String GetName() const
G4VSolid(const G4String &name)
Definition G4VSolid.cc:59
void DumpInfo() const
G4VoxelLimits represents limitation/restrictions of space, where restrictions are only made perpendic...
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