Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4BezierSurface.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// G4BezierSurface.cc
33//
34// ----------------------------------------------------------------------
35// History:
36// -------
37// - Replaced addition of coordinates by addition of 2 points
38// (L. Broglia, 10/10/98)
39// ----------------------------------------------------------------------
40
41#include "G4BezierSurface.hh"
42#include "G4ConvexHull.hh"
43
47
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}
58
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}
79
81 : G4Surface(), bezier_list(0), smin(0.), smax(0.),
82 average_u(0.), average_v(0.), dir(0), u_knots(0), v_knots(0),
83 ctl_points(0), new_knots(0), ord(0), oslo_m(0), lower(0), upper(0),
84 u_min(0.), u_max(0.), v_min(0.), v_max(0.), old_points(0)
85{
86 order[0]=0; order[1]=0;
87 u[0]=0.; u[1]=0.;
88 v[0]=0.; v[1]=0.;
89}
90
92{
93 return G4Vector3D(0,0,0);
94}
95
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}
140
141
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}
173
174
175void G4BezierSurface::CalcAverage()
176{
177 // Calculate the average point from the average clip-values.
178 average_u = (u_min + u_max)/2.0;
179 average_v = (v_min + v_max)/2.0;
180}
181
182
183void G4BezierSurface::CalcDistance(const G4Point3D& ray_start)
184{
185 // Calculate the distance between the average point and
186 // the ray starting point.
187 distance = ((((ray_start.x() - average_pt.x())*
188 (ray_start.x() - average_pt.x()))+
189 ((ray_start.y() - average_pt.y())*
190 (ray_start.y() - average_pt.y()))+
191 ((ray_start.z() - average_pt.z())*
192 (ray_start.z() - average_pt.z()))));
193}
194
195
196void G4BezierSurface::SetValues()
197{
198 if(dir)
199 {
200 v_min = smin;
201 v_max = smax;
202 }
203 else
204 {
205 u_min = smin;
206 u_max = smax;
207 }
208}
209
210
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}
299
300
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}
584
585
586void G4BezierSurface::GetClippedRegionFromSurface()
587{
588 // Returns the clipped part of the surface. First calculates the
589 // length of the new knotvector. Then uses the refinement function to
590 // get the new knotvector and controlmesh.
591
592 // G4cout << "\nBezier region clipped.";
593
594 delete new_knots;
595 if ( dir == ROW)
596 {
597 new_knots = new G4KnotVector(GetOrder(0) * 2);
598 for (register G4int i = 0; i < GetOrder(0); i++)
599 {
600 new_knots->PutKnot(i, smin);
601 new_knots->PutKnot(i+ GetOrder(0), smax);
602 }
603 }
604 else
605 {
606 new_knots = new G4KnotVector( GetOrder(1) * 2);
607 for ( register G4int i = 0; i < GetOrder(1); i++)
608 {
609 new_knots->PutKnot(i, smin);
610 new_knots->PutKnot(i+ GetOrder(1), smax);
611 }
612 }
613} // NURB_REGION_FROM_SURFACE
614
615
616void G4BezierSurface::RefineSurface()
617{
618 // Returns the new clipped surface. Calculates the new controlmesh
619 // and knotvectorvalues for the surface by using the Oslo-algorithm
620
621 delete old_points;
622 if (dir == ROW)
623 {
624 // Row (u) G4Vector3D
625 ord = GetOrder(0);
626 CalcOsloMatrix();
627 for(register G4int a=0;a<new_knots->GetSize();a++)
628 u_knots->PutKnot(a, new_knots->GetKnot(a));
629
630 lower = 0;
631 upper = new_knots->GetSize() - GetOrder(0);
632
633 // Copy of the old points.
634 old_points = new G4ControlPoints(*ctl_points);
635 MapSurface(this);
636 }
637 else
638 {
639 ord = GetOrder(1);
640 CalcOsloMatrix ();
641 for(register G4int a=0;a < new_knots->GetSize();a++)
642 v_knots->PutKnot(a, new_knots->GetKnot(a));
643
644 // Copy of the old points.
645 old_points = new G4ControlPoints(*ctl_points);
646
647 // Make new controlpoint matrix,
648 register G4int cols = ctl_points->GetCols();
649 delete ctl_points;
650
651 ctl_points = new G4ControlPoints(2,(new_knots->GetSize()-
652 GetOrder(1)),cols);
653 lower = 0;
654 upper = new_knots->GetSize() - GetOrder(1);
655 MapSurface(this);
656 }
657}// REFINE_SURFACE
658
659
660void G4BezierSurface::CalcOsloMatrix()
661{
662 // This algorithm is described in the paper "Making the Oslo-algorithm
663 // more efficient" in SIAM J.NUMER.ANAL. Vol.23, No. 3, June '86
664 // Calculates the oslo-matrix , which is used in mapping the new
665 // knotvector- and controlpoint-values.
666
667 register G4KnotVector *ah;
668 register G4KnotVector *newknots;
669 register G4int i;
670 register G4int j;
671 register G4int mu, muprim;
672 register G4int vv, p;
673 register G4int iu, il, ih, n1;
674 register G4int ahi;
675 register G4double beta1;
676 register G4double tj;
677
678 ah = new G4KnotVector(ord*(ord + 1)/2);
679 newknots = new G4KnotVector(ord * 2 );
680
681 n1 = new_knots->GetSize() - ord;
682 mu = 0;
683
684 if(oslo_m!=(G4OsloMatrix*)0)
685 {
686 G4OsloMatrix* tmp;
687
688 // while(oslo_m!=oslo_m->next)
689 while(oslo_m!=(G4OsloMatrix*)0)
690 {
691 tmp=oslo_m->GetNextNode();delete oslo_m; oslo_m=tmp;
692 }
693 }
694
695 delete oslo_m;
696 oslo_m = new G4OsloMatrix();
697
698 register G4OsloMatrix* o_ptr = oslo_m;
699
700 register G4KnotVector* old_knots;
701 if(dir)
702 old_knots = v_knots;
703 else
704 old_knots = u_knots;
705
706 for (j = 0; j < n1; j++)
707 {
708 if ( j != 0 )
709 {
710 oslo_m->SetNextNode(new G4OsloMatrix());
711 oslo_m = oslo_m->GetNextNode();
712 }
713
714 while (old_knots->GetKnot(mu + 1) <= new_knots->GetKnot(j))
715 mu = mu + 1; // find the bounding mu
716
717 i = j + 1;
718 muprim = mu;
719
720 while ((new_knots->GetKnot(i) == old_knots->GetKnot(muprim)) &&
721 i < (j + ord))
722 {
723 i++;
724 muprim--;
725 }
726
727 ih = muprim + 1;
728
729 for (vv = 0, p = 1; p < ord; p++)
730 {
731 if (new_knots->GetKnot(j + p) == old_knots->GetKnot(ih))
732 ih++;
733 else
734 newknots->PutKnot(++vv - 1,new_knots->GetKnot(j + p));
735 }
736
737 ahi = AhIndex(0, ord - 1,ord);
738 ah->PutKnot(ahi, 1.0);
739
740 for (p = 1; p <= vv; p++)
741 {
742 beta1 = 0.0;
743 tj = newknots->GetKnot(p-1);
744
745 if (p - 1 >= muprim)
746 {
747 beta1 = AhIndex(p - 1, ord - muprim,ord);
748 beta1 = ((tj - old_knots->GetKnot(0)) * beta1) /
749 (old_knots->GetKnot(p + ord - vv) - old_knots->GetKnot(0));
750 }
751
752 i = muprim - p + 1;
753 il = Amax (1, i);
754 i = n1 - 1 + vv - p;
755 iu = Amin (muprim, i);
756
757 for (i = il; i <= iu; i++)
758 {
759 register G4double d1, d2;
760 register G4double beta;
761
762 d1 = tj - old_knots->GetKnot(i);
763 d2 = old_knots->GetKnot(i + p + ord - vv - 1) - tj;
764
765 beta = ah->GetKnot(AhIndex(p - 1, i + ord - muprim - 1,ord)) /
766 (d1 + d2);
767
768
769 ah->PutKnot(AhIndex(p, i + ord - muprim - 2,ord), d2 * beta + beta1) ;
770 beta1 = d1 * beta;
771 }
772
773 ah->PutKnot(AhIndex(p, iu + ord - muprim - 1,ord), beta1);
774
775 if (iu < muprim)
776 {
777 register G4double kkk;
778 register G4double ahv;
779
780 kkk = old_knots->GetKnot(n1 - 1 + ord);
781 ahv = AhIndex (p - 1, iu + ord - muprim,ord);
782 ah->PutKnot(AhIndex(p, iu + ord - muprim - 1,ord),
783 beta1 + (kkk - tj) * ahv /
784 (kkk - old_knots->GetKnot(iu + 1)));
785 }
786 }
787
788 // Remove the oslo matrix List
789 G4OsloMatrix* temp_oslo = oslo_m;
790
791/*
792 if(oslo_m != (G4OsloMatrix*)0)
793 while(oslo_m->next != oslo_m)
794 {
795 oslo_m = oslo_m->next;
796 delete temp_oslo;
797 temp_oslo = oslo_m;
798 }
799
800 // Remove the last
801 delete oslo_m;
802*/
803
804 while(oslo_m != (G4OsloMatrix*)0)
805 {
806 oslo_m = oslo_m->GetNextNode();
807 delete temp_oslo;
808 temp_oslo = oslo_m;
809 }
810
811 delete oslo_m;
812
813 // Create a new oslo matrix
814 oslo_m = new G4OsloMatrix(vv+1, Amax(muprim - vv,0), vv);
815
816 for ( i = vv, p = 0; i >= 0; i--)
817 oslo_m->GetKnotVector()
818 ->PutKnot ( p++, ah->GetKnot(AhIndex (vv, (ord-1) - i,ord)));
819
820 }
821
822 delete ah;
823 delete newknots;
824 oslo_m->SetNextNode(0);
825 oslo_m = o_ptr;
826}
827
828
829void G4BezierSurface::MapSurface(G4Surface*)
830{
831 // This algorithm is described in the paper Making the Oslo-algorithm
832 // more efficient in SIAM J.NUMER.ANAL. Vol.23, No. 3, June '86
833 // Maps the new controlpoints into the new surface.
834
835 register G4ControlPoints *c_ptr;
836 register G4OsloMatrix *o_ptr;
837 register G4ControlPoints* new_pts;
838 register G4ControlPoints* old_pts;
839
840 new_pts = ctl_points;
841
842 // Copy the old points so they can be used in calculating the new ones.
843 // old_pts = new G4ControlPoints(*ctl_points);
844 old_pts = old_points;
845 register G4int j, // j loop
846 i; // oslo loop
847
848 c_ptr = new_pts;
849 register G4int size; // The number of rows or columns,
850 // depending on processing order
851
852 if(!dir)
853 size=new_pts->GetRows();
854 else
855 size=new_pts->GetCols();
856
857 for(G4int a=0; a<size;a++)
858 {
859 if ( lower != 0)
860 {
861 for ( i = 0, o_ptr = oslo_m;
862 i < lower;
863 i++, o_ptr = o_ptr->GetNextNode()){;}
864 }
865 else
866 {
867 o_ptr = oslo_m;
868 }
869
870 if(!dir)// Direction ROW
871 {
872 for ( j = lower; j < upper; j++, o_ptr = o_ptr->GetNextNode())
873 {
874 register G4double o_scale;
875 register G4int x;
876 x=a;
877
878/* L. Broglia
879 G4Point2d o_pts= (G4Point2d&)old_pts->Get2d(x, o_ptr->GetOffset());
880 G4Point2d tempc= (G4Point2d&)c_ptr->Get2d(j/upper,(j)%upper-lower);
881*/
882 G4Point3D o_pts = old_pts->Get3D(x, o_ptr->GetOffset());
883 G4Point3D tempc = c_ptr->Get3D(j/upper, (j)%upper-lower);
884
885 o_scale = o_ptr->GetKnotVector()->GetKnot(0);
886
887 tempc.setX(o_pts.x() * o_scale);
888 tempc.setY(o_pts.x() * o_scale);
889
890 for ( i = 1; i <= o_ptr->GetSize(); i++)
891 {
892 o_scale = o_ptr->GetKnotVector()->GetKnot(i);
893
894/* L. Broglia
895 o_pts = (G4Point2d&)old_pts->get(x, i+o_ptr->GetOffset());
896 tempc.X(tempc.X() + o_scale * o_pts.X());
897 tempc.Y(tempc.Y() + o_scale * o_pts.Y());
898*/
899 o_pts = old_pts->Get3D(x, i+o_ptr->GetOffset());
900 tempc.setX(tempc.x() + o_scale * o_pts.x());
901 tempc.setY(tempc.y() + o_scale * o_pts.y());
902
903 }
904
905 c_ptr->put(a,(j)%upper-lower,tempc);
906 }
907 }
908 else // dir = COL
909 {
910 for ( j = lower; j < upper; j++, o_ptr = o_ptr->GetNextNode())
911 {
912 register G4double o_scale;
913 register G4int x;
914 x=a;
915
916/* L. Broglia
917 G4Point2d o_pts= (G4Point2d&)old_pts->Get2d(o_ptr->GetOffset(), x);
918 G4Point2d tempc = (G4Point2d&)c_ptr->Get2d((j)%upper-lower,j/upper);
919*/
920 G4Point3D o_pts = old_pts->Get3D(o_ptr->GetOffset(), x);
921 G4Point3D tempc = c_ptr->Get3D((j)%upper-lower,j/upper);
922
923 o_scale = o_ptr->GetKnotVector()->GetKnot(0);
924
925 tempc.setX(o_pts.x() * o_scale);
926 tempc.setY(o_pts.y() * o_scale);
927
928 for ( i = 1; i <= o_ptr->GetSize(); i++)
929 {
930 o_scale = o_ptr->GetKnotVector()->GetKnot(i);
931/* L. Broglia
932 o_pts= (G4Point2d&)old_pts->get(i+o_ptr->GetOffset(),a);
933*/
934 o_pts= old_pts->Get3D(i+o_ptr->GetOffset(),a);
935 tempc.setX(tempc.x() + o_scale * o_pts.x());
936 tempc.setY(tempc.y() + o_scale * o_pts.y());
937 }
938
939 c_ptr->put((j)%upper-lower,a,tempc);
940 }
941 }
942 }
943}
944
945
946void G4BezierSurface::SplitNURBSurface()
947{
948 // Divides the surface in two parts. Uses the oslo-algorithm to calculate
949 // the new knotvectors and controlpoints for the subsurfaces.
950
951 // G4cout << "\nBezier splitted.";
952
953 register G4double value;
954 register G4int i;
955 register G4int k_index=0;
956 G4BezierSurface *srf1, *srf2;
957 G4int nr,nc;
958
959 if ( dir == ROW )
960 {
961 value = u_knots->GetKnot((u_knots->GetSize()-1)/2);
962
963 for( i = 0; i < u_knots->GetSize(); i++)
964 if( value == u_knots->GetKnot(i) )
965 {
966 k_index = i;
967 break;
968 }
969
970 if ( k_index == 0)
971 {
972 value = ( value + u_knots->GetKnot(u_knots->GetSize() -1))/2.0;
973 k_index = GetOrder(ROW);
974 }
975
976 new_knots = u_knots->MultiplyKnotVector(GetOrder(ROW), value);
977
978 ord = GetOrder(ROW);
979 CalcOsloMatrix();
980
981 srf1 = new G4BezierSurface(*this);
982 // srf1->dir=ROW;
983 srf1->dir=COL;
984
985 new_knots->ExtractKnotVector(srf1->u_knots, k_index +
986 srf1->GetOrder(ROW),0);
987
988 nr= srf1->v_knots->GetSize() - srf1->GetOrder(COL);
989 nc= srf1->u_knots->GetSize() - srf1->GetOrder(ROW);
990 delete srf1->ctl_points;
991
992 srf1->ctl_points= new G4ControlPoints(2, nr, nc);
993 srf2 = new G4BezierSurface(*this);
994
995 // srf2->dir = ROW;
996 srf2->dir = COL;
997
998 new_knots->ExtractKnotVector(srf2->u_knots,
999 new_knots->GetSize(), k_index);
1000
1001 nr= srf2->v_knots->GetSize() - srf2->GetOrder(COL);
1002 nc= srf2->u_knots->GetSize() - srf2->GetOrder(ROW);
1003
1004 delete srf2->ctl_points;
1005 srf2->ctl_points = new G4ControlPoints(2, nr, nc);
1006
1007 lower = 0;
1008 upper = k_index;
1009 MapSurface(srf1);
1010
1011 lower = k_index;
1012 upper = new_knots->GetSize() - srf2->GetOrder(ROW);
1013 MapSurface(srf2);
1014 }
1015 else // G4Vector3D = col
1016 {
1017 value = v_knots->GetKnot((v_knots->GetSize() -1)/2);
1018
1019 for( i = 0; i < v_knots->GetSize(); i++)
1020 if( value == v_knots->GetKnot(i))
1021 {
1022 k_index = i;
1023 break;
1024 }
1025 if ( k_index == 0)
1026 {
1027 value = ( value + v_knots->GetKnot(v_knots->GetSize() -1))/2.0;
1028 k_index = GetOrder(COL);
1029 }
1030
1031 new_knots = v_knots->MultiplyKnotVector( GetOrder(COL), value );
1032 ord = GetOrder(COL);
1033
1034 CalcOsloMatrix();
1035
1036 srf1 = new G4BezierSurface(*this);
1037 // srf1->dir = COL;
1038 srf1->dir = ROW;
1039
1040 new_knots->ExtractKnotVector(srf1->v_knots,
1041 k_index + srf1->GetOrder(COL), 0);
1042
1043 nr = srf1->v_knots->GetSize() - srf1->GetOrder(COL);
1044 nc = srf1->u_knots->GetSize() - srf1->GetOrder(ROW);
1045
1046 delete srf1->ctl_points;
1047 srf1->ctl_points = new G4ControlPoints(2, nr, nc);
1048
1049 srf2 = new G4BezierSurface(*this);
1050 // srf2->dir = COL;
1051 srf2->dir = ROW;
1052
1053 new_knots->ExtractKnotVector(srf2->v_knots, new_knots->GetSize(), k_index);
1054
1055 nr = srf2->v_knots->GetSize() - srf2->GetOrder(COL);
1056 nc = srf2->u_knots->GetSize() - srf2->GetOrder(ROW);
1057
1058 delete srf2->ctl_points;
1059 srf2->ctl_points = new G4ControlPoints(2,nr, nc);
1060
1061 lower = 0;
1062 upper = k_index;
1063 MapSurface(srf1);
1064
1065 // next->oslo_m = oslo_m;
1066 lower = k_index;
1067 upper = new_knots->GetSize() - srf2->GetOrder(COL);
1068 MapSurface(srf2);
1069 }
1070
1071 bezier_list->AddSurface(srf1);
1072 bezier_list->AddSurface(srf2);
1073 delete new_knots;
1074
1075 // Testing
1076 Splits++;
1077}
HepGeom::Point3D< G4double > G4Point3D
Definition: G4Point3D.hh:35
const G4Point3D PINFINITY(kInfinity, kInfinity, kInfinity)
#define COL
Definition: G4PointRat.hh:52
#define ROW
Definition: G4PointRat.hh:51
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
HepGeom::Vector3D< G4double > G4Vector3D
Definition: G4Vector3D.hh:35
G4DLLIMPORT std::ostream G4cout
static G4double Tolerance
static G4int Splits
virtual G4Vector3D SurfaceNormal(const G4Point3D &Pt) const
G4SurfaceList * bezier_list
virtual ~G4BezierSurface()
static G4int Clips
G4int GetOrder(G4int direction) const
G4int BIntersect(G4SurfaceList &)
G4int GetTestResult() const
G4int GetCols() const
G4int GetRows() const
G4Point3D Get3D(G4int i, G4int j) const
void put(G4int i, G4int j, const G4Point3D &tmp)
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
G4int GetSize() const
void ExtractKnotVector(G4KnotVector *kv, G4int upper, G4int lower)
G4double GetKnot(G4int knot_number) const
G4KnotVector * MultiplyKnotVector(G4int num, G4double value)
void PutKnot(G4int knot_number, G4double value)
G4OsloMatrix * GetNextNode()
void SetNextNode(G4OsloMatrix *)
G4KnotVector * GetKnotVector()
void AddSurface(G4Surface *srf)
void RemoveSurface(G4Surface *srf)
G4double distance
Definition: G4Surface.hh:203
G4BoundingBox3D * bbox
Definition: G4Surface.hh:185
G4double kCarTolerance
Definition: G4Surface.hh:192