Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VCSGfaceted.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// G4VCSGfaceted implementation; a virtual class of a CSG type shape
27// that is built entirely out of G4VCSGface faces.
28//
29// Author: David C. Williams (UCSC), 1998
30// --------------------------------------------------------------------
31
32#include "G4VCSGfaceted.hh"
33#include "G4VCSGface.hh"
34#include "G4SolidExtentList.hh"
35
36#include "G4VoxelLimits.hh"
37#include "G4AffineTransform.hh"
38
39#include "G4QuickRand.hh"
40
41#include "G4Polyhedron.hh"
42#include "G4VGraphicsScene.hh"
43#include "G4VisExtent.hh"
44
45#include "G4AutoLock.hh"
46
47namespace
48{
49 G4Mutex polyhedronMutex = G4MUTEX_INITIALIZER;
50 G4Mutex vcsgMutex = G4MUTEX_INITIALIZER;
51}
52
53//
54// Constructor
55//
57 : G4VSolid(name),
58 fStatistics(1000000), fCubVolEpsilon(0.001), fAreaAccuracy(-1.)
59{
60}
61
62
63//
64// Fake default constructor - sets only member data and allocates memory
65// for usage restricted to object persistency.
66//
68 : G4VSolid(a),
69 fStatistics(1000000), fCubVolEpsilon(0.001), fAreaAccuracy(-1.)
70{
71}
72
73//
74// Destructor
75//
81
82
83//
84// Copy constructor
85//
87 : G4VSolid( source )
88{
89 fStatistics = source.fStatistics;
90 fCubVolEpsilon = source.fCubVolEpsilon;
91 fAreaAccuracy = source.fAreaAccuracy;
92
93 CopyStuff( source );
94}
95
96
97//
98// Assignment operator
99//
101{
102 if (&source == this) { return *this; }
103
104 // Copy base class data
105 //
106 G4VSolid::operator=(source);
107
108 // Copy data
109 //
110 fStatistics = source.fStatistics;
111 fCubVolEpsilon = source.fCubVolEpsilon;
112 fAreaAccuracy = source.fAreaAccuracy;
113
114 DeleteStuff();
115 CopyStuff( source );
116
117 return *this;
118}
119
120
121//
122// CopyStuff (protected)
123//
124// Copy the contents of source
125//
127{
128 numFace = source.numFace;
129 if (numFace == 0) { return; } // odd, but permissable?
130
131 faces = new G4VCSGface*[numFace];
132
133 G4VCSGface **face = faces,
134 **sourceFace = source.faces;
135 do // Loop checking, 13.08.2015, G.Cosmo
136 {
137 *face = (*sourceFace)->Clone();
138 } while( ++sourceFace, ++face < faces+numFace );
139 fCubicVolume = source.fCubicVolume;
140 fSurfaceArea = source.fSurfaceArea;
141 fRebuildPolyhedron = false;
142 fpPolyhedron = nullptr;
143}
144
145
146//
147// DeleteStuff (protected)
148//
149// Delete all allocated objects
150//
152{
153 if (numFace != 0)
154 {
155 G4VCSGface **face = faces;
156 do // Loop checking, 13.08.2015, G.Cosmo
157 {
158 delete *face;
159 } while( ++face < faces + numFace );
160
161 delete [] faces;
162 }
163 delete fpPolyhedron; fpPolyhedron = nullptr;
164}
165
166
167//
168// CalculateExtent
169//
171 const G4VoxelLimits& voxelLimit,
172 const G4AffineTransform& transform,
173 G4double& min,
174 G4double& max ) const
175{
176 G4SolidExtentList extentList( axis, voxelLimit );
177
178 //
179 // Loop over all faces, checking min/max extent as we go.
180 //
181 G4VCSGface **face = faces;
182 do // Loop checking, 13.08.2015, G.Cosmo
183 {
184 (*face)->CalculateExtent( axis, voxelLimit, transform, extentList );
185 } while( ++face < faces + numFace );
186
187 //
188 // Return min/max value
189 //
190 return extentList.GetExtent( min, max );
191}
192
193
194//
195// Inside
196//
197// It could be a good idea to override this virtual
198// member to add first a simple test (such as spherical
199// test or whatnot) and to call this version only if
200// the simplier test fails.
201//
203{
204 EInside answer=kOutside;
205 G4VCSGface **face = faces;
206 G4double best = kInfinity;
207 do // Loop checking, 13.08.2015, G.Cosmo
208 {
209 G4double distance;
210 EInside result = (*face)->Inside( p, kCarTolerance/2, &distance );
211 if (result == kSurface) { return kSurface; }
212 if (distance < best)
213 {
214 best = distance;
215 answer = result;
216 }
217 } while( ++face < faces + numFace );
218
219 return answer;
220}
221
222
223//
224// SurfaceNormal
225//
227{
228 G4ThreeVector answer;
229 G4VCSGface **face = faces;
230 G4double best = kInfinity;
231 do // Loop checking, 13.08.2015, G.Cosmo
232 {
233 G4double distance = kInfinity;
234 G4ThreeVector normal = (*face)->Normal( p, &distance );
235 if (distance < best)
236 {
237 best = distance;
238 answer = normal;
239 }
240 } while( ++face < faces + numFace );
241
242 return answer;
243}
244
245
246//
247// DistanceToIn(p,v)
248//
250 const G4ThreeVector& v ) const
251{
252 G4double distance = kInfinity;
253 G4double distFromSurface = kInfinity;
254 G4VCSGface **face = faces;
255 G4VCSGface *bestFace = *face;
256 do // Loop checking, 13.08.2015, G.Cosmo
257 {
258 G4double faceDistance,
259 faceDistFromSurface;
260 G4ThreeVector faceNormal;
261 G4bool faceAllBehind;
262 if ((*face)->Intersect( p, v, false, kCarTolerance/2,
263 faceDistance, faceDistFromSurface,
264 faceNormal, faceAllBehind ) )
265 {
266 //
267 // Intersecting face
268 //
269 if (faceDistance < distance)
270 {
271 distance = faceDistance;
272 distFromSurface = faceDistFromSurface;
273 bestFace = *face;
274 if (distFromSurface <= 0) { return 0; }
275 }
276 }
277 } while( ++face < faces + numFace );
278
279 if (distance < kInfinity && distFromSurface<kCarTolerance/2)
280 {
281 if (bestFace->Distance(p,false) < kCarTolerance/2) { distance = 0; }
282 }
283
284 return distance;
285}
286
287
288//
289// DistanceToIn(p)
290//
292{
293 return DistanceTo( p, false );
294}
295
296
297//
298// DistanceToOut(p,v)
299//
301 const G4ThreeVector& v,
302 const G4bool calcNorm,
303 G4bool* validNorm,
304 G4ThreeVector* n ) const
305{
306 G4bool allBehind = true;
307 G4double distance = kInfinity;
308 G4double distFromSurface = kInfinity;
309 G4ThreeVector normal;
310
311 G4VCSGface **face = faces;
312 G4VCSGface *bestFace = *face;
313 do // Loop checking, 13.08.2015, G.Cosmo
314 {
315 G4double faceDistance,
316 faceDistFromSurface;
317 G4ThreeVector faceNormal;
318 G4bool faceAllBehind;
319 if ((*face)->Intersect( p, v, true, kCarTolerance/2,
320 faceDistance, faceDistFromSurface,
321 faceNormal, faceAllBehind ) )
322 {
323 //
324 // Intersecting face
325 //
326 if ( (distance < kInfinity) || (!faceAllBehind) ) { allBehind = false; }
327 if (faceDistance < distance)
328 {
329 distance = faceDistance;
330 distFromSurface = faceDistFromSurface;
331 normal = faceNormal;
332 bestFace = *face;
333 if (distFromSurface <= 0.) { break; }
334 }
335 }
336 } while( ++face < faces + numFace );
337
338 if (distance < kInfinity)
339 {
340 if (distFromSurface <= 0.)
341 {
342 distance = 0.;
343 }
344 else if (distFromSurface<kCarTolerance/2)
345 {
346 if (bestFace->Distance(p,true) < kCarTolerance/2) { distance = 0.; }
347 }
348
349 if (calcNorm)
350 {
351 *validNorm = allBehind;
352 *n = normal;
353 }
354 }
355 else
356 {
357 if (Inside(p) == kSurface) { distance = 0.; }
358 if (calcNorm) { *validNorm = false; }
359 }
360
361 return distance;
362}
363
364
365//
366// DistanceToOut(p)
367//
369{
370 return DistanceTo( p, true );
371}
372
373
374//
375// DistanceTo
376//
377// Protected routine called by DistanceToIn and DistanceToOut
378//
380 const G4bool outgoing ) const
381{
382 G4VCSGface **face = faces;
383 G4double best = kInfinity;
384 do // Loop checking, 13.08.2015, G.Cosmo
385 {
386 G4double distance = (*face)->Distance( p, outgoing );
387 if (distance < best) { best = distance; }
388 } while( ++face < faces + numFace );
389
390 return (best < 0.5*kCarTolerance) ? 0. : best;
391}
392
393
394//
395// DescribeYourselfTo
396//
398{
399 scene.AddSolid( *this );
400}
401
402
403//
404// GetExtent
405//
406// Define the sides of the box into which our solid instance would fit.
407//
409{
410 static const G4ThreeVector xMax(1,0,0), xMin(-1,0,0),
411 yMax(0,1,0), yMin(0,-1,0),
412 zMax(0,0,1), zMin(0,0,-1);
413 static const G4ThreeVector *axes[6] =
414 { &xMin, &xMax, &yMin, &yMax, &zMin, &zMax };
415
416 G4double answers[6] =
417 {-kInfinity, -kInfinity, -kInfinity, -kInfinity, -kInfinity, -kInfinity};
418
419 G4VCSGface **face = faces;
420 do // Loop checking, 13.08.2015, G.Cosmo
421 {
422 const G4ThreeVector **axis = axes+5 ;
423 G4double* answer = answers+5;
424 do // Loop checking, 13.08.2015, G.Cosmo
425 {
426 G4double testFace = (*face)->Extent( **axis );
427 if (testFace > *answer) { *answer = testFace; }
428 }
429 while( --axis, --answer >= answers );
430
431 } while( ++face < faces + numFace );
432
433 return { -answers[0], answers[1],
434 -answers[2], answers[3],
435 -answers[4], answers[5] };
436}
437
438
439//
440// GetEntityType
441//
443{
444 return {"G4CSGfaceted"};
445}
446
447
448//
449// Stream object contents to an output stream
450//
451std::ostream& G4VCSGfaceted::StreamInfo( std::ostream& os ) const
452{
453 os << "-----------------------------------------------------------\n"
454 << " *** Dump for solid - " << GetName() << " ***\n"
455 << " ===================================================\n"
456 << " Solid type: G4VCSGfaceted\n"
457 << " Parameters: \n"
458 << " number of faces: " << numFace << "\n"
459 << "-----------------------------------------------------------\n";
460
461 return os;
462}
463
464
465//
466// GetCubVolStatistics
467//
469{
470 return fStatistics;
471}
472
473
474//
475// GetCubVolEpsilon
476//
478{
479 return fCubVolEpsilon;
480}
481
482
483//
484// SetCubVolStatistics
485//
487{
488 fCubicVolume=0.;
489 fStatistics=st;
490}
491
492
493//
494// SetCubVolEpsilon
495//
497{
498 fCubicVolume=0.;
499 fCubVolEpsilon=ep;
500}
501
502
503//
504// GetAreaStatistics
505//
507{
508 return fStatistics;
509}
510
511
512//
513// GetAreaAccuracy
514//
516{
517 return fAreaAccuracy;
518}
519
520
521//
522// SetAreaStatistics
523//
525{
526 fSurfaceArea=0.;
527 fStatistics=st;
528}
529
530
531//
532// SetAreaAccuracy
533//
535{
536 fSurfaceArea=0.;
537 fAreaAccuracy=ep;
538}
539
540
541//
542// GetCubicVolume
543//
545{
546 if (fCubicVolume == 0)
547 {
548 G4AutoLock l(&vcsgMutex);
549 fCubicVolume = EstimateCubicVolume(fStatistics,fCubVolEpsilon);
550 l.unlock();
551 }
552 return fCubicVolume;
553}
554
555
556//
557// GetSurfaceArea
558//
560{
561 if (fSurfaceArea == 0)
562 {
563 G4AutoLock l(&vcsgMutex);
564 fSurfaceArea = EstimateSurfaceArea(fStatistics,fAreaAccuracy);
565 l.unlock();
566 }
567 return fSurfaceArea;
568}
569
570
571//
572// GetPolyhedron
573//
575{
576 if (fpPolyhedron == nullptr ||
578 fpPolyhedron->GetNumberOfRotationStepsAtTimeOfCreation() !=
579 fpPolyhedron->GetNumberOfRotationSteps())
580 {
581 G4AutoLock l(&polyhedronMutex);
582 delete fpPolyhedron;
584 fRebuildPolyhedron = false;
585 l.unlock();
586 }
587 return fpPolyhedron;
588}
589
590
591//
592// GetPointOnSurfaceGeneric proportional to Areas of faces
593// in case of GenericPolycone or GenericPolyhedra
594//
596{
597 // Preparing variables
598 //
599 G4ThreeVector answer=G4ThreeVector(0.,0.,0.);
600 G4VCSGface **face = faces;
601 G4double area = 0.;
602 G4int i;
603 std::vector<G4double> areas;
604
605 // First step: calculate surface areas
606 //
607 do // Loop checking, 13.08.2015, G.Cosmo
608 {
609 G4double result = (*face)->SurfaceArea( );
610 areas.push_back(result);
611 area=area+result;
612 } while( ++face < faces + numFace );
613
614 // Second Step: choose randomly one surface
615 //
616 G4VCSGface **face1 = faces;
617 G4double chose = area*G4QuickRand();
618 G4double Achose1, Achose2;
619 Achose1=0.; Achose2=0.;
620 i=0;
621
622 do
623 {
624 Achose2+=areas[i];
625 if(chose>=Achose1 && chose<Achose2)
626 {
627 G4ThreeVector point;
628 point= (*face1)->GetPointOnFace();
629 return point;
630 }
631 ++i;
632 Achose1=Achose2;
633 } while( ++face1 < faces + numFace );
634
635 return answer;
636}
G4TemplateAutoLock< G4Mutex > G4AutoLock
G4double G4QuickRand(uint32_t seed=0)
#define G4MUTEX_INITIALIZER
std::mutex G4Mutex
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4String G4GeometryType
Definition G4VSolid.hh:70
G4AffineTransform is a class for geometric affine transformations. It supports efficient arbitrary ro...
G4SolidExtentList is utility class designed for calculating the extent of a CSG-like solid for a voxe...
G4bool GetExtent(G4double &min, G4double &max) const
G4VCSGface is virtual base class, representing one side (or face) of a CSG-like solid....
virtual G4VCSGface * Clone()=0
virtual G4double Distance(const G4ThreeVector &p, G4bool outgoing)=0
void CopyStuff(const G4VCSGfaceted &source)
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const override
G4double GetCubVolEpsilon() const
G4Polyhedron * GetPolyhedron() const override
G4VCSGfaceted & operator=(const G4VCSGfaceted &source)
G4double fCubicVolume
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const override
G4int GetCubVolStatistics() const
void SetAreaStatistics(G4int st)
EInside Inside(const G4ThreeVector &p) const override
G4double GetSurfaceArea() override
G4VisExtent GetExtent() const override
G4VCSGface ** faces
void SetAreaAccuracy(G4double ep)
~G4VCSGfaceted() override
void DescribeYourselfTo(G4VGraphicsScene &scene) const override
G4VCSGfaceted(const G4String &name)
G4Polyhedron * fpPolyhedron
G4bool fRebuildPolyhedron
G4Polyhedron * CreatePolyhedron() const override=0
G4double GetCubicVolume() override
G4int GetAreaStatistics() const
std::ostream & StreamInfo(std::ostream &os) const override
G4double GetAreaAccuracy() const
G4double fSurfaceArea
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const override
G4ThreeVector GetPointOnSurfaceGeneric() const
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const override
virtual G4double DistanceTo(const G4ThreeVector &p, const G4bool outgoing) const
void SetCubVolEpsilon(G4double ep)
void SetCubVolStatistics(G4int st)
G4GeometryType GetEntityType() const override
virtual void AddSolid(const G4Box &)=0
G4String GetName() const
G4double EstimateCubicVolume(G4int nStat, G4double epsilon) const
Definition G4VSolid.cc:229
G4VSolid(const G4String &name)
Definition G4VSolid.cc:59
G4double kCarTolerance
Definition G4VSolid.hh:418
G4VSolid & operator=(const G4VSolid &rhs)
Definition G4VSolid.cc:108
G4double EstimateSurfaceArea(G4int nStat, G4double epsilon) const
Definition G4VSolid.cc:290
G4VoxelLimits represents limitation/restrictions of space, where restrictions are only made perpendic...
EAxis
Definition geomdefs.hh:54
EInside
Definition geomdefs.hh:67
@ kOutside
Definition geomdefs.hh:68
@ kSurface
Definition geomdefs.hh:69
const axis_t axis_to_type< N >::axis
Definition pugixml.cc:9668