Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4BSplineSurface.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 //
27 // $Id$
28 //
29 // ----------------------------------------------------------------------
30 // GEANT 4 class source file
31 //
32 // G4BSplineSurface.cc
33 //
34 // ----------------------------------------------------------------------
35 
36 #include "G4BSplineSurface.hh"
37 #include "G4BezierSurface.hh"
38 #include "G4ControlPoints.hh"
39 #include "G4BoundingBox3D.hh"
40 
42  : ord(0), k_index(0), param(0.), Rational(0)
43 {
44  distance = kInfinity;
45  dir=ROW;
46  first_hit = Hit = (G4UVHit*)0;
47  order[0] = 0; order[1] = 0;
48  ctl_points = (G4ControlPoints*)0;
49  u_knots = v_knots = tmp_knots = (G4KnotVector*)0;
50 }
51 
52 
54  : dir(0), ord(0), k_index(0), param(0.), Rational(0)
55 {
56  distance = kInfinity;
57  first_hit = Hit = (G4UVHit*)0;
58  order[0] = 0; order[1] = 0;
59  ctl_points = (G4ControlPoints*)0;
60  u_knots = v_knots = tmp_knots = (G4KnotVector*)0;
61 }
62 
63 
65  G4KnotVector& v_kv, G4ControlPoints& cp)
66  : dir(0), ord(0), k_index(0), param(0.), Rational(0)
67 {
68  first_hit = Hit = (G4UVHit*)0;
69 
70  order[0] = u+1; order[1] = v+1;
71 
72  u_knots = new G4KnotVector(u_kv);
73  v_knots = new G4KnotVector(v_kv);
74  tmp_knots = (G4KnotVector*)0;
75 
76  ctl_points = new G4ControlPoints(cp);
77 }
78 
79 
81 {
82  delete u_knots;
83  delete v_knots;
84  delete ctl_points;
85  G4UVHit* temphit=Hit;
86  Hit = first_hit;
87  while(Hit!=(G4UVHit*)0)
88  {
89  Hit=Hit->GetNext();
90  delete temphit;
91  temphit=Hit;
92  }
93  // delete temphit;// remove last
94 
95 }
96 
97 
99 {
100  Intersected = 1;
101  FindIntersections(rayref);
102  G4BezierSurface *bez_ptr;
103  bezier_list.MoveToFirst();
104  distance = kInfinity;
105 
106  while( bezier_list.GetSurface() != (G4Surface*)0)
107  {
108  bez_ptr = (G4BezierSurface*)bezier_list.GetSurface();
109 
110  if(bez_ptr->IsActive())
111  {
112  if(distance > bez_ptr->GetDistance())
113  {
114  // Put data from closest bezier to b-spline data struct
115  closest_hit = bez_ptr->AveragePoint();
116  distance = bez_ptr->GetDistance();
117  }
118  else
119  {
120  // Set other beziers as inactive
121  bez_ptr->SetActive(0);
122 
123  // Remove beziers that are not closest
124  // bezier_list.RemoveSurface(bez_ptr);
125  }
126  }
127 
128  bezier_list.Step();
129  }
130 
131  bezier_list.MoveToFirst();
132 
133  if(bezier_list.GetSize())
134  return 1;
135  else
136  {
137  active=0;
138  return 0;
139  }
140 }
141 
142 
143 G4Point3D G4BSplineSurface::FinalIntersection()
144 {
145  // Compute the real intersection point.
146  G4BezierSurface* bez_ptr;
147  while ( bezier_list.GetSize() > 0 &&
148  bezier_list.GetSurface() != (G4Surface*)0)
149  {
150  bez_ptr = (G4BezierSurface*)bezier_list.GetSurface();
151  G4int tmp = 0;
152 
153  // L. Broglia
154  // Modify G4BezierSurface intersection function name
155  // tmp = bez_ptr->Intersect( bezier_list);
156  tmp = bez_ptr->BIntersect( bezier_list);
157 
158  if(!tmp)
159  {
160  bezier_list.RemoveSurface(bez_ptr);
161  if(bezier_list.GetSurface() != (G4Surface*)0)
162  bezier_list.GetSurface()->SetActive(1);
163  }
164  else
165  if(tmp==1)
166  {
167  active=1;
168  // Hit found
169  AddHit(bez_ptr->GetU(), bez_ptr->GetV());
170 
171  // Delete beziers
172  bezier_list.EmptyList();
173  }
174  else
175  if(tmp==2)
176  {
177  // The bezier was split so the last
178  // two surfaces in the List should
179  // be bbox tested and if passed
180  // clipped in both dirs.
181 
182  // Move to first
183  bezier_list.MoveToFirst();
184  // Find the second last.
185 // What!? Casting a G4Surface* to a G4SurfaceList* !?!?!? - GC
186 //
187 // if(bezier_list.index != bezier_list.last)
188 // while ( ((G4SurfaceList*)bezier_list.index)->next !=
189 // bezier_list.last) bezier_list.Step();
190 //
191 // Try the following instead (if that's the wished behavior)...
192 //
193  if(bezier_list.GetSurface() != bezier_list.GetLastSurface())
194  while (bezier_list.GetNext() != bezier_list.GetLastSurface())
195  bezier_list.Step();
196 
197  G4BezierSurface* tmpbz = (G4BezierSurface*) bezier_list.GetSurface();
198  tmpbz->CalcBBox();
199 
200 // L. Broglia tmpbz->bbox->Test();
201 
202  G4int result=0;
203  if(tmpbz->bbox->GetTestResult())
204  {
205  // Clip
206  while(!result)
207  result = tmpbz->ClipBothDirs();
208  }
209  else
210  {
211  bezier_list.RemoveSurface(tmpbz);
212  }
213  // Second surface
214  tmpbz = (G4BezierSurface*) bezier_list.GetLastSurface();
215  tmpbz->CalcBBox();
216 
217 // L. Broglia tmpbz->bbox->Test();
218 
219  if(tmpbz->bbox->GetTestResult())
220  {
221  result = 0;
222  while(!result)
223  result = tmpbz->ClipBothDirs();
224  }
225  else
226  {
227  bezier_list.RemoveSurface(tmpbz);
228  }
229 
230  bezier_list.RemoveSurface(bez_ptr);
231  bezier_list.MoveToFirst();
232  }
233 
234  bezier_list.Step();
235  }//While....
236 
237  Hit = first_hit;
238  G4Point3D result;
239  if(Hit == (G4UVHit*)0)
240  active = 0;
241  else
242  {
243  while(Hit != (G4UVHit*)0)
244  {
245  // L. Broglia
246  // Modify function name
247  // result = Evaluate();
248  result = BSEvaluate();
249 
250  Hit = Hit->GetNext();
251  }
252 
253  Hit = first_hit;
254  }
255 
256  return result;
257 }
258 
259 
261 {
262 
263  // Finds the bounds of the b-spline surface iow
264  // calculates the bounds for a bounding box
265  // to the surface. The bounding box is used
266  // for a preliminary check of intersection.
267 
268  G4Point3D box_min = G4Point3D( PINFINITY);
269  G4Point3D box_max = G4Point3D(-PINFINITY);
270 
271  // Loop to search the whole control point mesh
272  // for the minimum and maximum values for x, y and z.
273 
274  for(register int a = ctl_points->GetRows()-1; a>=0;a--)
275  for(register int b = ctl_points->GetCols()-1; b>=0;b--)
276  {
277  G4Point3D tmp = ctl_points->Get3D(a,b);
278  if((box_min.x()) > (tmp.x())) box_min.setX(tmp.x());
279  if((box_min.y()) > (tmp.y())) box_min.setY(tmp.y());
280  if((box_min.z()) > (tmp.z())) box_min.setZ(tmp.z());
281  if((box_max.x()) < (tmp.x())) box_max.setX(tmp.x());
282  if((box_max.y()) < (tmp.y())) box_max.setY(tmp.y());
283  if((box_max.z()) < (tmp.z())) box_max.setZ(tmp.z());
284  }
285  bbox = new G4BoundingBox3D( box_min, box_max);
286 }
287 
288 
289 G4ProjectedSurface* G4BSplineSurface::CopyToProjectedSurface
290 (const G4Ray& rayref)
291 {
292  G4ProjectedSurface* proj_srf = new G4ProjectedSurface() ;
293  proj_srf->PutOrder(0,GetOrder(0));
294  proj_srf->PutOrder(1,GetOrder(1));
295  proj_srf->dir = dir;
296 
297  proj_srf->u_knots = new G4KnotVector(*u_knots);
298  proj_srf->v_knots = new G4KnotVector(*v_knots);
299  proj_srf->ctl_points = new G4ControlPoints
300  (2, ctl_points->GetRows(), ctl_points->GetCols());
301 
302  const G4Plane& plane1 = rayref.GetPlane(1);
303  const G4Plane& plane2 = rayref.GetPlane(2);
304  ProjectNURBSurfaceTo2D(plane1, plane2, proj_srf);
305 
306  return proj_srf;
307 }
308 
309 
310 void G4BSplineSurface::FindIntersections(const G4Ray& rayref)
311 {
312  // Do the projection to 2D
313  G4ProjectedSurface* proj_srf = CopyToProjectedSurface(rayref);
314 
315  // Put surface in projected List
316  projected_list.AddSurface(proj_srf);
317 
318  // Loop through List of projected surfaces
319  while(projected_list.GetSize() > 0)
320  {
321  // Get first in List
322  proj_srf = (G4ProjectedSurface*)projected_list.GetSurface();
323 
324  // Create the bounding box for the projected surface.
325  proj_srf->CalcBBox();
326 
327 // L. Broglia proj_srf->bbox->Test();
328 
329  // Check bbox test result is ok
330  if(proj_srf->bbox->GetTestResult())
331  // Convert the projected surface to a bezier. Split if necessary.
332  proj_srf->ConvertToBezier(projected_list, bezier_list);
333 
334  // Remove projected surface
335  projected_list.RemoveSurface(proj_srf);
336  }
337 
338  // Loop through the bezier List
339  G4BezierSurface* bez_ptr;
340  distance = kInfinity;
341 
342  while(bezier_list.GetSurface() != (G4Surface*)0)
343  {
344  bez_ptr = (G4BezierSurface*)bezier_list.GetSurface();
345 
346  // Add a temporary Hit
347  AddHit(bez_ptr->UAverage(), bez_ptr->VAverage());
348 
349  // Evaluate Hit
350 
351  // L. Broglia
352  // Modify function name
353  // bez_ptr->SetAveragePoint(Evaluate());
354  bez_ptr->SetAveragePoint(BSEvaluate());
355 
356  // Calculate distance to ray origin
357  bez_ptr->CalcDistance(rayref.GetStart());
358 
359  // Put closest to b_splines distance value
360  if(bez_ptr->GetDistance() < distance) distance = bez_ptr->GetDistance();
361 
362  // Remove the temporary Hit
363  if (first_hit == Hit) first_hit = (G4UVHit*)0;
364  delete Hit;
365  Hit = (G4UVHit*)0;
366 
367  // Move to next in the List
368  bezier_list.Step();
369  }
370 
371  bezier_list.MoveToFirst();
372  if(bezier_list.GetSize() == 0)
373  {
374  active=0;
375  return;
376  }
377 
378  // Check that approx Hit is in direction of ray
379  const G4Point3D& Pt = rayref.GetStart();
380  const G4Vector3D& Dir = rayref.GetDir();
381  G4Point3D TestPoint = G4Point3D( (0.00001*Dir) + Pt );
382  G4BezierSurface* Bsrf = (G4BezierSurface*)bezier_list.GetSurface(0);
383 
384  G4Point3D AveragePoint = Bsrf->AveragePoint();
385  G4double TestDistance = TestPoint.distance2(AveragePoint);
386 
387  if(TestDistance > distance)
388  // Hit behind ray starting point, no intersection.
389  active=0;
390 }
391 
392 
393 void G4BSplineSurface::AddHit(G4double u, G4double v)
394 {
395  if(Hit == (G4UVHit*)0)
396  {
397  first_hit = new G4UVHit(u,v);
398  first_hit->SetNext(0);
399  Hit = first_hit;
400  }
401  else
402  {
403  Hit->SetNext(new G4UVHit(u,v));
404  Hit = Hit->GetNext();
405  Hit->SetNext(0);
406  }
407 }
408 
409 
410 void G4BSplineSurface::ProjectNURBSurfaceTo2D
411  (const G4Plane& plane1, const G4Plane& plane2,
412  register G4ProjectedSurface* proj_srf)
413 {
414  // Projects the nurb surface so that the z-axis = ray.
415 
416  /* L. Broglia
417  G4Point* tmp = (G4Point*)&ctl_points->get(0,0);
418  */
419 
420  G4PointRat tmp = ctl_points->GetRat(0,0);
421  G4int rational = tmp.GetType();// Get the type of control point
422  G4Point3D psrfcoords;
423  register int rows = ctl_points->GetRows();
424  register int cols = ctl_points->GetCols();
425 
426  for (register int i=0; i< rows; i++)
427  for(register int j=0; j < cols;j++)
428  {
429  if ( rational==4 ) // 4 coordinates
430  {
431  G4PointRat& srfcoords = ctl_points->GetRat(i, j);
432 
433 // L. Broglia
434 // Changes for new G4PointRat
435 
436  // Calculate the x- and y-coordinates for the new
437  // 2-D surface.
438  psrfcoords.setX(( srfcoords.x() * plane1.a
439  +srfcoords.y() * plane1.b
440  +srfcoords.z() * plane1.c
441  -srfcoords.w() * plane1.d));
442  psrfcoords.setY(( srfcoords.x() * plane2.a
443  +srfcoords.y() * plane2.b
444  +srfcoords.z() * plane2.c
445  -srfcoords.w() * plane2.d));
446 
447  proj_srf->ctl_points->put(i,j,psrfcoords);
448  }
449  else // 3 coordinates
450  {
451  G4Point3D srfcoords = ctl_points->Get3D(i, j);
452 
453  psrfcoords.setX(( srfcoords.x() * plane1.a
454  +srfcoords.y() * plane1.b
455  +srfcoords.z() * plane1.c
456  - plane1.d));
457 
458  psrfcoords.setY(( srfcoords.x() * plane2.a
459  +srfcoords.y() * plane2.b
460  +srfcoords.z() * plane2.c
461  - plane2.d));
462 
463  proj_srf->ctl_points->put(i,j,psrfcoords);
464  }
465  }
466 }
467 
468 /* L. Broglia
469  Changes for new G4PointRat
470 G4Point& G4BSplineSurface::InternalEvalCrv(G4int i, G4ControlPoints* crv)*/
471 
472 G4PointRat& G4BSplineSurface::InternalEvalCrv(G4int i, G4ControlPoints* crv)
473 {
474  if ( ord <= 1 )
475  return crv->GetRat(i, k_index);
476 
477  register int j = k_index;
478 
479  while ( j > (k_index - ord + 1))
480  {
481  register G4double k1, k2;
482 
483  k1 = tmp_knots->GetKnot((j + ord - 1));
484  k2 = tmp_knots->GetKnot(j);
485 
486  if ((std::abs(k1 - k2)) > kCarTolerance )
487  {
488  /* L. Broglia
489  register G4PointRat* pts1 = &crv->get(i,j-1);
490  register G4PointRat* pts2 = &crv->get(i,j );
491  if(pts1->GetType()==3)
492  {
493  crv->CalcValues(k1, param, *(G4Point3D*)pts1, k2, *(G4Point3D*)pts2);
494  crv->put(0, j, *(G4Point3D*)pts2);
495  }
496  else
497  {
498  crv->CalcValues(k1, param, *(G4PointRat*)pts1, k2, *(G4PointRat*)pts2);
499  crv->put(0, j, *(G4PointRat*)pts2);
500  }
501  register G4PointRat* pts1 = &crv->GetRat(i,j-1);
502  register G4PointRat* pts2 = &crv->GetRat(i,j );
503  */
504  }
505 
506  j--;
507  }
508 
509  ord = ord-1;
510  return InternalEvalCrv(0, crv); // Recursion
511 }
512 
513 
514 G4Point3D G4BSplineSurface::BSEvaluate()
515 {
516  register int i;
517  register int row_size = ctl_points->GetRows();
518  register G4ControlPoints *diff_curve;
519  register G4ControlPoints* curves;
520  G4Point3D result;
521 
522  /* L. Broglia
523  G4Point* tmp = (G4Point*)&ctl_points->get(0,0);
524  */
525 
526  G4PointRat* tmp = &ctl_points->GetRat(0,0);
527 
528  register int point_type = tmp->GetType();
529  diff_curve = new G4ControlPoints(point_type, row_size, 1);
530  k_index = u_knots->GetKnotIndex(Hit->GetU(), GetOrder(ROW) );
531 
532  ord = GetOrder(ROW);
533  if(k_index==-1)
534  {
535  delete diff_curve;
536  active = 0;
537  return result;
538  }
539 
540  curves=new G4ControlPoints(*ctl_points);
541  tmp_knots = u_knots;
542  param = Hit->GetU();
543 
544  if(point_type == 4)
545  {
546  for ( i = 0; i < row_size; i++)
547  {
548  ord = GetOrder(ROW);
549  G4PointRat rtr_pt = InternalEvalCrv(i, curves);
550  diff_curve->put(0,i,rtr_pt);
551  }
552 
553  k_index = v_knots->GetKnotIndex( Hit->GetV(), GetOrder(COL) );
554  if(k_index==-1)
555  {
556  delete diff_curve;
557  delete curves;
558  active = 0;
559  return result;
560  }
561 
562  ord = GetOrder(COL);
563  tmp_knots = v_knots;
564  param = Hit->GetV();
565 
566  // Evaluate the diff_curve...
567  // G4PointRat rat_result = (G4PointRat&) InternalEvalCrv(0, diff_curve);
568  G4PointRat rat_result(InternalEvalCrv(0, diff_curve));
569 
570  // Calc the 3D values.
571  // L. Broglia
572  // Changes for new G4PointRat
573  result.setX(rat_result.x()/rat_result.w());
574  result.setY(rat_result.y()/rat_result.w());
575  result.setZ(rat_result.z()/rat_result.w());
576  }
577  else
578  if(point_type == 3)
579  {
580  for ( i = 0; i < row_size; i++)
581  {
582  ord = GetOrder(ROW);
583  // G4Point3D rtr_pt = (G4Point3D&) InternalEvalCrv(i, curves);
584  G4Point3D rtr_pt = (InternalEvalCrv(i, curves)).pt();
585  diff_curve->put(0,i,rtr_pt);
586  }
587 
588  k_index = v_knots->GetKnotIndex( Hit->GetV(), GetOrder(COL) );
589  if(k_index==-1)
590  {
591  delete diff_curve;
592  delete curves;
593  active = 0;
594  return result;
595  }
596 
597  ord = GetOrder(COL);
598  tmp_knots = v_knots;
599  param = Hit->GetV();
600 
601  // Evaluate the diff_curve...
602  result = (InternalEvalCrv(0, diff_curve)).pt();
603  }
604 
605  delete diff_curve;
606  delete curves;
607  closest_hit = result;
608  return result;
609 }
610 
611 
612 G4Point3D G4BSplineSurface::Evaluation(const G4Ray&)
613 {
614  // Delete old UVhits
615  G4UVHit* temphit=Hit;
616  while(Hit!=(G4UVHit*)0)
617  {
618  Hit=Hit->GetNext();
619  delete temphit;
620  temphit=Hit;
621  }
622 
623  // delete temphit;
624 
625  // Get the real Hit point
626  closest_hit = FinalIntersection();
627 
628  // The following part (commented out) is old bullshit
629  // Check that Hit is not in a void i.e. InnerBoundary.
630  // for(int a=0; a<NumberOfInnerBoundaries;a++)
631  // if(InnerBoundary[a]->Inside(closest_hit, rayref))
632  // {
633  // Active(0);
634  // Distance(kInfinity);
635  // return closest_hit;
636  // }
637  return closest_hit;
638 }
639 
640 
642 {
643  G4double PointDistance=0;
644  PointDistance = ctl_points->ClosestDistanceToPoint(Pt);
645  return PointDistance;
646 }