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

#include <G4BREPSolidPolyhedra.hh>

+ Inheritance diagram for G4BREPSolidPolyhedra:

Public Member Functions

 G4BREPSolidPolyhedra (const G4String &name, G4double phi1, G4double dphi, G4int sides, G4int num_z_planes, G4double z_start, G4double z_values[], G4double RMIN[], G4double RMAX[])
 
 ~G4BREPSolidPolyhedra ()
 
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
 
G4VSolidClone () const
 
std::ostream & StreamInfo (std::ostream &os) const
 
G4PolyhedronCreatePolyhedron () const
 
void Reset () const
 
 G4BREPSolidPolyhedra (__void__ &)
 
 G4BREPSolidPolyhedra (const G4BREPSolidPolyhedra &rhs)
 
G4BREPSolidPolyhedraoperator= (const G4BREPSolidPolyhedra &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 59 of file G4BREPSolidPolyhedra.hh.

Constructor & Destructor Documentation

◆ G4BREPSolidPolyhedra() [1/3]

G4BREPSolidPolyhedra::G4BREPSolidPolyhedra ( const G4String name,
G4double  phi1,
G4double  dphi,
G4int  sides,
G4int  num_z_planes,
G4double  z_start,
G4double  z_values[],
G4double  RMIN[],
G4double  RMAX[] 
)

Definition at line 71 of file G4BREPSolidPolyhedra.cc.

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}
@ JustWarning
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

◆ ~G4BREPSolidPolyhedra()

G4BREPSolidPolyhedra::~G4BREPSolidPolyhedra ( )

Definition at line 148 of file G4BREPSolidPolyhedra.cc.

149{
150 if( constructorParams.num_z_planes > 0 )
151 {
152 delete [] constructorParams.z_values;
153 delete [] constructorParams.RMIN;
154 delete [] constructorParams.RMAX;
155 }
156}

◆ G4BREPSolidPolyhedra() [2/3]

G4BREPSolidPolyhedra::G4BREPSolidPolyhedra ( __void__ &  a)

Definition at line 135 of file G4BREPSolidPolyhedra.cc.

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}

◆ G4BREPSolidPolyhedra() [3/3]

G4BREPSolidPolyhedra::G4BREPSolidPolyhedra ( const G4BREPSolidPolyhedra rhs)

Definition at line 158 of file G4BREPSolidPolyhedra.cc.

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}

Member Function Documentation

◆ Clone()

G4VSolid * G4BREPSolidPolyhedra::Clone ( ) const
virtual

Reimplemented from G4BREPSolid.

Definition at line 1340 of file G4BREPSolidPolyhedra.cc.

1341{
1342 return new G4BREPSolidPolyhedra(*this);
1343}

◆ CreatePolyhedron()

G4Polyhedron * G4BREPSolidPolyhedra::CreatePolyhedron ( ) const
virtual

Reimplemented from G4BREPSolid.

Definition at line 1517 of file G4BREPSolidPolyhedra.cc.

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}

◆ DistanceToIn() [1/2]

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

Reimplemented from G4BREPSolid.

Definition at line 1057 of file G4BREPSolidPolyhedra.cc.

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

◆ DistanceToIn() [2/2]

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

Reimplemented from G4BREPSolid.

Definition at line 1106 of file G4BREPSolidPolyhedra.cc.

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}
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
G4double GetDistance() const
G4double kCarTolerance
Definition: G4VSolid.hh:307

◆ DistanceToOut() [1/2]

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

Reimplemented from G4BREPSolid.

Definition at line 1285 of file G4BREPSolidPolyhedra.cc.

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}
virtual G4double HowNear(const G4Vector3D &x) const
Definition: G4Surface.cc:283

◆ DistanceToOut() [2/2]

G4double G4BREPSolidPolyhedra::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 1187 of file G4BREPSolidPolyhedra.cc.

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}
virtual G4int Intersect(const G4Ray &)
Definition: G4Surface.cc:170

◆ Initialize()

void G4BREPSolidPolyhedra::Initialize ( )
virtual

Reimplemented from G4BREPSolid.

Definition at line 906 of file G4BREPSolidPolyhedra.cc.

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

◆ Inside()

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

Reimplemented from G4BREPSolid.

Definition at line 929 of file G4BREPSolidPolyhedra.cc.

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}
EInside Inside(register const G4ThreeVector &Pt) const
const G4BoundingBox3D * GetBBox() const
G4int IsActive() const
@ kInside
Definition: geomdefs.hh:58
@ kOutside
Definition: geomdefs.hh:58
@ kSurface
Definition: geomdefs.hh:58

Referenced by Inside().

◆ operator=()

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

Definition at line 187 of file G4BREPSolidPolyhedra.cc.

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

◆ Reset()

void G4BREPSolidPolyhedra::Reset ( ) const
virtual

Reimplemented from G4BREPSolid.

Definition at line 919 of file G4BREPSolidPolyhedra.cc.

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

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

◆ StreamInfo()

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

Reimplemented from G4BREPSolid.

Definition at line 1345 of file G4BREPSolidPolyhedra.cc.

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

◆ SurfaceNormal()

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

Reimplemented from G4BREPSolid.

Definition at line 1008 of file G4BREPSolidPolyhedra.cc.

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}
bool G4bool
Definition: G4Types.hh:67
Hep3Vector unit() const
G4Point3D GetSrfPoint() const
G4Vector3D SurfaceNormal(const G4Point3D &Pt) const

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