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