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

#include <G4BezierSurface.hh>

+ Inheritance diagram for G4BezierSurface:

Public Member Functions

 G4BezierSurface ()
 
virtual ~G4BezierSurface ()
 
G4Point3D AveragePoint () const
 
void SetAveragePoint (G4Point3D p)
 
G4double UAverage () const
 
G4double VAverage () const
 
void Dir (G4int d)
 
void ChangeDir ()
 
G4double SMin () const
 
G4double SMax () const
 
G4int GetOrder (G4int direction) const
 
void PutOrder (G4int direction, G4int value)
 
G4double GetU () const
 
G4double GetV () const
 
void CalcBBox ()
 
G4int BIntersect (G4SurfaceList &)
 
G4int ClipBothDirs ()
 
void ClipSurface ()
 
virtual G4Vector3D SurfaceNormal (const G4Point3D &Pt) 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 Attributes

G4SurfaceListbezier_list
 
- 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
 

Static Protected Attributes

static G4int Clips =0
 
static G4int Splits =0
 
static G4double Tolerance =0
 

Friends

class G4BSplineSurface
 
class G4ProjectedSurface
 
void CopySurface (G4BezierSurface &bez)
 

Additional Inherited Members

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

Detailed Description

Definition at line 51 of file G4BezierSurface.hh.

Constructor & Destructor Documentation

◆ G4BezierSurface()

G4BezierSurface::G4BezierSurface ( )

Definition at line 48 of file G4BezierSurface.cc.

49 : G4Surface(), bezier_list(0), smin(0.), smax(0.),
50 average_u(0.), average_v(0.), dir(0), u_knots(0), v_knots(0),
51 ctl_points(0), new_knots(0), ord(0), oslo_m(0), lower(0), upper(0),
52 u_min(0.), u_max(0.), v_min(0.), v_max(0.), old_points(0)
53{
54 order[0]=0; order[1]=0;
55 u[0]=0.; u[1]=0.;
56 v[0]=0.; v[1]=0.;
57}
G4SurfaceList * bezier_list

◆ ~G4BezierSurface()

G4BezierSurface::~G4BezierSurface ( )
virtual

Definition at line 59 of file G4BezierSurface.cc.

60{
61 delete u_knots;
62 delete v_knots;
63 delete new_knots;
64 delete ctl_points;
65 delete old_points;
66
67 G4OsloMatrix* temp_oslo = oslo_m;
68
69 while(oslo_m != (G4OsloMatrix*)0)
70 {
71 oslo_m = oslo_m->GetNextNode();
72 delete temp_oslo;
73 temp_oslo = oslo_m;
74 }
75
76 delete oslo_m;
77 delete bbox;
78}
G4OsloMatrix * GetNextNode()
G4BoundingBox3D * bbox
Definition: G4Surface.hh:185

Member Function Documentation

◆ AveragePoint()

G4Point3D G4BezierSurface::AveragePoint ( ) const
inline

◆ BIntersect()

G4int G4BezierSurface::BIntersect ( G4SurfaceList bez_list)

Definition at line 211 of file G4BezierSurface.cc.

212{
213 bezier_list = &bez_list;
214 G4int clip_regions = 0; // Used for tolerance/efficiency-testing
215
216 do
217 {
218 // Calc bbox
219 CalcBBox();
220
221 // Test bbox
222/* L. Broglia
223 bbox->Test2dBBox();
224*/
225 // bbox->Test();
226
227 // Check result
228 if(!bbox->GetTestResult())
229 return 0;
230
231 // The first clipping has already been Done
232 // previously so we continue by doing the
233 // actual clip.
234
235 // Cut out the clipped region of the surface
236 GetClippedRegionFromSurface();
237 clip_regions++;
238
239 // Calculate the knot vectors and control points
240 // for the clipped surface
241 RefineSurface();
242
243 // Gets the u- and v-bounds for the clipped surface
244 u_min = u_knots->GetKnot(0);
245 u_max = u_knots->GetKnot(u_knots->GetSize() - 1);
246 v_min = v_knots->GetKnot(0);
247 v_max = v_knots->GetKnot(v_knots->GetSize() - 1);
248
249 // Choose the G4Vector3D for the next() clipping so that
250 // the larger side will be clipped.
251 if( (u_max - u_min) < (v_max - v_min) )
252 dir = 1;
253 else
254 dir = 0;
255
256 // Calculate the clip points
257 ClipSurface();
258 // G4cout << "\n SMINMAX : " << smin << " " << smax;
259
260 // The ray intersects with the bounding box
261 // but not with the surface itself.
262 if( smin > 1.0 || smax < 0.0 )
263 {
264 // G4cout << "\nG4BezierSurface::Intersect : bezier missed!";
265 // bezier_list->RemoveSurface(this);
266 return 0;
267 }
268
269 if( (smax - smin) > 0.8)
270 {
271 // Multiple intersections
272 // G4cout << "\nG4BezierSurface::Intersect : Bezier split.";
273 SplitNURBSurface();
274 // Now the two new surfaces should also be
275 // clipped in both G4Vector3Ds i.e the
276 // last and the second last surface
277 // in the List. This is Done after returning
278 // from this function.
279 // G4cout << "\n\n BEZ SPLIT in final Calc! \n\n";
280
281
282 return 2;
283 }
284
285 // Calculate the smin and smax values on the
286 // b_spline.
287 LocalizeClipValues();
288
289 // Check if the size of the remaining surface is within the
290 // Tolerance .
291 } while ((u_max - u_min > Tolerance) || (v_max - v_min) > Tolerance);
292
293 SetValues();
294 // G4cout << "\nG4BezierSurface::Intersect :Regions were cut "
295 // << clip_regions << " Times.\n";
296
297 return 1;
298}
int G4int
Definition: G4Types.hh:66
static G4double Tolerance
G4int GetTestResult() const
G4int GetSize() const
G4double GetKnot(G4int knot_number) const

◆ CalcBBox()

void G4BezierSurface::CalcBBox ( )
virtual

Reimplemented from G4Surface.

Definition at line 142 of file G4BezierSurface.cc.

143{
144 // Finds the bounds of the 2D-projected nurb iow
145 // calculates the bounds for a bounding rectangle
146 // to the surface. The bounding rectangle is used
147 // for a preliminary check of intersection.
148 G4Point3D box_min = G4Point3D(PINFINITY);
149 G4Point3D box_max = G4Point3D(-PINFINITY);
150
151
152 // Loop to search the whole control point mesh
153 // for the minimum and maximum values for.X() and y.
154 for(register G4int a = ctl_points->GetRows()-1; a>=0;a--)
155 for(register G4int b = ctl_points->GetCols()-1; b>=0;b--)
156 {
157/* L. Broglia
158 G4Point2d& tmp = (G4Point2d&)ctl_points->get(a,b);
159 if((box_min.X()) > (tmp.X())) box_min.X(tmp.X());
160 if((box_max.X()) < (tmp.X())) box_max.X(tmp.X());
161 if((box_min.Y()) > (tmp.Y())) box_min.Y(tmp.Y());
162 if((box_max.Y()) < (tmp.Y())) box_max.Y(tmp.Y());
163*/
164 G4Point3D tmp = ctl_points->Get3D(a,b);
165 if((box_min.x()) > (tmp.x())) box_min.setX(tmp.x());
166 if((box_max.x()) < (tmp.x())) box_max.setX(tmp.x());
167 if((box_min.y()) > (tmp.y())) box_min.setY(tmp.y());
168 if((box_max.y()) < (tmp.y())) box_max.setY(tmp.y());
169 }
170
171 bbox = new G4BoundingBox3D(box_min, box_max);
172}
HepGeom::Point3D< G4double > G4Point3D
Definition: G4Point3D.hh:35
const G4Point3D PINFINITY(kInfinity, kInfinity, kInfinity)
G4int GetCols() const
G4int GetRows() const
G4Point3D Get3D(G4int i, G4int j) const

Referenced by BIntersect().

◆ ChangeDir()

void G4BezierSurface::ChangeDir ( )
inline

◆ ClipBothDirs()

G4int G4BezierSurface::ClipBothDirs ( )

Definition at line 96 of file G4BezierSurface.cc.

97{
98 dir = ROW;
99 ClipSurface();
100
101 // G4cout << "\n CLIP BOTH DIRS 1: " << smin << " " << smax;
102
103 if(smin > 1.0 || smax < 0.0)
104 {
106 return 1;
107 }
108 else
109 if((smax - smin) > 0.8)
110 {
111 SplitNURBSurface();
112 return 0;
113 }
114
115 LocalizeClipValues();
116 SetValues();
117
118 // Other G4Vector3D clipping and testing.
119 dir = COL;
120 ClipSurface();
121 // G4cout << "\n CLIP BOTH DIRS 2: " << smin << " " << smax;
122
123 if(smin > 1.0 || smax < 0.0)
124 {
126 return 1;
127 }
128 else
129 if((smax - smin) > 0.8)
130 {
131 SplitNURBSurface();
132 return 0;
133 }
134
135 LocalizeClipValues();
136 SetValues();
137 CalcAverage();
138 return 1;
139}
#define COL
Definition: G4PointRat.hh:52
#define ROW
Definition: G4PointRat.hh:51
void RemoveSurface(G4Surface *srf)

◆ ClipSurface()

void G4BezierSurface::ClipSurface ( )

Definition at line 301 of file G4BezierSurface.cc.

302{
303 // This routine is described in Computer Graphics, Volume 24,
304 // Number 4, August 1990 under the title Ray Tracing Trimmed
305 // Rational Surface Patches.
306
307
308 // G4cout << "\nBezier clip.";
309
310 register G4int i,j;
311 register G4int col_size = ctl_points->GetCols();
312 register G4int row_size = ctl_points->GetRows();
313
314 G4ConvexHull *ch_tmp= new G4ConvexHull(0,1.0e8,-1.0e8);
315 G4ConvexHull *ch_ptr=0, *ch_first=0;
316
317 // The four cornerpoints of the controlpoint mesh.
318
319/* L. Broglia
320 G4Point2d pt1 = ctl_points->get(0,0);
321 G4Point2d pt2 = ctl_points->get(0,col_size-1);
322 G4Point2d pt3 = ctl_points->get(row_size-1,0);
323 G4Point2d pt4 = ctl_points->get(row_size-1,col_size-1);
324 G4Point2d v1,v2,v3;
325*/
326 G4Point3D pt1 = ctl_points->Get3D(0,0);
327 G4Point3D pt2 = ctl_points->Get3D(0,col_size-1);
328 G4Point3D pt3 = ctl_points->Get3D(row_size-1,0);
329 G4Point3D pt4 = ctl_points->Get3D(row_size-1,col_size-1);
330 G4Point3D v1,v2,v3;
331
332 if ( dir == ROW)
333 {
334 // Vectors from cornerpoints
335 v1 = (pt1 - pt3);
336 // v1.X() = pt1.X() - pt3.X();
337 // v1.Y() = pt1.Y() - pt3.Y();
338 v2 = (pt2 - pt4);
339 // v2.X() = pt2.X() - pt4.X();
340 // v2.Y() = pt2.Y() - pt4.Y();
341 }
342 else
343 {
344 v1 = pt1 - pt2;
345 v2 = pt3 - pt4;
346 // v1.X() = pt1.X() - pt2.X();
347 // v1.Y() = pt1.Y() - pt2.Y();
348 // v2.X() = pt3.X() - pt4.X();
349 // v2.Y() = pt3.Y() - pt4.Y();
350 }
351/* L. Broglia
352 v3.X(v1.X() + v2.X());
353 v3.Y(v1.Y() + v1.Y());
354*/
355 v3 = v1 + v2 ;
356
357 smin = 1.0e8;
358 smax = -1.0e8;
359
360 G4double norm = std::sqrt(v3.x() * v3.x() + v3.y() * v3.y());
361 if(!norm)
362 {
363 G4cout << "\nNormal zero!";
364 G4cout << "\nLINE & DIR: " << line.x() << " " << line.y() << " " << dir;
365 G4cout << "\n";
366
367 if((std::abs(line.x())) > kCarTolerance)
368 line.setX(-line.x());
369 else
370 if((std::abs(line.y())) > kCarTolerance)
371 line.setY(-line.y());
372 else
373 {
374 G4cout << "\n RETURNING FROm CLIP..";
375 smin = 0; smax = 1; delete ch_tmp;
376 return;
377 }
378
379 G4cout << "\nCHANGED LINE & DIR: " << line.x() << " "
380 << line.y() << " " << dir;
381 }
382 else
383 {
384 line.setX( v3.y() / norm);
385 line.setY(-v3.x() / norm);
386 }
387
388 // smin = 1.0e8;
389 // smax = -1.0e8;
390 // G4cout << "\n FINAL LINE & DIR: " << line.X() << " "
391 // << line.Y() << " " << dir;
392
393 if( dir == ROW)
394 {
395 // Create a Convex() hull List
396 for(G4int a = 0; a < col_size; a++)
397 {
398 ch_ptr = new G4ConvexHull(a/(col_size - 1.0),1.0e8,-1.0e8);
399 if(! a)
400 {
401 ch_first=ch_ptr; delete ch_tmp;
402 }
403 else ch_tmp->SetNextHull(ch_ptr);
404
405 ch_tmp=ch_ptr;
406 }
407
408 ch_ptr=ch_first;
409 register G4double value;
410
411 // Loops through the control point mesh and calculates
412 // the nvex() hull for the surface.
413
414 for( G4int h = 0; h < row_size; h++)
415 {
416 for(G4int k = 0; k < col_size; k++)
417 {
418/* L. Broglia
419 G4Point2d& coordstmp = (G4Point2d&)ctl_points->get(h,k);
420 value = - ((coordstmp.X() * line.X() + coordstmp.Y() * line.Y()));
421*/
422 G4Point3D coordstmp = ctl_points->Get3D(h,k);
423 value = - ((coordstmp.x() * line.x() + coordstmp.y() * line.y()));
424
425 if( value <= (ch_ptr->GetMin()+kCarTolerance)) ch_ptr->SetMin(value);
426 if( value >= (ch_ptr->GetMax()-kCarTolerance)) ch_ptr->SetMax(value);
427
428 ch_ptr=ch_ptr->GetNextHull();
429 }
430
431 ch_ptr=ch_first;
432 }
433
434 ch_ptr=ch_first;
435 // Finds the points where the nvex() hull intersects
436 // with the coordinate .X()is. These points are the
437 // minimum and maximum values to where to clip the
438 // surface.
439
440 for(G4int l = 0; l < col_size - 1; l++)
441 {
442 ch_tmp=ch_ptr->GetNextHull();
443 for(G4int q = l+1; q < col_size; q++)
444 {
445 register G4double d, param1, param2;
446 param1 = ch_ptr->GetParam();
447 param2 = ch_tmp->GetParam();
448
449 if(ch_tmp->GetMax() - ch_ptr->GetMax())
450 {
451 d = Findzero( param1, param2, ch_ptr->GetMax(), ch_tmp->GetMax());
452 if( d <= (smin + kCarTolerance) ) smin = d * .99;
453 if( d >= (smax - kCarTolerance) ) smax = d * .99 + .01;
454 }
455
456 if(ch_tmp->GetMin() - ch_ptr->GetMin())
457 {
458 d = Findzero( param1, param2, ch_ptr->GetMin(), ch_tmp->GetMin());
459 if( d <= (smin + kCarTolerance)) smin = d * .99;
460 if( d >= (smax - kCarTolerance)) smax = d * .99 + .01;
461 }
462
463 ch_tmp=ch_tmp->GetNextHull();
464 }
465
466 ch_ptr=ch_ptr->GetNextHull();
467 }
468
469 ch_ptr=ch_first;
470
471 if (smin <= 0.0) smin = 0.0;
472 if (smax >= 1.0) smax = 1.0;
473
474 if ( (ch_ptr)
475 && (Sign(ch_ptr->GetMin()) != Sign(ch_ptr->GetMax()))) smin = 0.0;
476
477 i = Sign(ch_tmp->GetMin()); // ch_tmp points to last nvex()_hull in List
478 j = Sign(ch_tmp->GetMax());
479
480 if ( std::abs(i-j) > kCarTolerance ) smax = 1.0;
481 // if ( i != j) smax = 1.0;
482
483 }
484 else // Other G4Vector3D
485 {
486 for(G4int n = 0; n < row_size; n++)
487 {
488 ch_ptr = new G4ConvexHull(n/(row_size - 1.0),1.0e8,-1.0e8);
489 if(!n)
490 {
491 ch_first=ch_ptr; delete ch_tmp;
492 }
493 else ch_tmp->SetNextHull(ch_ptr);
494
495 ch_tmp=ch_ptr;
496 }
497
498 ch_ptr=ch_first;
499
500 for( G4int o = 0; o < col_size; o++)
501 {
502 for(G4int p = 0; p < row_size; p++)
503 {
504 register G4double value;
505
506/* L. Broglia
507 G4Point2d& coordstmp =(G4Point2d&) ctl_points->get(p,o);
508 value = - ((coordstmp.X() * line.X() + coordstmp.Y() * line.Y()));
509*/
510 G4Point3D coordstmp = ctl_points->Get3D(p,o);
511 value = - ((coordstmp.x() * line.x() + coordstmp.y() * line.y()));
512
513 if( value <= (ch_ptr->GetMin()+kCarTolerance)) ch_ptr->SetMin(value);
514 if( value >= (ch_ptr->GetMax()-kCarTolerance)) ch_ptr->SetMax(value);
515
516 ch_ptr=ch_ptr->GetNextHull();
517 }
518
519 ch_ptr=ch_first;
520 }
521
522 ch_ptr=ch_first;
523 delete ch_tmp;
524
525 for(G4int q = 0; q < row_size - 1; q++)
526 {
527 ch_tmp=ch_ptr->GetNextHull();
528 for(G4int r = q+1; r < row_size; r++)
529 {
530 register G4double param1 = ch_ptr->GetParam();
531 register G4double param2 = ch_tmp->GetParam();
532 register G4double d;
533
534 if(ch_tmp->GetMax() - ch_ptr->GetMax())
535 {
536 d = Findzero( param1, param2, ch_ptr->GetMax(), ch_tmp->GetMax());
537 if( d <= (smin + kCarTolerance) ) smin = d * .99;
538 if( d >= (smax - kCarTolerance) ) smax = d * .99 + .01;
539 }
540
541 if(ch_tmp->GetMin()-ch_ptr->GetMin())
542 {
543 d = Findzero( param1, param2, ch_ptr->GetMin(), ch_tmp->GetMin());
544 if( d <= (smin + kCarTolerance) ) smin = d * .99;
545 if( d >= (smax - kCarTolerance) ) smax = d * .99 + .01;
546 }
547
548 ch_tmp=ch_tmp->GetNextHull();
549 }
550
551 ch_ptr=ch_ptr->GetNextHull();
552 }
553
554 ch_tmp=ch_ptr;
555 ch_ptr=ch_first;
556
557 if (smin <= 0.0) smin = 0.0;
558 if (smax >= 1.0) smax = 1.0;
559
560 if ( (ch_ptr)
561 && (Sign(ch_ptr->GetMin()) != Sign(ch_ptr->GetMax()))) smin = 0.0;
562
563 if ( ch_tmp )
564 {
565 i = Sign(ch_tmp->GetMin()); // ch_tmp points to last nvex()_hull in List
566 j = Sign(ch_tmp->GetMax());
567 if ( (std::abs(i-j) > kCarTolerance)) smax = 1.0;
568 }
569 }
570
571 ch_ptr=ch_first;
572 while(ch_ptr && (ch_ptr!=ch_ptr->GetNextHull()))
573 {
574 ch_tmp=ch_ptr;
575 ch_ptr=ch_ptr->GetNextHull();
576 delete ch_tmp;
577 }
578
579 delete ch_ptr;
580
581 // Testing...
582 Clips++;
583}
double G4double
Definition: G4Types.hh:64
G4DLLIMPORT std::ostream G4cout
static G4int Clips
G4double GetMin() const
Definition: G4ConvexHull.hh:57
G4ConvexHull * GetNextHull()
Definition: G4ConvexHull.hh:55
G4double GetMax() const
Definition: G4ConvexHull.hh:58
G4double GetParam() const
Definition: G4ConvexHull.hh:56
void SetMin(G4double x)
Definition: G4ConvexHull.hh:62
void SetNextHull(G4ConvexHull *n)
Definition: G4ConvexHull.hh:60
void SetMax(G4double y)
Definition: G4ConvexHull.hh:63
G4double kCarTolerance
Definition: G4Surface.hh:192

Referenced by BIntersect(), and ClipBothDirs().

◆ Dir()

void G4BezierSurface::Dir ( G4int  d)
inline

◆ GetOrder()

G4int G4BezierSurface::GetOrder ( G4int  direction) const
inline

◆ GetU()

G4double G4BezierSurface::GetU ( ) const
inline

◆ GetV()

G4double G4BezierSurface::GetV ( ) const
inline

◆ PutOrder()

void G4BezierSurface::PutOrder ( G4int  direction,
G4int  value 
)
inline

◆ SetAveragePoint()

void G4BezierSurface::SetAveragePoint ( G4Point3D  p)
inline

◆ SMax()

G4double G4BezierSurface::SMax ( ) const
inline

◆ SMin()

G4double G4BezierSurface::SMin ( ) const
inline

◆ SurfaceNormal()

G4Vector3D G4BezierSurface::SurfaceNormal ( const G4Point3D Pt) const
virtual

Implements G4Surface.

Definition at line 91 of file G4BezierSurface.cc.

92{
93 return G4Vector3D(0,0,0);
94}
HepGeom::Vector3D< G4double > G4Vector3D
Definition: G4Vector3D.hh:35

◆ UAverage()

G4double G4BezierSurface::UAverage ( ) const
inline

◆ VAverage()

G4double G4BezierSurface::VAverage ( ) const
inline

Friends And Related Function Documentation

◆ CopySurface

void CopySurface ( G4BezierSurface bez)
friend

◆ G4BSplineSurface

friend class G4BSplineSurface
friend

Definition at line 53 of file G4BezierSurface.hh.

◆ G4ProjectedSurface

friend class G4ProjectedSurface
friend

Definition at line 54 of file G4BezierSurface.hh.

Member Data Documentation

◆ bezier_list

G4SurfaceList* G4BezierSurface::bezier_list
protected

Definition at line 93 of file G4BezierSurface.hh.

Referenced by BIntersect(), and ClipBothDirs().

◆ Clips

G4int G4BezierSurface::Clips =0
staticprotected

Definition at line 95 of file G4BezierSurface.hh.

Referenced by ClipSurface().

◆ Splits

G4int G4BezierSurface::Splits =0
staticprotected

Definition at line 96 of file G4BezierSurface.hh.

◆ Tolerance

G4double G4BezierSurface::Tolerance =0
staticprotected

Definition at line 97 of file G4BezierSurface.hh.

Referenced by BIntersect().


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