Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BREPSolidPolyhedra.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// $Id$
27//
28// ----------------------------------------------------------------------
29// GEANT 4 class source file
30//
31// G4BREPSolidPolyhedra.cc
32//
33// ----------------------------------------------------------------------
34// The polygonal solid G4BREPSolidPolyhedra is a shape defined by an inner
35// and outer polygonal surface and two planes perpendicular to the Z axis.
36// Each polygonal surface is created by linking a series of polygons created
37// at different planes perpendicular to the Z-axis. All these polygons all
38// have the same number of sides (sides) and are defined at the same Z planes
39// for both inner and outer polygonal surfaces.
40// ----------------------------------------------------------------------
41//
42// History
43// -------
44//
45// Bug-fix #266 by R.Chytracek:
46// - The situation when phi1 = 0 dphi1 = 2*pi and all RMINs = 0.0 is handled
47// now. In this case the inner planes are not created. The fix goes even
48// further this means it considers more than 2 z-planes and inner planes
49// are not created whenever two consecutive RMINs are = 0.0 .
50//
51// Corrections by S.Giani:
52// - Xaxis now corresponds to phi=0
53// - partial angle = phiTotal / Nsides
54// - end planes exact boundary calculation for phiTotal < 2pi
55// (also including case with RMIN=RMAX)
56// - Xaxis now properly rotated to compute correct scope of vertixes
57// - corrected surface orientation for outer faces parallel to Z
58// - completed explicit setting of the orientation for all faces
59// - some comparison between doubles avoided by using tolerances
60// - visualisation parameters made consistent with the use made by
61// constructor of the input arguments (i.e. circumscribed radius).
62// ----------------------------------------------------------------------
63
64#include <sstream>
65
68#include "G4SystemOfUnits.hh"
69#include "G4FPlane.hh"
70
72 G4double start_angle,
73 G4double opening_angle,
74 G4int sides,
75 G4int num_z_planes,
76 G4double z_start,
77 G4double z_values[],
78 G4double RMIN[],
79 G4double RMAX[] )
80 : G4BREPSolid(name)
81{
82 // Store the original parameters, to be used in visualisation
83 // Note radii are not scaled because this BREP uses the radius of the
84 // circumscribed circle and also graphics_reps/G4Polyhedron uses the
85 // radius of the circumscribed circle.
86
87 // Save contructor parameters
88 //
89 constructorParams.start_angle = start_angle;
90 constructorParams.opening_angle = opening_angle;
91 constructorParams.sides = sides;
92 constructorParams.num_z_planes = num_z_planes;
93 constructorParams.z_start = z_start;
94 constructorParams.z_values = 0;
95 constructorParams.RMIN = 0;
96 constructorParams.RMAX = 0;
97
98 if( num_z_planes > 0 )
99 {
100 constructorParams.z_values = new G4double[num_z_planes];
101 constructorParams.RMIN = new G4double[num_z_planes];
102 constructorParams.RMAX = new G4double[num_z_planes];
103 for( G4int idx = 0; idx < num_z_planes; ++idx )
104 {
105 constructorParams.z_values[idx] = z_values[idx];
106 constructorParams.RMIN[idx] = RMIN[idx];
107 constructorParams.RMAX[idx] = RMAX[idx];
108 }
109 }
110
111 // z_values[0] should be equal to z_start, for consistency
112 // with what the constructor does.
113 // Otherwise the z_values that are shifted by (z_values[0] - z_start) ,
114 // because z_values are only used in the form
115 // length = z_values[d+1] - z_values[d]; // JA Apr 2, 97
116
117 if( z_values[0] != z_start )
118 {
119 std::ostringstream message;
120 message << "Construction Error. z_values[0] must be equal to z_start!"
121 << G4endl
122 << " Wrong solid parameters: "
123 << " z_values[0]= " << z_values[0] << " is not equal to "
124 << " z_start= " << z_start << ".";
125 G4Exception( "G4BREPSolidPolyhedra::G4BREPSolidPolyhedra()",
126 "GeomSolids1002", JustWarning, message );
127 if( num_z_planes <= 0 ) { constructorParams.z_values = new G4double[1]; }
128 constructorParams.z_values[0]= z_start;
129 }
130
131 active=1;
132 InitializePolyhedra();
133}
134
136 : G4BREPSolid(a)
137{
138 constructorParams.start_angle = 0.;
139 constructorParams.opening_angle = 0.;
140 constructorParams.sides = 0;
141 constructorParams.num_z_planes = 0;
142 constructorParams.z_start = 0.;
143 constructorParams.z_values = 0;
144 constructorParams.RMIN = 0;
145 constructorParams.RMAX = 0;
146}
147
149{
150 if( constructorParams.num_z_planes > 0 )
151 {
152 delete [] constructorParams.z_values;
153 delete [] constructorParams.RMIN;
154 delete [] constructorParams.RMAX;
155 }
156}
157
159 : G4BREPSolid(rhs)
160{
161 constructorParams.start_angle = rhs.constructorParams.start_angle;
162 constructorParams.opening_angle = rhs.constructorParams.opening_angle;
163 constructorParams.sides = rhs.constructorParams.sides;
164 constructorParams.num_z_planes = rhs.constructorParams.num_z_planes;
165 constructorParams.z_start = rhs.constructorParams.z_start;
166 constructorParams.z_values = 0;
167 constructorParams.RMIN = 0;
168 constructorParams.RMAX = 0;
169 G4int num_z_planes = constructorParams.num_z_planes;
170 if( num_z_planes > 0 )
171 {
172 constructorParams.z_values = new G4double[num_z_planes];
173 constructorParams.RMIN = new G4double[num_z_planes];
174 constructorParams.RMAX = new G4double[num_z_planes];
175 for( G4int idx = 0; idx < num_z_planes; ++idx )
176 {
177 constructorParams.z_values[idx] = rhs.constructorParams.z_values[idx];
178 constructorParams.RMIN[idx] = rhs.constructorParams.RMIN[idx];
179 constructorParams.RMAX[idx] = rhs.constructorParams.RMAX[idx];
180 }
181 }
182
183 InitializePolyhedra();
184}
185
188{
189 // Check assignment to self
190 //
191 if (this == &rhs) { return *this; }
192
193 // Copy base class data
194 //
196
197 // Copy data
198 //
199 constructorParams.start_angle = rhs.constructorParams.start_angle;
200 constructorParams.opening_angle = rhs.constructorParams.opening_angle;
201 constructorParams.sides = rhs.constructorParams.sides;
202 constructorParams.num_z_planes = rhs.constructorParams.num_z_planes;
203 constructorParams.z_start = rhs.constructorParams.z_start;
204 G4int num_z_planes = constructorParams.num_z_planes;
205 if( num_z_planes > 0 )
206 {
207 delete [] constructorParams.z_values;
208 delete [] constructorParams.RMIN;
209 delete [] constructorParams.RMAX;
210 constructorParams.z_values = new G4double[num_z_planes];
211 constructorParams.RMIN = new G4double[num_z_planes];
212 constructorParams.RMAX = new G4double[num_z_planes];
213 for( G4int idx = 0; idx < num_z_planes; ++idx )
214 {
215 constructorParams.z_values[idx] = rhs.constructorParams.z_values[idx];
216 constructorParams.RMIN[idx] = rhs.constructorParams.RMIN[idx];
217 constructorParams.RMAX[idx] = rhs.constructorParams.RMAX[idx];
218 }
219 }
220
221 InitializePolyhedra();
222
223 return *this;
224}
225
226void G4BREPSolidPolyhedra::InitializePolyhedra()
227{
228 G4double start_angle = constructorParams.start_angle;
229 G4double opening_angle = constructorParams.opening_angle;
230 G4int sides = constructorParams.sides;
231 G4int num_z_planes = constructorParams.num_z_planes;
232 G4double z_start = constructorParams.z_start;
233 G4double* z_values = constructorParams.z_values;
234 G4double* RMIN = constructorParams.RMIN;
235 G4double* RMAX = constructorParams.RMAX;
236 G4int sections = num_z_planes - 1;
237
238 if( opening_angle >= 2*pi-perMillion )
239 {
240 nb_of_surfaces = 2*(sections * sides) + 2;
241 }
242 else
243 {
244 nb_of_surfaces = 2*(sections * sides) + 4;
245 }
246
247 G4int MaxNbOfSurfaces = nb_of_surfaces;
248 G4Surface** MaxSurfaceVec = new G4Surface*[MaxNbOfSurfaces];
249
250 G4Vector3D Axis(0,0,1);
251 G4Vector3D XAxis(1,0,0);
252 G4Vector3D TmpAxis;
253 G4Point3D Origin(0,0,z_start);
254 G4Point3D LocalOrigin(0,0,z_start);
255 G4double Length;
256 G4int Count = 0 ;
257 G4double PartAngle = (opening_angle)/sides;
258
259
260 ///////////////////////////////////////////////////
261 // Preconditions check
262
263 // Detecting minimal required number of sides
264 //
265 if( sides < 3 )
266 {
267 G4Exception( "G4BREPSolidPolyhedra::G4BREPSolidPolyhedra()",
268 "InvalidSetup", FatalException,
269 "The solid must have at least 3 sides!" );
270 }
271
272 // Detecting minimal required number of z-sections
273 //
274 if( num_z_planes < 2 )
275 {
276 G4Exception( "G4BREPSolidPolyhedra::G4BREPSolidPolyhedra()",
277 "GeomSolids0002", FatalException,
278 "The solid must have at least 2 z-sections!" );
279 }
280
281 // Detect invalid configurations at the ends of polyhedra which
282 // would not lead to a valid solid creation and likely to a crash
283 //
284 if( z_values[0] == z_values[1]
285 || z_values[sections-1] == z_values[sections] )
286 {
287 G4Exception( "G4BREPSolidPolyhedra::G4BREPSolidPolyhedra()",
288 "GeomSolids0002", FatalException,
289 "The solid must have the first 2 and the last 2 z-values different!" );
290 }
291
292 // Find out how the z-values sequence is ordered
293 //
294 G4bool increasing;
295 if( z_values[0] < z_values[1] )
296 {
297 increasing = true;
298 }
299 else
300 {
301 increasing = false;
302 }
303
304 // Detecting polyhedra teeth.
305 // It's forbidden to specify unordered, e.g. non-increasing or
306 // non-decreasing sequence of z-values. It may be provided by a
307 // specific solid in a future.
308 //
309 for( G4int idx = 0; idx < sections; idx++ )
310 {
311 if( ( z_values[idx] > z_values[idx+1] && increasing ) ||
312 ( z_values[idx] < z_values[idx+1] && !increasing ) )
313 {
314 // ERROR! Invalid sequence of z-values
315 //
316 std::ostringstream message;
317 message << "Unordered, non-increasing or non-decreasing sequence."
318 << G4endl
319 << " Unordered z_values sequence detected !" << G4endl
320 << " Check z_values with indexes: "
321 << idx << " " << (idx+1) << ".";
322 G4Exception( "G4BREPSolidPolyhedra::G4BREPSolidPolyhedra()",
323 "GeomSolids0002", FatalException, message );
324 }
325 }
326
327///////////////////////////////////////////////////
328#ifdef G4_EXPERIMENTAL_CODE
329 // There is one problem when sequence of z values is not increasing in a
330 // regular way, in other words, it's not purely increasing or decreasing
331 // Irregular sequence can be provided in order to define a polyhedra having
332 // teeth as shown on the picture bellow
333 // In this sequence can happen the following:
334 // z[a-1] > z[a] < z[a+1] && z[a+1] >= z[a-1]
335 // One has to check the RMAX and RMIN values due to the possible
336 // intersections.
337 //
338 // 1 2 3
339 // ___ ___ ____
340 // 00/ 00/ _ 000/
341 // 0/ 0/ |0 00|
342 // V___ V__+0 00+--
343 // 0000 00000 00000
344 // ---- ----- -----
345 // ------------------------------------ z-axis
346 //
347 //
348 // NOTE: This picture doesn't show all the possible configurations of
349 // a polyhedra having teeth when looking at its profile.
350 // The picture shows only one half of the polyhedra's profile
351 /////////////////////////////////////////////////////////////////////////
352
353 // Experimental code! Not recommended for production, it's incomplete!
354 // The task is to identify invalid combination of z, RMIN and RMAX values
355 // in the case of toothydra :-)
356 //
357 G4int toothIdx;
358
359 for( G4int idx = 1; idx < sections+1; idx++ )
360 {
361 if( z_values[idx-1] > z_values[idx] )
362 {
363 G4double toothdist = std::fabs( z_values[idx-1] - z_values[idx] );
364 G4double aftertoothdist = std::fabs( z_values[idx+1] - z_values[idx] );
365 if( toothdist > aftertoothdist )
366 {
367 // Check for possible intersection
368 //
369 if( RMAX[idx-1] < RMAX[idx+1] || RMIN[idx-1] > RMIN[idx+1] )
370 {
371 // ERROR! The surface conflict!
372 //
373 std::ostringstream message;
374 message << "Unordered sequence of z_values detected." << G4endl
375 << " Conflicting RMAX or RMIN values!" << G4endl
376 << " Check z_values with indexes: "
377 << (idx-1) << " " << idx << " " << (idx+1) << ".";
378 G4Exception( "G4BREPSolidPolyhedra::G4BREPSolidPolyhedra()",
379 "GeomSolids0002", FatalException, message );
380 }
381 }
382 }
383 }
384#endif // G4_EXPERIMENTAL_CODE
385///////////////////////////////////////////////////
386
387 for(G4int a=0;a<sections;a++)
388 {
389 Length = z_values[a+1] - z_values[a];
390
391 if( Length != 0.0 )
392 {
393 TmpAxis= XAxis;
394 TmpAxis.rotateZ(start_angle);
395
396 // L. Broglia: Be careful in the construction of the planes,
397 // see G4FPlane
398 //
399 for( G4int b = 0; b < sides; b++ )
400 {
401 // Create inner side by calculation of points for the planar surface
402 // boundary. The order of the points gives the surface sense -> changed
403 // to explicit sense set-up by R. Chytracek, 12/02/2002
404 // We must check if a pair of two consecutive RMINs is not = 0.0,
405 // this means no inner plane exists!
406 //
407 if( RMIN[a] != 0.0 )
408 {
409 if( RMIN[a+1] != 0.0 )
410 {
411 // Standard case
412 //
413 MaxSurfaceVec[Count] =
414 CreateTrapezoidalSurface( RMIN[a], RMIN[a+1],
415 LocalOrigin, Length,
416 TmpAxis, PartAngle, EInverse );
417 }
418 else
419 {
420 // The special case of r1 > r2 where we end at the
421 // point (0,0,z[a+1])
422 //
423 MaxSurfaceVec[Count] =
424 CreateTriangularSurface( RMIN[a], RMIN[a+1],
425 LocalOrigin, Length,
426 TmpAxis, PartAngle, EInverse );
427 }
428 }
429 else if( RMIN[a+1] != 0.0 )
430 {
431 // The special case of r1 < r2 where we start at the
432 // point ( 0,0,z[a])
433 //
434 MaxSurfaceVec[Count] =
435 CreateTriangularSurface( RMIN[a], RMIN[a+1], LocalOrigin, Length,
436 TmpAxis, PartAngle, EInverse );
437 }
438 else
439 {
440 // Insert nothing into the vector of sufaces, we'll replicate
441 // the vector anyway later
442 //
443 MaxSurfaceVec[Count] = 0;
444
445 // We need to reduce the number of planes by 1,
446 // one we have just skipped
447 //
449 }
450
451 if( MaxSurfaceVec[Count] != 0 )
452 {
453 // Rotate axis back for the other surface point calculation
454 // only in the case any of the Create* methods above have been
455 // called because they modify the passed in TmpAxis
456 //
457 TmpAxis.rotateZ(-PartAngle);
458 }
459
460 Count++;
461
462 // Create outer side
463
464 if( RMAX[a] != 0.0 )
465 {
466 if( RMAX[a+1] != 0.0 )
467 {
468 // Standard case
469 //
470 MaxSurfaceVec[Count] =
471 CreateTrapezoidalSurface( RMAX[a], RMAX[a+1],
472 LocalOrigin, Length,
473 TmpAxis, PartAngle, ENormal );
474 }
475 else
476 {
477 // The special case of r1 > r2 where we end
478 // at the point (0,0,z[a+1])
479 //
480 MaxSurfaceVec[Count] =
481 CreateTriangularSurface( RMAX[a], RMAX[a+1],
482 LocalOrigin, Length,
483 TmpAxis, PartAngle, ENormal );
484 }
485 }
486 else if( RMAX[a+1] != 0.0 )
487 {
488 // The special case of r1 < r2 where we start
489 // at the point ( 0,0,z[a])
490 //
491 MaxSurfaceVec[Count] =
492 CreateTriangularSurface( RMAX[a], RMAX[a+1], LocalOrigin, Length,
493 TmpAxis, PartAngle, ENormal );
494 }
495 else
496 {
497 // Two consecutive RMAX values can't be zero as
498 // it's against the definition of BREP polyhedra
499 //
500 G4Exception( "G4BREPSolidPolyhedra::G4BREPSolidPolyhedra()",
501 "GeomSolids0002", FatalException,
502 "Two consecutive RMAX values cannot be zero!" );
503 }
504
505 Count++;
506 } // End of for loop over sides
507 }
508 else
509 {
510 // Create planar surfaces perpendicular to z-axis
511
512 ESurfaceSense OuterSurfSense, InnerSurfSense;
513
514 if( RMAX[a] != RMAX[a+1] && RMIN[a] != RMIN[a+1] )
515 {
516 // We're about to create a planar surface perpendicular to z-axis
517 // We can have the 8 following configurations here:
518 //
519 // 1. 2. 3. 4.
520 // --+ +-- --+ +--
521 // xx|-> <-|xx xx| |xx
522 // xx+-- --+xx --+ +--
523 // xxxxx xxxxx | |
524 // xxxxx xxxxx +-- --+
525 // xx+-- --+xx |xx xx|
526 // xx|-> <-|xx +-- --+
527 // --+ +--
528 // -------------------------- Z axis
529 //
530 //////////////////////////////////////////////////////////////
531 //////////////////////////////////////////////////////////////
532 //
533 // 5. 6. 7. 8.
534 // --+ +-- --+ +--
535 // xx|-> <-|xx xx|-> <-|xx
536 // --+-- --+-- xx+-- --+xx
537 // <-|xx xx|-> xxxxx xxxxx
538 // +-- --+ --+xx xx+--
539 // <-|xx xx|->
540 // +-- --+
541 // -------------------------- Z axis
542 //
543 // NOTE: The pictures shows only one half of polyhedra!
544 // The arrows show the expected surface normal direction.
545 // The configuration No. 3 and 4 are not valid solids!
546
547 // Eliminate the invalid cases 3 and 4.
548 // At this point is guaranteed that each RMIN[i] < RMAX[i]
549 // where i in in interval 0 < i < num_z_planes-1. So:
550 //
551 if( RMIN[a] > RMAX[a+1] || RMAX[a] < RMIN[a+1] )
552 {
553 std::ostringstream message;
554 message << "The values of RMIN[" << a << "] & RMAX[" << a+1
555 << "] or RMAX[" << a << "] & RMIN[" << a+1 << "]" << G4endl
556 << "generate an invalid configuration of solid: "
557 << GetName() << "!";
558 G4Exception( "G4BREPSolidPolyhedra::G4BREPSolidPolyhedra()",
559 "GeomSolids0002", FatalException, message );
560 }
561
562 // We need to clasify all the cases in order to figure out
563 // the planar surface sense
564 //
565 if( RMAX[a] > RMAX[a+1] )
566 {
567 // Cases 1, 5, 7
568 //
569 if( RMIN[a] < RMIN[a+1] )
570 {
571 // Case 1
572 OuterSurfSense = EInverse;
573 InnerSurfSense = EInverse;
574 }
575 else if( RMAX[a+1] != RMIN[a])
576 {
577 // Case 7
578 OuterSurfSense = EInverse;
579 InnerSurfSense = ENormal;
580 }
581 else
582 {
583 // Case 5
584 OuterSurfSense = EInverse;
585 InnerSurfSense = ENormal;
586 }
587 }
588 else
589 {
590 // Cases 2, 6, 8
591 if( RMIN[a] > RMIN[a+1] )
592 {
593 // Case 2
594 OuterSurfSense = ENormal;
595 InnerSurfSense = ENormal;
596 }
597 else if( RMIN[a+1] != RMAX[a] )
598 {
599 // Case 8
600 OuterSurfSense = ENormal;
601 InnerSurfSense = EInverse;
602 }
603 else
604 {
605 // Case 6
606 OuterSurfSense = ENormal;
607 InnerSurfSense = EInverse;
608 }
609 }
610
611 TmpAxis= XAxis;
612 TmpAxis.rotateZ(start_angle);
613
614 // Compute the outer planar surface
615 //
616 MaxSurfaceVec[Count] =
617 ComputePlanarSurface( RMAX[a], RMAX[a+1], LocalOrigin, TmpAxis,
618 sides, PartAngle, OuterSurfSense );
619 if( MaxSurfaceVec[Count] == 0 )
620 {
621 // No surface was created
622 //
624 }
625 Count++;
626
627 TmpAxis= XAxis;
628 TmpAxis.rotateZ(start_angle);
629
630 // Compute the inner planar surface
631 //
632 MaxSurfaceVec[Count] =
633 ComputePlanarSurface( RMIN[a], RMIN[a+1], LocalOrigin, TmpAxis,
634 sides, PartAngle, InnerSurfSense );
635 if( MaxSurfaceVec[Count] == 0 )
636 {
637 // No surface was created
638 //
640 }
641 Count++;
642
643 // Since we can create here at maximum 2 surfaces
644 // we need to reflect this in the total
645 //
646 nb_of_surfaces -= (2*(sides-1));
647 }
648 else
649 {
650 // The case where only one of the radius values has changed
651 //
652 // RMAX RMIN
653 // change change
654 //
655 // 1 2 3 4
656 // --+ +-- ----- -----
657 // 00|-> <-|00 00000 00000
658 // 00+-- --+00 --+00 00+--
659 // 00000 00000 <-|00 00|->
660 // +-- --+
661 // --------------------------- Z axis
662 //
663 // NOTE: The picture shows only one half of polyhedra!
664
665 G4double R1, R2;
666 ESurfaceSense SurfSense;
667
668 // The case by case classification
669 //
670 if( RMAX[a] != RMAX[a+1] )
671 {
672 // Cases 1, 2
673 //
674 R1 = RMAX[a];
675 R2 = RMAX[a+1];
676 if( R1 > R2 )
677 {
678 // Case 1
679 //
680 SurfSense = EInverse;
681 }
682 else
683 {
684 // Case 2
685 //
686 SurfSense = ENormal;
687 }
688 }
689 else if(RMIN[a] != RMIN[a+1])
690 {
691 // Cases 3, 4
692 //
693 R1 = RMIN[a];
694 R2 = RMIN[a+1];
695 if( R1 > R2 )
696 {
697 // Case 3
698 //
699 SurfSense = ENormal;
700 }
701 else
702 {
703 // Case 4
704 //
705 SurfSense = EInverse;
706 }
707 }
708 else
709 {
710 std::ostringstream message;
711 message << "Error in construction." << G4endl
712 << " Exactly the same z, rmin and rmax given for"
713 << G4endl
714 << " consecutive indices, " << a << " and " << a+1;
715 G4Exception( "G4BREPSolidPolyhedra::G4BREPSolidPolyhedra()",
716 "GeomSolids1001", JustWarning, message );
717 continue;
718 }
719 TmpAxis= XAxis;
720 TmpAxis.rotateZ(start_angle);
721
722 MaxSurfaceVec[Count] =
723 ComputePlanarSurface( R1, R2, LocalOrigin, TmpAxis,
724 sides, PartAngle, SurfSense );
725 if( MaxSurfaceVec[Count] == 0 )
726 {
727 // No surface was created
728 //
730 }
731 Count++;
732
733 // Since we can create here at maximum 1 surface
734 // we need to reflect this in the total
735 //
736 nb_of_surfaces -= ((2*sides) - 1);
737 }
738 } // End of if( Length != 0.0 )
739
740 LocalOrigin = LocalOrigin + (Length*Axis);
741
742 } // End of for loop over z sections
743
744 if(opening_angle >= 2*pi-perMillion)
745 {
746 // Create the end planes for the configuration where delta phi >= 2*PI
747
748 TmpAxis = XAxis;
749 TmpAxis.rotateZ(start_angle);
750
751 MaxSurfaceVec[Count] =
752 ComputePlanarSurface( RMIN[0], RMAX[0], Origin, TmpAxis,
753 sides, PartAngle, ENormal );
754
755 if( MaxSurfaceVec[Count] == 0 )
756 {
757 // No surface was created
758 //
760 }
761 Count++;
762
763 // Reset plane axis
764 //
765 TmpAxis = XAxis;
766 TmpAxis.rotateZ(start_angle);
767
768 MaxSurfaceVec[Count] =
769 ComputePlanarSurface( RMIN[sections], RMAX[sections],
770 LocalOrigin, TmpAxis,
771 sides, PartAngle, EInverse );
772
773 if( MaxSurfaceVec[Count] == 0 )
774 {
775 // No surface was created
776 //
778 }
779 Count++;
780 }
781 else
782 {
783 // If delta phi < 2*PI then create a single boundary
784 // (case with RMIN=0 included)
785
786 // Create the lateral planars
787 //
788 TmpAxis = XAxis;
789 G4Vector3D TmpAxis2 = XAxis;
790 TmpAxis.rotateZ(start_angle);
791 TmpAxis2.rotateZ(start_angle);
792 TmpAxis2.rotateZ(start_angle);
793
794 LocalOrigin = Origin;
795 G4int points = sections*2+2;
796 G4int PointCount = 0;
797
798 G4Point3DVector GapPointList(points);
799 G4Point3DVector GapPointList2(points);
800
801 for(G4int d=0;d<sections+1;d++)
802 {
803 GapPointList[PointCount] = LocalOrigin + (RMAX[d]*TmpAxis);
804 GapPointList[points-1-PointCount] = LocalOrigin + (RMIN[d]*TmpAxis);
805
806 GapPointList2[PointCount] = LocalOrigin + (RMAX[d]*TmpAxis2);
807 GapPointList2[points-1-PointCount] = LocalOrigin + (RMIN[d]*TmpAxis2);
808
809 PointCount++;
810
811 Length = z_values[d+1] - z_values[d];
812 LocalOrigin = LocalOrigin+(Length*Axis);
813 }
814
815 // Add the lateral planars to the surfaces list and set/reverse sense
816 //
817 MaxSurfaceVec[Count++] = new G4FPlane( &GapPointList, 0, ENormal );
818 MaxSurfaceVec[Count++] = new G4FPlane( &GapPointList2, 0, EInverse );
819
820 TmpAxis = XAxis;
821 TmpAxis.rotateZ(start_angle);
822 TmpAxis.rotateZ(opening_angle);
823
824 // Create end planes
825 //
826 G4Point3DVector EndPointList ((sides+1)*2);
827 G4Point3DVector EndPointList2((sides+1)*2);
828
829 for(G4int c=0;c<sides+1;c++)
830 {
831 // outer polylines for origin end and opposite side
832 EndPointList[c] = Origin + (RMAX[0] * TmpAxis);
833 EndPointList[(sides+1)*2-1-c] = Origin + (RMIN[0] * TmpAxis);
834 EndPointList2[c] = LocalOrigin + (RMAX[sections] * TmpAxis);
835 EndPointList2[(sides+1)*2-1-c] = LocalOrigin + (RMIN[sections] * TmpAxis);
836 TmpAxis.rotateZ(-PartAngle);
837 }
838
839 // Add the end planes to the surfaces list
840 // Note the surface sense in this case is reversed
841 // It's because here we have created the end planes in reversed order
842 // than it's done by ComputePlanarSurface() method
843 //
844 if(RMAX[0]-RMIN[0] >= perMillion)
845 {
846 MaxSurfaceVec[Count] = new G4FPlane( &EndPointList, 0, EInverse );
847 }
848 else
849 {
850 MaxSurfaceVec[Count] = 0;
852 }
853
854 Count++;
855
856 if(RMAX[sections]-RMIN[sections] >= perMillion)
857 {
858 MaxSurfaceVec[Count] = new G4FPlane( &EndPointList2, 0, ENormal );
859 }
860 else
861 {
862 MaxSurfaceVec[Count] = 0;
864 }
865 }
866
867 // Now let's replicate the relevant surfaces into
868 // G4BREPSolid's vector of surfaces
869 //
871 G4int sf = 0; G4int zeroCount = 0;
872 for( G4int srf = 0; srf < MaxNbOfSurfaces; srf++ )
873 {
874 if( MaxSurfaceVec[srf] != 0 )
875 {
876 if( sf < nb_of_surfaces )
877 {
878 SurfaceVec[sf] = MaxSurfaceVec[srf];
879 }
880 sf++;
881 }
882 else
883 {
884 zeroCount++;
885 }
886 }
887
888 if( sf != nb_of_surfaces )
889 {
890 std::ostringstream message;
891 message << "Bad number of surfaces!" << G4endl
892 << " sf : " << sf
893 << " nb_of_surfaces: " << nb_of_surfaces
894 << " Count : " << Count;
895 G4Exception( "G4BREPSolidPolyhedra::G4BREPSolidPolyhedra()",
896 "GeomSolids0002", FatalException, message);
897 }
898
899 // Clean up the temporary vector of surfaces
900 //
901 delete [] MaxSurfaceVec;
902
903 Initialize();
904}
905
907{
908 // Calc bounding box for solids and surfaces
909 // Convert concave planes to convex
910 //
911 ShortestDistance=1000000;
913 if(!Box || !AxisBox)
914 IsConvex();
915
916 CalcBBoxes();
917}
918
920{
921 Active(1);
922 ((G4BREPSolidPolyhedra*)this)->intersectionDistance=kInfinity;
923 StartInside(0);
924 for(register G4int a=0;a<nb_of_surfaces;a++)
925 SurfaceVec[a]->Reset();
926 ShortestDistance = kInfinity;
927}
928
930{
931 // This function find if the point Pt is inside,
932 // outside or on the surface of the solid
933
934 const G4double sqrHalfTolerance = kCarTolerance*kCarTolerance*0.25;
935
936 G4Vector3D v(1, 0, 0.01);
937 G4Vector3D Pttmp(Pt);
938 G4Vector3D Vtmp(v);
939 G4Ray r(Pttmp, Vtmp);
940
941 // Check if point is inside the Polyhedra bounding box
942 //
943 if( !GetBBox()->Inside(Pttmp) )
944 {
945 return kOutside;
946 }
947
948 // Set the surfaces to active again
949 //
950 Reset();
951
952 // Test if the bounding box of each surface is intersected
953 // by the ray. If not, the surface is deactivated.
954 //
956
957 G4int hits=0, samehit=0;
958
959 for(G4int a=0; a < nb_of_surfaces; a++)
960 {
961 G4Surface* surface = SurfaceVec[a];
962
963 if(surface->IsActive())
964 {
965 // count the number of intersections.
966 // if this number is odd, the start of the ray is
967 // inside the volume bounded by the surfaces, so
968 // increment the number of intersection by 1 if the
969 // point is not on the surface and if this intersection
970 // was not found before
971 //
972 if( (surface->Intersect(r)) & 1 )
973 {
974 // test if the point is on the surface
975 //
976 if(surface->GetDistance() < sqrHalfTolerance)
977 {
978 return kSurface;
979 }
980 // test if this intersection was found before
981 //
982 for(G4int i=0; i<a; i++)
983 {
984 if(surface->GetDistance() == SurfaceVec[i]->GetDistance())
985 {
986 samehit++;
987 break;
988 }
989 }
990
991 // count the number of surfaces intersected by the ray
992 //
993 if(!samehit)
994 {
995 hits++;
996 }
997 }
998 }
999 }
1000
1001 // if the number of surfaces intersected is odd,
1002 // the point is inside the solid
1003 //
1004 return ( (hits&1) ? kInside : kOutside );
1005}
1006
1009{
1010 // This function calculates the normal of the closest surface
1011 // to the given point
1012 // Note : the sense of the normal depends on the sense of the surface
1013
1014 G4int iplane;
1015 G4bool normflag = false;
1016 const G4double sqrHalfTolerance = kCarTolerance*kCarTolerance*0.25;
1017
1018 // Determine if the point is on the surface
1019 //
1020 G4double minDist = kInfinity;
1021 G4int normPlane = 0;
1022 for(iplane = 0; iplane < nb_of_surfaces; iplane++)
1023 {
1024 G4double dist = std::fabs(SurfaceVec[iplane]->HowNear(Pt));
1025 if( minDist > dist )
1026 {
1027 minDist = dist;
1028 normPlane = iplane;
1029 }
1030 if( dist < sqrHalfTolerance)
1031 {
1032 // the point is on this surface, so take this as the
1033 // the surface to consider for computing the normal
1034 //
1035 normflag = true;
1036 break;
1037 }
1038 }
1039
1040 // Calculate the normal at this point, if the point is on the
1041 // surface, otherwise compute the normal to the closest surface
1042 //
1043 if ( normflag ) // point on surface
1044 {
1045 G4ThreeVector norm = SurfaceVec[iplane]->SurfaceNormal(Pt);
1046 return norm.unit();
1047 }
1048 else // point not on surface
1049 {
1050 G4FPlane* nPlane = (G4FPlane*)(SurfaceVec[normPlane]);
1051 G4ThreeVector hitPt = nPlane->GetSrfPoint();
1052 G4ThreeVector hitNorm = nPlane->SurfaceNormal(hitPt);
1053 return hitNorm.unit();
1054 }
1055}
1056
1058{
1059 // Calculates the shortest distance ("safety") from a point
1060 // outside the solid to any boundary of this solid.
1061 // Return 0 if the point is already inside.
1062
1063 G4double *dists = new G4double[nb_of_surfaces];
1064 G4int a;
1065
1066 // Set the surfaces to active again
1067 //
1068 Reset();
1069
1070 // compute the shortest distance of the point to each surfaces
1071 // Be careful : it's a signed value
1072 //
1073 for(a=0; a< nb_of_surfaces; a++)
1074 dists[a] = SurfaceVec[a]->HowNear(Pt);
1075
1076 G4double Dist = kInfinity;
1077
1078 // if dists[] is positive, the point is outside
1079 // so take the shortest of the shortest positive distances
1080 // dists[] can be equal to 0 : point on a surface
1081 // ( Problem with the G4FPlane : there is no inside and no outside...
1082 // So, to test if the point is inside to return 0, utilize the Inside
1083 // function. But I don`t know if it is really needed because dToIn is
1084 // called only if the point is outside )
1085 //
1086 for(a = 0; a < nb_of_surfaces; a++)
1087 if( std::fabs(Dist) > std::fabs(dists[a]) )
1088 //if( dists[a] >= 0)
1089 Dist = dists[a];
1090
1091 delete[] dists;
1092
1093 if(Dist == kInfinity)
1094 {
1095 // the point is inside the solid or on a surface
1096 //
1097 return 0;
1098 }
1099 else
1100 {
1101 return std::fabs(Dist);
1102 }
1103}
1104
1107 register const G4ThreeVector& V) const
1108{
1109 // Calculates the distance from a point outside the solid
1110 // to the solid`s boundary along a specified direction vector.
1111 //
1112 // Note : Intersections with boundaries less than the
1113 // tolerance must be ignored if the direction
1114 // is away from the boundary
1115
1116 G4int a;
1117
1118 // Set the surfaces to active again
1119 //
1120 Reset();
1121
1122 const G4double sqrHalfTolerance = kCarTolerance*kCarTolerance*0.25;
1123 G4Vector3D Pttmp(Pt);
1124 G4Vector3D Vtmp(V);
1125 G4Ray r(Pttmp, Vtmp);
1126
1127 // Test if the bounding box of each surface is intersected
1128 // by the ray. If not, the surface become deactive.
1129 //
1131
1132 ShortestDistance = kInfinity;
1133
1134 for(a=0; a< nb_of_surfaces; a++)
1135 {
1136 if( SurfaceVec[a]->IsActive() )
1137 {
1138 // test if the ray intersect the surface
1139 //
1140 if( SurfaceVec[a]->Intersect(r) )
1141 {
1142 G4double surfDistance = SurfaceVec[a]->GetDistance();
1143
1144 // if more than 1 surface is intersected,
1145 // take the nearest one
1146 //
1147 if( surfDistance < ShortestDistance )
1148 {
1149 if( surfDistance > sqrHalfTolerance )
1150 {
1151 ShortestDistance = surfDistance;
1152 }
1153 else
1154 {
1155 // the point is within the boundary
1156 // ignore it if the direction is away from the boundary
1157 //
1158 G4Vector3D Norm = SurfaceVec[a]->SurfaceNormal(Pttmp);
1159
1160 if( (Norm * Vtmp) < 0 )
1161 {
1162 ShortestDistance = surfDistance;
1163// ShortestDistance = surfDistance==0
1164// ? sqrHalfTolerance
1165// : surfDistance;
1166 }
1167 }
1168 }
1169 }
1170 }
1171 }
1172
1173 // Be careful !
1174 // SurfaceVec->Distance is in fact the squared distance
1175 //
1176 if(ShortestDistance != kInfinity)
1177 {
1178 return std::sqrt(ShortestDistance);
1179 }
1180 else // no intersection
1181 {
1182 return kInfinity;
1183 }
1184}
1185
1188 register const G4ThreeVector& V,
1189 const G4bool,
1190 G4bool *validNorm,
1191 G4ThreeVector * ) const
1192{
1193 // Calculates the distance from a point inside the solid
1194 // to the solid`s boundary along a specified direction vector.
1195 // Return 0 if the point is already outside (even number of
1196 // intersections greater than the tolerance).
1197 //
1198 // Note : If the shortest distance to a boundary is less
1199 // than the tolerance, it is ignored. This allows
1200 // for a point within a tolerant boundary to leave
1201 // immediately
1202
1203 G4int parity = 0;
1204
1205 // Set the surfaces to active again
1206 //
1207 Reset();
1208
1209 const G4double sqrHalfTolerance = kCarTolerance*kCarTolerance*0.25;
1210 G4Vector3D Ptv = Pt;
1211 G4int a;
1212
1213 // I don`t understand this line
1214 //
1215 if(validNorm)
1216 *validNorm=false;
1217
1218 G4Vector3D Pttmp(Pt);
1219 G4Vector3D Vtmp(V);
1220
1221 G4Ray r(Pttmp, Vtmp);
1222
1223 // Test if the bounding box of each surface is intersected
1224 // by the ray. If not, the surface become deactive.
1225 //
1227
1228 ShortestDistance = kInfinity; // this is actually the square of the distance
1229
1230 for(a=0; a< nb_of_surfaces; a++)
1231 {
1232 G4double surfDistance = SurfaceVec[a]->GetDistance();
1233
1234 if(SurfaceVec[a]->IsActive())
1235 {
1236 G4int intersects = SurfaceVec[a]->Intersect(r);
1237
1238 // test if the ray intersects the surface
1239 //
1240 if( intersects != 0 )
1241 {
1242 parity += 1;
1243
1244 // if more than 1 surface is intersected, take the nearest one
1245 //
1246 if( surfDistance < ShortestDistance )
1247 {
1248 if( surfDistance > sqrHalfTolerance )
1249 {
1250 ShortestDistance = surfDistance;
1251 }
1252 else
1253 {
1254 // the point is within the boundary: ignore it
1255 //
1256 parity -= 1;
1257 }
1258 }
1259 }
1260 }
1261 }
1262
1263 G4double distance = 0.;
1264
1265 // Be careful !
1266 // SurfaceVec->Distance is in fact the squared distance
1267 //
1268 // This condition was changed in order to give not zero answer
1269 // when particle is passing the border of two Touching Surfaces
1270 // and the distance to this surfaces is not zero.
1271 // parity is for the points on the boundary,
1272 // parity is counting only surfDistance<kCarTolerance/2.
1273 //
1274 // if((ShortestDistance != kInfinity) && (parity&1))
1275 //
1276 //
1277 if((ShortestDistance != kInfinity) || (parity&1))
1278 {
1279 distance = std::sqrt(ShortestDistance);
1280 }
1281
1282 return distance;
1283}
1284
1286{
1287 // Calculates the shortest distance ("safety") from a point
1288 // inside the solid to any boundary of this solid.
1289 // Return 0 if the point is already outside.
1290
1291 G4double *dists = new G4double[nb_of_surfaces];
1292 G4int a;
1293
1294 // Set the surfaces to active again
1295 //
1296 Reset();
1297
1298 // calculate the shortest distance of the point to each surfaces
1299 // Be careful : it's a signed value
1300 //
1301 for(a=0; a< nb_of_surfaces; a++)
1302 {
1303 dists[a] = SurfaceVec[a]->HowNear(Pt);
1304 }
1305
1306 G4double Dist = kInfinity;
1307
1308 // if dists[] is negative, the point is inside
1309 // so take the shortest of the shortest negative distances
1310 // dists[] can be equal to 0 : point on a surface
1311 // ( Problem with the G4FPlane : there is no inside and no outside...
1312 // So, to test if the point is outside to return 0, utilize the Inside
1313 // function. But I don`t know if it is really needed because dToOut is
1314 // called only if the point is inside )
1315
1316 for(a = 0; a < nb_of_surfaces; a++)
1317 {
1318 if( std::fabs(Dist) > std::fabs(dists[a]) )
1319 {
1320 //if( dists[a] <= 0)
1321 Dist = dists[a];
1322 }
1323 }
1324
1325 delete[] dists;
1326
1327 if(Dist == kInfinity)
1328 {
1329 // the point is ouside the solid or on a surface
1330 //
1331 return 0;
1332 }
1333 else
1334 {
1335 // return Dist;
1336 return std::fabs(Dist);
1337 }
1338}
1339
1341{
1342 return new G4BREPSolidPolyhedra(*this);
1343}
1344
1345std::ostream& G4BREPSolidPolyhedra::StreamInfo(std::ostream& os) const
1346{
1347 // Streams solid contents to output stream.
1348
1350 << "\n start_angle: " << constructorParams.start_angle
1351 << "\n opening_angle: " << constructorParams.opening_angle
1352 << "\n sides: " << constructorParams.sides
1353 << "\n num_z_planes: " << constructorParams.num_z_planes
1354 << "\n z_start: " << constructorParams.z_start
1355 << "\n z_values: ";
1356 G4int idx;
1357 for( idx = 0; idx < constructorParams.num_z_planes; idx++ )
1358 {
1359 os << constructorParams.z_values[idx] << " ";
1360 }
1361 os << "\n RMIN: ";
1362 for( idx = 0; idx < constructorParams.num_z_planes; idx++ )
1363 {
1364 os << constructorParams.RMIN[idx] << " ";
1365 }
1366 os << "\n RMAX: ";
1367 for( idx = 0; idx < constructorParams.num_z_planes; idx++ )
1368 {
1369 os << constructorParams.RMAX[idx] << " ";
1370 }
1371 os << "\n-----------------------------------------------------------\n";
1372
1373 return os;
1374}
1375
1376G4Surface*
1377G4BREPSolidPolyhedra::CreateTrapezoidalSurface( G4double r1,
1378 G4double r2,
1379 const G4Point3D& origin,
1380 G4double distance,
1381 G4Vector3D& xAxis,
1382 G4double partAngle,
1383 ESurfaceSense sense )
1384{
1385 // The surface to be returned
1386 //
1387 G4Surface* trapsrf = 0;
1388 G4Point3DVector PointList(4);
1389 G4Vector3D zAxis(0,0,1);
1390
1391 PointList[0] = origin + ( r1 * xAxis);
1392 PointList[3] = origin + ( distance * zAxis) + (r2 * xAxis);
1393
1394 xAxis.rotateZ( partAngle );
1395
1396 PointList[2] = origin + ( distance * zAxis) + (r2 * xAxis);
1397 PointList[1] = origin + ( r1 * xAxis);
1398
1399 // Return the planar trapezoidal surface
1400 //
1401 trapsrf = new G4FPlane( &PointList, 0, sense );
1402
1403 return trapsrf;
1404}
1405
1406G4Surface*
1407G4BREPSolidPolyhedra::CreateTriangularSurface( G4double r1,
1408 G4double r2,
1409 const G4Point3D& origin,
1410 G4double distance,
1411 G4Vector3D& xAxis,
1412 G4double partAngle,
1413 ESurfaceSense sense )
1414{
1415 // The surface to be returned
1416 //
1417 G4Surface* trapsrf = 0;
1418 G4Point3DVector PointList(3);
1419 G4Vector3D zAxis(0,0,1);
1420
1421 PointList[0] = origin + ( r1 * xAxis);
1422 PointList[2] = origin + ( distance * zAxis) + (r2 * xAxis);
1423
1424 xAxis.rotateZ( partAngle );
1425
1426 if( r1 < r2 )
1427 {
1428 PointList[1] = origin + ( distance * zAxis) + (r2 * xAxis);
1429 }
1430 else
1431 {
1432 PointList[1] = origin + ( r1 * xAxis);
1433 }
1434
1435 // Return the planar trapezoidal surface
1436 //
1437 trapsrf = new G4FPlane( &PointList, 0, sense );
1438
1439 return trapsrf;
1440}
1441
1442G4Surface*
1443G4BREPSolidPolyhedra::ComputePlanarSurface( G4double r1,
1444 G4double r2,
1445 const G4Point3D& origin,
1446 G4Vector3D& xAxis,
1447 G4int sides,
1448 G4double partAngle,
1449 ESurfaceSense sense )
1450{
1451 // This method can be called only when r1 != r2,
1452 // otherwise it returns 0 which means that no surface can be
1453 // created out of the given radius pair.
1454 // This method requires the xAxis to be pre-rotated properly.
1455
1456 G4Point3DVector OuterPointList( sides );
1457 G4Point3DVector InnerPointList( sides );
1458
1459 G4double rIn, rOut;
1460 G4Surface* planarSrf = 0;
1461
1462 if( r1 < r2 )
1463 {
1464 rIn = r1;
1465 rOut = r2;
1466 }
1467 else if( r1 > r2 )
1468 {
1469 rIn = r2;
1470 rOut = r1;
1471 }
1472 else
1473 {
1474 // Invalid precondition, the radius values are r1 == r2,
1475 // which means we can create only polyline but no surface
1476 //
1477 return 0;
1478 }
1479
1480 for( G4int pidx = 0; pidx < sides; pidx++ )
1481 {
1482 // Outer polyline
1483 //
1484 OuterPointList[pidx] = origin + ( rOut * xAxis);
1485 // Inner polyline
1486 //
1487 InnerPointList[pidx] = origin + ( rIn * xAxis);
1488 xAxis.rotateZ( partAngle );
1489 }
1490
1491 if( rIn != 0.0 && rOut != 0.0 )
1492 {
1493 // Standard case
1494 //
1495 planarSrf = new G4FPlane( &OuterPointList, &InnerPointList, sense );
1496 }
1497 else if( rOut != 0.0 )
1498 {
1499 // Special case where inner radius is zero so no polyline
1500 // is actually created
1501 //
1502 planarSrf = new G4FPlane( &OuterPointList, 0, sense );
1503 }
1504 else
1505 {
1506 // No surface being created
1507 // This should not happen as filtered out by precondition check above
1508 }
1509
1510 return planarSrf;
1511}
1512
1513// In graphics_reps:
1514
1515#include "G4Polyhedron.hh"
1516
1518{
1519 return new G4PolyhedronPgon( constructorParams.start_angle,
1520 constructorParams.opening_angle,
1521 constructorParams.sides,
1522 constructorParams.num_z_planes,
1523 constructorParams.z_values,
1524 constructorParams.RMIN,
1525 constructorParams.RMAX);
1526}
ESurfaceSense
@ JustWarning
@ FatalException
std::vector< G4Point3D > G4Point3DVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
Hep3Vector unit() const
G4ThreeVector SurfaceNormal(const G4ThreeVector &) const
G4Polyhedron * CreatePolyhedron() const
std::ostream & StreamInfo(std::ostream &os) const
G4BREPSolidPolyhedra & operator=(const G4BREPSolidPolyhedra &rhs)
EInside Inside(register const G4ThreeVector &Pt) const
G4BREPSolidPolyhedra(const G4String &name, G4double phi1, G4double dphi, G4int sides, G4int num_z_planes, G4double z_start, G4double z_values[], G4double RMIN[], G4double RMAX[])
G4double DistanceToIn(const G4ThreeVector &) const
G4double DistanceToOut(register const G4ThreeVector &Pt, register const G4ThreeVector &V, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const
G4int Active() const
G4Surface ** SurfaceVec
Definition: G4BREPSolid.hh:231
G4BREPSolid & operator=(const G4BREPSolid &rhs)
Definition: G4BREPSolid.cc:138
G4int StartInside() const
G4int Intersect(register const G4Ray &) const
G4bool IsConvex()
Definition: G4BREPSolid.cc:460
G4int nb_of_surfaces
Definition: G4BREPSolid.hh:229
virtual std::ostream & StreamInfo(std::ostream &os) const
void CheckSurfaceNormals()
Definition: G4BREPSolid.cc:213
virtual void CalcBBoxes()
const G4BoundingBox3D * GetBBox() const
static G4double ShortestDistance
Definition: G4BREPSolid.hh:221
void TestSurfaceBBoxes(register const G4Ray &) const
const G4String & GetName() const
G4Point3D GetSrfPoint() const
G4Vector3D SurfaceNormal(const G4Point3D &Pt) const
Definition: G4Ray.hh:49
virtual G4int Intersect(const G4Ray &)
Definition: G4Surface.cc:170
G4int IsActive() const
virtual G4Vector3D SurfaceNormal(const G4Point3D &Pt) const =0
virtual G4double HowNear(const G4Vector3D &x) const
Definition: G4Surface.cc:283
G4double GetDistance() const
G4double kCarTolerance
Definition: G4VSolid.hh:307
BasicVector3D< T > & rotateZ(T a)
EInside
Definition: geomdefs.hh:58
@ kInside
Definition: geomdefs.hh:58
@ kOutside
Definition: geomdefs.hh:58
@ kSurface
Definition: geomdefs.hh:58
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41