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

#include <G4FPlane.hh>

+ Inheritance diagram for G4FPlane:

Public Member Functions

 G4FPlane ()
 
virtual ~G4FPlane ()
 
 G4FPlane (const G4Vector3D &direction, const G4Vector3D &axis, const G4Point3D &Pt0, G4int sense=1)
 
 G4FPlane (const G4Point3DVector *pVec, const G4Point3DVector *iVec=0, G4int sense=1)
 
G4int Intersect (const G4Ray &G4Rayref)
 
void CalcBBox ()
 
void Project ()
 
G4int GetConvex () const
 
G4int GetNumberOfPoints () const
 
G4Point3D GetSrfPoint () const
 
const G4Point3DGetPoint (G4int Count) const
 
void CalcNormal ()
 
G4Vector3D SurfaceNormal (const G4Point3D &Pt) const
 
const char * Name () const
 
G4double ClosestDistanceToPoint (const G4Point3D &Pt)
 
G4double HowNear (const G4Vector3D &x) const
 
G4Axis2Placement3D GetPplace () const
 
G4Plane GetPplane () const
 
G4int MyType () const
 
G4int IsConvex () const
 
void Deactivate ()
 
G4RayNorm ()
 
const G4Point3DGetHitPoint () const
 
- Public Member Functions inherited from G4Surface
 G4Surface ()
 
virtual ~G4Surface ()
 
G4int operator== (const G4Surface &s)
 
virtual G4String GetEntityType () const
 
virtual const char * Name () const
 
virtual G4int MyType () const
 
void SetBoundaries (G4CurveVector *)
 
virtual G4double HowNear (const G4Vector3D &x) const
 
virtual G4double ClosestDistanceToPoint (const G4Point3D &Pt)
 
G4Vector3D GetOrigin () const
 
G4double GetDistance () const
 
void SetDistance (G4double Dist)
 
G4int IsActive () const
 
void SetActive (G4int act)
 
void Deactivate ()
 
void SetSameSense (G4int sameSense0)
 
G4int GetSameSense () const
 
G4BoundingBox3DGetBBox ()
 
const G4Point3DGetClosestHit () const
 
void SetNextNode (G4Surface *)
 
G4SurfaceGetNextNode ()
 
virtual void Reset ()
 
virtual G4int Intersect (const G4Ray &)
 
virtual G4Vector3D Normal (const G4Vector3D &p) const
 
virtual void CalcBBox ()
 
virtual G4double GetUHit () const
 
virtual G4double GetVHit () const
 
virtual G4Point3D Evaluation (const G4Ray &G4Rayref)
 
virtual G4int Evaluate (register const G4Ray &Rayref)
 
virtual void Project ()
 
virtual void CalcNormal ()
 
virtual G4int IsConvex () const
 
virtual G4int GetConvex () const
 
virtual G4int GetNumberOfPoints () const
 
virtual const G4Point3DGetPoint (G4int Count) const
 
virtual G4RayNorm ()
 
virtual G4Vector3D SurfaceNormal (const G4Point3D &Pt) const =0
 
- Public Member Functions inherited from G4STEPEntity
 G4STEPEntity ()
 
virtual ~G4STEPEntity ()
 
virtual G4String GetEntityType () const =0
 

Protected Member Functions

void InitBounded ()
 
virtual void InitBounded ()
 

Protected Attributes

G4Point3D hitpoint
 
- Protected Attributes inherited from G4Surface
G4BoundingBox3Dbbox
 
G4Point3D closest_hit
 
G4Surfacenext
 
G4SurfaceBoundary surfaceBoundary
 
G4double kCarTolerance
 
G4double kAngTolerance
 
G4int Intersected
 
G4Vector3D origin
 
G4int Type
 
G4int AdvancedFace
 
G4int active
 
G4double distance
 
G4double uhit
 
G4double vhit
 
G4int sameSense
 

Additional Inherited Members

- Static Public Member Functions inherited from G4Surface
static void Project (G4double &Coord, const G4Point3D &Pt, const G4Plane &Pl)
 

Detailed Description

Definition at line 67 of file G4FPlane.hh.

Constructor & Destructor Documentation

◆ G4FPlane() [1/3]

G4FPlane::G4FPlane ( )

◆ ~G4FPlane()

G4FPlane::~G4FPlane ( )
virtual

Definition at line 124 of file G4FPlane.cc.

125{
126 delete NormalX;
127}

◆ G4FPlane() [2/3]

G4FPlane::G4FPlane ( const G4Vector3D direction,
const G4Vector3D axis,
const G4Point3D Pt0,
G4int  sense = 1 
)

Definition at line 49 of file G4FPlane.cc.

52 : pplace(direction, axis, Pt0), Convex(0), projectedBoundary(0)
53{
54 G4Point3D Pt1 = G4Point3D( Pt0 + direction );
55
56 // The plane include direction and axis is the normal,
57 // so axis^direction is included in the plane
58 G4Point3D Pt2 = G4Point3D( Pt0 + axis.cross(direction) );
59
60 G4Ray::CalcPlane3Pts( Pl, Pt0, Pt1, Pt2 );
61
62 active = 1;
63 sameSense = sense;
64 CalcNormal();
65 distance = kInfinity;
66 Type = 1;
67}
HepGeom::Point3D< G4double > G4Point3D
Definition: G4Point3D.hh:35
void CalcNormal()
Definition: G4FPlane.cc:147
static G4int CalcPlane3Pts(G4Plane &plane, const G4Point3D &a, const G4Point3D &b, const G4Point3D &c)
Definition: G4Ray.cc:212
G4int sameSense
Definition: G4Surface.hh:207
G4int Type
Definition: G4Surface.hh:200
G4double distance
Definition: G4Surface.hh:203
G4int active
Definition: G4Surface.hh:202
BasicVector3D< T > cross(const BasicVector3D< T > &v) const

◆ G4FPlane() [3/3]

G4FPlane::G4FPlane ( const G4Point3DVector pVec,
const G4Point3DVector iVec = 0,
G4int  sense = 1 
)

Definition at line 70 of file G4FPlane.cc.

73 : pplace( (*pVec)[0]-(*pVec)[1], // direction
74 ((*pVec)[pVec->size()-1]-(*pVec)[0])
75 .cross((*pVec)[0]-(*pVec)[1]), // axis
76 (*pVec)[0] ) // location
77
78{
79 G4Ray::CalcPlane3Pts( Pl, (*pVec)[0], (*pVec)[1], (*pVec)[2] );
80
81 G4CurveVector bounds;
82 G4CompositeCurve* polygon;
83
84 projectedBoundary = new G4SurfaceBoundary;
85
86 sameSense = sense;
87
88 // Outer boundary
89
90 polygon= new G4CompositeCurve(*pVec);
91
92 for (size_t i=0; i< polygon->GetSegments().size(); i++)
93 polygon->GetSegments()[i]->SetSameSense(sameSense);
94
95 bounds.push_back(polygon);
96
97 // Eventual inner boundary
98
99 if (iVec)
100 {
101 polygon= new G4CompositeCurve(*iVec);
102
103 for (size_t i=0; i< polygon->GetSegments().size(); i++)
104 polygon->GetSegments()[i]->SetSameSense(sameSense);
105
106 bounds.push_back(polygon);
107 }
108
109 // Set sense for boundaries
110
111 for (size_t j=0; j< bounds.size(); j++)
112 bounds[j]->SetSameSense(sameSense);
113
114
115 SetBoundaries(&bounds);
116
117 CalcNormal();
118 IsConvex();
119 distance = kInfinity;
120 Type=1;
121}
std::vector< G4Curve * > G4CurveVector
const G4CurveVector & GetSegments() const
G4int IsConvex() const
Definition: G4FPlane.cc:204
void SetBoundaries(G4CurveVector *)
Definition: G4Surface.cc:140
void SetSameSense(G4int sameSense0)

Member Function Documentation

◆ CalcBBox()

void G4FPlane::CalcBBox ( )
virtual

Reimplemented from G4Surface.

Definition at line 130 of file G4FPlane.cc.

131{
132 // This is needed since the bounds are used for the Solid
133 // bbox calculation. The bbox test is NOT performed for
134 // planar surfaces.
135
136 // Finds the bounds of the G4Plane surface iow
137 // calculates the bounds for a bounding box
138 // to the surface. The bounding box is used
139 // for a preliminary check of intersection.
140
143
144}
G4Point3D GetBoxMin() const
G4Point3D GetBoxMax() const
const G4BoundingBox3D & BBox() const
G4BoundingBox3D * bbox
Definition: G4Surface.hh:185
G4SurfaceBoundary surfaceBoundary
Definition: G4Surface.hh:189

◆ CalcNormal()

void G4FPlane::CalcNormal ( )
virtual

Reimplemented from G4Surface.

Definition at line 147 of file G4FPlane.cc.

148{
149/*
150 // Calc Normal for surface which is used for the projection
151 // Make planes
152 G4Vector3D norm;
153
154 G4Vector3D RefDirection = pplace.GetRefDirection();
155 G4Vector3D Axis = pplace.GetAxis();
156
157 // L. Broglia : before in G4Placement
158 if( RefDirection == Axis )
159 norm = RefDirection;
160 else
161 {
162 // L. Broglia : error on setY, and it`s better to use cross function
163 // norm.setX( RefDirection.y() * Axis.z() - RefDirection.z() * Axis.y() );
164 // norm.setY( RefDirection.x() * Axis.z() - RefDirection.z() * Axis.x() );
165 // norm.setZ( RefDirection.x() * Axis.y() - RefDirection.y() * Axis.x() );
166
167 norm = RefDirection.cross(Axis);
168 }
169
170 // const G4Point3D& tmp = pplace.GetSrfPoint();
171 const G4Point3D tmp = pplace.GetLocation();
172*/
173
174 // L. Broglia
175 // The direction of the normal is the axis of his location
176 // Its sense depend on the orientation of the bounded curve
177 const G4Point3D tmp = pplace.GetLocation();
178 G4Vector3D norm;
179 G4int sense = GetSameSense();
180
181 if (sense)
182 norm = pplace.GetAxis();
183 else
184 norm = - pplace.GetAxis();
185
186 NormalX = new G4Ray(tmp, norm);
187 NormalX->RayCheck();
188 NormalX->CreatePlanes();
189}
int G4int
Definition: G4Types.hh:66
G4Point3D GetLocation() const
G4Vector3D GetAxis() const
Definition: G4Ray.hh:49
void RayCheck()
Definition: G4Ray.cc:251
void CreatePlanes()
Definition: G4Ray.cc:63
G4int GetSameSense() const

Referenced by G4FPlane().

◆ ClosestDistanceToPoint()

G4double G4FPlane::ClosestDistanceToPoint ( const G4Point3D Pt)
virtual

Reimplemented from G4Surface.

Definition at line 341 of file G4FPlane.cc.

342{
343 // Calculates signed distance of point Pt to G4Plane Pl
344 // Be careful, the equation of the plane is :
345 // ax + by + cz = d
346 G4double dist = Pt.x()*Pl.a + Pt.y()*Pl.b + Pt.z()*Pl.c - Pl.d;
347
348 return dist;
349}
double G4double
Definition: G4Types.hh:64
G4double d
Definition: G4Plane.hh:52
G4double a
Definition: G4Plane.hh:52
G4double c
Definition: G4Plane.hh:52
G4double b
Definition: G4Plane.hh:52

◆ Deactivate()

void G4FPlane::Deactivate ( )
inline

◆ GetConvex()

G4int G4FPlane::GetConvex ( ) const
inlinevirtual

Reimplemented from G4Surface.

◆ GetHitPoint()

const G4Point3D & G4FPlane::GetHitPoint ( ) const
inline

◆ GetNumberOfPoints()

G4int G4FPlane::GetNumberOfPoints ( ) const
inlinevirtual

Reimplemented from G4Surface.

◆ GetPoint()

const G4Point3D & G4FPlane::GetPoint ( G4int  Count) const
inlinevirtual

Reimplemented from G4Surface.

◆ GetPplace()

G4Axis2Placement3D G4FPlane::GetPplace ( ) const
inline

◆ GetPplane()

G4Plane G4FPlane::GetPplane ( ) const
inline

◆ GetSrfPoint()

G4Point3D G4FPlane::GetSrfPoint ( ) const
inline

◆ HowNear()

G4double G4FPlane::HowNear ( const G4Vector3D x) const
virtual

Reimplemented from G4Surface.

Definition at line 360 of file G4FPlane.cc.

361{
362 G4double hownear = Pt.x()*Pl.a + Pt.y()*Pl.b + Pt.z()*Pl.c - Pl.d;
363
364 return hownear;
365}

◆ InitBounded()

void G4FPlane::InitBounded ( )
protectedvirtual

Reimplemented from G4Surface.

Definition at line 352 of file G4FPlane.cc.

353{
354 // L. Broglia
355
356 projectedBoundary =
358}
const G4Transform3D & GetToPlacementCoordinates() const
G4SurfaceBoundary * Project(const G4Transform3D &tr=G4Transform3D::Identity)

◆ Intersect()

G4int G4FPlane::Intersect ( const G4Ray G4Rayref)
virtual

Reimplemented from G4Surface.

Definition at line 210 of file G4FPlane.cc.

211{
212 // This function count the number of intersections of a
213 // bounded surface by a ray.
214
215
216 // Find the intersection with the infinite plane
217 Intersected =1;
218
219 // s is solution, line is p + tq, n is G4Plane Normal, r is point on G4Plane
220 // all parameters are pointers to arrays of three elements
221
223 register G4double a, b, t;
224
225 register const G4Vector3D& RayDir = rayref.GetDir();
226 register const G4Point3D& RayStart = rayref.GetStart();
227
228 G4double dirx = RayDir.x();
229 G4double diry = RayDir.y();
230 G4double dirz = RayDir.z();
231
232 G4Vector3D norm = (*NormalX).GetDir();
233 G4Point3D srf_point = pplace.GetLocation();
234
235 b = norm.x() * dirx + norm.y() * diry + norm.z() * dirz;
236
237 if ( std::fabs(b) < perMillion )
238 {
239 // G4cout << "\nLine is parallel to G4Plane.No Hit.";
240 }
241 else
242 {
243 G4double startx = RayStart.x();
244 G4double starty = RayStart.y();
245 G4double startz = RayStart.z();
246
247 a = norm.x() * (srf_point.x() - startx) +
248 norm.y() * (srf_point.y() - starty) +
249 norm.z() * (srf_point.z() - startz) ;
250
251 t = a/b;
252
253 // substitute t into line equation
254 // to calculate final solution
255 G4double solx,soly,solz;
256 solx = startx + t * dirx;
257 soly = starty + t * diry;
258 solz = startz + t * dirz;
259
260 // solve tolerance problem
261 if( (t*dirx >= -kCarTolerance/2) && (t*dirx <= kCarTolerance/2) )
262 solx = startx;
263
264 if( (t*diry >= -kCarTolerance/2) && (t*diry <= kCarTolerance/2) )
265 soly = starty;
266
267 if( (t*dirz >= -kCarTolerance/2) && (t*dirz <= kCarTolerance/2) )
268 solz = startz;
269
270 G4bool xhit = (dirx < 0 && solx <= startx) || (dirx >= 0 && solx >= startx);
271 G4bool yhit = (diry < 0 && soly <= starty) || (diry >= 0 && soly >= starty);
272 G4bool zhit = (dirz < 0 && solz <= startz) || (dirz >= 0 && solz >= startz);
273
274 if( xhit && yhit && zhit ) {
275 hitpoint= G4Point3D(solx, soly, solz);
276 }
277 }
278
279 // closest_hit is a public Point3D in G4Surface
281
282 if(closest_hit.x() == kInfinity)
283 {
284 // no hit
285 active=0;
286 SetDistance(kInfinity);
287 return 0;
288 }
289 else
290 {
291 // calculate the squared distance from the point to the intersection
292 // and set it in the distance data member (all clients know they have
293 // to take the sqrt)
294 SetDistance( RayStart.distance2(closest_hit) );
295
296 // now, we have to verify that the hit point founded
297 // is included into the G4FPlane boundaries
298
299 // project the hit to the xy plane,
300 // with the same projection that took the boundary
301 // into projectedBoundary
302 G4Point3D projectedHit= pplace.GetToPlacementCoordinates() * closest_hit;
303
304 // test ray from the hit on the xy plane
305 G4Ray testRay( projectedHit, G4Vector3D(1, 0.01, 0) );
306
307 // check if it intersects the boundary
308 G4int nbinter = projectedBoundary->IntersectRay2D(testRay);
309
310 // If this number is par, it`s signify that the projected point
311 // is outside the projected surface, so the hit point is outside
312 // the bounded surface
313 if(nbinter&1)
314 {
315 // the intersection point is into the boundaries
316 // check if the intersection point is on the surface
318 {
319 // the point is on the surface, set the distance to 0
320 SetDistance(0);
321 }
322 else
323 {
324 // the point is outside the surface
325 }
326
327 return 1 ;
328 }
329 else
330 {
331 // the intersection point is out the boundaries
332 // it is not a real intersection
333 active=0;
334 SetDistance(kInfinity);
335 return 0;
336 }
337 }
338}
const G4Point3D PINFINITY(kInfinity, kInfinity, kInfinity)
bool G4bool
Definition: G4Types.hh:67
HepGeom::Vector3D< G4double > G4Vector3D
Definition: G4Vector3D.hh:35
G4Point3D hitpoint
Definition: G4FPlane.hh:151
G4int IntersectRay2D(const G4Ray &ray)
void SetDistance(G4double Dist)
G4int Intersected
Definition: G4Surface.hh:194
G4Point3D closest_hit
Definition: G4Surface.hh:186
G4double kCarTolerance
Definition: G4Surface.hh:192

◆ IsConvex()

G4int G4FPlane::IsConvex ( ) const
virtual

Reimplemented from G4Surface.

Definition at line 204 of file G4FPlane.cc.

205{
206 return -1;
207}

Referenced by G4FPlane().

◆ MyType()

G4int G4FPlane::MyType ( ) const
inlinevirtual

Reimplemented from G4Surface.

◆ Name()

const char * G4FPlane::Name ( ) const
inlinevirtual

Reimplemented from G4Surface.

◆ Norm()

G4Ray * G4FPlane::Norm ( )
inlinevirtual

Reimplemented from G4Surface.

◆ Project()

void G4FPlane::Project ( )
virtual

Reimplemented from G4Surface.

Definition at line 192 of file G4FPlane.cc.

193{
194 // Project
195 // const G4Plane& Plane1 = NormalX->GetPlane(1);
196 // const G4Plane& Plane2 = NormalX->GetPlane(2);
197
198 // probably not necessary
199 // projections of the boundary should be handled by the intersection
200 // OuterBoundary->ProjectBoundaryTo2D(Plane1, Plane2, 0);
201}

◆ SurfaceNormal()

G4Vector3D G4FPlane::SurfaceNormal ( const G4Point3D Pt) const
inlinevirtual

Member Data Documentation

◆ hitpoint

G4Point3D G4FPlane::hitpoint
protected

Definition at line 151 of file G4FPlane.hh.

Referenced by Intersect().


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