updated to modern VTK
[engrid-github.git] / src / libengrid / geometrytools.cpp
blobce3ae2a30ed9a1838111e1689ad0c13efce542b4
1 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 // + +
3 // + This file is part of enGrid. +
4 // + +
5 // + Copyright 2008-2014 enGits GmbH +
6 // + +
7 // + enGrid is free software: you can redistribute it and/or modify +
8 // + it under the terms of the GNU General Public License as published by +
9 // + the Free Software Foundation, either version 3 of the License, or +
10 // + (at your option) any later version. +
11 // + +
12 // + enGrid is distributed in the hope that it will be useful, +
13 // + but WITHOUT ANY WARRANTY; without even the implied warranty of +
14 // + MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the +
15 // + GNU General Public License for more details. +
16 // + +
17 // + You should have received a copy of the GNU General Public License +
18 // + along with enGrid. If not, see <http://www.gnu.org/licenses/>. +
19 // + +
20 // ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
21 #include "geometrytools.h"
22 #include "containertricks.h"
23 #include "engrid.h"
24 #include "utilities.h"
25 #include "math/linsolve.h"
27 #include <vtkCellType.h>
28 #include <cmath>
30 namespace GeometryTools {
32 double rad2deg( double rad )
34 return rad/M_PI*180;
37 double deg2rad( double deg )
39 return deg/180*M_PI;
42 void rotate(vec3_t g1, vec3_t g2,vec3_t g3, vec3_t &b, double theta)
44 //cout << b << endl;
45 mat3_t g_t;
46 g_t[0] = g1;
47 g_t[1] = g2;
48 g_t[2] = g3;
49 mat3_t g = g_t.transp();
50 //cout << g << endl;
51 vec3_t bb = g_t*b;
52 //cout << bb << endl;
53 mat3_t rot = mat3_t::identity();
54 rot[0][0] = cos(theta);
55 rot[0][1] = -sin(theta);
56 rot[1][0] = sin(theta);
57 rot[1][1] = cos(theta);
58 //cout << rot << endl;
59 vec3_t bbb = rot*bb;
60 //cout << bbb << endl;
61 b = g*bbb;
62 //cout << b << endl;
65 vec3_t rotate(vec3_t v, vec3_t axis, double theta)
67 axis.normalise();
69 // transposed base of rotate system
70 mat3_t g_t;
72 // compute projection of v in axis direction
73 vec3_t v_axis = (axis*v)*axis;
75 // compute first orthogonal vector (first base vector)
76 g_t[0] = v-v_axis;
78 //In case of points on the rotation axis, do nothing
79 if(g_t[0].abs()==0) return v;
81 g_t[0].normalise();
83 // second base vector is the normalised axis
84 g_t[1] = axis;
86 // compute second orthogonal vector (third base vector)
87 g_t[2] = g_t[0].cross(g_t[1]);
89 // base of rotate system
90 mat3_t g = g_t.transp();
92 // matrix for rotation around g_t[1];
93 mat3_t rot = mat3_t::identity();
94 rot[0][0] = cos(theta);
95 rot[0][2] = sin(theta);
96 rot[2][0] = -sin(theta);
97 rot[2][2] = cos(theta);
99 // transfer v to rotate system
100 vec3_t v_r = g_t*v;
102 // rotate the vector and transfer it back
103 v_r = rot*v_r;
104 v = g*v_r;
106 return v;
109 vec3_t orthogonalVector(vec3_t v)
111 // get absolute values
112 double xx = v[0] < 0.0 ? -v[0] : v[0];
113 double yy = v[1] < 0.0 ? -v[1] : v[1];
114 double zz = v[2] < 0.0 ? -v[2] : v[2];
115 // switch both biggest values and set the other one to zero
116 vec3_t u;
117 if (xx < yy) {
118 u = xx < zz ? vec3_t(0,v[2],-v[1]) : vec3_t(v[1],-v[0],0);
119 } else {
120 u = yy < zz ? vec3_t(-v[2],0,v[0]) : vec3_t(v[1],-v[0],0);
122 u.normalise();
123 return u;
127 double intersection(vec3_t x_straight, vec3_t v_straight, vec3_t x_plane, vec3_t u_plane, vec3_t v_plane)
129 vec3_t n = u_plane.cross(v_plane);
130 return intersection(x_straight,v_straight,x_plane,n);
133 double intersection(vec3_t x_straight, vec3_t v_straight, vec3_t x_plane, vec3_t n_plane)
135 double k = (x_plane*n_plane - x_straight*n_plane)/(v_straight*n_plane);
136 return k;
139 bool intersection (double &k1, double &k2, vec2_t r1, vec2_t u1, vec2_t r2, vec2_t u2)
141 double ave_length = .5*(sqrt(u1[0]*u1[0]+u1[1]*u1[1]) + sqrt(u2[0]*u2[0]+u2[1]*u2[1]));
142 double DET = (u1[0]*u2[1]-u1[1]*u2[0]);
143 if (fabs(DET) > 1e-6*ave_length) {
144 k1 = -(u2[0]*r2[1]-u2[0]*r1[1]-r2[0]*u2[1]+r1[0]*u2[1])/DET;
145 k2 = -(-u1[1]*r2[0]+u1[0]*r2[1]-u1[0]*r1[1]+u1[1]*r1[0])/DET;
146 return true;
147 } else {
148 return false;
152 void sliceTriangle(const vector<vec3_t> &Tin, vec3_t x, vec3_t n, vector<vector<vec3_t> > &Tout)
154 vec3_t a = Tin[0];
155 vec3_t b = Tin[1];
156 vec3_t c = Tin[2];
157 double kab = intersection(a,b-a,x,n);
158 double kbc = intersection(b,c-b,x,n);
159 double kca = intersection(c,a-c,x,n);
160 bool ab_cut = ((kab >= 0) && (kab <= 1));
161 bool bc_cut = ((kbc >= 0) && (kbc <= 1));
162 bool ca_cut = ((kca >= 0) && (kca <= 1));
163 if (ab_cut && bc_cut && ca_cut) {
164 //cerr << "invalid triangle (SliceTriangle) A" << endl;
165 //exit(EXIT_FAILURE);
166 if ((kab <= kbc) && (kab <= kca)) ab_cut = false;
167 else if ((kbc <= kab) && (kbc <= kca)) bc_cut = false;
168 else ca_cut = false;
170 if (ab_cut && bc_cut) {
171 vec3_t ab = a + kab*(b-a);
172 vec3_t bc = b + kbc*(c-b);
173 Tout.resize(3,vector<vec3_t>(3));
174 clinit(Tout[0]) = a,ab,bc;
175 clinit(Tout[1]) = ab,b,bc;
176 clinit(Tout[2]) = bc,c,a;
177 } else if (bc_cut && ca_cut) {
178 vec3_t bc = b + kbc*(c-b);
179 vec3_t ca = c + kca*(a-c);
180 Tout.resize(3,vector<vec3_t>(3));
181 clinit(Tout[0]) = a,bc,ca;
182 clinit(Tout[1]) = a,b,bc;
183 clinit(Tout[2]) = bc,c,ca;
184 } else if (ca_cut && ab_cut) {
185 vec3_t ca = c + kca*(a-c);
186 vec3_t ab = a + kab*(b-a);
187 Tout.resize(3,vector<vec3_t>(3));
188 clinit(Tout[0]) = a,ab,ca;
189 clinit(Tout[1]) = ab,b,ca;
190 clinit(Tout[2]) = b,c,ca;
191 } else {
192 Tout.resize(1,vector<vec3_t>(3));
193 clinit(Tout[0]) = a,b,c;
197 double tetraVol(const vec3_t& x0, const vec3_t& x1, const vec3_t& x2, const vec3_t& x3, bool neg)
199 static double f16 = 1.0/6.0;
200 vec3_t v1(x1[0]-x0[0], x1[1]-x0[1], x1[2]-x0[2]);
201 vec3_t v2(x2[0]-x0[0], x2[1]-x0[1], x2[2]-x0[2]);
202 vec3_t v3(x3[0]-x0[0], x3[1]-x0[1], x3[2]-x0[2]);
203 double V = v1[0]*(v2[1]*v3[2]-v2[2]*v3[1]) + v1[1]*(v2[2]*v3[0]-v2[0]*v3[2]) + v1[2]*(v2[0]*v3[1]-v2[1]*v3[0]);
204 V *= f16;
205 if (!neg && (V < 0)) {
206 V = -1e99;
208 return V; //fabs(V);
211 double pyraVol(vec3_t x0, vec3_t x1, vec3_t x2, vec3_t x3, vec3_t x4, bool neg)
213 double V1 = tetraVol(x0, x1, x3, x4, neg) + tetraVol(x1, x2, x3, x4, neg);
214 double V2 = tetraVol(x0, x1, x2, x4, neg) + tetraVol(x2, x3, x0, x4, neg);
215 return min(V1,V2);
217 double V = 0;
218 vec3_t m0 = .25*(x0+x1+x2+x3);
219 V += tetraVol(x0, x1, m0, x4, neg);
220 V += tetraVol(x1, x2, m0, x4, neg);
221 V += tetraVol(x2, x3, m0, x4, neg);
222 V += tetraVol(x3, x0, m0, x4, neg);
223 return V;
227 double prismVol(vec3_t x0, vec3_t x1, vec3_t x2, vec3_t x3, vec3_t x4, vec3_t x5, bool neg)
229 double V = 0;
230 vec3_t p = 1.0/6.0*(x0+x1+x2+x3+x4+x5);
231 V += tetraVol(x0, x2, x1, p, neg);
232 V += tetraVol(x3, x4, x5, p, neg);
233 V += pyraVol (x0, x1, x4, x3, p, neg);
234 V += pyraVol (x1, x2, x5, x4, p, neg);
235 V += pyraVol (x0, x3, x5, x2, p, neg);
236 return V;
239 double hexaVol(vec3_t x0, vec3_t x1, vec3_t x2, vec3_t x3, vec3_t x4, vec3_t x5, vec3_t x6, vec3_t x7, bool neg)
241 double V = 0;
242 vec3_t p = 1.0/8.0*(x0+x1+x2+x3+x4+x5+x6+x7);
243 V += pyraVol(x0, x1, x3, x2, p, neg);
244 V += pyraVol(x0, x4, x5, x1, p, neg);
245 V += pyraVol(x4, x6, x7, x5, p, neg);
246 V += pyraVol(x2, x3, x7, x6, p, neg);
247 V += pyraVol(x1, x5, x7, x3, p, neg);
248 V += pyraVol(x0, x2, x6, x4, p, neg);
249 return V;
252 double triArea(vec3_t x0, vec3_t x1, vec3_t x2)
254 vec3_t a = x1-x0;
255 vec3_t b = x2-x0;
256 double A = 0.5*((a.cross(b)).abs());
257 return A;
260 double quadArea(vec3_t x0, vec3_t x1, vec3_t x2, vec3_t x3)
262 double A = 0;
263 vec3_t p = .25*(x0+x1+x2+x3);
264 A += triArea(x0,x1,p);
265 A += triArea(x1,x2,p);
266 A += triArea(x2,x3,p);
267 A += triArea(x3,x0,p);
268 return A;
271 vec3_t triNormal(vec3_t x0, vec3_t x1, vec3_t x2)
273 vec3_t a = x1-x0;
274 vec3_t b = x2-x0;
275 vec3_t n = 0.5*(a.cross(b));
276 return n;
279 vec3_t quadNormal(vec3_t x0, vec3_t x1, vec3_t x2, vec3_t x3)
281 vec3_t n;
282 clinit(n) = 0,0,0;
283 vec3_t p = .25*(x0+x1+x2+x3);
284 n += triNormal(x0,x1,p);
285 n += triNormal(x1,x2,p);
286 n += triNormal(x2,x3,p);
287 n += triNormal(x3,x0,p);
288 return n;
291 vec3_t triNormal(vtkUnstructuredGrid *grid, vtkIdType p1, vtkIdType p2, vtkIdType p3)
293 vec3_t x1, x2, x3;
294 grid->GetPoint(p1,x1.data());
295 grid->GetPoint(p2,x2.data());
296 grid->GetPoint(p3,x3.data());
297 return triNormal(x1,x2,x3);
300 vec3_t quadNormal(vtkUnstructuredGrid *grid, vtkIdType p1, vtkIdType p2, vtkIdType p3, vtkIdType p4)
302 vec3_t x1, x2, x3, x4;
303 grid->GetPoint(p1,x1.data());
304 grid->GetPoint(p2,x2.data());
305 grid->GetPoint(p3,x3.data());
306 grid->GetPoint(p4,x4.data());
307 return quadNormal(x1,x2,x3,x4);
310 vec3_t cellNormal(vtkUnstructuredGrid *grid, vtkIdType i)
312 vec3_t n(0,0,0);
313 EG_GET_CELL(i, grid);
314 if (num_pts < 3) {
315 EG_BUG;
316 } else if (num_pts == 3) {
317 return triNormal(grid,pts[0],pts[1],pts[2]);
318 } else if (num_pts == 4) {
319 return quadNormal(grid,pts[0],pts[1],pts[2],pts[3]);
320 } else {
321 QList<vec3_t> x;
322 for (int i = 0; i < num_pts; ++i) {
323 vec3_t xp;
324 grid->GetPoint(pts[i], xp.data());
325 x << xp;
327 x.append(x.first());
328 return polyNormal(x, true);
330 return n;
333 double cellVA(vtkUnstructuredGrid *grid, vtkIdType cellId, bool neg)
335 vec3_t p[8];
336 EG_GET_CELL(cellId, grid);
337 for (int i_pts = 0; i_pts < num_pts; ++i_pts) {
338 grid->GetPoints()->GetPoint(pts[i_pts], p[i_pts].data());
340 vtkIdType cellType = grid->GetCellType(cellId);
341 if (cellType == VTK_TRIANGLE) return triArea (p[0], p[1], p[2]);
342 else if (cellType == VTK_QUAD) return quadArea(p[0], p[1], p[2], p[3]);
343 else if (cellType == VTK_TETRA) return tetraVol(p[0], p[1], p[2], p[3], neg);
344 else if (cellType == VTK_PYRAMID) return pyraVol (p[0], p[1], p[2], p[3], p[4], neg);
345 else if (cellType == VTK_WEDGE) return prismVol(p[0], p[1], p[2], p[3], p[4], p[5], neg);
346 else if (cellType == VTK_HEXAHEDRON) return hexaVol (p[0], p[1], p[2], p[3], p[4], p[5], p[6], p[7], neg);
347 return 0.0;
350 double angle(const vec3_t & u, const vec3_t & v)
352 // return the angle w.r.t. another 3-vector
353 double ptot2 = u.abs2()*v.abs2();
354 if(ptot2 <= 0) {
355 return 0.0;
356 } else {
357 double arg = (u*v)/sqrt(ptot2);
358 if(arg > 1.0) arg = 1.0;
359 if(arg < -1.0) arg = -1.0;
360 return acos(arg);
364 double deviation(vtkUnstructuredGrid *grid, vtkIdType p1, vtkIdType p2, vtkIdType p3)
366 vec3_t x1, x2, x3;
367 grid->GetPoint(p1,x1.data());
368 grid->GetPoint(p2,x2.data());
369 grid->GetPoint(p3,x3.data());
370 vec3_t u=x2-x1;
371 vec3_t v=x3-x2;
372 return angle(u,v);
375 double angle(vtkUnstructuredGrid *grid, vtkIdType p1, vtkIdType p2, vtkIdType p3)
377 vec3_t x1, x2, x3;
378 grid->GetPoint(p1,x1.data());
379 grid->GetPoint(p2,x2.data());
380 grid->GetPoint(p3,x3.data());
381 vec3_t u=x1-x2;
382 vec3_t v=x3-x2;
383 return angle(u,v);
386 double cosAngle(vtkUnstructuredGrid *grid, vtkIdType cell1, vtkIdType cell2)
388 vec3_t u1 = cellNormal(grid, cell1);
389 vec3_t u2 = cellNormal(grid, cell2);
390 u1.normalise();
391 u2.normalise();
392 return(u1*u2);
395 vec3_t getCenter(vtkUnstructuredGrid *grid, vtkIdType cellId, double& Rmin, double& Rmax)
397 EG_GET_CELL(cellId, grid);
398 if (num_pts <= 0) {
399 cout<<"FATAL ERROR: num_pts <= 0"<<endl;
400 abort();
403 //calculate center position
404 vec3_t xc(0,0,0);
405 for (vtkIdType i = 0; i < num_pts; ++i) {
406 vec3_t xp;
407 grid->GetPoints()->GetPoint(pts[i], xp.data());
408 xc += xp;
410 xc = 1.0/(double)num_pts * xc;
412 //calculate Rmin+Rmax
413 vec3_t xp;
414 grid->GetPoints()->GetPoint(pts[0], xp.data());
415 Rmin = 0.25*(xp-xc).abs();
416 Rmax = 0.25*(xp-xc).abs();
417 for (vtkIdType i = 1; i < num_pts; ++i) {
418 grid->GetPoints()->GetPoint(pts[i], xp.data());
419 Rmin = min(Rmin, 0.25*(xp-xc).abs());
420 Rmax = max(Rmax, 0.25*(xp-xc).abs());
423 return(xc);
426 bool isInsideTriangle(vec2_t r, double tol)
428 if (r[0] < 0-tol || 1+tol < r[0] || r[1] < 0-tol || 1+tol < r[1] || r[0]+r[1] > 1+tol) {
429 return false;
431 return true;
434 bool intersectEdgeAndTriangle(const vec3_t& a, const vec3_t& b, const vec3_t& c,
435 const vec3_t& x1, const vec3_t& x2, vec3_t& xi, vec3_t& ri, double tol)
437 // triangle base
438 vec3_t g1 = b - a;
439 vec3_t g2 = c - a;
440 vec3_t g3 = g1.cross(g2);
441 g3.normalise();
443 // direction of the edge
444 vec3_t v = x2 - x1;
446 // parallel?
447 if (fabs(g3*v) < 1e-6) {
448 return false;
451 // compute intersection between straight and triangular plane
452 double k = intersection(x1, v, a, g3);
453 xi = x1 + k*v;
455 // transform xi to triangular base
456 mat3_t G;
457 G.column(0, g1);
458 G.column(1, g2);
459 G.column(2, g3);
461 mat3_t GI = G.inverse();
462 ri = xi - a;
463 ri = GI*ri;
466 vec3_t b = xi-a;
467 linsolve(G, b, ri);
471 // intersection outside of edge range?
472 if (k < 0 - tol) {
473 return false;
475 if (k > 1 + tol) {
476 return false;
479 // intersection outside of triangle?
480 if (!isInsideTriangle(vec2_t(ri[0],ri[1]),tol)) {
481 return false;
483 return true;
486 double distance(vtkUnstructuredGrid *grid, vtkIdType id_node1, vtkIdType id_node2) {
487 vec3_t A;
488 vec3_t B;
489 grid->GetPoints()->GetPoint(id_node1, A.data());
490 grid->GetPoints()->GetPoint(id_node2, B.data());
491 return((B-A).abs());
494 double distance2(vtkUnstructuredGrid *grid, vtkIdType id_node1, vtkIdType id_node2) {
495 vec3_t A;
496 vec3_t B;
497 grid->GetPoints()->GetPoint(id_node1, A.data());
498 grid->GetPoints()->GetPoint(id_node2, B.data());
499 return((B-A).abs2());
502 double areaOfCircumscribedCircle(vtkUnstructuredGrid *grid, vtkIdType id_cell) {
503 EG_GET_CELL(id_cell, grid);
504 vec3_t A,B,C;
505 grid->GetPoints()->GetPoint(pts[0], A.data());
506 grid->GetPoints()->GetPoint(pts[1], B.data());
507 grid->GetPoints()->GetPoint(pts[2], C.data());
508 double a=(C-B).abs();
509 double alpha=angle((B-A),(C-A));
510 double R=a/(2*sin(alpha));
511 return(M_PI*R*R);
514 void computeCircumscribedCircle(vec3_t a, vec3_t b, vec3_t c, vec3_t &x, double &radius)
516 double la = (b-c).abs();
517 double lb = (a-c).abs();
518 double lc = (a-b).abs();
519 double bca = sqr(la)*(sqr(lb) + sqr(lc) - sqr(la));
520 double bcb = sqr(lb)*(sqr(la) + sqr(lc) - sqr(lb));
521 double bcc = sqr(lc)*(sqr(la) + sqr(lb) - sqr(lc));
522 double sum = bca + bcb + bcc;
523 if (fabs(sum) < 1e-6) {
524 x = (1.0/3.0)*(a + b + c);
525 radius = 1e99;
526 } else {
527 x = bca*a + bcb*b + bcc*c;
528 x *= 1.0/sum;
529 radius = (x-a).abs();
533 vec3_t getBarycentricCoordinates(double x, double y)
535 if(isnan(x) || isinf(x) || isnan(y) || isinf(y)) {
536 qWarning()<<"x="<<x;
537 qWarning()<<"y="<<y;
538 EG_BUG;
541 double x_1=0;
542 double y_1=0;
543 double x_2=1;
544 double y_2=0;
545 double x_3=0;
546 double y_3=1;
548 mat2_t T;
549 T[0][0]=x_1-x_3; T[0][1]=x_2-x_3;
550 T[1][0]=y_1-y_3; T[1][1]=y_2-y_3;
552 if(T.det()==0) {
553 qWarning()<<"T.det()="<<T.det();
554 qWarning()<<T[0][0]<<T[0][1];
555 qWarning()<<T[1][0]<<T[1][1];
556 qWarning()<<"T[0][0]*T[1][1]-T[1][0]*T[0][1]="<<T[0][0]*T[1][1]-T[1][0]*T[0][1];
557 EG_BUG;
560 double lambda_1 = ((y_2-y_3)*(x-x_3)-(x_2-x_3)*(y-y_3))/(T.det());
561 double lambda_2 = (-(y_1-y_3)*(x-x_3)+(x_1-x_3)*(y-y_3))/(T.det());
562 double lambda_3 = 1-lambda_1-lambda_2;
564 vec3_t bary_coords(lambda_1,lambda_2,lambda_3);
565 return bary_coords;
567 // initialize
568 /* double t1=0;
569 double t2=0;
570 double t3=0;*/
572 /* if(x==0) {
573 t3=y;
574 t1=1-y;
575 t2=0;
577 else if(y==0) {
578 t2=x;
579 t1=1-x;
580 t3=0;
582 else if((x+y)==1) {
585 else {
588 double k1,k2;
589 if(!intersection (k1, k2, t_A, t_M-t_A, t_B, t_C-t_B)) EG_BUG;
590 vec2_t t_I1 = t_A+k1*(t_M-t_A);
591 vec3_t g_nI1 = (1-k2)*g_nB + k2*g_nC;
592 vec2_t pm1_M(1.0/k1,0);
594 // normalize
595 double total = t1+t2+t3;
596 t1=t1/total;
597 t2=t2/total;
598 t3=t3/total;*/
600 /* t2 = x;
601 t3 = y;
602 t1 = 1-t2-t3;*/
604 // return value
605 // vec3_t bary_coords(t1,t2,t3);
606 // return bary_coords;
609 vec3_t intersectionOnPlane(vec3_t v, vec3_t A, vec3_t nA, vec3_t B, vec3_t nB)
611 vec3_t u = B-A;
612 // u.normalise();
613 v.normalise();
614 v = u.abs()*v;
616 //cout<<"u="<<u<<" v="<<v<<endl;
618 vec2_t p_A(0,0);
619 vec2_t p_B(1,0);
620 vec2_t p_nA = projectVectorOnPlane(nA,u,v);
621 vec2_t p_nB = projectVectorOnPlane(nB,u,v);
623 vec2_t p_tA = turnRight(p_nA);
624 vec2_t p_tB = turnRight(p_nB);
626 double k1, k2;
627 vec2_t p_K;
628 if(!intersection(k1, k2, p_A, p_tA, p_B, p_tB)) {
629 //qDebug()<<"WARNING: No intersection found!!!";
630 p_K = 0.5*(p_A + p_B);
632 else {
633 p_K = p_A + k1*p_tA;
636 //cout<<"nA="<<nA<<endl;
637 //cout<<"p_nA="<<p_nA<<endl;
638 //cout<<"p_tA="<<p_tA<<endl;
639 //cout<<"p_K="<<p_K<<endl;
640 if(p_K[0]<0) p_K[0] = 0;
641 if(p_K[0]>1) p_K[0] = 1;
642 vec3_t K = A + p_K[0]*u + p_K[1]*v;
643 //cout<<"K="<<K<<endl;
644 return K;
647 vec2_t projectVectorOnPlane(vec3_t V,vec3_t i,vec3_t j)
649 if(i.abs2()==0) EG_BUG;
650 if(j.abs2()==0) EG_BUG;
651 double x = V*i/i.abs2();
652 double y = V*j/j.abs2();
653 return vec2_t(x,y);
656 vec3_t projectPointOnPlane(const vec3_t& M, const vec3_t& A, const vec3_t& N)
658 double k = ((M-A)*N)/N.abs2();
659 return( M - k*N );
662 vec3_t projectPointOnEdge(const vec3_t& M,const vec3_t& A, const vec3_t& u)
664 checkVector(u);
665 if(u.abs2()==0) EG_BUG;
666 double k = ((M-A)*u)/u.abs2();
667 return A + k*u;
670 void cart2spherical(vec3_t x, double &alpha, double &beta, double &r)
672 r = x.abs();
673 static const vec3_t ex(1,0,0);
674 vec3_t xy(x[0],x[1],0);
675 if (x[1] >= 0) {
676 alpha = angle(ex, xy);
677 } else {
678 alpha = 2*M_PI - angle(xy, ex);
680 if (xy.abs2() > 0) {
681 if (x[2] >= 0) {
682 beta = angle(xy, x);
683 } else {
684 beta = -angle(xy, x);
686 } else {
687 beta = 0.5*M_PI;
691 } // namespace