Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BREPSolidPCone.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// G4BREPSolidPCone.cc
32//
33// ----------------------------------------------------------------------
34// The polyconical solid G4BREPSolidPCone is a shape defined by a set of
35// inner and outer conical or cylindrical surface sections and two planes
36// perpendicular to the Z axis. Each conical surface is defined by its
37// radius at two different planes perpendicular to the Z-axis. Inner and
38// outer conical surfaces are defined using common Z planes.
39// ----------------------------------------------------------------------
40
41#include "G4Types.hh"
42#include <sstream>
43
44#include "G4BREPSolidPCone.hh"
46#include "G4SystemOfUnits.hh"
48#include "G4FConicalSurface.hh"
49#include "G4CircularCurve.hh"
50#include "G4FPlane.hh"
51
52typedef enum
53{
55 ENormal = 1
57
59 G4double start_angle,
60 G4double opening_angle,
61 G4int num_z_planes, // sections,
62 G4double z_start,
63 G4double z_values[],
64 G4double RMIN[],
65 G4double RMAX[] )
66 : G4BREPSolid(name)
67{
68 // Save contructor parameters
69 //
70 constructorParams.start_angle = start_angle;
71 constructorParams.opening_angle = opening_angle;
72 constructorParams.num_z_planes = num_z_planes;
73 constructorParams.z_start = z_start;
74 constructorParams.z_values = 0;
75 constructorParams.RMIN = 0;
76 constructorParams.RMAX = 0;
77
78 if( num_z_planes > 0 )
79 {
80 constructorParams.z_values = new G4double[num_z_planes];
81 constructorParams.RMIN = new G4double[num_z_planes];
82 constructorParams.RMAX = new G4double[num_z_planes];
83 for( G4int idx = 0; idx < num_z_planes; ++idx )
84 {
85 constructorParams.z_values[idx] = z_values[idx];
86 constructorParams.RMIN[idx] = RMIN[idx];
87 constructorParams.RMAX[idx] = RMAX[idx];
88 }
89 }
90
91 active=1;
92 InitializePCone();
93}
94
96 : G4BREPSolid(a)
97{
98 constructorParams.start_angle = 0.;
99 constructorParams.opening_angle = 0.;
100 constructorParams.num_z_planes = 0;
101 constructorParams.z_start = 0.;
102 constructorParams.z_values = 0;
103 constructorParams.RMIN = 0;
104 constructorParams.RMAX = 0;
105}
106
108{
109 if( constructorParams.num_z_planes > 0 )
110 {
111 delete [] constructorParams.z_values;
112 delete [] constructorParams.RMIN;
113 delete [] constructorParams.RMAX;
114 }
115}
116
118 : G4BREPSolid(rhs)
119{
120 constructorParams.start_angle = rhs.constructorParams.start_angle;
121 constructorParams.opening_angle = rhs.constructorParams.opening_angle;
122 constructorParams.num_z_planes = rhs.constructorParams.num_z_planes;
123 constructorParams.z_start = rhs.constructorParams.z_start;
124 constructorParams.z_values = 0;
125 constructorParams.RMIN = 0;
126 constructorParams.RMAX = 0;
127 G4int nplanes = constructorParams.num_z_planes;
128 if( nplanes > 0 )
129 {
130 constructorParams.z_values = new G4double[nplanes];
131 constructorParams.RMIN = new G4double[nplanes];
132 constructorParams.RMAX = new G4double[nplanes];
133 for( G4int idx = 0; idx < nplanes; ++idx )
134 {
135 constructorParams.z_values[idx] = rhs.constructorParams.z_values[idx];
136 constructorParams.RMIN[idx] = rhs.constructorParams.RMIN[idx];
137 constructorParams.RMAX[idx] = rhs.constructorParams.RMAX[idx];
138 }
139 }
140
141 InitializePCone();
142}
143
146{
147 // Check assignment to self
148 //
149 if (this == &rhs) { return *this; }
150
151 // Copy base class data
152 //
154
155 // Copy data
156 //
157 constructorParams.start_angle = rhs.constructorParams.start_angle;
158 constructorParams.opening_angle = rhs.constructorParams.opening_angle;
159 constructorParams.num_z_planes = rhs.constructorParams.num_z_planes;
160 constructorParams.z_start = rhs.constructorParams.z_start;
161 G4int nplanes = constructorParams.num_z_planes;
162 if( nplanes > 0 )
163 {
164 delete [] constructorParams.z_values;
165 delete [] constructorParams.RMIN;
166 delete [] constructorParams.RMAX;
167 constructorParams.z_values = new G4double[nplanes];
168 constructorParams.RMIN = new G4double[nplanes];
169 constructorParams.RMAX = new G4double[nplanes];
170 for( G4int idx = 0; idx < nplanes; ++idx )
171 {
172 constructorParams.z_values[idx] = rhs.constructorParams.z_values[idx];
173 constructorParams.RMIN[idx] = rhs.constructorParams.RMIN[idx];
174 constructorParams.RMAX[idx] = rhs.constructorParams.RMAX[idx];
175 }
176 }
177
178 InitializePCone();
179
180 return *this;
181}
182
183void G4BREPSolidPCone::InitializePCone()
184{
185 G4double opening_angle = constructorParams.opening_angle;
186 G4int num_z_planes = constructorParams.num_z_planes;
187 G4double z_start = constructorParams.z_start;
188 G4double* z_values = constructorParams.z_values;
189 G4double* RMIN = constructorParams.RMIN;
190 G4double* RMAX = constructorParams.RMAX;
191
192 G4int sections= constructorParams.num_z_planes-1;
193 nb_of_surfaces = 2*sections+2;
194
196 G4ThreeVector Axis(0,0,1);
197 G4ThreeVector Origin(0,0,z_start);
198 G4double Length;
199 G4ThreeVector LocalOrigin(0,0,z_start);
200 G4int a, b = 0;
201
202 G4ThreeVector PlaneAxis(0, 0, 1);
203 G4ThreeVector PlaneDir (0, 1, 0);
204
205 ///////////////////////////////////////////////////
206 // Test delta phi
207
208 // At the moment (11/03/2002) the phi section is not implemented
209 // so we take a G4 application down if there is a request for such
210 // a configuration
211 //
212 if( opening_angle < 2*pi-perMillion )
213 {
214 G4Exception("G4BREPSolidPCone::G4BREPSolidPCone()",
215 "GeomSolids0001", FatalException,
216 "Sorry, phi-section not supported yet, try to use G4Polycone instead!");
217 }
218
219 ///////////////////////////////////////////////////
220 // Test the validity of the R values
221
222 // RMIN[0] and RMIN[num_z_planes-1] cannot be = 0
223 // when RMIN[0] or RMIN[num_z_planes-1] are = 0
224 //
225 if( ((RMIN[0] == 0) && (RMAX[0] == 0)) ||
226 ((RMIN[num_z_planes-1] == 0) && (RMAX[num_z_planes-1] == 0)) )
227 G4Exception("G4BREPSolidPCone::G4BREPSolidPCone()",
228 "GeomSolids0002", FatalException,
229 "RMIN at the extremities cannot be 0 when RMAX=0 !");
230
231 // only RMAX[0] and RMAX[num_z_planes-1] can be = 0
232 //
233 for(a = 1; a < num_z_planes-1; a++)
234 if (RMAX[a] == 0)
235 G4Exception("G4BREPSolidPCone::G4BREPSolidPCone()",
236 "GeomSolids0002", FatalException,
237 "RMAX inside the solid cannot be 0 !");
238
239 // RMAX[a] must be greater than RMIN[a]
240 //
241 for(a = 1; a < num_z_planes-1; a++)
242 if (RMIN[a] >= RMAX[a])
243 G4Exception("G4BREPSolidPCone::G4BREPSolidPCone()",
244 "GeomSolids0002", FatalException,
245 "RMAX must be greater than RMIN in the middle Z planes !");
246
247 if( (RMIN[num_z_planes-1] > RMAX[num_z_planes-1] )
248 || (RMIN[0] > RMAX[0]) )
249 G4Exception("G4BREPSolidPCone::G4BREPSolidPCone()",
250 "GeomSolids0002", FatalException,
251 "RMAX must be greater or equal than RMIN at the ends !");
252
253 // Create surfaces
254
255 for( a = 0; a < sections; a++)
256 {
257 // Surface length
258 //
259 Length = z_values[a+1] - z_values[a];
260
261 if( Length == 0 )
262 {
263 // We need to create planar surface(s)
264 //
265 if( RMAX[a] != RMAX[a+1] && RMIN[a] != RMIN[a+1] )
266 {
267 // We can have the 8 following configurations here:
268 //
269 // 1. 2. 3. 4.
270 // --+ +-- --+ +--
271 // xx|-> <-|xx xx| |xx
272 // xx+-- --+xx --+ +--
273 // xxxxx xxxxx | |
274 // xxxxx xxxxx +-- --+
275 // xx+-- --+xx |xx xx|
276 // xx|-> <-|xx +-- --+
277 // --+ +--
278 // -------------------------- Z axis
279 //
280 //////////////////////////////////////////////////////////////
281 //////////////////////////////////////////////////////////////
282 //
283 // 5. 6. 7. 8.
284 // --+ +-- --+ +--
285 // xx|-> <-|xx xx|-> <-|xx
286 // --+-- --+-- xx+-- --+xx
287 // <-|xx xx|-> xxxxx xxxxx
288 // +-- --+ --+xx xx+--
289 // <-|xx xx|->
290 // +-- --+
291 // -------------------------- Z axis
292 //
293 // NOTE: The pictures show only one half of polycone!
294 // The arrows show the expected surface normal direction.
295 // The configuration No. 3 and 4 are not valid solids!
296
297 // Eliminate the invalid cases 3 and 4.
298 // At this point is guaranteed that each RMIN[i] < RMAX[i]
299 // where i in in interval 0 < i < num_z_planes-1. So:
300 //
301 if( RMIN[a] > RMAX[a+1] || RMAX[a] < RMIN[a+1] )
302 {
303 std::ostringstream message;
304 message << "The values of RMIN[" << a
305 << "] & RMAX[" << a+1 << "] or RMAX[" << a
306 << "] & RMIN[" << a+1 << "]" << G4endl
307 << "generate an invalid configuration for solid: "
308 << GetName() << "!";
309 G4Exception("G4BREPSolidPCone::G4BREPSolidPCone()",
310 "GeomSolids0002", FatalException, message );
311 }
312
313 ESurfaceSense UpSurfSense, LowSurfSense;
314
315 // We need to classify all the cases in order to figure out
316 // the planar surface sense
317 //
318 if( RMAX[a] > RMAX[a+1] )
319 {
320 // Cases 1, 5, 7
321 //
322 if( RMIN[a] < RMIN[a+1] )
323 {
324 // Case 1
325 //
326 UpSurfSense = ENormal;
327 LowSurfSense = ENormal;
328 }
329 else if( RMAX[a+1] != RMIN[a])
330 {
331 // Case 7
332 //
333 UpSurfSense = ENormal;
334 LowSurfSense = EInverse;
335 }
336 else
337 {
338 // Case 5
339 //
340 UpSurfSense = ENormal;
341 LowSurfSense = EInverse;
342 }
343 }
344 else
345 {
346 // Cases 2, 6, 8
347 //
348 if( RMIN[a] > RMIN[a+1] )
349 {
350 // Case 2
351 //
352 UpSurfSense = EInverse;
353 LowSurfSense = EInverse;
354 }
355 else if( RMIN[a+1] != RMAX[a] )
356 {
357 // Case 8
358 //
359 UpSurfSense = EInverse;
360 LowSurfSense = ENormal;
361 }
362 else
363 {
364 // Case 6
365 UpSurfSense = EInverse;
366 LowSurfSense = ENormal;
367 }
368 }
369
370 SurfaceVec[b] =
371 ComputePlanarSurface( RMAX[a], RMAX[a+1], LocalOrigin, PlaneAxis,
372 PlaneDir, UpSurfSense );
373 //SurfaceVec[b]->SetSameSense( UpSurfSense );
374 b++;
375
376 SurfaceVec[b] =
377 ComputePlanarSurface( RMIN[a], RMIN[a+1], LocalOrigin, PlaneAxis,
378 PlaneDir, LowSurfSense );
379 //SurfaceVec[b]->SetSameSense( LowSurfSense );
380 b++;
381 }
382 else
383 {
384 // The original code creating single planar surface
385 // in case where only either RMAX or RMIN have changed
386 // at the same Z value; e.g.:
387 //
388 // RMAX RMIN
389 // change change
390 //
391 // 1 2 3 4
392 // --+ +-- ----- -----
393 // 00|-> <-|00 00000 00000
394 // 00+-- --+00 --+00 00+--
395 // 00000 00000 <-|00 00|->
396 // +-- --+
397 // --------------------------- Z axis
398 //
399 // NOTE: The picture shows only one half of polycone!
400
401 G4double R1, R2;
402 ESurfaceSense SurfSense;
403
404 // test where is the plane surface
405 // if( RMAX[a] != RMAX[a+1] )
406 // {
407 // R1 = RMAX[a];
408 // R2 = RMAX[a+1];
409 // }
410 // else if(RMIN[a] != RMIN[a+1])
411 // {
412 // R1 = RMIN[a];
413 // R2 = RMIN[a+1];
414 // }
415 // else
416 // {
417 // std::ostringstream message;
418 // message << "Error in construction." << G4endl
419 // << "Exactly the same z, rmin and rmax given for "
420 // << "consecutive indices, " << a << " and " << a+1;
421 // G4Exception("G4BREPSolidPCone::InitializePCone()",
422 // "GeomSolids1002", JustWarning, message);
423 // continue;
424 // }
425 if( RMAX[a] != RMAX[a+1] )
426 {
427 // Cases 1, 2
428 //
429 R1 = RMAX[a];
430 R2 = RMAX[a+1];
431 if( R1 > R2 )
432 {
433 // Case 1
434 //
435 SurfSense = ENormal;
436 }
437 else
438 {
439 // Case 2
440 //
441 SurfSense = EInverse;
442 }
443 }
444 else if(RMIN[a] != RMIN[a+1])
445 {
446 // Cases 1, 2
447 //
448 R1 = RMIN[a];
449 R2 = RMIN[a+1];
450 if( R1 > R2 )
451 {
452 // Case 3
453 //
454 SurfSense = EInverse;
455 }
456 else
457 {
458 // Case 4
459 //
460 SurfSense = ENormal;
461 }
462 }
463 else
464 {
465 std::ostringstream message;
466 message << "Error in construction." << G4endl
467 << "Exactly the same z, rmin and rmax given for"
468 << "consecutive indices, " << a << " and " << a+1;
469 G4Exception("G4BREPSolidPCone::G4BREPSolidPCone()",
470 "GeomSolids1001", JustWarning, message);
471 continue;
472 }
473
474 SurfaceVec[b] =
475 ComputePlanarSurface( R1, R2, LocalOrigin, PlaneAxis,
476 PlaneDir, SurfSense );
477 //SurfaceVec[b]->SetSameSense( SurfSense );
478 b++;
479
480 // SurfaceVec[b]->SetSameSense(1);
482 }
483 }
484 else
485 {
486 // The surface to create is conical or cylindrical
487
488 // Inner PCone
489 //
490 if(RMIN[a] != RMIN[a+1])
491 {
492 // Create cone
493 //
494 if(RMIN[a] > RMIN[a+1])
495 {
496 G4Vector3D ConeOrigin = G4Vector3D(LocalOrigin) ;
497 SurfaceVec[b] = new G4FConicalSurface(ConeOrigin, Axis, Length,
498 RMIN[a+1], RMIN[a]);
500 }
501 else
502 {
503 G4Vector3D Axis2 = G4Vector3D( -1*Axis );
504 G4Vector3D LocalOrigin2 = G4Vector3D( LocalOrigin + (Length*Axis) );
505 G4Vector3D ConeOrigin = LocalOrigin2;
506 SurfaceVec[b] = new G4FConicalSurface(ConeOrigin, Axis2, Length,
507 RMIN[a], RMIN[a+1]);
509 }
510 b++;
511 }
512 else
513 {
514 if( RMIN[a] == 0 )
515 {
516 // Do not create any surface and decrease nb_of_surfaces
517 //
519 }
520 else
521 {
522 // Create cylinder
523 //
524 G4Vector3D CylOrigin = G4Vector3D( LocalOrigin );
525 SurfaceVec[b] = new G4FCylindricalSurface(CylOrigin, Axis,
526 RMIN[a], Length );
528 b++;
529 }
530 }
531
532 // Outer PCone
533 //
534 if(RMAX[a] != RMAX[a+1])
535 {
536 // Create cone
537 //
538 if(RMAX[a] > RMAX[a+1])
539 {
540 G4Vector3D ConeOrigin = G4Vector3D( LocalOrigin );
541 SurfaceVec[b] = new G4FConicalSurface(ConeOrigin, Axis, Length, RMAX[a+1], RMAX[a]);
543 }
544 else
545 {
546 G4Vector3D Axis2 = G4Vector3D( -1*Axis );
547 G4Vector3D LocalOrigin2 = G4Vector3D( LocalOrigin + (Length*Axis) );
548 G4Vector3D ConeOrigin = LocalOrigin2;
549 SurfaceVec[b] = new G4FConicalSurface(ConeOrigin, Axis2, Length,
550 RMAX[a], RMAX[a+1]);
552 }
553 b++;
554 }
555 else
556 {
557 // Create cylinder
558 //
559 G4Vector3D CylOrigin = G4Vector3D( LocalOrigin );
560
561 if (RMAX[a] == 0)
562 {
563 // Do not create any surface and decrease nb_of_surfaces
564 //
566 }
567 else
568 {
569 SurfaceVec[b] = new G4FCylindricalSurface(CylOrigin, Axis,
570 RMAX[a], Length );
572 b++;
573 }
574 }
575 }
576
577 // Move surface origin to next section
578 //
579 LocalOrigin = LocalOrigin + (Length*Axis);
580 }
581
582 ///////////////////////////////////////////////////
583 // Create two end planes
584
585 G4CurveVector cv;
586 G4CircularCurve* tmp;
587
588 if(RMIN[0] < RMAX[0])
589 {
590 // Create start G4Plane & boundaries
591 //
592 G4Point3D ArcStart1a = G4Point3D( Origin + (RMIN[0]*PlaneDir) );
593 G4Point3D ArcStart1b = G4Point3D( Origin + (RMAX[0]*PlaneDir) );
594
595 if( RMIN[0] > 0.0 )
596 {
597 tmp = new G4CircularCurve;
598 tmp->Init(G4Axis2Placement3D(PlaneDir, PlaneAxis, Origin), RMIN[0]);
599 tmp->SetBounds(ArcStart1a, ArcStart1a);
600 tmp->SetSameSense(0);
601 cv.push_back(tmp);
602 }
603
604 tmp = new G4CircularCurve;
605 tmp->Init(G4Axis2Placement3D(PlaneDir, PlaneAxis, Origin), RMAX[0]);
606 tmp->SetBounds(ArcStart1b, ArcStart1b);
607 tmp->SetSameSense(1);
608 cv.push_back(tmp);
609
610 SurfaceVec[nb_of_surfaces-2] = new G4FPlane(PlaneDir, -PlaneAxis, Origin);
613 }
614 else
615 {
616 // RMIN[0] == RMAX[0] so no surface is needed, it is a line!
617 //
619 }
620
621 if(RMIN[sections] < RMAX[sections])
622 {
623 // Create end G4Plane & boundaries
624 //
625 G4Point3D ArcStart2a = G4Point3D( LocalOrigin+(RMIN[sections]*PlaneDir) );
626 G4Point3D ArcStart2b = G4Point3D( LocalOrigin+(RMAX[sections]*PlaneDir) );
627
628 cv.clear();
629
630 if( RMIN[sections] > 0.0 )
631 {
632 tmp = new G4CircularCurve;
633 tmp->Init(G4Axis2Placement3D(PlaneDir, PlaneAxis, LocalOrigin),
634 RMIN[sections]);
635 tmp->SetBounds(ArcStart2a, ArcStart2a);
636 tmp->SetSameSense(0);
637 cv.push_back(tmp);
638 }
639
640 tmp = new G4CircularCurve;
641 tmp->Init(G4Axis2Placement3D(PlaneDir, PlaneAxis, LocalOrigin),
642 RMAX[sections]);
643 tmp->SetBounds(ArcStart2b, ArcStart2b);
644 tmp->SetSameSense(1);
645 cv.push_back(tmp);
646
647 SurfaceVec[nb_of_surfaces-1] = new G4FPlane(PlaneDir, PlaneAxis,
648 LocalOrigin);
650
651 // set sense of the surface
652 //
654 }
655 else
656 {
657 // RMIN[0] == RMAX[0] so no surface is needed, it is a line!
658 //
660 }
661
662 Initialize();
663}
664
666{
667 // Computes the bounding box for solids and surfaces
668 // Converts concave planes to convex
669 //
670 ShortestDistance=1000000;
672
673 if(!Box || !AxisBox)
674 IsConvex();
675
676 CalcBBoxes();
677}
678
680{
681 // This function find if the point Pt is inside,
682 // outside or on the surface of the solid
683
684 G4Vector3D v(1, 0, 0.01);
685 //G4Vector3D v(1, 0, 0); // This will miss the planar surface perp. to Z axis
686 //G4Vector3D v(0, 0, 1); // This works, however considered as hack not a fix
687 G4Vector3D Pttmp(Pt);
688 G4Vector3D Vtmp(v);
689 G4Ray r(Pttmp, Vtmp);
690
691 // Check if point is inside the PCone bounding box
692 //
693 if( !GetBBox()->Inside(Pttmp) )
694 {
695 return kOutside;
696 }
697
698 // Set the surfaces to active again
699 //
700 Reset();
701
702 // Test if the bounding box of each surface is intersected
703 // by the ray. If not, the surface is deactivated.
704 //
706
707 G4double dist = kInfinity;
708 G4bool isIntersected = false;
709 G4int WhichSurface = 0;
710
711 const G4double sqrHalfTolerance = kCarTolerance*kCarTolerance*0.25;
712
713 // Chech if the point is on the surface, otherwise
714 // find the nearest intersected suface. If there are not intersections the
715 // point is outside
716
717 for(G4int a=0; a < nb_of_surfaces; a++)
718 {
719 G4Surface* surface = SurfaceVec[a];
720
721 if( surface->IsActive() )
722 {
723 G4double hownear = surface->HowNear(Pt);
724
725 if( std::fabs( hownear ) < sqrHalfTolerance )
726 {
727 return kSurface;
728 }
729
730 if( surface->Intersect(r) )
731 {
732 isIntersected = true;
733 hownear = surface->GetDistance();
734
735 if ( std::fabs( hownear ) < dist )
736 {
737 dist = hownear;
738 WhichSurface = a;
739 }
740 }
741 }
742 }
743
744 if ( !isIntersected )
745 {
746 return kOutside;
747 }
748
749 // Find the point of intersection on the surface and the normal
750 // !!!! be carefull the distance is std::sqrt(dist) !!!!
751
752 dist = std::sqrt( dist );
753 G4Vector3D IntersectionPoint = Pttmp + dist*Vtmp;
754 G4Vector3D Normal =
755 SurfaceVec[WhichSurface]->SurfaceNormal( IntersectionPoint );
756
757 G4double dot = Normal*Vtmp;
758
759 return ( (dot > 0) ? kInside : kOutside );
760}
761
764{
765 // This function calculates the normal of the closest surface
766 // at a point on the surface
767 // Note : the sense of the normal depends on the sense of the surface
768
769 G4int isurf;
770 G4bool normflag = false;
771
772 const G4double sqrHalfTolerance = kCarTolerance*kCarTolerance*0.25;
773
774 // Determine if the point is on the surface
775 //
776 G4double minDist = kInfinity;
777 G4int normSurface = 0;
778 for(isurf = 0; isurf < nb_of_surfaces; isurf++)
779 {
780 G4double dist = std::fabs(SurfaceVec[isurf]->HowNear(Pt));
781 if( minDist > dist )
782 {
783 minDist = dist;
784 normSurface = isurf;
785 }
786 if( dist < sqrHalfTolerance)
787 {
788 // the point is on this surface
789 //
790 normflag = true;
791 break;
792 }
793 }
794
795 // Calculate the normal at this point, if the point is on the
796 // surface, otherwise compute the normal to the closest surface
797 //
798 if ( normflag ) // point on surface
799 {
800 G4ThreeVector norm = SurfaceVec[isurf]->SurfaceNormal(Pt);
801 return norm.unit();
802 }
803 else // point not on surface
804 {
805 G4Surface* nSurface = SurfaceVec[normSurface];
806 G4ThreeVector hitPt = nSurface->GetClosestHit();
807 G4ThreeVector hitNorm = nSurface->SurfaceNormal(hitPt);
808 return hitNorm.unit();
809 }
810}
811
813{
814 // Calculates the shortest distance ("safety") from a point
815 // outside the solid to any boundary of this solid.
816 // Return 0 if the point is already inside.
817
818 G4double *dists = new G4double[nb_of_surfaces];
819 G4int a;
820
821 // Set the surfaces to active again
822 //
823 Reset();
824
825 // Compute the shortest distance of the point to each surfaces
826 // Be careful : it's a signed value
827 //
828 for(a=0; a< nb_of_surfaces; a++)
829 dists[a] = SurfaceVec[a]->HowNear(Pt);
830
831 G4double Dist = kInfinity;
832
833 // if dists[] is positive, the point is outside
834 // so take the shortest of the shortest positive distances
835 // dists[] can be equal to 0 : point on a surface
836 // ( Problem with the G4FPlane : there is no inside and no outside...
837 // So, to test if the point is inside to return 0, utilize the Inside
838 // function. But I don`t know if it is really needed because dToIn is
839 // called only if the point is outside )
840 //
841 for(a = 0; a < nb_of_surfaces; a++)
842 if( std::fabs(Dist) > std::fabs(dists[a]) )
843 //if( dists[a] >= 0)
844 Dist = dists[a];
845
846 delete[] dists;
847
848 if(Dist == kInfinity)
849 // the point is inside the solid or on a surface
850 //
851 return 0;
852 else
853 return std::fabs(Dist);
854}
855
857 register const G4ThreeVector& V) const
858{
859 // Calculates the distance from a point outside the solid
860 // to the solid`s boundary along a specified direction vector.
861 //
862 // Note : Intersections with boundaries less than the
863 // tolerance must be ignored if the direction
864 // is away from the boundary
865
866 G4int a;
867
868 // Set the surfaces to active again
869 //
870 Reset();
871
872 G4double sqrHalfTolerance = kCarTolerance*kCarTolerance*0.25;
873 G4Vector3D Pttmp(Pt);
874 G4Vector3D Vtmp(V);
875 G4Ray r(Pttmp, Vtmp);
876
877 // Test if the bounding box of each surface is intersected
878 // by the ray. If not, the surface become deactive.
879 //
881
882 ShortestDistance = kInfinity;
883
884 for(a=0; a< nb_of_surfaces; a++)
885 {
886 if(SurfaceVec[a]->IsActive())
887 {
888 // test if the ray intersect the surface
889 //
890 G4Vector3D Norm = SurfaceVec[a]->SurfaceNormal(Pttmp);
891 G4double hownear = SurfaceVec[a]->HowNear(Pt);
892
893 if( (Norm * Vtmp) < 0 && std::fabs(hownear) < sqrHalfTolerance )
894 return 0;
895
896 if( (SurfaceVec[a]->Intersect(r)) )
897 {
898 // if more than 1 surface is intersected,
899 // take the nearest one
900 G4double distance = SurfaceVec[a]->GetDistance();
901
902 if( distance < ShortestDistance )
903 {
904 if( distance > sqrHalfTolerance )
905 {
906 ShortestDistance = distance;
907 }
908 }
909 }
910 }
911 }
912
913 // Be careful !
914 // SurfaceVec->Distance is in fact the squared distance
915 //
916 if(ShortestDistance != kInfinity)
917 return( std::sqrt(ShortestDistance) );
918 else
919 // no intersection
920 //
921 return kInfinity;
922}
923
925 register const G4ThreeVector& V,
926 const G4bool,
927 G4bool *validNorm,
928 G4ThreeVector *) const
929{
930 // Calculates the distance from a point inside the solid
931 // to the solid`s boundary along a specified direction vector.
932 // Return 0 if the point is already outside.
933 //
934 // Note : If the shortest distance to a boundary is less
935 // than the tolerance, it is ignored. This allows
936 // for a point within a tolerant boundary to leave
937 // immediately
938
939 // Set the surfaces to active again
940 //
941 Reset();
942
943 const G4double sqrHalfTolerance = kCarTolerance*kCarTolerance*0.25;
944 G4Vector3D Ptv = G4Vector3D( Pt );
945 G4int a;
946
947 // I don`t understand this line
948 //
949 if(validNorm)
950 *validNorm=false;
951
952 G4Vector3D Pttmp(Pt);
953 G4Vector3D Vtmp(V);
954
955 G4Ray r(Pttmp, Vtmp);
956
957 // Test if the bounding box of each surface is intersected
958 // by the ray. If not, the surface become deactive.
959 //
961
962 ShortestDistance = kInfinity;
963
964 for(a=0; a< nb_of_surfaces; a++)
965 {
966 if( SurfaceVec[a]->IsActive() )
967 {
968 G4Vector3D Norm = SurfaceVec[a]->SurfaceNormal(Pttmp);
969 G4double hownear = SurfaceVec[a]->HowNear(Pt);
970
971 if( (Norm * Vtmp) > 0 && std::fabs( hownear ) < sqrHalfTolerance )
972 return 0;
973
974 // test if the ray intersect the surface
975 //
976 if( SurfaceVec[a]->Intersect(r) )
977 {
978 // if more than 1 surface is intersected,
979 // take the nearest one
980 //
981 G4double distance = SurfaceVec[a]->GetDistance();
982
983 if( distance < ShortestDistance )
984 {
985 if( distance > sqrHalfTolerance )
986 {
987 ShortestDistance = distance;
988 }
989 }
990 }
991 }
992 }
993
994 // Be careful !
995 // SurfaceVec->Distance is in fact the squared distance
996 //
997 if(ShortestDistance != kInfinity)
998 return std::sqrt(ShortestDistance);
999 else
1000 // if no intersection is founded, the point is outside
1001 //
1002 return 0;
1003}
1004
1006{
1007 // Calculates the shortest distance ("safety") from a point
1008 // inside the solid to any boundary of this solid.
1009 // Return 0 if the point is already outside.
1010
1011 G4double *dists = new G4double[nb_of_surfaces];
1012 G4int a;
1013
1014 // Set the surfaces to active again
1015 Reset();
1016
1017 // calcul of the shortest distance of the point to each surfaces
1018 // Be careful : it's a signed value
1019 //
1020 for(a=0; a< nb_of_surfaces; a++)
1021 dists[a] = SurfaceVec[a]->HowNear(Pt);
1022
1023 G4double Dist = kInfinity;
1024
1025 // if dists[] is negative, the point is inside
1026 // so take the shortest of the shortest negative distances
1027 // dists[] can be equal to 0 : point on a surface
1028 // ( Problem with the G4FPlane : there is no inside and no outside...
1029 // So, to test if the point is outside to return 0, utilize the Inside
1030 // function. But I don`t know if it is really needed because dToOut is
1031 // called only if the point is inside )
1032
1033 for(a = 0; a < nb_of_surfaces; a++)
1034 if( std::fabs(Dist) > std::fabs(dists[a]) )
1035 //if( dists[a] <= 0)
1036 Dist = dists[a];
1037
1038 delete[] dists;
1039
1040 if(Dist == kInfinity)
1041 // the point is ouside the solid or on a surface
1042 //
1043 return 0;
1044 else
1045 return std::fabs(Dist);
1046}
1047
1049{
1050 return new G4BREPSolidPCone(*this);
1051}
1052
1053std::ostream& G4BREPSolidPCone::StreamInfo(std::ostream& os) const
1054{
1055 // Streams solid contents to output stream.
1056
1058 << "\n start_angle: " << constructorParams.start_angle
1059 << "\n opening_angle: " << constructorParams.opening_angle
1060 << "\n num_z_planes: " << constructorParams.num_z_planes
1061 << "\n z_start: " << constructorParams.z_start
1062 << "\n z_values: ";
1063 G4int idx;
1064 for( idx = 0; idx < constructorParams.num_z_planes; idx++ )
1065 {
1066 os << constructorParams.z_values[idx] << " ";
1067 }
1068 os << "\n RMIN: ";
1069 for( idx = 0; idx < constructorParams.num_z_planes; idx++ )
1070 {
1071 os << constructorParams.RMIN[idx] << " ";
1072 }
1073 os << "\n RMAX: ";
1074 for( idx = 0; idx < constructorParams.num_z_planes; idx++ )
1075 {
1076 os << constructorParams.RMAX[idx] << " ";
1077 }
1078 os << "\n-----------------------------------------------------------\n";
1079
1080 return os;
1081}
1082
1084{
1085 Active(1);
1086 ((G4BREPSolidPCone*)this)->intersectionDistance=kInfinity;
1087 StartInside(0);
1088 for(register int a=0;a<nb_of_surfaces;a++)
1089 SurfaceVec[a]->Reset();
1090 ShortestDistance = kInfinity;
1091}
1092
1093G4Surface*
1094G4BREPSolidPCone::ComputePlanarSurface( G4double r1,
1095 G4double r2,
1096 G4ThreeVector& origin,
1097 G4ThreeVector& planeAxis,
1098 G4ThreeVector& planeDirection,
1099 G4int surfSense )
1100{
1101 // The planar surface to return
1102 G4Surface* planarFace = 0;
1103
1104 G4CurveVector cv1;
1105 G4CircularCurve *tmp1, *tmp2;
1106
1107 // Create plane surface
1108 G4Point3D ArcStart1 = G4Point3D( origin + (r1 * planeDirection) );
1109 G4Point3D ArcStart2 = G4Point3D( origin + (r2 * planeDirection) );
1110
1111 if(r1 != 0)
1112 {
1113 tmp1 = new G4CircularCurve;
1114 tmp1->Init(G4Axis2Placement3D( planeDirection, planeAxis, origin), r1);
1115 tmp1->SetBounds(ArcStart1, ArcStart1);
1116
1117 if( surfSense )
1118 tmp1->SetSameSense(1);
1119 else
1120 tmp1->SetSameSense(0);
1121
1122 cv1.push_back(tmp1);
1123 }
1124
1125 if(r2 != 0)
1126 {
1127 tmp2 = new G4CircularCurve;
1128 tmp2->Init(G4Axis2Placement3D( planeDirection, planeAxis, origin), r2);
1129 tmp2->SetBounds(ArcStart2, ArcStart2);
1130
1131 if( surfSense )
1132 tmp2->SetSameSense(0);
1133 else
1134 tmp2->SetSameSense(1);
1135
1136 cv1.push_back(tmp2);
1137 }
1138
1139 planarFace = new G4FPlane( planeDirection, planeAxis, origin, surfSense );
1140
1141 planarFace->SetBoundaries(&cv1);
1142
1143 return planarFace;
1144}
1145
1146// In graphics_reps:
1147
1148#include "G4Polyhedron.hh"
1149
1151{
1152 return new G4PolyhedronPcon( constructorParams.start_angle,
1153 constructorParams.opening_angle,
1154 constructorParams.num_z_planes,
1155 constructorParams.z_values,
1156 constructorParams.RMIN,
1157 constructorParams.RMAX );
1158}
ESurfaceSense
@ ENormal
@ EInverse
std::vector< G4Curve * > G4CurveVector
@ JustWarning
@ FatalException
HepGeom::Point3D< G4double > G4Point3D
Definition: G4Point3D.hh:35
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
HepGeom::Vector3D< G4double > G4Vector3D
Definition: G4Vector3D.hh:35
#define G4endl
Definition: G4ios.hh:52
Hep3Vector unit() const
std::ostream & StreamInfo(std::ostream &os) const
G4BREPSolidPCone(const G4String &name, G4double start_angle, G4double opening_angle, G4int num_z_planes, G4double z_start, G4double z_values[], G4double RMIN[], G4double RMAX[])
G4VSolid * Clone() const
G4double DistanceToIn(const G4ThreeVector &) const
EInside Inside(register const G4ThreeVector &Pt) const
G4ThreeVector SurfaceNormal(const G4ThreeVector &) const
G4double DistanceToOut(register const G4ThreeVector &Pt, register const G4ThreeVector &V, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const
G4Polyhedron * CreatePolyhedron() const
G4BREPSolidPCone & operator=(const G4BREPSolidPCone &rhs)
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
void Init(const G4Axis2Placement3D &position0, G4double radius0)
void SetSameSense(G4int sameSense0)
void SetBounds(G4double p1, G4double p2)
Definition: G4Ray.hh:49
virtual G4int Intersect(const G4Ray &)
Definition: G4Surface.cc:170
void SetBoundaries(G4CurveVector *)
Definition: G4Surface.cc:140
G4int IsActive() const
virtual G4Vector3D SurfaceNormal(const G4Point3D &Pt) const =0
void SetSameSense(G4int sameSense0)
const G4Point3D & GetClosestHit() const
virtual G4double HowNear(const G4Vector3D &x) const
Definition: G4Surface.cc:283
G4double GetDistance() const
G4double kCarTolerance
Definition: G4VSolid.hh:307
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