Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BooleanSolid.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 the base class for solids created by Boolean
27// operations between other solids
28//
29// 1998.09.10 V.Grichine - created
30// --------------------------------------------------------------------
31
32#include "G4BooleanSolid.hh"
33#include "G4VSolid.hh"
34#include "G4DisplacedSolid.hh"
35#include "G4ReflectedSolid.hh"
36#include "G4ScaledSolid.hh"
37#include "G4Polyhedron.hh"
39#include "G4QuickRand.hh"
40
41#include "G4AutoLock.hh"
42
43namespace
44{
45 G4RecursiveMutex polyhedronMutex = G4MUTEX_INITIALIZER;
46 G4Mutex boolMutex = G4MUTEX_INITIALIZER;
47}
48
50
51//////////////////////////////////////////////////////////////////
52//
53// Constructor
54
56 G4VSolid* pSolidA ,
57 G4VSolid* pSolidB )
58 : G4VSolid(pName), fPtrSolidA(pSolidA), fPtrSolidB(pSolidB)
59{
60}
61
62//////////////////////////////////////////////////////////////////
63//
64// Constructor
65
67 G4VSolid* pSolidA ,
68 G4VSolid* pSolidB ,
69 G4RotationMatrix* rotMatrix,
70 const G4ThreeVector& transVector )
71 : G4VSolid(pName), createdDisplacedSolid(true)
72{
73 fPtrSolidA = pSolidA ;
74 fPtrSolidB = new G4DisplacedSolid("placedB",pSolidB,rotMatrix,transVector) ;
75}
76
77//////////////////////////////////////////////////////////////////
78//
79// Constructor
80
82 G4VSolid* pSolidA ,
83 G4VSolid* pSolidB ,
84 const G4Transform3D& transform )
85 : G4VSolid(pName), createdDisplacedSolid(true)
86{
87 fPtrSolidA = pSolidA ;
88 fPtrSolidB = new G4DisplacedSolid("placedB",pSolidB,transform) ;
89}
90
91///////////////////////////////////////////////////////////////
92//
93// Fake default constructor - sets only member data and allocates memory
94// for usage restricted to object persistency.
95
97 : G4VSolid(a)
98{
99}
100
101///////////////////////////////////////////////////////////////
102//
103// Destructor deletes transformation contents of the created displaced solid
104
106{
107 if(createdDisplacedSolid)
108 {
109 ((G4DisplacedSolid*)fPtrSolidB)->CleanTransformations();
110 }
111 delete fpPolyhedron; fpPolyhedron = nullptr;
112}
113
114///////////////////////////////////////////////////////////////
115//
116// Copy constructor
117
121 fCubVolStatistics(rhs.fCubVolStatistics),
122 fAreaStatistics(rhs.fAreaStatistics),
123 fCubVolEpsilon(rhs.fCubVolEpsilon),
124 fAreaAccuracy(rhs.fAreaAccuracy),
125 createdDisplacedSolid(rhs.createdDisplacedSolid)
126{
127 fPrimitives.resize(0); fPrimitivesSurfaceArea = 0.;
128}
129
130///////////////////////////////////////////////////////////////
131//
132// Assignment operator
133
135{
136 // Check assignment to self
137 //
138 if (this == &rhs) { return *this; }
139
140 // Copy base class data
141 //
143
144 // Copy data
145 //
148 fCubVolStatistics = rhs.fCubVolStatistics; fCubVolEpsilon = rhs.fCubVolEpsilon;
149 fAreaStatistics = rhs.fAreaStatistics; fAreaAccuracy = rhs.fAreaAccuracy;
150 createdDisplacedSolid= rhs.createdDisplacedSolid;
151
152 fRebuildPolyhedron = false;
153 delete fpPolyhedron; fpPolyhedron = nullptr;
154 fPrimitives.resize(0); fPrimitivesSurfaceArea = 0.;
155
156 return *this;
157}
158
159///////////////////////////////////////////////////////////////
160//
161// If solid is made up from a Boolean operation of two solids,
162// return the corresponding solid (for no=0 and 1)
163// If the solid is not a "Boolean", return 0
164
166{
167 const G4VSolid* subSolid = nullptr;
168 if( no == 0 )
169 {
170 subSolid = fPtrSolidA;
171 }
172 else if( no == 1 )
173 {
174 subSolid = fPtrSolidB;
175 }
176 else
177 {
178 DumpInfo();
179 G4Exception("G4BooleanSolid::GetConstituentSolid()",
180 "GeomSolids0002", FatalException, "Invalid solid index.");
181 }
182 return subSolid;
183}
184
185///////////////////////////////////////////////////////////////
186//
187// If solid is made up from a Boolean operation of two solids,
188// return the corresponding solid (for no=0 and 1)
189// If the solid is not a "Boolean", return 0
190
192{
193 G4VSolid* subSolid = nullptr;
194 if( no == 0 )
195 {
196 subSolid = fPtrSolidA;
197 }
198 else if( no == 1 )
199 {
200 subSolid = fPtrSolidB;
201 }
202 else
203 {
204 DumpInfo();
205 G4Exception("G4BooleanSolid::GetConstituentSolid()",
206 "GeomSolids0002", FatalException, "Invalid solid index.");
207 }
208 return subSolid;
209}
210
211//////////////////////////////////////////////////////////////////////////
212//
213// Returns entity type
214
216{
217 return {"G4BooleanSolid"};
218}
219
220//////////////////////////////////////////////////////////////////////////
221//
222// Set number of random points to be used for computing cubic volume
223
225{
226 if (st != fCubVolStatistics) { fCubicVolume = -1.; }
227 fCubVolStatistics = st;
228
229 // Propagate st to all components of the 1st solid
230 if (fPtrSolidA->GetNumOfConstituents() > 1)
231 {
232 G4VSolid* ptr = fPtrSolidA;
233 while(true)
234 {
235 G4String type = ptr->GetEntityType();
236 if (type == "G4DisplacedSolid")
237 {
238 ptr = ((G4DisplacedSolid*)ptr)->GetConstituentMovedSolid();
239 continue;
240 }
241 if (type == "G4ReflectedSolid")
242 {
243 ptr = ((G4ReflectedSolid*)ptr)->GetConstituentMovedSolid();
244 continue;
245 }
246 if (type == "G4ScaledSolid")
247 {
248 ptr = ((G4ScaledSolid*)ptr)->GetUnscaledSolid();
249 continue;
250 }
251 if (type != "G4MultiUnion") // G4MultiUnion doesn't have SetCubVolStatistics()
252 {
253 ((G4BooleanSolid*)ptr)->SetCubVolStatistics(st);
254 }
255 break;
256 }
257 }
258
259 // Propagate st to all components of the 2nd solid
260 if (fPtrSolidB->GetNumOfConstituents() > 1)
261 {
262 G4VSolid* ptr = fPtrSolidB;
263 while(true)
264 {
265 G4String type = ptr->GetEntityType();
266 if (type == "G4DisplacedSolid")
267 {
268 ptr = ((G4DisplacedSolid*)ptr)->GetConstituentMovedSolid();
269 continue;
270 }
271 if (type == "G4ReflectedSolid")
272 {
273 ptr = ((G4ReflectedSolid*)ptr)->GetConstituentMovedSolid();
274 continue;
275 }
276 if (type == "G4ScaledSolid")
277 {
278 ptr = ((G4ScaledSolid*)ptr)->GetUnscaledSolid();
279 continue;
280 }
281 if (type != "G4MultiUnion") // G4MultiUnion doesn't have SetCubVolStatistics()
282 {
283 ((G4BooleanSolid*)ptr)->SetCubVolStatistics(st);
284 }
285 break;
286 }
287 }
288}
289
290//////////////////////////////////////////////////////////////////////////
291//
292// Set epsilon for computing cubic volume
293
295{
296 if (ep != fCubVolEpsilon) { fCubicVolume = -1.; }
297 fCubVolEpsilon = ep;
298
299 // Propagate ep to all components of the 1st solid
300 if (fPtrSolidA->GetNumOfConstituents() > 1)
301 {
302 G4VSolid* ptr = fPtrSolidA;
303 while(true)
304 {
305 G4String type = ptr->GetEntityType();
306 if (type == "G4DisplacedSolid")
307 {
308 ptr = ((G4DisplacedSolid*)ptr)->GetConstituentMovedSolid();
309 continue;
310 }
311 if (type == "G4ReflectedSolid")
312 {
313 ptr = ((G4ReflectedSolid*)ptr)->GetConstituentMovedSolid();
314 continue;
315 }
316 if (type == "G4ScaledSolid")
317 {
318 ptr = ((G4ScaledSolid*)ptr)->GetUnscaledSolid();
319 continue;
320 }
321 if (type != "G4MultiUnion") // G4MultiUnion doesn't have SetCubVolEpsilon()
322 {
323 ((G4BooleanSolid*)ptr)->SetCubVolEpsilon(ep);
324 }
325 break;
326 }
327 }
328
329 // Propagate ep to all components of the 2nd solid
330 if (fPtrSolidB->GetNumOfConstituents() > 1)
331 {
332 G4VSolid* ptr = fPtrSolidB;
333 while(true)
334 {
335 G4String type = ptr->GetEntityType();
336 if (type == "G4DisplacedSolid")
337 {
338 ptr = ((G4DisplacedSolid*)ptr)->GetConstituentMovedSolid();
339 continue;
340 }
341 if (type == "G4ReflectedSolid")
342 {
343 ptr = ((G4ReflectedSolid*)ptr)->GetConstituentMovedSolid();
344 continue;
345 }
346 if (type == "G4ScaledSolid")
347 {
348 ptr = ((G4ScaledSolid*)ptr)->GetUnscaledSolid();
349 continue;
350 }
351 if (type != "G4MultiUnion") // G4MultiUnion doesn't have SetCubVolEpsilon()
352 {
353 ((G4BooleanSolid*)ptr)->SetCubVolEpsilon(ep);
354 }
355 break;
356 }
357 }
358}
359
360//////////////////////////////////////////////////////////////////////////
361//
362// Stream object contents to an output stream
363
364std::ostream& G4BooleanSolid::StreamInfo(std::ostream& os) const
365{
366 os << "-----------------------------------------------------------\n"
367 << " *** Dump for Boolean solid - " << GetName() << " ***\n"
368 << " ===================================================\n"
369 << " Solid type: " << GetEntityType() << "\n"
370 << " Parameters of constituent solids: \n"
371 << "===========================================================\n";
372 fPtrSolidA->StreamInfo(os);
373 fPtrSolidB->StreamInfo(os);
374 os << "===========================================================\n";
375
376 return os;
377}
378
379//////////////////////////////////////////////////////////////////////////
380//
381// Creates list of constituent primitives of and their placements
382
384 std::vector<std::pair<G4VSolid*,G4Transform3D>>& primitives,
385 const G4Transform3D& curPlacement) const
386{
387 G4Transform3D transform;
388 G4VSolid* solid;
389 G4String type;
390
391 // Repeat two times, first time for fPtrSolidA and then for fPtrSolidB
392 //
393 for (auto i=0; i<2; ++i)
394 {
395 transform = curPlacement;
396 solid = (i == 0) ? fPtrSolidA : fPtrSolidB;
397 type = solid->GetEntityType();
398
399 // While current solid is a trasformed solid just modify transform
400 //
401 while (type == "G4DisplacedSolid" ||
402 type == "G4ReflectedSolid" ||
403 type == "G4ScaledSolid")
404 {
405 if (type == "G4DisplacedSolid")
406 {
407 transform = transform * G4Transform3D(
408 ((G4DisplacedSolid*)solid)->GetObjectRotation(),
409 ((G4DisplacedSolid*)solid)->GetObjectTranslation());
410 solid = ((G4DisplacedSolid*)solid)->GetConstituentMovedSolid();
411 }
412 else if (type == "G4ReflectedSolid")
413 {
414 transform= transform*((G4ReflectedSolid*)solid)->GetDirectTransform3D();
415 solid = ((G4ReflectedSolid*)solid)->GetConstituentMovedSolid();
416 }
417 else if (type == "G4ScaledSolid")
418 {
419 transform = transform * ((G4ScaledSolid*)solid)->GetScaleTransform();
420 solid = ((G4ScaledSolid*)solid)->GetUnscaledSolid();
421 }
422 type = solid->GetEntityType();
423 }
424
425 // If current solid is a Boolean solid then continue recursion,
426 // otherwise add it to the list of primitives
427 //
428 if (type == "G4UnionSolid" ||
429 type == "G4SubtractionSolid" ||
430 type == "G4IntersectionSolid" ||
431 type == "G4BooleanSolid")
432 {
433 ((G4BooleanSolid *)solid)->GetListOfPrimitives(primitives,transform);
434 }
435 else
436 {
437 primitives.emplace_back(solid,transform);
438 }
439 }
440}
441
442//////////////////////////////////////////////////////////////////////////
443//
444// Returns a point (G4ThreeVector) randomly and uniformly selected
445// on the surface of the solid
446
448{
449 std::size_t nprims = fPrimitives.size();
450 std::pair<G4VSolid *, G4Transform3D> prim;
451
452 // Get list of primitives and find the total area of their surfaces
453 //
454 if (nprims == 0)
455 {
456 GetListOfPrimitives(fPrimitives, G4Transform3D());
457 nprims = fPrimitives.size();
458 fPrimitivesSurfaceArea = 0.;
459 for (std::size_t i=0; i<nprims; ++i)
460 {
461 fPrimitivesSurfaceArea += fPrimitives[i].first->GetSurfaceArea();
462 }
463 }
464
465 // Select random primitive, get random point on its surface and
466 // check that the point belongs to the surface of the solid
467 //
469 for (std::size_t k=0; k<100000; ++k) // try 100k times
470 {
471 G4double rand = fPrimitivesSurfaceArea * G4QuickRand();
472 G4double area = 0.;
473 for (std::size_t i=0; i<nprims; ++i)
474 {
475 prim = fPrimitives[i];
476 area += prim.first->GetSurfaceArea();
477 if (rand < area) { break; }
478 }
479 p = prim.first->GetPointOnSurface();
480 p = prim.second * G4Point3D(p);
481 if (Inside(p) == kSurface) { return p; }
482 }
483 std::ostringstream message;
484 message << "Solid - " << GetName() << "\n"
485 << "All 100k attempts to generate a point on the surface have failed!\n"
486 << "The solid created may be an invalid Boolean construct!";
487 G4Exception("G4BooleanSolid::GetPointOnSurface()",
488 "GeomSolids1001", JustWarning, message);
489 return p;
490}
491
492//////////////////////////////////////////////////////////////////////////
493//
494// Return total number of constituents used for construction of the solid
495
497{
498 return (fPtrSolidA->GetNumOfConstituents() + fPtrSolidB->GetNumOfConstituents());
499}
500
501//////////////////////////////////////////////////////////////////////////
502//
503// Return true if the resulting solid has only planar faces
504
506{
507 return (fPtrSolidA->IsFaceted() && fPtrSolidB->IsFaceted());
508}
509
510//////////////////////////////////////////////////////////////////////////
511//
512// Returns polyhedron for visualization
513
515{
516 if (fpPolyhedron == nullptr ||
517 fRebuildPolyhedron ||
518 fpPolyhedron->GetNumberOfRotationStepsAtTimeOfCreation() !=
519 fpPolyhedron->GetNumberOfRotationSteps())
520 {
521 G4RecursiveAutoLock l(&polyhedronMutex);
522 delete fpPolyhedron;
523 fpPolyhedron = CreatePolyhedron();
524 fRebuildPolyhedron = false;
525 l.unlock();
526 }
527 return fpPolyhedron;
528}
529
530//////////////////////////////////////////////////////////////////////////
531//
532// Stacks polyhedra for processing. Returns top polyhedron.
533
536 const G4VSolid* solid) const
537{
539 const G4String& type = solid->GetEntityType();
540 if (type == "G4UnionSolid")
541 { operation = HepPolyhedronProcessor::UNION; }
542 else if (type == "G4IntersectionSolid")
544 else if (type == "G4SubtractionSolid")
546 else
547 {
548 std::ostringstream message;
549 message << "Solid - " << solid->GetName()
550 << " - Unrecognised composite solid" << G4endl
551 << " Returning NULL !";
552 G4Exception("StackPolyhedron()", "GeomSolids1001", JustWarning, message);
553 return nullptr;
554 }
555
556 G4Polyhedron* top = nullptr;
557 const G4VSolid* solidA = solid->GetConstituentSolid(0);
558 const G4VSolid* solidB = solid->GetConstituentSolid(1);
559
560 if (solidA->GetConstituentSolid(0) != nullptr)
561 {
562 top = StackPolyhedron(processor, solidA);
563 }
564 else
565 {
566 top = solidA->GetPolyhedron();
567 }
568 G4Polyhedron* operand = solidB->GetPolyhedron();
569 if (operand != nullptr)
570 {
571 processor.push_back (operation, *operand);
572 }
573 else
574 {
575 std::ostringstream message;
576 message << "Solid - " << solid->GetName()
577 << " - No G4Polyhedron for Boolean component";
578 G4Exception("G4BooleanSolid::StackPolyhedron()",
579 "GeomSolids2001", JustWarning, message);
580 }
581
582 return top;
583}
584
585
586//////////////////////////////////////////////////////////////////////////
587//
588// Estimate Cubic Volume (capacity) and cache it for reuse.
589
591{
592 if(fCubicVolume < 0.)
593 {
594 G4AutoLock l(&boolMutex);
595 fCubicVolume = EstimateCubicVolume(fCubVolStatistics, fCubVolEpsilon);
596 l.unlock();
597 }
598 return fCubicVolume;
599}
600
601//////////////////////////////////////////////////////////////////////////
602//
603// Estimate Surface Area and cache it for reuse.
604
606{
607 if(fSurfaceArea < 0.)
608 {
609 G4AutoLock l(&boolMutex);
610 fSurfaceArea = EstimateSurfaceArea(fAreaStatistics, fAreaAccuracy);
611 l.unlock();
612 }
613 return fSurfaceArea;
614}
615
616//////////////////////////////////////////////////////////////////////////
617//
618// Set external Boolean processor.
619
620void
625
626//////////////////////////////////////////////////////////////////////////
627//
628// Get external Boolean processor.
629
G4TemplateAutoLock< G4RecursiveMutex > G4RecursiveAutoLock
G4TemplateAutoLock< G4Mutex > G4AutoLock
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
HepGeom::Point3D< G4double > G4Point3D
Definition G4Point3D.hh:34
G4double G4QuickRand(uint32_t seed=0)
CLHEP::HepRotation G4RotationMatrix
#define G4MUTEX_INITIALIZER
std::recursive_mutex G4RecursiveMutex
std::mutex G4Mutex
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
#define G4endl
Definition G4ios.hh:67
G4double GetSurfaceArea() override
G4VSolid * fPtrSolidA
G4BooleanSolid & operator=(const G4BooleanSolid &rhs)
void SetCubVolEpsilon(G4double ep)
G4Polyhedron * GetPolyhedron() const override
G4VSolid * fPtrSolidB
const G4VSolid * GetConstituentSolid(G4int no) const override
~G4BooleanSolid() override
void GetListOfPrimitives(std::vector< std::pair< G4VSolid *, G4Transform3D > > &, const G4Transform3D &) const
G4double GetCubicVolume() override
G4GeometryType GetEntityType() const override
G4BooleanSolid(const G4String &pName, G4VSolid *pSolidA, G4VSolid *pSolidB)
static G4VBooleanProcessor * GetExternalBooleanProcessor()
G4int GetNumOfConstituents() const override
static G4VBooleanProcessor * fExternalBoolProcessor
G4Polyhedron * StackPolyhedron(HepPolyhedronProcessor &, const G4VSolid *) const
std::ostream & StreamInfo(std::ostream &os) const override
void SetCubVolStatistics(G4int st)
G4bool IsFaceted() const override
static void SetExternalBooleanProcessor(G4VBooleanProcessor *extProcessor)
G4ThreeVector GetPointOnSurface() const override
G4DisplacedSolid is a solid that has been shifted from its original frame of reference to a new one....
G4VSolid * GetConstituentMovedSolid() const
G4ReflectedSolid is a solid that has been shifted from its original frame of reference to a new refle...
G4ScaledSolid is a solid that has been scaled in dimensions in X, Y or Z, from its original descripti...
G4VBooleanProcessor is a virtual base class for Boolean solid processing.
G4String GetName() const
virtual const G4VSolid * GetConstituentSolid(G4int no) const
Definition G4VSolid.cc:185
G4double EstimateCubicVolume(G4int nStat, G4double epsilon) const
Definition G4VSolid.cc:229
G4VSolid(const G4String &name)
Definition G4VSolid.cc:59
virtual EInside Inside(const G4ThreeVector &p) const =0
void DumpInfo() const
virtual G4Polyhedron * GetPolyhedron() const
Definition G4VSolid.cc:731
virtual G4Polyhedron * CreatePolyhedron() const
Definition G4VSolid.cc:726
G4VSolid & operator=(const G4VSolid &rhs)
Definition G4VSolid.cc:108
G4double EstimateSurfaceArea(G4int nStat, G4double epsilon) const
Definition G4VSolid.cc:290
virtual G4GeometryType GetEntityType() const =0
void push_back(Operation, const HepPolyhedron &)
@ kSurface
Definition geomdefs.hh:69