Geant4 11.2.2
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}
47
49
50//////////////////////////////////////////////////////////////////
51//
52// Constructor
53
55 G4VSolid* pSolidA ,
56 G4VSolid* pSolidB )
57 : G4VSolid(pName), fPtrSolidA(pSolidA), fPtrSolidB(pSolidB)
58{
59}
60
61//////////////////////////////////////////////////////////////////
62//
63// Constructor
64
66 G4VSolid* pSolidA ,
67 G4VSolid* pSolidB ,
68 G4RotationMatrix* rotMatrix,
69 const G4ThreeVector& transVector )
70 : G4VSolid(pName), createdDisplacedSolid(true)
71{
72 fPtrSolidA = pSolidA ;
73 fPtrSolidB = new G4DisplacedSolid("placedB",pSolidB,rotMatrix,transVector) ;
74}
75
76//////////////////////////////////////////////////////////////////
77//
78// Constructor
79
81 G4VSolid* pSolidA ,
82 G4VSolid* pSolidB ,
83 const G4Transform3D& transform )
84 : G4VSolid(pName), createdDisplacedSolid(true)
85{
86 fPtrSolidA = pSolidA ;
87 fPtrSolidB = new G4DisplacedSolid("placedB",pSolidB,transform) ;
88}
89
90///////////////////////////////////////////////////////////////
91//
92// Fake default constructor - sets only member data and allocates memory
93// for usage restricted to object persistency.
94
96 : G4VSolid(a)
97{
98}
99
100///////////////////////////////////////////////////////////////
101//
102// Destructor deletes transformation contents of the created displaced solid
103
105{
106 if(createdDisplacedSolid)
107 {
108 ((G4DisplacedSolid*)fPtrSolidB)->CleanTransformations();
109 }
110 delete fpPolyhedron; fpPolyhedron = nullptr;
111}
112
113///////////////////////////////////////////////////////////////
114//
115// Copy constructor
116
118 : G4VSolid (rhs), fPtrSolidA(rhs.fPtrSolidA), fPtrSolidB(rhs.fPtrSolidB),
119 fCubicVolume(rhs.fCubicVolume), fSurfaceArea(rhs.fSurfaceArea),
120 fCubVolStatistics(rhs.fCubVolStatistics),
121 fAreaStatistics(rhs.fAreaStatistics),
122 fCubVolEpsilon(rhs.fCubVolEpsilon),
123 fAreaAccuracy(rhs.fAreaAccuracy),
124 createdDisplacedSolid(rhs.createdDisplacedSolid)
125{
126 fPrimitives.resize(0); fPrimitivesSurfaceArea = 0.;
127}
128
129///////////////////////////////////////////////////////////////
130//
131// Assignment operator
132
134{
135 // Check assignment to self
136 //
137 if (this == &rhs) { return *this; }
138
139 // Copy base class data
140 //
142
143 // Copy data
144 //
147 fCubVolStatistics = rhs.fCubVolStatistics; fCubVolEpsilon = rhs.fCubVolEpsilon;
148 fAreaStatistics = rhs.fAreaStatistics; fAreaAccuracy = rhs.fAreaAccuracy;
149 createdDisplacedSolid= rhs.createdDisplacedSolid;
150
151 fRebuildPolyhedron = false;
152 delete fpPolyhedron; fpPolyhedron = nullptr;
153 fPrimitives.resize(0); fPrimitivesSurfaceArea = 0.;
154
155 return *this;
156}
157
158///////////////////////////////////////////////////////////////
159//
160// If solid is made up from a Boolean operation of two solids,
161// return the corresponding solid (for no=0 and 1)
162// If the solid is not a "Boolean", return 0
163
165{
166 const G4VSolid* subSolid = nullptr;
167 if( no == 0 )
168 subSolid = fPtrSolidA;
169 else if( no == 1 )
170 subSolid = fPtrSolidB;
171 else
172 {
173 DumpInfo();
174 G4Exception("G4BooleanSolid::GetConstituentSolid()",
175 "GeomSolids0002", FatalException, "Invalid solid index.");
176 }
177 return subSolid;
178}
179
180///////////////////////////////////////////////////////////////
181//
182// If solid is made up from a Boolean operation of two solids,
183// return the corresponding solid (for no=0 and 1)
184// If the solid is not a "Boolean", return 0
185
187{
188 G4VSolid* subSolid = nullptr;
189 if( no == 0 )
190 subSolid = fPtrSolidA;
191 else if( no == 1 )
192 subSolid = fPtrSolidB;
193 else
194 {
195 DumpInfo();
196 G4Exception("G4BooleanSolid::GetConstituentSolid()",
197 "GeomSolids0002", FatalException, "Invalid solid index.");
198 }
199 return subSolid;
200}
201
202//////////////////////////////////////////////////////////////////////////
203//
204// Returns entity type
205
207{
208 return {"G4BooleanSolid"};
209}
210
211//////////////////////////////////////////////////////////////////////////
212//
213// Stream object contents to an output stream
214
215std::ostream& G4BooleanSolid::StreamInfo(std::ostream& os) const
216{
217 os << "-----------------------------------------------------------\n"
218 << " *** Dump for Boolean solid - " << GetName() << " ***\n"
219 << " ===================================================\n"
220 << " Solid type: " << GetEntityType() << "\n"
221 << " Parameters of constituent solids: \n"
222 << "===========================================================\n";
225 os << "===========================================================\n";
226
227 return os;
228}
229
230//////////////////////////////////////////////////////////////////////////
231//
232// Creates list of constituent primitives of and their placements
233
235 std::vector<std::pair<G4VSolid*,G4Transform3D>>& primitives,
236 const G4Transform3D& curPlacement) const
237{
238 G4Transform3D transform;
239 G4VSolid* solid;
240 G4String type;
241
242 // Repeat two times, first time for fPtrSolidA and then for fPtrSolidB
243 //
244 for (auto i=0; i<2; ++i)
245 {
246 transform = curPlacement;
247 solid = (i == 0) ? fPtrSolidA : fPtrSolidB;
248 type = solid->GetEntityType();
249
250 // While current solid is a trasformed solid just modify transform
251 //
252 while (type == "G4DisplacedSolid" ||
253 type == "G4ReflectedSolid" ||
254 type == "G4ScaledSolid")
255 {
256 if (type == "G4DisplacedSolid")
257 {
258 transform = transform * G4Transform3D(
259 ((G4DisplacedSolid*)solid)->GetObjectRotation(),
260 ((G4DisplacedSolid*)solid)->GetObjectTranslation());
261 solid = ((G4DisplacedSolid*)solid)->GetConstituentMovedSolid();
262 }
263 else if (type == "G4ReflectedSolid")
264 {
265 transform= transform*((G4ReflectedSolid*)solid)->GetDirectTransform3D();
266 solid = ((G4ReflectedSolid*)solid)->GetConstituentMovedSolid();
267 }
268 else if (type == "G4ScaledSolid")
269 {
270 transform = transform * ((G4ScaledSolid*)solid)->GetScaleTransform();
271 solid = ((G4ScaledSolid*)solid)->GetUnscaledSolid();
272 }
273 type = solid->GetEntityType();
274 }
275
276 // If current solid is a Boolean solid then continue recursion,
277 // otherwise add it to the list of primitives
278 //
279 if (type == "G4UnionSolid" ||
280 type == "G4SubtractionSolid" ||
281 type == "G4IntersectionSolid" ||
282 type == "G4BooleanSolid")
283 {
284 ((G4BooleanSolid *)solid)->GetListOfPrimitives(primitives,transform);
285 }
286 else
287 {
288 primitives.emplace_back(solid,transform);
289 }
290 }
291}
292
293//////////////////////////////////////////////////////////////////////////
294//
295// Returns a point (G4ThreeVector) randomly and uniformly selected
296// on the surface of the solid
297
299{
300 std::size_t nprims = fPrimitives.size();
301 std::pair<G4VSolid *, G4Transform3D> prim;
302
303 // Get list of primitives and find the total area of their surfaces
304 //
305 if (nprims == 0)
306 {
307 GetListOfPrimitives(fPrimitives, G4Transform3D());
308 nprims = fPrimitives.size();
309 fPrimitivesSurfaceArea = 0.;
310 for (std::size_t i=0; i<nprims; ++i)
311 {
312 fPrimitivesSurfaceArea += fPrimitives[i].first->GetSurfaceArea();
313 }
314 }
315
316 // Select random primitive, get random point on its surface and
317 // check that the point belongs to the surface of the solid
318 //
320 for (std::size_t k=0; k<100000; ++k) // try 100k times
321 {
322 G4double rand = fPrimitivesSurfaceArea * G4QuickRand();
323 G4double area = 0.;
324 for (std::size_t i=0; i<nprims; ++i)
325 {
326 prim = fPrimitives[i];
327 area += prim.first->GetSurfaceArea();
328 if (rand < area) break;
329 }
330 p = prim.first->GetPointOnSurface();
331 p = prim.second * G4Point3D(p);
332 if (Inside(p) == kSurface) return p;
333 }
334 std::ostringstream message;
335 message << "Solid - " << GetName() << "\n"
336 << "All 100k attempts to generate a point on the surface have failed!\n"
337 << "The solid created may be an invalid Boolean construct!";
338 G4Exception("G4BooleanSolid::GetPointOnSurface()",
339 "GeomSolids1001", JustWarning, message);
340 return p;
341}
342
343//////////////////////////////////////////////////////////////////////////
344//
345// Returns polyhedron for visualization
346
348{
349 if (fpPolyhedron == nullptr ||
350 fRebuildPolyhedron ||
352 fpPolyhedron->GetNumberOfRotationSteps())
353 {
354 G4RecursiveAutoLock l(&polyhedronMutex);
355 delete fpPolyhedron;
356 fpPolyhedron = CreatePolyhedron();
357 fRebuildPolyhedron = false;
358 l.unlock();
359 }
360 return fpPolyhedron;
361}
362
363//////////////////////////////////////////////////////////////////////////
364//
365// Stacks polyhedra for processing. Returns top polyhedron.
366
369 const G4VSolid* solid) const
370{
372 const G4String& type = solid->GetEntityType();
373 if (type == "G4UnionSolid")
374 { operation = HepPolyhedronProcessor::UNION; }
375 else if (type == "G4IntersectionSolid")
377 else if (type == "G4SubtractionSolid")
379 else
380 {
381 std::ostringstream message;
382 message << "Solid - " << solid->GetName()
383 << " - Unrecognised composite solid" << G4endl
384 << " Returning NULL !";
385 G4Exception("StackPolyhedron()", "GeomSolids1001", JustWarning, message);
386 return nullptr;
387 }
388
389 G4Polyhedron* top = nullptr;
390 const G4VSolid* solidA = solid->GetConstituentSolid(0);
391 const G4VSolid* solidB = solid->GetConstituentSolid(1);
392
393 if (solidA->GetConstituentSolid(0) != nullptr)
394 {
395 top = StackPolyhedron(processor, solidA);
396 }
397 else
398 {
399 top = solidA->GetPolyhedron();
400 }
401 G4Polyhedron* operand = solidB->GetPolyhedron();
402 if (operand != nullptr)
403 {
404 processor.push_back (operation, *operand);
405 }
406 else
407 {
408 std::ostringstream message;
409 message << "Solid - " << solid->GetName()
410 << " - No G4Polyhedron for Boolean component";
411 G4Exception("G4BooleanSolid::StackPolyhedron()",
412 "GeomSolids2001", JustWarning, message);
413 }
414
415 return top;
416}
417
418
419//////////////////////////////////////////////////////////////////////////
420//
421// Estimate Cubic Volume (capacity) and cache it for reuse.
422
424{
425 if(fCubicVolume < 0.)
426 {
427 fCubicVolume = EstimateCubicVolume(fCubVolStatistics, fCubVolEpsilon);
428 }
429 return fCubicVolume;
430}
431
432//////////////////////////////////////////////////////////////////////////
433//
434// Set external Boolean processor.
435
436void
441
442//////////////////////////////////////////////////////////////////////////
443//
444// Get external Boolean processor.
445
@ 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()
#define G4MUTEX_INITIALIZER
std::recursive_mutex G4RecursiveMutex
HepGeom::Transform3D G4Transform3D
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4VSolid * fPtrSolidA
G4BooleanSolid & operator=(const G4BooleanSolid &rhs)
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()
static G4VBooleanProcessor * fExternalBoolProcessor
G4Polyhedron * StackPolyhedron(HepPolyhedronProcessor &, const G4VSolid *) const
std::ostream & StreamInfo(std::ostream &os) const override
static void SetExternalBooleanProcessor(G4VBooleanProcessor *extProcessor)
G4ThreeVector GetPointOnSurface() const override
G4VSolid * GetConstituentMovedSolid() const
G4int GetNumberOfRotationStepsAtTimeOfCreation() const
G4String GetName() const
virtual const G4VSolid * GetConstituentSolid(G4int no) const
Definition G4VSolid.cc:167
virtual std::ostream & StreamInfo(std::ostream &os) const =0
G4double EstimateCubicVolume(G4int nStat, G4double epsilon) const
Definition G4VSolid.cc:203
virtual EInside Inside(const G4ThreeVector &p) const =0
void DumpInfo() const
virtual G4Polyhedron * GetPolyhedron() const
Definition G4VSolid.cc:705
virtual G4Polyhedron * CreatePolyhedron() const
Definition G4VSolid.cc:700
G4VSolid & operator=(const G4VSolid &rhs)
Definition G4VSolid.cc:107
virtual G4GeometryType GetEntityType() const =0
void push_back(Operation, const HepPolyhedron &)
static G4int GetNumberOfRotationSteps()
@ kSurface
Definition geomdefs.hh:69