Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BREPSolidPCone Class Reference

#include <G4BREPSolidPCone.hh>

+ Inheritance diagram for G4BREPSolidPCone:

Public Member Functions

 G4BREPSolidPCone (const G4String &name, G4double start_angle, G4double opening_angle, G4int num_z_planes, G4double z_start, G4double z_values[], G4double RMIN[], G4double RMAX[])
 
 ~G4BREPSolidPCone ()
 
void Initialize ()
 
EInside Inside (register const G4ThreeVector &Pt) const
 
G4ThreeVector SurfaceNormal (const G4ThreeVector &) const
 
G4double DistanceToIn (const G4ThreeVector &) const
 
G4double DistanceToIn (register const G4ThreeVector &Pt, register const G4ThreeVector &V) const
 
G4double DistanceToOut (register const G4ThreeVector &Pt, register const G4ThreeVector &V, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const
 
G4double DistanceToOut (const G4ThreeVector &) const
 
void Reset () const
 
G4PolyhedronCreatePolyhedron () const
 
G4VSolidClone () const
 
std::ostream & StreamInfo (std::ostream &os) const
 
 G4BREPSolidPCone (__void__ &)
 
 G4BREPSolidPCone (const G4BREPSolidPCone &rhs)
 
G4BREPSolidPConeoperator= (const G4BREPSolidPCone &rhs)
 
- Public Member Functions inherited from G4BREPSolid
 G4BREPSolid (const G4String &name)
 
 G4BREPSolid (const G4String &, G4Surface **, G4int)
 
virtual ~G4BREPSolid ()
 
virtual void Initialize ()
 
G4bool CalculateExtent (const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
 
virtual EInside Inside (register const G4ThreeVector &Pt) const
 
virtual G4ThreeVector SurfaceNormal (const G4ThreeVector &) const
 
virtual G4double DistanceToIn (const G4ThreeVector &) const
 
virtual G4double DistanceToIn (register const G4ThreeVector &Pt, register const G4ThreeVector &V) const
 
virtual G4double DistanceToOut (const G4ThreeVector &) const
 
virtual G4double DistanceToOut (register const G4ThreeVector &Pt, register const G4ThreeVector &V, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const
 
G4Point3D Scope () const
 
virtual G4String GetEntityType () const
 
virtual G4VSolidClone () const
 
virtual std::ostream & StreamInfo (std::ostream &os) const
 
void DescribeYourselfTo (G4VGraphicsScene &scene) const
 
G4PolyhedronCreatePolyhedron () const
 
G4NURBSCreateNURBS () const
 
virtual G4PolyhedronGetPolyhedron () const
 
G4int Intersect (register const G4Ray &) const
 
G4SurfaceGetSurface (G4int) const
 
void Active (G4int) const
 
G4int Active () const
 
G4double GetShortestDistance () const
 
G4int GetId () const
 
void SetId (G4int)
 
const G4StringGetName () const
 
void SetName (const G4String &name)
 
G4int GetNumberOfFaces () const
 
G4int GetNumberOfSolids () const
 
const G4Axis2Placement3DGetPlace () const
 
const G4BoundingBox3DGetBBox () const
 
G4int GetCubVolStatistics () const
 
G4double GetCubVolEpsilon () const
 
void SetCubVolStatistics (G4int st)
 
void SetCubVolEpsilon (G4double ep)
 
G4int GetAreaStatistics () const
 
G4double GetAreaAccuracy () const
 
void SetAreaStatistics (G4int st)
 
void SetAreaAccuracy (G4double ep)
 
G4double GetCubicVolume ()
 
G4double GetSurfaceArea ()
 
G4double IntersectionDistance () const
 
void IntersectionDistance (G4double) const
 
virtual void Reset () const
 
 G4BREPSolid (__void__ &)
 
 G4BREPSolid (const G4BREPSolid &rhs)
 
G4BREPSolidoperator= (const G4BREPSolid &rhs)
 
- Public Member Functions inherited from G4VSolid
 G4VSolid (const G4String &name)
 
virtual ~G4VSolid ()
 
G4bool operator== (const G4VSolid &s) const
 
G4String GetName () const
 
void SetName (const G4String &name)
 
G4double GetTolerance () const
 
virtual G4bool CalculateExtent (const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const =0
 
virtual EInside Inside (const G4ThreeVector &p) const =0
 
virtual G4ThreeVector SurfaceNormal (const G4ThreeVector &p) const =0
 
virtual G4double DistanceToIn (const G4ThreeVector &p, const G4ThreeVector &v) const =0
 
virtual G4double DistanceToIn (const G4ThreeVector &p) const =0
 
virtual G4double DistanceToOut (const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const =0
 
virtual G4double DistanceToOut (const G4ThreeVector &p) const =0
 
virtual void ComputeDimensions (G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
 
virtual G4double GetCubicVolume ()
 
virtual G4double GetSurfaceArea ()
 
virtual G4GeometryType GetEntityType () const =0
 
virtual G4ThreeVector GetPointOnSurface () const
 
virtual G4VSolidClone () const
 
virtual std::ostream & StreamInfo (std::ostream &os) const =0
 
void DumpInfo () const
 
virtual void DescribeYourselfTo (G4VGraphicsScene &scene) const =0
 
virtual G4VisExtent GetExtent () const
 
virtual G4PolyhedronCreatePolyhedron () const
 
virtual G4NURBSCreateNURBS () const
 
virtual G4PolyhedronGetPolyhedron () const
 
virtual const G4VSolidGetConstituentSolid (G4int no) const
 
virtual G4VSolidGetConstituentSolid (G4int no)
 
virtual const G4DisplacedSolidGetDisplacedSolidPtr () const
 
virtual G4DisplacedSolidGetDisplacedSolidPtr ()
 
 G4VSolid (__void__ &)
 
 G4VSolid (const G4VSolid &rhs)
 
G4VSolidoperator= (const G4VSolid &rhs)
 

Additional Inherited Members

- Protected Member Functions inherited from G4BREPSolid
G4ThreeVectorListCreateRotatedVertices (const G4AffineTransform &) const
 
G4bool IsConvex ()
 
virtual void CalcBBoxes ()
 
void CheckSurfaceNormals ()
 
void RemoveHiddenFaces (register const G4Ray &G4Rayref, G4int) const
 
void TestSurfaceBBoxes (register const G4Ray &) const
 
G4int StartInside () const
 
void StartInside (G4int si) const
 
void QuickSort (register G4Surface **SrfVec, register G4int left, register G4int right) const
 
- Protected Member Functions inherited from G4VSolid
void CalculateClippedPolygonExtent (G4ThreeVectorList &pPolygon, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
 
void ClipCrossSection (G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
 
void ClipBetweenSections (G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
 
void ClipPolygon (G4ThreeVectorList &pPolygon, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis) const
 
G4double EstimateCubicVolume (G4int nStat, G4double epsilon) const
 
G4double EstimateSurfaceArea (G4int nStat, G4double ell) const
 
- Protected Attributes inherited from G4BREPSolid
G4int Box
 
G4int Convex
 
G4int AxisBox
 
G4int PlaneSolid
 
G4Axis2Placement3Dplace
 
G4BoundingBox3Dbbox
 
G4double intersectionDistance
 
G4int active
 
G4int startInside
 
G4int nb_of_surfaces
 
G4Point3D intersection_point
 
G4Surface ** SurfaceVec
 
G4double RealDist
 
G4String solidname
 
G4int Id
 
- Protected Attributes inherited from G4VSolid
G4double kCarTolerance
 
- Static Protected Attributes inherited from G4BREPSolid
static G4int NumberOfSolids =0
 
static G4Ray Track
 
static G4double ShortestDistance = kInfinity
 

Detailed Description

Definition at line 57 of file G4BREPSolidPCone.hh.

Constructor & Destructor Documentation

◆ G4BREPSolidPCone() [1/3]

G4BREPSolidPCone::G4BREPSolidPCone ( const G4String name,
G4double  start_angle,
G4double  opening_angle,
G4int  num_z_planes,
G4double  z_start,
G4double  z_values[],
G4double  RMIN[],
G4double  RMAX[] 
)

Definition at line 58 of file G4BREPSolidPCone.cc.

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}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66

◆ ~G4BREPSolidPCone()

G4BREPSolidPCone::~G4BREPSolidPCone ( )

Definition at line 107 of file G4BREPSolidPCone.cc.

108{
109 if( constructorParams.num_z_planes > 0 )
110 {
111 delete [] constructorParams.z_values;
112 delete [] constructorParams.RMIN;
113 delete [] constructorParams.RMAX;
114 }
115}

◆ G4BREPSolidPCone() [2/3]

G4BREPSolidPCone::G4BREPSolidPCone ( __void__ &  a)

Definition at line 95 of file G4BREPSolidPCone.cc.

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}

◆ G4BREPSolidPCone() [3/3]

G4BREPSolidPCone::G4BREPSolidPCone ( const G4BREPSolidPCone rhs)

Definition at line 117 of file G4BREPSolidPCone.cc.

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}

Member Function Documentation

◆ Clone()

G4VSolid * G4BREPSolidPCone::Clone ( ) const
virtual

Reimplemented from G4BREPSolid.

Definition at line 1048 of file G4BREPSolidPCone.cc.

1049{
1050 return new G4BREPSolidPCone(*this);
1051}

◆ CreatePolyhedron()

G4Polyhedron * G4BREPSolidPCone::CreatePolyhedron ( ) const
virtual

Reimplemented from G4BREPSolid.

Definition at line 1150 of file G4BREPSolidPCone.cc.

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}

◆ DistanceToIn() [1/2]

G4double G4BREPSolidPCone::DistanceToIn ( const G4ThreeVector Pt) const
virtual

Reimplemented from G4BREPSolid.

Definition at line 812 of file G4BREPSolidPCone.cc.

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}
G4Surface ** SurfaceVec
Definition: G4BREPSolid.hh:231
G4int nb_of_surfaces
Definition: G4BREPSolid.hh:229

◆ DistanceToIn() [2/2]

G4double G4BREPSolidPCone::DistanceToIn ( register const G4ThreeVector Pt,
register const G4ThreeVector V 
) const
virtual

Reimplemented from G4BREPSolid.

Definition at line 856 of file G4BREPSolidPCone.cc.

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}
G4int Intersect(register const G4Ray &) const
static G4double ShortestDistance
Definition: G4BREPSolid.hh:221
void TestSurfaceBBoxes(register const G4Ray &) const
Definition: G4Ray.hh:49
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

◆ DistanceToOut() [1/2]

G4double G4BREPSolidPCone::DistanceToOut ( const G4ThreeVector Pt) const
virtual

Reimplemented from G4BREPSolid.

Definition at line 1005 of file G4BREPSolidPCone.cc.

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}

◆ DistanceToOut() [2/2]

G4double G4BREPSolidPCone::DistanceToOut ( register const G4ThreeVector Pt,
register const G4ThreeVector V,
const G4bool  calcNorm = false,
G4bool validNorm = 0,
G4ThreeVector n = 0 
) const
virtual

Reimplemented from G4BREPSolid.

Definition at line 924 of file G4BREPSolidPCone.cc.

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}
HepGeom::Vector3D< G4double > G4Vector3D
Definition: G4Vector3D.hh:35

◆ Initialize()

void G4BREPSolidPCone::Initialize ( )
virtual

Reimplemented from G4BREPSolid.

Definition at line 665 of file G4BREPSolidPCone.cc.

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}
G4bool IsConvex()
Definition: G4BREPSolid.cc:460
void CheckSurfaceNormals()
Definition: G4BREPSolid.cc:213
virtual void CalcBBoxes()

◆ Inside()

EInside G4BREPSolidPCone::Inside ( register const G4ThreeVector Pt) const
virtual

Reimplemented from G4BREPSolid.

Definition at line 679 of file G4BREPSolidPCone.cc.

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}
bool G4bool
Definition: G4Types.hh:67
EInside Inside(register const G4ThreeVector &Pt) const
const G4BoundingBox3D * GetBBox() const
virtual G4int Intersect(const G4Ray &)
Definition: G4Surface.cc:170
G4int IsActive() const
@ kInside
Definition: geomdefs.hh:58
@ kOutside
Definition: geomdefs.hh:58
@ kSurface
Definition: geomdefs.hh:58

Referenced by Inside().

◆ operator=()

G4BREPSolidPCone & G4BREPSolidPCone::operator= ( const G4BREPSolidPCone rhs)

Definition at line 145 of file G4BREPSolidPCone.cc.

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}
G4BREPSolid & operator=(const G4BREPSolid &rhs)
Definition: G4BREPSolid.cc:138

◆ Reset()

void G4BREPSolidPCone::Reset ( ) const
virtual

Reimplemented from G4BREPSolid.

Definition at line 1083 of file G4BREPSolidPCone.cc.

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}
G4int Active() const
G4int StartInside() const

Referenced by DistanceToIn(), DistanceToOut(), Inside(), and Reset().

◆ StreamInfo()

std::ostream & G4BREPSolidPCone::StreamInfo ( std::ostream &  os) const
virtual

Reimplemented from G4BREPSolid.

Definition at line 1053 of file G4BREPSolidPCone.cc.

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}
virtual std::ostream & StreamInfo(std::ostream &os) const

◆ SurfaceNormal()

G4ThreeVector G4BREPSolidPCone::SurfaceNormal ( const G4ThreeVector Pt) const
virtual

Reimplemented from G4BREPSolid.

Definition at line 763 of file G4BREPSolidPCone.cc.

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}
Hep3Vector unit() const
const G4Point3D & GetClosestHit() const

The documentation for this class was generated from the following files: