source: trunk/source/geometry/solids/BREPS/src/G4ProjectedSurface.cc@ 1224

Last change on this file since 1224 was 1058, checked in by garnier, 17 years ago

file release beta

File size: 19.4 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: G4ProjectedSurface.cc,v 1.12 2008/03/13 14:18:57 gcosmo Exp $
28// GEANT4 tag $Name: geant4-09-02-ref-02 $
29//
30// ----------------------------------------------------------------------
31// GEANT 4 class source file
32//
33// G4ProjectedSurface.cc
34//
35// ----------------------------------------------------------------------
36
37#include "G4ProjectedSurface.hh"
38
39G4int G4ProjectedSurface::Splits=0;
40
41G4ProjectedSurface::G4ProjectedSurface()
42{
43 distance = 0;
44 oslo_m =(G4OsloMatrix*)0;
45}
46
47
48G4ProjectedSurface::~G4ProjectedSurface()
49{
50 delete u_knots;
51 delete v_knots;
52 delete ctl_points;
53
54 G4OsloMatrix* temp_oslo;
55 if(oslo_m!=(G4OsloMatrix*)0)
56 {
57 while(oslo_m->GetNextNode() != oslo_m)
58 {
59 temp_oslo = oslo_m;
60 oslo_m = oslo_m->GetNextNode();
61
62 delete temp_oslo;
63 }
64
65 delete oslo_m;
66 }
67
68 delete bbox;
69}
70
71G4ProjectedSurface::G4ProjectedSurface(const G4ProjectedSurface&)
72 : G4Surface()
73{
74}
75
76void G4ProjectedSurface::CopySurface()
77 // Copies the projected surface into a bezier surface
78 // and adds it to the List of bezier surfaces.
79{
80 G4BezierSurface *bez = new G4BezierSurface();
81
82 bez->SetDistance(distance);
83 bez->PutOrder(0, order[0]);
84 bez->PutOrder(1, order[1]);
85 bez->Dir(dir);
86 bez->u_knots = new G4KnotVector(*u_knots);
87 bez->v_knots = new G4KnotVector(*v_knots);
88 bez->ctl_points = new G4ControlPoints(*ctl_points);
89
90 bezier_list->AddSurface(bez);
91}
92
93
94void G4ProjectedSurface::CalcBBox()
95{
96 // Finds the bounds of the 2D-projected nurb iow
97 // calculates the bounds for a bounding rectangle
98 // to the surface. The bounding rectangle is used
99 // for a preliminary check of intersection.
100
101 // Loop to search the whole control point mesh
102 // for the minimum and maximum values for x and y.
103 G4double box_minx,box_miny,box_maxx,box_maxy;
104 box_minx = kInfinity;
105 box_miny = kInfinity;
106 box_maxx = -kInfinity;
107 box_maxy = -kInfinity;
108
109 G4double bminx,bminy,bmaxx,bmaxy,tmpx,tmpy;
110 bminx = box_minx; bminy = box_miny;
111 bmaxx = box_maxx; bmaxy = box_maxy;
112
113 for(register G4int a = ctl_points->GetRows()-1; a>=0;a--)
114 for(register G4int b = ctl_points->GetCols()-1; b>=0;b--)
115 {
116/* L. Broglia
117 G4Point2d& tmp = (G4Point2d&)ctl_points->get(a,b);
118*/
119 G4Point3D tmp = ctl_points->Get3D(a,b);
120
121 tmpx = tmp.x(); tmpy = tmp.y();
122 if(bminx > tmpx) box_minx=tmpx;
123 if(bmaxx < tmpx) box_maxx=tmpx;
124 if(bminy > tmpy) box_miny=tmpy;
125 if(bmaxy < tmpy) box_maxy=tmpy;
126 }
127
128 G4Point3D box_min(box_minx,box_miny,0.);
129 G4Point3D box_max(box_maxx,box_maxy,0.);
130
131 delete bbox;
132 bbox = new G4BoundingBox3D(box_min, box_max);
133}
134
135
136void G4ProjectedSurface::ConvertToBezier(G4SurfaceList& proj_list,
137 G4SurfaceList& bez_list)
138{
139 projected_list = &proj_list;
140 bezier_list = &bez_list;
141
142 // Check wether the surface is a bezier surface by checking
143 // if internal knots exist.
144 if(CheckBezier())
145 {
146 // Make it a G4BezierSurface -object and add it to the bezier
147 // surface List
148 CopySurface();
149
150 // Retrieve a pointer to the newly added surface iow the
151 // last in the List
152 G4BezierSurface* bez_ptr = (G4BezierSurface*)bezier_list->GetLastSurface();
153
154 // Do the first clip to the bezier.
155 bez_ptr->ClipSurface();
156 G4double dMin = bez_ptr->SMin();
157 G4double dMax = bez_ptr->SMax();
158 G4double dMaxMinusdMin = dMax - dMin;
159
160 if(( dMaxMinusdMin > kCarTolerance ))
161 {
162 if( dMaxMinusdMin > 0.8 )
163 {
164 // The clipping routine selected a larger Area than one
165 // knot interval which indicates that we have a case of
166 // multiple intersections. The projected surface has to
167 // be split again in order to separate the intersections
168 // to different surfaces.
169
170 // Check tolerance of clipping
171 // G4cout << "\nClip Area too big -> Split";
172 dir = bez_ptr->dir;
173 bezier_list->RemoveSurface(bez_ptr);
174
175 SplitNURBSurface();
176 return;
177 //}
178 }
179 else
180 if( dMin > 0.0 || dMax < 0.0 )
181 {
182 // The ray intersects with the bounding box
183 // but not with the surface itself.
184 // G4cout << "\nConvex hull missed.";
185 bezier_list->RemoveSurface(bez_ptr);
186 return;
187 }
188 }
189 else
190 if(dMaxMinusdMin < kCarTolerance && dMaxMinusdMin > -kCarTolerance)
191 {
192 bezier_list->RemoveSurface(bez_ptr);
193 return;
194 }
195
196 bez_ptr->LocalizeClipValues();
197 bez_ptr->SetValues();
198
199 // Other G4ThreeVec clipping and testing.
200 bez_ptr->ChangeDir();//bez->dir = !bez_ptr->dir;
201 bez_ptr->ClipSurface();
202 // G4cout<<"\nSMIN: " << bez_ptr->smin << " SMAX: "
203 // << bez_ptr->smax << " DIR: " << bez_ptr->dir;
204
205 dMin = bez_ptr->SMin();
206 dMax = bez_ptr->SMax();
207 dMaxMinusdMin = dMax-dMin;
208
209 if((dMaxMinusdMin > kCarTolerance ))// ||
210 // (dMaxMinusdMin < -kCarTolerance))
211 {
212 if( (dMaxMinusdMin) > 0.8 )
213 {
214 // G4cout << "\nClip Area too big -> Split";
215 dir = bez_ptr->dir;//1.2 klo 18.30
216
217 // dir=!dir;
218 bezier_list->RemoveSurface(bez_ptr);
219 SplitNURBSurface();
220 return;
221 //}
222 }
223 else
224 if( dMin > 1.0 || dMax < 0.0 )
225 {
226 // G4cout << "\nConvex hull missed.";
227 bezier_list->RemoveSurface(bez_ptr);
228 return;
229 }
230 }
231 else
232 if(dMaxMinusdMin < kCarTolerance && dMaxMinusdMin > -kCarTolerance)
233 {
234 bezier_list->RemoveSurface(bez_ptr);
235 return;
236 }
237
238 bez_ptr->LocalizeClipValues();
239 bez_ptr->SetValues();
240 bez_ptr->CalcAverage();
241 }
242 else
243 {
244 // Split the surface into two new surfaces. The G4ThreeVec
245 // is set in the CheckBezier function.
246 // G4cout << "\nNot a bezier surface -> Split";
247 SplitNURBSurface();
248 }
249}
250
251
252G4int G4ProjectedSurface::CheckBezier()
253{
254 // Checks if the surface is a bezier surface by
255 // checking wether internal knots exist. If no internal
256 // knots exist the quantity of knots is 2*order of the
257 // surface. Returns 1 if the surface
258 // is a bezier.
259
260 if( u_knots->GetSize() > (2.0 * GetOrder(ROW)))
261 {dir=0;return 0;}
262
263 if( v_knots->GetSize() > (2.0 * GetOrder(COL)))
264 {dir=1;return 0;}
265
266 return 1;
267}
268
269
270void G4ProjectedSurface::SplitNURBSurface()
271{
272 // Divides the surface in two parts. Uses the oslo-algorithm to calculate
273 // the new knotvectors and controlpoints for the subsurfaces.
274
275 // G4cout << "\nProjected splitted.";
276
277 register G4double value;
278 register G4int i;
279 register G4int k_index=0;
280 register G4ProjectedSurface *srf1, *srf2;
281 register G4int nr,nc;
282
283 if ( dir == ROW )
284 {
285 value = u_knots->GetKnot((u_knots->GetSize()-1)/2);
286
287 for( i = 0; i < u_knots->GetSize(); i++)
288 if( (std::abs(value - u_knots->GetKnot(i))) < kCarTolerance )
289 {
290 k_index = i;
291 break;
292 }
293
294 if ( k_index == 0)
295 {
296 value = ( value + u_knots->GetKnot(u_knots->GetSize() -1))/2.0;
297 k_index = GetOrder(ROW);
298 }
299
300 new_knots = u_knots->MultiplyKnotVector(GetOrder(ROW), value);
301
302 ord = GetOrder(ROW);
303 CalcOsloMatrix();
304
305 srf1 = new G4ProjectedSurface(*this);
306
307 //srf1->dir=ROW;
308 srf1->dir=COL;
309
310 new_knots->ExtractKnotVector(srf1->u_knots,
311 k_index + srf1->GetOrder(ROW),0);
312
313 nr= srf1->v_knots->GetSize() - srf1->GetOrder(COL);
314 nc= srf1->u_knots->GetSize() - srf1->GetOrder(ROW);
315
316 delete srf1->ctl_points;
317 srf1->ctl_points= new G4ControlPoints(2, nr, nc);
318
319 srf2 = new G4ProjectedSurface(*this);
320
321 //srf2->dir = ROW;
322 srf2->dir = COL;
323
324 new_knots->ExtractKnotVector(srf2->u_knots,
325 new_knots->GetSize(), k_index);
326
327 nr= srf2->v_knots->GetSize() - srf2->GetOrder(COL);
328 nc= srf2->u_knots->GetSize() - srf2->GetOrder(ROW);
329
330 delete srf2->ctl_points;
331 srf2->ctl_points = new G4ControlPoints(2, nr, nc);
332
333 lower = 0;
334 upper = k_index;
335 MapSurface(srf1);
336
337 lower = k_index;
338 upper = new_knots->GetSize() - srf2->GetOrder(ROW);
339 MapSurface(srf2);
340 }
341 else // G4ThreeVec = col
342 {
343 value = v_knots->GetKnot((v_knots->GetSize() -1)/2);
344
345 for( i = 0; i < v_knots->GetSize(); i++)
346 if( (std::abs(value - v_knots->GetKnot(i))) < kCarTolerance )
347 {
348 k_index = i;
349 break;
350 }
351
352 if ( k_index == 0)
353 {
354 value = ( value + v_knots->GetKnot(v_knots->GetSize() -1))/2.0;
355 k_index = GetOrder(COL);
356 }
357
358 new_knots = v_knots->MultiplyKnotVector( GetOrder(COL), value );
359 ord = GetOrder(COL);
360
361 CalcOsloMatrix();
362
363 srf1 = new G4ProjectedSurface(*this);
364
365 //srf1->dir = COL;
366 srf1->dir = ROW;
367
368 new_knots->ExtractKnotVector(srf1->v_knots,
369 k_index + srf1->GetOrder(COL), 0);
370
371 nr = srf1->v_knots->GetSize() - srf1->GetOrder(COL);
372 nc = srf1->u_knots->GetSize() - srf1->GetOrder(ROW);
373
374 delete srf1->ctl_points;
375 srf1->ctl_points = new G4ControlPoints(2, nr, nc);
376
377 srf2 = new G4ProjectedSurface(*this);
378
379 //srf2->dir = COL;
380 srf2->dir = ROW;
381
382 new_knots->ExtractKnotVector(srf2->v_knots, new_knots->GetSize(), k_index);
383
384 nr = srf2->v_knots->GetSize() - srf2->GetOrder(COL);
385 nc = srf2->u_knots->GetSize() - srf2->GetOrder(ROW);
386
387 delete srf2->ctl_points;
388 srf2->ctl_points = new G4ControlPoints(2,nr, nc);
389
390 lower = 0;
391 upper = k_index;
392 MapSurface(srf1);
393
394 lower = k_index;
395 upper = new_knots->GetSize() - srf2->GetOrder(COL);
396 MapSurface(srf2);
397 }
398
399 // Check that surfaces are ok.
400 G4int col_size = srf1->ctl_points->GetCols();
401 G4int row_size = srf1->ctl_points->GetRows();
402
403/* L. Broglia
404 // get three cornerpoints of the controlpoint mesh.
405 G4Point2d pt1 = srf1->ctl_points->get(0,0);
406 G4Point2d pt2 = srf1->ctl_points->get(0,col_size-1);
407 G4Point2d pt3 = srf1->ctl_points->get(row_size-1,0);
408
409 // Calc distance between points
410 G4double pointDist1 = pt1.Distance(pt2);
411 G4double pointDist2 = pt1.Distance(pt3);
412*/
413
414 // get three cornerpoints of the controlpoint mesh.
415 G4Point3D pt1 = srf1->ctl_points->Get3D(0,0);
416 G4Point3D pt2 = srf1->ctl_points->Get3D(0,col_size-1);
417 G4Point3D pt3 = srf1->ctl_points->Get3D(row_size-1,0);
418
419 // Calc distance squared between points
420 G4double pointDist1 = pt1.distance2(pt2);
421 G4double pointDist2 = pt1.distance2(pt3);
422
423
424 // Add surfaces to List of projected surfaces
425 if(pointDist1 > kCarTolerance && pointDist2 > kCarTolerance)
426 projected_list->AddSurface(srf1);
427 else
428 delete srf1;
429
430 col_size = srf2->ctl_points->GetCols();
431 row_size = srf2->ctl_points->GetRows();
432
433/* L. Broglia
434 // get three cornerpoints of the controlpoint mesh.
435 pt1 = srf2->ctl_points->get(0,0);
436 pt2 = srf2->ctl_points->get(0,col_size-1);
437 pt3 = srf2->ctl_points->get(row_size-1,0);
438
439 // Calc distance between points
440 pointDist1 = pt1.Distance(pt2);
441 pointDist2 = pt1.Distance(pt3);
442*/
443
444 // get three cornerpoints of the controlpoint mesh.
445 pt1 = srf2->ctl_points->Get3D(0,0);
446 pt2 = srf2->ctl_points->Get3D(0,col_size-1);
447 pt3 = srf2->ctl_points->Get3D(row_size-1,0);
448
449 // Calc distance squared between points
450 pointDist1 = pt1.distance2(pt2);
451 pointDist2 = pt1.distance2(pt3);
452
453 // Add surfaces to List of projected surfaces
454 if(pointDist1 > kCarTolerance && pointDist2 > kCarTolerance)
455 projected_list->AddSurface(srf2);
456 else
457 delete srf2;
458
459 delete new_knots;
460
461 Splits++;
462}
463
464void G4ProjectedSurface::CalcOsloMatrix()
465{
466 // This algorithm is described in the paper "Making the Oslo-algorithm
467 // more efficient" in SIAM J.NUMER.ANAL. Vol.23, No. 3, June '86
468 // Calculates the oslo-matrix , which is used in mapping the new
469 // knotvector- and controlpoint-values.
470
471 register G4KnotVector *ah;
472 static G4KnotVector *newknots;
473 register G4int i;
474 register G4int j;
475 register G4int mu, muprim;
476 register G4int v, p;
477 register G4int iu, il, ih, n1;
478 register G4int ahi;
479 register G4double beta1;
480 register G4double tj;
481
482 ah = new G4KnotVector(ord*(ord + 1)/2);
483
484 newknots = new G4KnotVector(ord * 2 );
485
486 n1 = new_knots->GetSize() - ord;
487 mu = 0;
488
489 if(oslo_m!=(G4OsloMatrix*)0)
490 {
491 G4OsloMatrix* tmp;
492 while(oslo_m!=oslo_m->GetNextNode())
493 {
494 tmp=oslo_m->GetNextNode();
495 delete oslo_m;
496 oslo_m=tmp;
497 }
498 }
499
500 delete oslo_m;
501
502 oslo_m = new G4OsloMatrix();
503 register G4OsloMatrix* o_ptr = oslo_m;
504
505 register G4KnotVector* old_knots;
506
507 if(dir)
508 old_knots = v_knots;
509 else
510 old_knots = u_knots;
511
512 for (j = 0; j < n1; j++)
513 {
514 if ( j != 0 )
515 {
516 oslo_m->SetNextNode(new G4OsloMatrix());
517 oslo_m = oslo_m->GetNextNode();
518 }
519
520 //while (old_knots->GetKnot(mu + 1) <= new_knots->GetKnot(j))
521 while ( (new_knots->GetKnot(j) - old_knots->GetKnot(mu + 1)) >
522 kCarTolerance )
523 mu = mu + 1; // find the bounding mu
524
525 i = j + 1;
526 muprim = mu;
527
528 while ( ((std::abs(new_knots->GetKnot(i) - old_knots->GetKnot(muprim))) <
529 kCarTolerance) && i < (j + ord) )
530 {
531 i++;
532 muprim--;
533 }
534
535 ih = muprim + 1;
536
537 for (v = 0, p = 1; p < ord; p++)
538 {
539 // if (new_knots->GetKnot(j + p) == old_knots->GetKnot(ih))
540 if ( (std::abs((new_knots->GetKnot(j + p)) - (old_knots->GetKnot(ih)))) <
541 kCarTolerance )
542 ih++;
543 else
544 newknots->PutKnot(++v - 1,new_knots->GetKnot(j + p));
545 }
546
547 ahi = AhIndex(0, ord - 1,ord);
548 ah->PutKnot(ahi, 1.0);
549
550 for (p = 1; p <= v; p++)
551 {
552 beta1 = 0.0;
553 tj = newknots->GetKnot(p-1);
554
555 if (p - 1 >= muprim)
556 {
557 beta1 = AhIndex(p - 1, ord - muprim,ord);
558 beta1 = ((tj - old_knots->GetKnot(0)) * beta1) /
559 (old_knots->GetKnot(p + ord - v) - old_knots->GetKnot(0));
560 }
561
562 i = muprim - p + 1;
563 il = Amax (1, i);
564 i = n1 - 1 + v - p;
565 iu = Amin (muprim, i);
566
567 for (i = il; i <= iu; i++)
568 {
569 register G4double d1, d2;
570 register G4double beta;
571
572 d1 = tj - old_knots->GetKnot(i);
573 d2 = old_knots->GetKnot(i + p + ord - v - 1) - tj;
574
575 beta = ah->GetKnot(AhIndex(p - 1, i + ord - muprim - 1,ord)) /
576 (d1 + d2);
577
578 ah->PutKnot(AhIndex(p, i + ord - muprim - 2,ord), d2 * beta + beta1) ;
579 beta1 = d1 * beta;
580 }
581
582 ah->PutKnot(AhIndex(p, iu + ord - muprim - 1,ord), beta1);
583
584 if (iu < muprim)
585 {
586 register G4double kkk;
587 register G4double ahv;
588
589 kkk = old_knots->GetKnot(n1 - 1 + ord);
590 ahv = AhIndex (p - 1, iu + ord - muprim,ord);
591 ah->PutKnot(AhIndex(p, iu + ord - muprim - 1,ord),
592 beta1 + (kkk - tj) * ahv /
593 (kkk - old_knots->GetKnot(iu + 1)));
594 }
595 }
596
597 delete oslo_m->GetKnotVector();
598 oslo_m->SetKnotVector(new G4KnotVector(v+1));
599 oslo_m->SetOffset(Amax(muprim - v, 0));
600 oslo_m->SetSize(v);
601
602 for ( i = v, p = 0; i >= 0; i--)
603 oslo_m->GetKnotVector()
604 ->PutKnot( p++, ah->GetKnot(AhIndex (v, (ord-1) - i,ord)) );
605
606 }
607
608 delete ah;
609 delete newknots;
610 oslo_m->SetNextNode(oslo_m);
611 oslo_m = o_ptr;
612}
613
614void G4ProjectedSurface::MapSurface(G4ProjectedSurface* srf)
615{
616 // This algorithm is described in the paper "Making the Oslo-algorithm
617 // more efficient" in SIAM J.NUMER.ANAL. Vol.23, No. 3, June '86
618 // Maps the new controlpoints into the new surface.
619
620 register G4ControlPoints *c_ptr;
621 register G4OsloMatrix *o_ptr;
622 register G4ControlPoints* new_pts;
623 register G4ControlPoints* old_pts;
624
625 new_pts = srf->ctl_points;
626
627 // Copy the old points so they can be used in calculating the new ones.
628 // In this version, where the splitted surfaces are given
629 // as parameters the copying is not necessary.
630
631 old_pts = new G4ControlPoints(*ctl_points);
632 register G4int j, // j loop
633 i; // oslo loop
634 c_ptr = new_pts;
635
636 register G4int size; // The number of rows or columns,
637 // depending on processing order
638
639 if(!dir)
640 size=new_pts->GetRows();
641 else
642 size=new_pts->GetCols();
643
644 for( register G4int a=0; a<size;a++)
645 {
646 if ( lower != 0)
647 for ( i = 0, o_ptr = oslo_m; i < lower; i++, o_ptr = o_ptr->GetNextNode()){;}
648 else
649 o_ptr = oslo_m;
650
651 if(!dir)// Direction ROW
652 {
653 for ( j = lower; j < upper; j++, o_ptr = o_ptr->GetNextNode())
654 {
655 register G4double o_scale;
656 register G4int x;
657 x=a;
658
659/* L. Broglia
660 G4Point2d o_pts = (G4Point2d&)old_pts->get(x,o_ptr->GetOffset());
661 G4Point2d tempc = (G4Point2d&)c_ptr->get(j/upper,(j)%upper-lower);
662*/
663 G4Point3D o_pts = old_pts->Get3D(x, o_ptr->GetOffset());
664 G4Point3D tempc = c_ptr->Get3D(j/upper, (j)%upper-lower);
665 o_scale = o_ptr->GetKnotVector()->GetKnot(0);
666
667 tempc.setX(o_pts.x() * o_scale);
668 tempc.setY(o_pts.y() * o_scale);
669
670 for ( i = 1; i <= o_ptr->GetSize(); i++)
671 {
672 o_scale = o_ptr->GetKnotVector()->GetKnot(i);
673
674/* L. Broglia
675 o_pts = (G4Point2d&)old_pts->get(x,i+o_ptr->GetOffset());
676 tempc.X(tempc.X() + o_scale * o_pts.X());
677 tempc.Y(tempc.Y() + o_scale * o_pts.Y());
678*/
679
680 o_pts = old_pts->Get3D(x,i+o_ptr->GetOffset());
681 tempc.setX(tempc.x() + o_scale * o_pts.x());
682 tempc.setY(tempc.y() + o_scale * o_pts.y());
683 }
684
685 c_ptr->put(a,(j)%upper-lower,tempc);
686 }
687 }
688 else // dir = COL
689 {
690 for ( j = lower; j < upper; j++, o_ptr = o_ptr->GetNextNode())
691 {
692 register G4double o_scale;
693 register G4int x;
694 x=a;
695
696/* L.Broglia
697 G4Point2d o_pts = (G4Point2d&)old_pts->get(o_ptr->GetOffset(),x);
698 G4Point2d tempc = (G4Point2d&)c_ptr->get((j)%upper-lower,j/upper);
699*/
700 G4Point3D o_pts = old_pts->Get3D(o_ptr->GetOffset(),x);
701 G4Point3D tempc = c_ptr->Get3D((j)%upper-lower, j/upper);
702
703 o_scale = o_ptr->GetKnotVector()->GetKnot(0);
704
705 tempc.setX(o_pts.x() * o_scale);
706 tempc.setY(o_pts.y() * o_scale);
707
708 for ( i = 1; i <= o_ptr->GetSize(); i++)
709 {
710 o_scale = o_ptr->GetKnotVector()->GetKnot(i);
711 o_pts= old_pts->Get3D(i+o_ptr->GetOffset(),a);
712
713 tempc.setX(tempc.x() + o_scale * o_pts.x());
714 tempc.setY(tempc.y() + o_scale * o_pts.y());
715 }
716
717 c_ptr->put((j)%upper-lower,a,tempc);
718 }
719 }
720 }
721
722 delete old_pts;
723}
Note: See TracBrowser for help on using the repository browser.