source: trunk/source/geometry/solids/BREPS/src/G4BezierSurface.cc@ 1330

Last change on this file since 1330 was 1228, checked in by garnier, 16 years ago

update geant4.9.3 tag

File size: 26.8 KB
Line 
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: G4BezierSurface.cc,v 1.10 2008/03/13 14:18:57 gcosmo Exp $
28// GEANT4 tag $Name: geant4-09-03 $
29//
30// ----------------------------------------------------------------------
31// GEANT 4 class source file
32//
33// G4BezierSurface.cc
34//
35// ----------------------------------------------------------------------
36// History:
37// -------
38// - Replaced addition of coordinates by addition of 2 points
39// (L. Broglia, 10/10/98)
40// ----------------------------------------------------------------------
41
42#include "G4BezierSurface.hh"
43#include "G4ConvexHull.hh"
44
45G4double G4BezierSurface::Tolerance=0;
46G4int G4BezierSurface::Clips=0;
47G4int G4BezierSurface::Splits=0;
48
49
50G4BezierSurface::G4BezierSurface()
51{
52 oslo_m = (G4OsloMatrix*)0;
53 new_knots = (G4KnotVector*)0;
54 old_points = (G4ControlPoints*)0;
55
56 u[0]=0; u[1]=0;
57 v[0]=0; v[1]=0;
58}
59
60G4BezierSurface::~G4BezierSurface()
61{
62 delete u_knots;
63 delete v_knots;
64 delete new_knots;
65 delete ctl_points;
66 delete old_points;
67
68 G4OsloMatrix* temp_oslo = oslo_m;
69
70 while(oslo_m != (G4OsloMatrix*)0)
71 {
72 oslo_m = oslo_m->GetNextNode();
73 delete temp_oslo;
74 temp_oslo = oslo_m;
75 }
76
77 delete oslo_m;
78 delete bbox;
79}
80
81G4BezierSurface::G4BezierSurface(const G4BezierSurface&)
82 : G4Surface()
83{
84}
85
86G4Vector3D G4BezierSurface::SurfaceNormal(const G4Point3D&) const
87{
88 return G4Vector3D(0,0,0);
89}
90
91G4int G4BezierSurface::ClipBothDirs()
92{
93 dir = ROW;
94 ClipSurface();
95
96 // G4cout << "\n CLIP BOTH DIRS 1: " << smin << " " << smax;
97
98 if(smin > 1.0 || smax < 0.0)
99 {
100 bezier_list->RemoveSurface(this);
101 return 1;
102 }
103 else
104 if((smax - smin) > 0.8)
105 {
106 SplitNURBSurface();
107 return 0;
108 }
109
110 LocalizeClipValues();
111 SetValues();
112
113 // Other G4Vector3D clipping and testing.
114 dir = COL;
115 ClipSurface();
116 // G4cout << "\n CLIP BOTH DIRS 2: " << smin << " " << smax;
117
118 if(smin > 1.0 || smax < 0.0)
119 {
120 bezier_list->RemoveSurface(this);
121 return 1;
122 }
123 else
124 if((smax - smin) > 0.8)
125 {
126 SplitNURBSurface();
127 return 0;
128 }
129
130 LocalizeClipValues();
131 SetValues();
132 CalcAverage();
133 return 1;
134}
135
136
137void G4BezierSurface::CalcBBox()
138{
139 // Finds the bounds of the 2D-projected nurb iow
140 // calculates the bounds for a bounding rectangle
141 // to the surface. The bounding rectangle is used
142 // for a preliminary check of intersection.
143 G4Point3D box_min = G4Point3D(PINFINITY);
144 G4Point3D box_max = G4Point3D(-PINFINITY);
145
146
147 // Loop to search the whole control point mesh
148 // for the minimum and maximum values for.X() and y.
149 for(register G4int a = ctl_points->GetRows()-1; a>=0;a--)
150 for(register G4int b = ctl_points->GetCols()-1; b>=0;b--)
151 {
152/* L. Broglia
153 G4Point2d& tmp = (G4Point2d&)ctl_points->get(a,b);
154 if((box_min.X()) > (tmp.X())) box_min.X(tmp.X());
155 if((box_max.X()) < (tmp.X())) box_max.X(tmp.X());
156 if((box_min.Y()) > (tmp.Y())) box_min.Y(tmp.Y());
157 if((box_max.Y()) < (tmp.Y())) box_max.Y(tmp.Y());
158*/
159 G4Point3D tmp = ctl_points->Get3D(a,b);
160 if((box_min.x()) > (tmp.x())) box_min.setX(tmp.x());
161 if((box_max.x()) < (tmp.x())) box_max.setX(tmp.x());
162 if((box_min.y()) > (tmp.y())) box_min.setY(tmp.y());
163 if((box_max.y()) < (tmp.y())) box_max.setY(tmp.y());
164 }
165
166 bbox = new G4BoundingBox3D(box_min, box_max);
167}
168
169
170void G4BezierSurface::CalcAverage()
171{
172 // Calculate the average point from the average clip-values.
173 average_u = (u_min + u_max)/2.0;
174 average_v = (v_min + v_max)/2.0;
175}
176
177
178void G4BezierSurface::CalcDistance(const G4Point3D& ray_start)
179{
180 // Calculate the distance between the average point and
181 // the ray starting point.
182 distance = ((((ray_start.x() - average_pt.x())*
183 (ray_start.x() - average_pt.x()))+
184 ((ray_start.y() - average_pt.y())*
185 (ray_start.y() - average_pt.y()))+
186 ((ray_start.z() - average_pt.z())*
187 (ray_start.z() - average_pt.z()))));
188}
189
190
191void G4BezierSurface::SetValues()
192{
193 if(dir)
194 {
195 v_min = smin;
196 v_max = smax;
197 }
198 else
199 {
200 u_min = smin;
201 u_max = smax;
202 }
203}
204
205
206G4int G4BezierSurface::BIntersect(G4SurfaceList& bez_list)
207{
208 bezier_list = &bez_list;
209 G4int clip_regions = 0; // Used for tolerance/efficiency-testing
210
211 do
212 {
213 // Calc bbox
214 CalcBBox();
215
216 // Test bbox
217/* L. Broglia
218 bbox->Test2dBBox();
219*/
220 // bbox->Test();
221
222 // Check result
223 if(!bbox->GetTestResult())
224 return 0;
225
226 // The first clipping has already been Done
227 // previously so we continue by doing the
228 // actual clip.
229
230 // Cut out the clipped region of the surface
231 GetClippedRegionFromSurface();
232 clip_regions++;
233
234 // Calculate the knot vectors and control points
235 // for the clipped surface
236 RefineSurface();
237
238 // Gets the u- and v-bounds for the clipped surface
239 u_min = u_knots->GetKnot(0);
240 u_max = u_knots->GetKnot(u_knots->GetSize() - 1);
241 v_min = v_knots->GetKnot(0);
242 v_max = v_knots->GetKnot(v_knots->GetSize() - 1);
243
244 // Choose the G4Vector3D for the next() clipping so that
245 // the larger side will be clipped.
246 if( (u_max - u_min) < (v_max - v_min) )
247 dir = 1;
248 else
249 dir = 0;
250
251 // Calculate the clip points
252 ClipSurface();
253 // G4cout << "\n SMINMAX : " << smin << " " << smax;
254
255 // The ray intersects with the bounding box
256 // but not with the surface itself.
257 if( smin > 1.0 || smax < 0.0 )
258 {
259 // G4cout << "\nG4BezierSurface::Intersect : bezier missed!";
260 // bezier_list->RemoveSurface(this);
261 return 0;
262 }
263
264 if( (smax - smin) > 0.8)
265 {
266 // Multiple intersections
267 // G4cout << "\nG4BezierSurface::Intersect : Bezier split.";
268 SplitNURBSurface();
269 // Now the two new surfaces should also be
270 // clipped in both G4Vector3Ds i.e the
271 // last and the second last surface
272 // in the List. This is Done after returning
273 // from this function.
274 // G4cout << "\n\n BEZ SPLIT in final Calc! \n\n";
275
276
277 return 2;
278 }
279
280 // Calculate the smin and smax values on the
281 // b_spline.
282 LocalizeClipValues();
283
284 // Check if the size of the remaining surface is within the
285 // Tolerance .
286 } while ((u_max - u_min > Tolerance) || (v_max - v_min) > Tolerance);
287
288 SetValues();
289 // G4cout << "\nG4BezierSurface::Intersect :Regions were cut "
290 // << clip_regions << " Times.\n";
291
292 return 1;
293}
294
295
296void G4BezierSurface::ClipSurface()
297{
298 // This routine is described in Computer Graphics, Volume 24,
299 // Number 4, August 1990 under the title Ray Tracing Trimmed
300 // Rational Surface Patches.
301
302
303 // G4cout << "\nBezier clip.";
304
305 register G4int i,j;
306 register G4ConvexHull *ch_ptr=0, *ch_tmp=0, *ch_first=0;
307 register G4int col_size = ctl_points->GetCols();
308 register G4int row_size = ctl_points->GetRows();
309
310 // The four cornerpoints of the controlpoint mesh.
311
312/* L. Broglia
313 G4Point2d pt1 = ctl_points->get(0,0);
314 G4Point2d pt2 = ctl_points->get(0,col_size-1);
315 G4Point2d pt3 = ctl_points->get(row_size-1,0);
316 G4Point2d pt4 = ctl_points->get(row_size-1,col_size-1);
317 G4Point2d v1,v2,v3;
318*/
319 G4Point3D pt1 = ctl_points->Get3D(0,0);
320 G4Point3D pt2 = ctl_points->Get3D(0,col_size-1);
321 G4Point3D pt3 = ctl_points->Get3D(row_size-1,0);
322 G4Point3D pt4 = ctl_points->Get3D(row_size-1,col_size-1);
323 G4Point3D v1,v2,v3;
324
325 if ( dir == ROW)
326 {
327 // Vectors from cornerpoints
328 v1 = (pt1 - pt3);
329 // v1.X() = pt1.X() - pt3.X();
330 // v1.Y() = pt1.Y() - pt3.Y();
331 v2 = (pt2 - pt4);
332 // v2.X() = pt2.X() - pt4.X();
333 // v2.Y() = pt2.Y() - pt4.Y();
334 }
335 else
336 {
337 v1 = pt1 - pt2;
338 v2 = pt3 - pt4;
339 // v1.X() = pt1.X() - pt2.X();
340 // v1.Y() = pt1.Y() - pt2.Y();
341 // v2.X() = pt3.X() - pt4.X();
342 // v2.Y() = pt3.Y() - pt4.Y();
343 }
344/* L. Broglia
345 v3.X(v1.X() + v2.X());
346 v3.Y(v1.Y() + v1.Y());
347*/
348 v3 = v1 + v2 ;
349
350 smin = 1.0e8;
351 smax = -1.0e8;
352
353 G4double norm = std::sqrt(v3.x() * v3.x() + v3.y() * v3.y());
354 if(!norm)
355 {
356 G4cout << "\nNormal zero!";
357 G4cout << "\nLINE & DIR: " << line.x() << " " << line.y() << " " << dir;
358 G4cout << "\n";
359
360 if((std::abs(line.x())) > kCarTolerance)
361 line.setX(-line.x());
362 else
363 if((std::abs(line.y())) > kCarTolerance)
364 line.setY(-line.y());
365 else
366 {
367 G4cout << "\n RETURNING FROm CLIP..";
368 smin = 0; smax = 1;
369 return;
370 }
371
372 G4cout << "\nCHANGED LINE & DIR: " << line.x() << " "
373 << line.y() << " " << dir;
374 }
375 else
376 {
377 line.setX( v3.y() / norm);
378 line.setY(-v3.x() / norm);
379 }
380
381 // smin = 1.0e8;
382 // smax = -1.0e8;
383 // G4cout << "\n FINAL LINE & DIR: " << line.X() << " "
384 // << line.Y() << " " << dir;
385
386 if( dir == ROW)
387 {
388 // Create a Convex() hull List
389 for(G4int a = 0; a < col_size; a++)
390 {
391 ch_ptr = new G4ConvexHull(a/(col_size - 1.0),1.0e8,-1.0e8);
392 if(! a)
393 {
394 ch_first=ch_ptr;ch_tmp=ch_ptr;
395 }
396 else ch_tmp->SetNextHull(ch_ptr);
397
398 ch_tmp=ch_ptr;
399 }
400
401 ch_ptr=ch_first;
402 register G4double value;
403
404 // Loops through the control point mesh and calculates
405 // the nvex() hull for the surface.
406
407 for( G4int h = 0; h < row_size; h++)
408 {
409 for(G4int k = 0; k < col_size; k++)
410 {
411/* L. Broglia
412 G4Point2d& coordstmp = (G4Point2d&)ctl_points->get(h,k);
413 value = - ((coordstmp.X() * line.X() + coordstmp.Y() * line.Y()));
414*/
415 G4Point3D coordstmp = ctl_points->Get3D(h,k);
416 value = - ((coordstmp.x() * line.x() + coordstmp.y() * line.y()));
417
418 if( value <= (ch_ptr->GetMin()+kCarTolerance)) ch_ptr->SetMin(value);
419 if( value >= (ch_ptr->GetMax()-kCarTolerance)) ch_ptr->SetMax(value);
420
421 ch_ptr=ch_ptr->GetNextHull();
422 }
423
424 ch_ptr=ch_first;
425 }
426
427 ch_ptr=ch_first;
428 // Finds the points where the nvex() hull intersects
429 // with the coordinate .X()is. These points are the
430 // minimum and maximum values to where to clip the
431 // surface.
432
433 for(G4int l = 0; l < col_size - 1; l++)
434 {
435 ch_tmp=ch_ptr->GetNextHull();
436 for(G4int m = l+1; m < col_size; m++)
437 {
438 register G4double d;
439 register G4double param1, param2;
440 param1 = ch_ptr->GetParam();
441 param2 = ch_tmp->GetParam();
442
443 if(ch_tmp->GetMax() - ch_ptr->GetMax())
444 {
445 d = Findzero( param1, param2, ch_ptr->GetMax(), ch_tmp->GetMax());
446 if( d <= (smin + kCarTolerance) ) smin = d * .99;
447 if( d >= (smax - kCarTolerance) ) smax = d * .99 + .01;
448 }
449
450 if(ch_tmp->GetMin() - ch_ptr->GetMin())
451 {
452 d = Findzero( param1, param2, ch_ptr->GetMin(), ch_tmp->GetMin());
453 if( d <= (smin + kCarTolerance)) smin = d * .99;
454 if( d >= (smax - kCarTolerance)) smax = d * .99 + .01;
455 }
456
457 ch_tmp=ch_tmp->GetNextHull();
458 }
459
460 ch_ptr=ch_ptr->GetNextHull();
461 }
462
463 ch_ptr=ch_first;
464
465 if (smin <= 0.0) smin = 0.0;
466 if (smax >= 1.0) smax = 1.0;
467
468 if ( Sign(ch_ptr->GetMin()) != Sign(ch_ptr->GetMax())) smin = 0.0;
469
470 i = Sign(ch_tmp->GetMin()); // ch_tmp points to last nvex()_hull in List
471 j = Sign(ch_tmp->GetMax());
472
473 if ( std::abs(i-j) > kCarTolerance ) smax = 1.0;
474 // if ( i != j) smax = 1.0;
475
476 }
477 else // Other G4Vector3D
478 {
479 for(G4int n = 0; n < row_size; n++)
480 {
481 ch_ptr = new G4ConvexHull(n/(row_size - 1.0),1.0e8,-1.0e8);
482 if(!n)
483 {
484 ch_first=ch_ptr;
485 ch_tmp=ch_ptr;
486 }
487 else ch_tmp->SetNextHull(ch_ptr);
488
489 ch_tmp=ch_ptr;
490 }
491
492 ch_ptr=ch_first;
493
494 for( G4int o = 0; o < col_size; o++)
495 {
496 for(G4int p = 0; p < row_size; p++)
497 {
498 register G4double value;
499
500/* L. Broglia
501 G4Point2d& coordstmp =(G4Point2d&) ctl_points->get(p,o);
502 value = - ((coordstmp.X() * line.X() + coordstmp.Y() * line.Y()));
503*/
504 G4Point3D coordstmp = ctl_points->Get3D(p,o);
505 value = - ((coordstmp.x() * line.x() + coordstmp.y() * line.y()));
506
507 if( value <= (ch_ptr->GetMin()+kCarTolerance)) ch_ptr->SetMin(value);
508 if( value >= (ch_ptr->GetMax()-kCarTolerance)) ch_ptr->SetMax(value);
509
510 ch_ptr=ch_ptr->GetNextHull();
511 }
512
513 ch_ptr=ch_first;
514 }
515
516 ch_ptr=ch_first;
517 ch_tmp=ch_first;
518
519 for(G4int q = 0; q < row_size - 1; q++)
520 {
521 ch_tmp=ch_ptr->GetNextHull();
522 for(G4int r = q+1; r < row_size; r++)
523 {
524 register G4double param1 = ch_ptr->GetParam();
525 register G4double param2 = ch_tmp->GetParam();
526 register G4double d;
527
528 if(ch_tmp->GetMax() - ch_ptr->GetMax())
529 {
530 d = Findzero( param1, param2, ch_ptr->GetMax(), ch_tmp->GetMax());
531 if( d <= (smin + kCarTolerance) ) smin = d * .99;
532 if( d >= (smax - kCarTolerance) ) smax = d * .99 + .01;
533 }
534
535 if(ch_tmp->GetMin()-ch_ptr->GetMin())
536 {
537 d = Findzero( param1, param2, ch_ptr->GetMin(), ch_tmp->GetMin());
538 if( d <= (smin + kCarTolerance) ) smin = d * .99;
539 if( d >= (smax - kCarTolerance) ) smax = d * .99 + .01;
540 }
541
542 ch_tmp=ch_tmp->GetNextHull();
543 }
544
545 ch_ptr=ch_ptr->GetNextHull();
546 }
547
548 ch_tmp=ch_ptr;
549 ch_ptr=ch_first;
550
551 if (smin <= 0.0) smin = 0.0;
552 if (smax >= 1.0) smax = 1.0;
553
554 if ( Sign(ch_ptr->GetMin()) != Sign(ch_ptr->GetMax())) smin = 0.0;
555
556 i = Sign(ch_tmp->GetMin()); // ch_tmp points to last nvex()_hull in List
557 j = Sign(ch_tmp->GetMax());
558
559 //
560 if ( (std::abs(i-j) > kCarTolerance)) smax = 1.0;
561 }
562
563 ch_ptr=ch_first;
564 while(ch_ptr!=ch_ptr->GetNextHull())
565 {
566 ch_tmp=ch_ptr;
567 ch_ptr=ch_ptr->GetNextHull();
568 delete ch_tmp;
569 }
570
571 delete ch_ptr;
572
573 // Testing...
574 Clips++;
575}
576
577
578void G4BezierSurface::GetClippedRegionFromSurface()
579{
580 // Returns the clipped part of the surface. First calculates the
581 // length of the new knotvector. Then uses the refinement function to
582 // get the new knotvector and controlmesh.
583
584 // G4cout << "\nBezier region clipped.";
585
586 delete new_knots;
587 if ( dir == ROW)
588 {
589 new_knots = new G4KnotVector(GetOrder(0) * 2);
590 for (register G4int i = 0; i < GetOrder(0); i++)
591 {
592 new_knots->PutKnot(i, smin);
593 new_knots->PutKnot(i+ GetOrder(0), smax);
594 }
595 }
596 else
597 {
598 new_knots = new G4KnotVector( GetOrder(1) * 2);
599 for ( register G4int i = 0; i < GetOrder(1); i++)
600 {
601 new_knots->PutKnot(i, smin);
602 new_knots->PutKnot(i+ GetOrder(1), smax);
603 }
604 }
605} // NURB_REGION_FROM_SURFACE
606
607
608void G4BezierSurface::RefineSurface()
609{
610 // Returns the new clipped surface. Calculates the new controlmesh
611 // and knotvectorvalues for the surface by using the Oslo-algorithm
612
613 delete old_points;
614 if (dir == ROW)
615 {
616 // Row (u) G4Vector3D
617 ord = GetOrder(0);
618 CalcOsloMatrix();
619 for(register G4int a=0;a<new_knots->GetSize();a++)
620 u_knots->PutKnot(a, new_knots->GetKnot(a));
621
622 lower = 0;
623 upper = new_knots->GetSize() - GetOrder(0);
624
625 // Copy of the old points.
626 old_points = new G4ControlPoints(*ctl_points);
627 MapSurface(this);
628 }
629 else
630 {
631 ord = GetOrder(1);
632 CalcOsloMatrix ();
633 for(register G4int a=0;a < new_knots->GetSize();a++)
634 v_knots->PutKnot(a, new_knots->GetKnot(a));
635
636 // Copy of the old points.
637 old_points = new G4ControlPoints(*ctl_points);
638
639 // Make new controlpoint matrix,
640 register G4int cols = ctl_points->GetCols();
641 delete ctl_points;
642
643 ctl_points = new G4ControlPoints(2,(new_knots->GetSize()-
644 GetOrder(1)),cols);
645 lower = 0;
646 upper = new_knots->GetSize() - GetOrder(1);
647 MapSurface(this);
648 }
649}// REFINE_SURFACE
650
651
652void G4BezierSurface::CalcOsloMatrix()
653{
654 // This algorithm is described in the paper "Making the Oslo-algorithm
655 // more efficient" in SIAM J.NUMER.ANAL. Vol.23, No. 3, June '86
656 // Calculates the oslo-matrix , which is used in mapping the new
657 // knotvector- and controlpoint-values.
658
659 register G4KnotVector *ah;
660 register G4KnotVector *newknots;
661 register G4int i;
662 register G4int j;
663 register G4int mu, muprim;
664 register G4int vv, p;
665 register G4int iu, il, ih, n1;
666 register G4int ahi;
667 register G4double beta1;
668 register G4double tj;
669
670 ah = new G4KnotVector(ord*(ord + 1)/2);
671 newknots = new G4KnotVector(ord * 2 );
672
673 n1 = new_knots->GetSize() - ord;
674 mu = 0;
675
676 if(oslo_m!=(G4OsloMatrix*)0)
677 {
678 G4OsloMatrix* tmp;
679
680 // while(oslo_m!=oslo_m->next)
681 while(oslo_m!=(G4OsloMatrix*)0)
682 {
683 tmp=oslo_m->GetNextNode();delete oslo_m; oslo_m=tmp;
684 }
685 }
686
687 delete oslo_m;
688 oslo_m = new G4OsloMatrix();
689
690 register G4OsloMatrix* o_ptr = oslo_m;
691
692 register G4KnotVector* old_knots;
693 if(dir)
694 old_knots = v_knots;
695 else
696 old_knots = u_knots;
697
698 for (j = 0; j < n1; j++)
699 {
700 if ( j != 0 )
701 {
702 oslo_m->SetNextNode(new G4OsloMatrix());
703 oslo_m = oslo_m->GetNextNode();
704 }
705
706 while (old_knots->GetKnot(mu + 1) <= new_knots->GetKnot(j))
707 mu = mu + 1; // find the bounding mu
708
709 i = j + 1;
710 muprim = mu;
711
712 while ((new_knots->GetKnot(i) == old_knots->GetKnot(muprim)) &&
713 i < (j + ord))
714 {
715 i++;
716 muprim--;
717 }
718
719 ih = muprim + 1;
720
721 for (vv = 0, p = 1; p < ord; p++)
722 {
723 if (new_knots->GetKnot(j + p) == old_knots->GetKnot(ih))
724 ih++;
725 else
726 newknots->PutKnot(++vv - 1,new_knots->GetKnot(j + p));
727 }
728
729 ahi = AhIndex(0, ord - 1,ord);
730 ah->PutKnot(ahi, 1.0);
731
732 for (p = 1; p <= vv; p++)
733 {
734 beta1 = 0.0;
735 tj = newknots->GetKnot(p-1);
736
737 if (p - 1 >= muprim)
738 {
739 beta1 = AhIndex(p - 1, ord - muprim,ord);
740 beta1 = ((tj - old_knots->GetKnot(0)) * beta1) /
741 (old_knots->GetKnot(p + ord - vv) - old_knots->GetKnot(0));
742 }
743
744 i = muprim - p + 1;
745 il = Amax (1, i);
746 i = n1 - 1 + vv - p;
747 iu = Amin (muprim, i);
748
749 for (i = il; i <= iu; i++)
750 {
751 register G4double d1, d2;
752 register G4double beta;
753
754 d1 = tj - old_knots->GetKnot(i);
755 d2 = old_knots->GetKnot(i + p + ord - vv - 1) - tj;
756
757 beta = ah->GetKnot(AhIndex(p - 1, i + ord - muprim - 1,ord)) /
758 (d1 + d2);
759
760
761 ah->PutKnot(AhIndex(p, i + ord - muprim - 2,ord), d2 * beta + beta1) ;
762 beta1 = d1 * beta;
763 }
764
765 ah->PutKnot(AhIndex(p, iu + ord - muprim - 1,ord), beta1);
766
767 if (iu < muprim)
768 {
769 register G4double kkk;
770 register G4double ahv;
771
772 kkk = old_knots->GetKnot(n1 - 1 + ord);
773 ahv = AhIndex (p - 1, iu + ord - muprim,ord);
774 ah->PutKnot(AhIndex(p, iu + ord - muprim - 1,ord),
775 beta1 + (kkk - tj) * ahv /
776 (kkk - old_knots->GetKnot(iu + 1)));
777 }
778 }
779
780 // Remove the oslo matrix List
781 G4OsloMatrix* temp_oslo = oslo_m;
782
783/*
784 if(oslo_m != (G4OsloMatrix*)0)
785 while(oslo_m->next != oslo_m)
786 {
787 oslo_m = oslo_m->next;
788 delete temp_oslo;
789 temp_oslo = oslo_m;
790 }
791
792 // Remove the last
793 delete oslo_m;
794*/
795
796 while(oslo_m != (G4OsloMatrix*)0)
797 {
798 oslo_m = oslo_m->GetNextNode();
799 delete temp_oslo;
800 temp_oslo = oslo_m;
801 }
802
803 delete oslo_m;
804
805 // Create a new oslo matrix
806 oslo_m = new G4OsloMatrix(vv+1, Amax(muprim - vv,0), vv);
807
808 for ( i = vv, p = 0; i >= 0; i--)
809 oslo_m->GetKnotVector()
810 ->PutKnot ( p++, ah->GetKnot(AhIndex (vv, (ord-1) - i,ord)));
811
812 }
813
814 delete ah;
815 delete newknots;
816 oslo_m->SetNextNode(0);
817 oslo_m = o_ptr;
818}
819
820
821void G4BezierSurface::MapSurface(G4Surface*)
822{
823 // This algorithm is described in the paper Making the Oslo-algorithm
824 // more efficient in SIAM J.NUMER.ANAL. Vol.23, No. 3, June '86
825 // Maps the new controlpoints into the new surface.
826
827 register G4ControlPoints *c_ptr;
828 register G4OsloMatrix *o_ptr;
829 register G4ControlPoints* new_pts;
830 register G4ControlPoints* old_pts;
831
832 new_pts = ctl_points;
833
834 // Copy the old points so they can be used in calculating the new ones.
835 // old_pts = new G4ControlPoints(*ctl_points);
836 old_pts = old_points;
837 register G4int j, // j loop
838 i; // oslo loop
839
840 c_ptr = new_pts;
841 register G4int size; // The number of rows or columns,
842 // depending on processing order
843
844 if(!dir)
845 size=new_pts->GetRows();
846 else
847 size=new_pts->GetCols();
848
849 for(G4int a=0; a<size;a++)
850 {
851 if ( lower != 0)
852 {
853 for ( i = 0, o_ptr = oslo_m;
854 i < lower;
855 i++, o_ptr = o_ptr->GetNextNode()){;}
856 }
857 else
858 {
859 o_ptr = oslo_m;
860 }
861
862 if(!dir)// Direction ROW
863 {
864 for ( j = lower; j < upper; j++, o_ptr = o_ptr->GetNextNode())
865 {
866 register G4double o_scale;
867 register G4int x;
868 x=a;
869
870/* L. Broglia
871 G4Point2d o_pts= (G4Point2d&)old_pts->Get2d(x, o_ptr->GetOffset());
872 G4Point2d tempc= (G4Point2d&)c_ptr->Get2d(j/upper,(j)%upper-lower);
873*/
874 G4Point3D o_pts = old_pts->Get3D(x, o_ptr->GetOffset());
875 G4Point3D tempc = c_ptr->Get3D(j/upper, (j)%upper-lower);
876
877 o_scale = o_ptr->GetKnotVector()->GetKnot(0);
878
879 tempc.setX(o_pts.x() * o_scale);
880 tempc.setY(o_pts.x() * o_scale);
881
882 for ( i = 1; i <= o_ptr->GetSize(); i++)
883 {
884 o_scale = o_ptr->GetKnotVector()->GetKnot(i);
885
886/* L. Broglia
887 o_pts = (G4Point2d&)old_pts->get(x, i+o_ptr->GetOffset());
888 tempc.X(tempc.X() + o_scale * o_pts.X());
889 tempc.Y(tempc.Y() + o_scale * o_pts.Y());
890*/
891 o_pts = old_pts->Get3D(x, i+o_ptr->GetOffset());
892 tempc.setX(tempc.x() + o_scale * o_pts.x());
893 tempc.setY(tempc.y() + o_scale * o_pts.y());
894
895 }
896
897 c_ptr->put(a,(j)%upper-lower,tempc);
898 }
899 }
900 else // dir = COL
901 {
902 for ( j = lower; j < upper; j++, o_ptr = o_ptr->GetNextNode())
903 {
904 register G4double o_scale;
905 register G4int x;
906 x=a;
907
908/* L. Broglia
909 G4Point2d o_pts= (G4Point2d&)old_pts->Get2d(o_ptr->GetOffset(), x);
910 G4Point2d tempc = (G4Point2d&)c_ptr->Get2d((j)%upper-lower,j/upper);
911*/
912 G4Point3D o_pts = old_pts->Get3D(o_ptr->GetOffset(), x);
913 G4Point3D tempc = c_ptr->Get3D((j)%upper-lower,j/upper);
914
915 o_scale = o_ptr->GetKnotVector()->GetKnot(0);
916
917 tempc.setX(o_pts.x() * o_scale);
918 tempc.setY(o_pts.y() * o_scale);
919
920 for ( i = 1; i <= o_ptr->GetSize(); i++)
921 {
922 o_scale = o_ptr->GetKnotVector()->GetKnot(i);
923/* L. Broglia
924 o_pts= (G4Point2d&)old_pts->get(i+o_ptr->GetOffset(),a);
925*/
926 o_pts= old_pts->Get3D(i+o_ptr->GetOffset(),a);
927 tempc.setX(tempc.x() + o_scale * o_pts.x());
928 tempc.setY(tempc.y() + o_scale * o_pts.y());
929 }
930
931 c_ptr->put((j)%upper-lower,a,tempc);
932 }
933 }
934 }
935}
936
937
938void G4BezierSurface::SplitNURBSurface()
939{
940 // Divides the surface in two parts. Uses the oslo-algorithm to calculate
941 // the new knotvectors and controlpoints for the subsurfaces.
942
943 // G4cout << "\nBezier splitted.";
944
945 register G4double value;
946 register G4int i;
947 register G4int k_index=0;
948 G4BezierSurface *srf1, *srf2;
949 G4int nr,nc;
950
951 if ( dir == ROW )
952 {
953 value = u_knots->GetKnot((u_knots->GetSize()-1)/2);
954
955 for( i = 0; i < u_knots->GetSize(); i++)
956 if( value == u_knots->GetKnot(i) )
957 {
958 k_index = i;
959 break;
960 }
961
962 if ( k_index == 0)
963 {
964 value = ( value + u_knots->GetKnot(u_knots->GetSize() -1))/2.0;
965 k_index = GetOrder(ROW);
966 }
967
968 new_knots = u_knots->MultiplyKnotVector(GetOrder(ROW), value);
969
970 ord = GetOrder(ROW);
971 CalcOsloMatrix();
972
973 srf1 = new G4BezierSurface(*this);
974 // srf1->dir=ROW;
975 srf1->dir=COL;
976
977 new_knots->ExtractKnotVector(srf1->u_knots, k_index +
978 srf1->GetOrder(ROW),0);
979
980 nr= srf1->v_knots->GetSize() - srf1->GetOrder(COL);
981 nc= srf1->u_knots->GetSize() - srf1->GetOrder(ROW);
982 delete srf1->ctl_points;
983
984 srf1->ctl_points= new G4ControlPoints(2, nr, nc);
985 srf2 = new G4BezierSurface(*this);
986
987 // srf2->dir = ROW;
988 srf2->dir = COL;
989
990 new_knots->ExtractKnotVector(srf2->u_knots,
991 new_knots->GetSize(), k_index);
992
993 nr= srf2->v_knots->GetSize() - srf2->GetOrder(COL);
994 nc= srf2->u_knots->GetSize() - srf2->GetOrder(ROW);
995
996 delete srf2->ctl_points;
997 srf2->ctl_points = new G4ControlPoints(2, nr, nc);
998
999 lower = 0;
1000 upper = k_index;
1001 MapSurface(srf1);
1002
1003 lower = k_index;
1004 upper = new_knots->GetSize() - srf2->GetOrder(ROW);
1005 MapSurface(srf2);
1006 }
1007 else // G4Vector3D = col
1008 {
1009 value = v_knots->GetKnot((v_knots->GetSize() -1)/2);
1010
1011 for( i = 0; i < v_knots->GetSize(); i++)
1012 if( value == v_knots->GetKnot(i))
1013 {
1014 k_index = i;
1015 break;
1016 }
1017 if ( k_index == 0)
1018 {
1019 value = ( value + v_knots->GetKnot(v_knots->GetSize() -1))/2.0;
1020 k_index = GetOrder(COL);
1021 }
1022
1023 new_knots = v_knots->MultiplyKnotVector( GetOrder(COL), value );
1024 ord = GetOrder(COL);
1025
1026 CalcOsloMatrix();
1027
1028 srf1 = new G4BezierSurface(*this);
1029 // srf1->dir = COL;
1030 srf1->dir = ROW;
1031
1032 new_knots->ExtractKnotVector(srf1->v_knots,
1033 k_index + srf1->GetOrder(COL), 0);
1034
1035 nr = srf1->v_knots->GetSize() - srf1->GetOrder(COL);
1036 nc = srf1->u_knots->GetSize() - srf1->GetOrder(ROW);
1037
1038 delete srf1->ctl_points;
1039 srf1->ctl_points = new G4ControlPoints(2, nr, nc);
1040
1041 srf2 = new G4BezierSurface(*this);
1042 // srf2->dir = COL;
1043 srf2->dir = ROW;
1044
1045 new_knots->ExtractKnotVector(srf2->v_knots, new_knots->GetSize(), k_index);
1046
1047 nr = srf2->v_knots->GetSize() - srf2->GetOrder(COL);
1048 nc = srf2->u_knots->GetSize() - srf2->GetOrder(ROW);
1049
1050 delete srf2->ctl_points;
1051 srf2->ctl_points = new G4ControlPoints(2,nr, nc);
1052
1053 lower = 0;
1054 upper = k_index;
1055 MapSurface(srf1);
1056
1057 // next->oslo_m = oslo_m;
1058 lower = k_index;
1059 upper = new_knots->GetSize() - srf2->GetOrder(COL);
1060 MapSurface(srf2);
1061 }
1062
1063 bezier_list->AddSurface(srf1);
1064 bezier_list->AddSurface(srf2);
1065 delete new_knots;
1066
1067 // Testing
1068 Splits++;
1069}
Note: See TracBrowser for help on using the repository browser.