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