gmsh-TingyuanDoc  0.1
An Open-Source Timing-driven Analytical Mixed-size FPGA Placer
Numeric.cpp
Go to the documentation of this file.
1 // Gmsh - Copyright (C) 1997-2022 C. Geuzaine, J.-F. Remacle
2 //
3 // See the LICENSE.txt file in the Gmsh root directory for license information.
4 // Please report all issues on https://gitlab.onelab.info/gmsh/gmsh/issues.
5 
6 #include <algorithm>
7 #include "GmshConfig.h"
8 #include "GmshMessage.h"
9 #include "fullMatrix.h"
10 #include "Numeric.h"
11 
12 double myatan2(double a, double b)
13 {
14  if(a == 0.0 && b == 0) return 0.0;
15  return atan2(a, b);
16 }
17 
18 double myasin(double a)
19 {
20  if(a <= -1.)
21  return -M_PI / 2.;
22  else if(a >= 1.)
23  return M_PI / 2.;
24  else
25  return asin(a);
26 }
27 
28 double myacos(double a)
29 {
30  if(a <= -1.)
31  return M_PI;
32  else if(a >= 1.)
33  return 0.;
34  else
35  return acos(a);
36 }
37 
38 double norm2(double a[3][3])
39 {
40  double norm2sq =
41  std::pow(a[0][0], 2) + std::pow(a[0][1], 2) + std::pow(a[0][2], 2) +
42  std::pow(a[1][0], 2) + std::pow(a[1][1], 2) + std::pow(a[1][2], 2) +
43  std::pow(a[2][0], 2) + std::pow(a[2][1], 2) + std::pow(a[2][2], 2);
44  return std::sqrt(norm2sq);
45 }
46 
47 void matvec(double mat[3][3], double vec[3], double res[3])
48 {
49  res[0] = mat[0][0] * vec[0] + mat[0][1] * vec[1] + mat[0][2] * vec[2];
50  res[1] = mat[1][0] * vec[0] + mat[1][1] * vec[1] + mat[1][2] * vec[2];
51  res[2] = mat[2][0] * vec[0] + mat[2][1] * vec[1] + mat[2][2] * vec[2];
52 }
53 
54 void matmat(double mat1[3][3], double mat2[3][3], double res[3][3])
55 {
56  res[0][0] =
57  mat1[0][0] * mat2[0][0] + mat1[0][1] * mat2[1][0] + mat1[0][2] * mat2[2][0];
58  res[0][1] =
59  mat1[0][0] * mat2[0][1] + mat1[0][1] * mat2[1][1] + mat1[0][2] * mat2[2][1];
60  res[0][2] =
61  mat1[0][0] * mat2[0][2] + mat1[0][1] * mat2[1][2] + mat1[0][2] * mat2[2][2];
62  res[1][0] =
63  mat1[1][0] * mat2[0][0] + mat1[1][1] * mat2[1][0] + mat1[1][2] * mat2[2][0];
64  res[1][1] =
65  mat1[1][0] * mat2[0][1] + mat1[1][1] * mat2[1][1] + mat1[1][2] * mat2[2][1];
66  res[1][2] =
67  mat1[1][0] * mat2[0][2] + mat1[1][1] * mat2[1][2] + mat1[1][2] * mat2[2][2];
68  res[2][0] =
69  mat1[2][0] * mat2[0][0] + mat1[2][1] * mat2[1][0] + mat1[2][2] * mat2[2][0];
70  res[2][1] =
71  mat1[2][0] * mat2[0][1] + mat1[2][1] * mat2[1][1] + mat1[2][2] * mat2[2][1];
72  res[2][2] =
73  mat1[2][0] * mat2[0][2] + mat1[2][1] * mat2[1][2] + mat1[2][2] * mat2[2][2];
74 }
75 
76 void normal3points(double x0, double y0, double z0, double x1, double y1,
77  double z1, double x2, double y2, double z2, double n[3])
78 {
79  double t1[3] = {x1 - x0, y1 - y0, z1 - z0};
80  double t2[3] = {x2 - x0, y2 - y0, z2 - z0};
81  prodve(t1, t2, n);
82  norme(n);
83 }
84 
85 void normal2points(double x0, double y0, double z0, double x1, double y1,
86  double z1, double n[3])
87 {
88  // this computes one of the normals to the edge
89  double t[3] = {x1 - x0, y1 - y0, z1 - z0};
90  double ex[3] = {0., 0., 0.};
91  if(t[0] == 0.)
92  ex[0] = 1.;
93  else if(t[1] == 0.)
94  ex[1] = 1.;
95  else
96  ex[2] = 1.;
97  prodve(t, ex, n);
98  norme(n);
99 }
100 
101 int sys2x2(double mat[2][2], double b[2], double res[2])
102 {
103  double const norm = std::pow(mat[0][0], 2) + std::pow(mat[1][1], 2) +
104  std::pow(mat[0][1], 2) + std::pow(mat[1][0], 2);
105 
106  double const det = mat[0][0] * mat[1][1] - mat[1][0] * mat[0][1];
107 
108  // TOLERANCE ! WARNING WARNING
109  if(norm == 0.0 || std::abs(det) / norm < 1.e-16) {
110  // if(norm)
111  // Msg::Debug("Assuming 2x2 matrix is singular (det/norm == %.16g)",
112  // std::abs(det) / norm);
113  res[0] = res[1] = 0.0;
114  return 0;
115  }
116  double const ud = 1.0 / det;
117 
118  res[0] = b[0] * mat[1][1] - mat[0][1] * b[1];
119  res[1] = mat[0][0] * b[1] - mat[1][0] * b[0];
120 
121  for(int i = 0; i < 2; i++) res[i] *= ud;
122 
123  return 1;
124 }
125 
126 double det3x3(double mat[3][3])
127 {
128  return (mat[0][0] * (mat[1][1] * mat[2][2] - mat[1][2] * mat[2][1]) -
129  mat[0][1] * (mat[1][0] * mat[2][2] - mat[1][2] * mat[2][0]) +
130  mat[0][2] * (mat[1][0] * mat[2][1] - mat[1][1] * mat[2][0]));
131 }
132 
133 double trace3x3(double mat[3][3]) { return mat[0][0] + mat[1][1] + mat[2][2]; }
134 
135 double trace2(double mat[3][3])
136 {
137  double a00 =
138  mat[0][0] * mat[0][0] + mat[1][0] * mat[0][1] + mat[2][0] * mat[0][2];
139  double a11 =
140  mat[1][0] * mat[0][1] + mat[1][1] * mat[1][1] + mat[1][2] * mat[2][1];
141  double a22 =
142  mat[2][0] * mat[0][2] + mat[2][1] * mat[1][2] + mat[2][2] * mat[2][2];
143 
144  return a00 + a11 + a22;
145 }
146 
147 int sys3x3(double mat[3][3], double b[3], double res[3], double *det)
148 {
149  double ud;
150  int i;
151 
152  *det = det3x3(mat);
153 
154  if(*det == 0.0) {
155  res[0] = res[1] = res[2] = 0.0;
156  return (0);
157  }
158 
159  ud = 1. / (*det);
160 
161  res[0] = b[0] * (mat[1][1] * mat[2][2] - mat[1][2] * mat[2][1]) -
162  mat[0][1] * (b[1] * mat[2][2] - mat[1][2] * b[2]) +
163  mat[0][2] * (b[1] * mat[2][1] - mat[1][1] * b[2]);
164 
165  res[1] = mat[0][0] * (b[1] * mat[2][2] - mat[1][2] * b[2]) -
166  b[0] * (mat[1][0] * mat[2][2] - mat[1][2] * mat[2][0]) +
167  mat[0][2] * (mat[1][0] * b[2] - b[1] * mat[2][0]);
168 
169  res[2] = mat[0][0] * (mat[1][1] * b[2] - b[1] * mat[2][1]) -
170  mat[0][1] * (mat[1][0] * b[2] - b[1] * mat[2][0]) +
171  b[0] * (mat[1][0] * mat[2][1] - mat[1][1] * mat[2][0]);
172 
173  for(i = 0; i < 3; i++) res[i] *= ud;
174  return (1);
175 }
176 
177 int sys3x3_with_tol(double mat[3][3], double b[3], double res[3], double *det)
178 {
179  int out;
180  double norm;
181 
182  out = sys3x3(mat, b, res, det);
183  norm =
184  std::pow(mat[0][0], 2) + std::pow(mat[0][1], 2) + std::pow(mat[0][2], 2) +
185  std::pow(mat[1][0], 2) + std::pow(mat[1][1], 2) + std::pow(mat[1][2], 2) +
186  std::pow(mat[2][0], 2) + std::pow(mat[2][1], 2) + std::pow(mat[2][2], 2);
187 
188  // TOLERANCE ! WARNING WARNING
189  if(norm == 0.0 || std::abs(*det) / norm < 1.e-12) {
190  if(norm)
191  Msg::Debug("Assuming 3x3 matrix is singular (det/norm == %.16g)",
192  std::abs(*det) / norm);
193  res[0] = res[1] = res[2] = 0.0;
194  return 0;
195  }
196  return out;
197 }
198 
199 double det2x2(double mat[2][2])
200 {
201  return mat[0][0] * mat[1][1] - mat[1][0] * mat[0][1];
202 }
203 
204 double det2x3(double mat[2][3])
205 {
206  double v1[3] = {mat[0][0], mat[0][1], mat[0][2]};
207  double v2[3] = {mat[1][0], mat[1][1], mat[1][2]};
208  double n[3];
209 
210  prodve(v1, v2, n);
211  return norm3(n);
212 }
213 
214 double inv2x2(double mat[2][2], double inv[2][2])
215 {
216  const double det = det2x2(mat);
217  if(det) {
218  double ud = 1. / det;
219  inv[0][0] = mat[1][1] * ud;
220  inv[1][0] = -mat[1][0] * ud;
221  inv[0][1] = -mat[0][1] * ud;
222  inv[1][1] = mat[0][0] * ud;
223  }
224  else {
225  Msg::Error("Singular matrix 2x2");
226  for(int i = 0; i < 2; i++)
227  for(int j = 0; j < 2; j++) inv[i][j] = 0.;
228  }
229  return det;
230 }
231 
232 double inv3x3(double mat[3][3], double inv[3][3])
233 {
234  double det = det3x3(mat);
235  if(det) {
236  double ud = 1. / det;
237  inv[0][0] = (mat[1][1] * mat[2][2] - mat[1][2] * mat[2][1]) * ud;
238  inv[1][0] = -(mat[1][0] * mat[2][2] - mat[1][2] * mat[2][0]) * ud;
239  inv[2][0] = (mat[1][0] * mat[2][1] - mat[1][1] * mat[2][0]) * ud;
240  inv[0][1] = -(mat[0][1] * mat[2][2] - mat[0][2] * mat[2][1]) * ud;
241  inv[1][1] = (mat[0][0] * mat[2][2] - mat[0][2] * mat[2][0]) * ud;
242  inv[2][1] = -(mat[0][0] * mat[2][1] - mat[0][1] * mat[2][0]) * ud;
243  inv[0][2] = (mat[0][1] * mat[1][2] - mat[0][2] * mat[1][1]) * ud;
244  inv[1][2] = -(mat[0][0] * mat[1][2] - mat[0][2] * mat[1][0]) * ud;
245  inv[2][2] = (mat[0][0] * mat[1][1] - mat[0][1] * mat[1][0]) * ud;
246  }
247  else {
248  Msg::Error("Singular matrix 3x3");
249  for(int i = 0; i < 3; i++)
250  for(int j = 0; j < 3; j++) inv[i][j] = 0.;
251  }
252  return det;
253 }
254 
255 double angle_02pi(double A3)
256 {
257  double DP = 2 * M_PI;
258  while(A3 > DP || A3 < 0.) {
259  if(A3 > 0)
260  A3 -= DP;
261  else
262  A3 += DP;
263  }
264  return A3;
265 }
266 
267 double angle_plan(double v[3], double p1[3], double p2[3], double n[3])
268 {
269  double PA[3], PB[3], angplan;
270 
271  PA[0] = p1[0] - v[0];
272  PA[1] = p1[1] - v[1];
273  PA[2] = p1[2] - v[2];
274 
275  PB[0] = p2[0] - v[0];
276  PB[1] = p2[1] - v[1];
277  PB[2] = p2[2] - v[2];
278 
279  norme(PA);
280  norme(PB);
281 
282  double c[3];
283  prodve(PA, PB, c);
284 
285  double const cosc = prosca(PA, PB);
286  double const sinc = prosca(c, n);
287 
288  angplan = myatan2(sinc, cosc);
289 
290  return angplan;
291 }
292 
293 double triangle_area(double p0[3], double p1[3], double p2[3])
294 {
295  double a[3], b[3], c[3];
296 
297  a[0] = p2[0] - p1[0];
298  a[1] = p2[1] - p1[1];
299  a[2] = p2[2] - p1[2];
300 
301  b[0] = p0[0] - p1[0];
302  b[1] = p0[1] - p1[1];
303  b[2] = p0[2] - p1[2];
304 
305  prodve(a, b, c);
306  return 0.5 * std::sqrt(c[0] * c[0] + c[1] * c[1] + c[2] * c[2]);
307 }
308 
309 double triangle_area2d(double p0[2], double p1[2], double p2[2])
310 {
311  const double c =
312  (p2[0] - p1[0]) * (p0[1] - p1[1]) - (p2[1] - p1[1]) * (p0[0] - p1[0]);
313 
314  return 0.5 * std::sqrt(c * c);
315 }
316 
317 void circumCenterXY(double *p1, double *p2, double *p3, double *res)
318 {
319  double d, a1, a2, a3;
320 
321  const double x1 = p1[0];
322  const double x2 = p2[0];
323  const double x3 = p3[0];
324  const double y1 = p1[1];
325  const double y2 = p2[1];
326  const double y3 = p3[1];
327 
328  d = 2. * (double)(y1 * (x2 - x3) + y2 * (x3 - x1) + y3 * (x1 - x2));
329  if(d == 0.0) {
330  // Msg::Warning("Colinear points in circum circle computation");
331  res[0] = res[1] = -99999.;
332  return;
333  }
334 
335  a1 = x1 * x1 + y1 * y1;
336  a2 = x2 * x2 + y2 * y2;
337  a3 = x3 * x3 + y3 * y3;
338  res[0] = (double)((a1 * (y3 - y2) + a2 * (y1 - y3) + a3 * (y2 - y1)) / d);
339  res[1] = (double)((a1 * (x2 - x3) + a2 * (x3 - x1) + a3 * (x1 - x2)) / d);
340 }
341 
342 void circumCenterXYZ(double *p1, double *p2, double *p3, double *res,
343  double *uv)
344 {
345  double v1[3] = {p2[0] - p1[0], p2[1] - p1[1], p2[2] - p1[2]};
346  double v2[3] = {p3[0] - p1[0], p3[1] - p1[1], p3[2] - p1[2]};
347  double vx[3] = {p2[0] - p1[0], p2[1] - p1[1], p2[2] - p1[2]};
348  double vy[3] = {p3[0] - p1[0], p3[1] - p1[1], p3[2] - p1[2]};
349  double vz[3];
350  prodve(vx, vy, vz);
351  prodve(vz, vx, vy);
352  norme(vx);
353  norme(vy);
354  norme(vz);
355 
356  double p1P[2] = {0.0, 0.0};
357  double p2P[2] = {prosca(v1, vx), prosca(v1, vy)};
358  double p3P[2] = {prosca(v2, vx), prosca(v2, vy)};
359 
360  double resP[2];
361 
362  circumCenterXY(p1P, p2P, p3P, resP);
363 
364  if(uv) {
365  double mat[2][2] = {{p2P[0] - p1P[0], p3P[0] - p1P[0]},
366  {p2P[1] - p1P[1], p3P[1] - p1P[1]}};
367  double rhs[2] = {resP[0] - p1P[0], resP[1] - p1P[1]};
368  sys2x2(mat, rhs, uv);
369  }
370 
371  res[0] = p1[0] + resP[0] * vx[0] + resP[1] * vy[0];
372  res[1] = p1[1] + resP[0] * vx[1] + resP[1] * vy[1];
373  res[2] = p1[2] + resP[0] * vx[2] + resP[1] * vy[2];
374 }
375 
376 void planarQuad_xyz2xy(double *x, double *y, double *z, double *xn, double *yn)
377 {
378  double v1[3] = {x[1] - x[0], y[1] - y[0], z[1] - z[0]};
379  double v2[3] = {x[2] - x[0], y[2] - y[0], z[2] - z[0]};
380  double v3[3] = {x[3] - x[0], y[3] - y[0], z[3] - z[0]};
381 
382  double vx[3] = {x[1] - x[0], y[1] - y[0], z[1] - z[0]};
383  double vy[3] = {x[2] - x[0], y[2] - y[0], z[2] - z[0]};
384  double vz[3];
385  prodve(vx, vy, vz);
386  prodve(vz, vx, vy);
387 
388  norme(vx);
389  norme(vy);
390  norme(vz);
391 
392  double p1P[2] = {0.0, 0.0};
393  double p2P[2] = {prosca(v1, vx), prosca(v1, vy)};
394  double p3P[2] = {prosca(v2, vx), prosca(v2, vy)};
395  double p4P[2] = {prosca(v3, vx), prosca(v3, vy)};
396 
397  xn[0] = p1P[0];
398  xn[1] = p2P[0];
399  xn[2] = p3P[0];
400  xn[3] = p4P[0];
401  yn[0] = p1P[1];
402  yn[1] = p2P[1];
403  yn[2] = p3P[1];
404  yn[3] = p4P[1];
405 }
406 
407 double computeInnerRadiusForQuad(double *x, double *y, int i)
408 {
409  // parameters of the equations of the 3 edges
410  double a1 = y[(4 + i) % 4] - y[(5 + i) % 4];
411  double a2 = y[(5 + i) % 4] - y[(6 + i) % 4];
412  double a3 = y[(6 + i) % 4] - y[(7 + i) % 4];
413 
414  double b1 = x[(5 + i) % 4] - x[(4 + i) % 4];
415  double b2 = x[(6 + i) % 4] - x[(5 + i) % 4];
416  double b3 = x[(7 + i) % 4] - x[(6 + i) % 4];
417 
418  double c1 = y[(5 + i) % 4] * x[(4 + i) % 4] - y[(4 + i) % 4] * x[(5 + i) % 4];
419  double c2 = y[(6 + i) % 4] * x[(5 + i) % 4] - y[(5 + i) % 4] * x[(6 + i) % 4];
420  double c3 = y[(7 + i) % 4] * x[(6 + i) % 4] - y[(6 + i) % 4] * x[(7 + i) % 4];
421 
422  // length of the 3 edges
423  double l1 = std::sqrt(a1 * a1 + b1 * b1);
424  double l2 = std::sqrt(a2 * a2 + b2 * b2);
425  double l3 = std::sqrt(a3 * a3 + b3 * b3);
426 
427  // parameters of the 2 bisectors
428  double a12 = a1 / l1 - a2 / l2;
429  double a23 = a2 / l2 - a3 / l3;
430 
431  double b12 = b1 / l1 - b2 / l2;
432  double b23 = b2 / l2 - b3 / l3;
433 
434  double c12 = c1 / l1 - c2 / l2;
435  double c23 = c2 / l2 - c3 / l3;
436 
437  // compute the coordinates of the center of the incircle,
438  // that is the point where the 2 bisectors meet
439  double x_s = (c12 * b23 - c23 * b12) / (a23 * b12 - a12 * b23);
440  double y_s = 0.;
441  if(b12 != 0) { y_s = -a12 / b12 * x_s - c12 / b12; }
442  else {
443  y_s = -a23 / b23 * x_s - c23 / b23;
444  }
445 
446  // finally get the radius of the circle
447  double r = (a1 * x_s + b1 * y_s + c1) / l1;
448 
449  return r;
450 }
451 
452 char float2char(float f)
453 {
454  // float normalized in [-1, 1], char in [-127, 127]
455  float c = f * 127.;
456  if(c > 127.)
457  return 127;
458  else if(c < -127.)
459  return -127;
460  else
461  return (int)c;
462 }
463 
464 float char2float(char c)
465 {
466  float f = c;
467  if(f > 127.)
468  return 1.;
469  else if(f < -127)
470  return -1.;
471  else
472  return f / 127.;
473 }
474 
475 void gradSimplex(double *x, double *y, double *z, double *v, double *grad)
476 {
477  // p = p1 * (1-u-v-w) + p2 u + p3 v + p4 w
478 
479  double mat[3][3];
480  double det, b[3];
481  mat[0][0] = x[1] - x[0];
482  mat[1][0] = x[2] - x[0];
483  mat[2][0] = x[3] - x[0];
484  mat[0][1] = y[1] - y[0];
485  mat[1][1] = y[2] - y[0];
486  mat[2][1] = y[3] - y[0];
487  mat[0][2] = z[1] - z[0];
488  mat[1][2] = z[2] - z[0];
489  mat[2][2] = z[3] - z[0];
490  b[0] = v[1] - v[0];
491  b[1] = v[2] - v[0];
492  b[2] = v[3] - v[0];
493  sys3x3(mat, b, grad, &det);
494 }
495 
496 double ComputeVonMises(double *V)
497 {
498  double tr = (V[0] + V[4] + V[8]) / 3.;
499  double v11 = V[0] - tr, v12 = V[1], v13 = V[2];
500  double v21 = V[3], v22 = V[4] - tr, v23 = V[5];
501  double v31 = V[6], v32 = V[7], v33 = V[8] - tr;
502  return sqrt(1.5 * (v11 * v11 + v12 * v12 + v13 * v13 + v21 * v21 + v22 * v22 +
503  v23 * v23 + v31 * v31 + v32 * v32 + v33 * v33));
504 }
505 
506 double ComputeScalarRep(int numComp, double *val, int tensorRep)
507 {
508  if(numComp == 1)
509  return val[0];
510  else if(numComp == 3)
511  return sqrt(val[0] * val[0] + val[1] * val[1] + val[2] * val[2]);
512  else if(numComp == 9) {
513  if(tensorRep == 0) { // Von-Mises
514  return ComputeVonMises(val);
515  }
516  else {
517  fullMatrix<double> tensor(3, 3);
518  fullVector<double> S(3), imS(3);
519  fullMatrix<double> V(3, 3);
520  fullMatrix<double> rightV(3, 3);
521  for(int j = 0; j < 3; j++) {
522  tensor(j, 0) = val[0 + j * 3];
523  tensor(j, 1) = val[1 + j * 3];
524  tensor(j, 2) = val[2 + j * 3];
525  }
526  tensor.eig(S, imS, V, rightV, true);
527  if(tensorRep == 1) { // max eigenvalue
528  return S(2);
529  }
530  else { // min eigenvalue
531  return S(0);
532  }
533  }
534  }
535  return 0.;
536 }
537 
538 void eigenvalue2x2(double mat[2][2], double v[2])
539 {
540  double a = 1;
541  double b = -(mat[0][0] + mat[1][1]);
542  double c = (mat[0][0] * mat[1][1]) - (mat[0][1] * mat[1][0]);
543 
544  double det = b * b - 4. * a * c;
545 
546  v[0] = (-b + sqrt(det)) / (2 * a);
547  v[1] = (-b - sqrt(det)) / (2 * a);
548 }
549 
550 void eigenvalue(double mat[3][3], double v[3])
551 {
552  // characteristic polynomial of T : find v root of
553  // v^3 - I1 v^2 + I2 T - I3 = 0
554  // I1 : first invariant , trace(T)
555  // I2 : second invariant , 1/2 (I1^2 -trace(T^2))
556  // I3 : third invariant , det T
557 
558  double c[4];
559  c[3] = 1.0;
560  c[2] = -trace3x3(mat);
561  c[1] = 0.5 * (c[2] * c[2] - trace2(mat));
562  c[0] = -det3x3(mat);
563 
564  // printf("%g %g %g\n", mat[0][0], mat[0][1], mat[0][2]);
565  // printf("%g %g %g\n", mat[1][0], mat[1][1], mat[1][2]);
566  // printf("%g %g %g\n", mat[2][0], mat[2][1], mat[2][2]);
567  // printf("%g x^3 + %g x^2 + %g x + %g = 0\n", c[3], c[2], c[1], c[0]);
568 
569  double imag[3];
570  FindCubicRoots(c, v, imag);
571  eigsort(v);
572 }
573 
574 void FindCubicRoots(const double coef[4], double real[3], double imag[3])
575 {
576  double a = coef[3];
577  double b = coef[2];
578  double c = coef[1];
579  double d = coef[0];
580 
581  if(!a || !d) {
582  // Msg::Error("Degenerate cubic: use a second degree solver!");
583  return;
584  }
585 
586  b /= a;
587  c /= a;
588  d /= a;
589 
590  double q = (3.0 * c - (b * b)) / 9.0;
591  double r = -(27.0 * d) + b * (9.0 * c - 2.0 * (b * b));
592  r /= 54.0;
593 
594  double discrim = q * q * q + r * r;
595  imag[0] = 0.0; // The first root is always real.
596  double term1 = (b / 3.0);
597 
598  if(discrim > 0) { // one root is real, two are complex
599  double s = r + sqrt(discrim);
600  s = ((s < 0) ? -pow(-s, (1.0 / 3.0)) : pow(s, (1.0 / 3.0)));
601  double t = r - sqrt(discrim);
602  t = ((t < 0) ? -pow(-t, (1.0 / 3.0)) : pow(t, (1.0 / 3.0)));
603  real[0] = -term1 + s + t;
604  term1 += (s + t) / 2.0;
605  real[1] = real[2] = -term1;
606  term1 = sqrt(3.0) * (-t + s) / 2;
607  imag[1] = term1;
608  imag[2] = -term1;
609  return;
610  }
611 
612  // The remaining options are all real
613  imag[1] = imag[2] = 0.0;
614 
615  double r13;
616  if(discrim == 0) { // All roots real, at least two are equal.
617  r13 = ((r < 0) ? -pow(-r, (1.0 / 3.0)) : pow(r, (1.0 / 3.0)));
618  real[0] = -term1 + 2.0 * r13;
619  real[1] = real[2] = -(r13 + term1);
620  return;
621  }
622 
623  // Only option left is that all roots are real and unequal (to get
624  // here, q < 0)
625  q = -q;
626  double dum1 = q * q * q;
627  dum1 = acos(r / sqrt(dum1));
628  r13 = 2.0 * sqrt(q);
629  real[0] = -term1 + r13 * cos(dum1 / 3.0);
630  real[1] = -term1 + r13 * cos((dum1 + 2.0 * M_PI) / 3.0);
631  real[2] = -term1 + r13 * cos((dum1 + 4.0 * M_PI) / 3.0);
632 }
633 
634 void eigsort(double d[3])
635 {
636  int k, j, i;
637  double p;
638 
639  for(i = 0; i < 3; i++) {
640  p = d[k = i];
641  for(j = i + 1; j < 3; j++)
642  if(d[j] >= p) p = d[k = j];
643  if(k != i) {
644  d[k] = d[i];
645  d[i] = p;
646  }
647  }
648 }
649 
650 void invert_singular_matrix3x3(double MM[3][3], double II[3][3])
651 {
652  int i, j, k, n = 3;
653  double TT[3][3];
654 
655  for(i = 1; i <= n; i++) {
656  for(j = 1; j <= n; j++) {
657  II[i - 1][j - 1] = 0.0;
658  TT[i - 1][j - 1] = 0.0;
659  }
660  }
661 
662  fullMatrix<double> M(3, 3), V(3, 3);
663  fullVector<double> W(3);
664  for(i = 1; i <= n; i++) {
665  for(j = 1; j <= n; j++) { M(i - 1, j - 1) = MM[i - 1][j - 1]; }
666  }
667  M.svd(V, W);
668  for(i = 1; i <= n; i++) {
669  for(j = 1; j <= n; j++) {
670  double ww = W(i - 1);
671  if(fabs(ww) > 1.e-16) { // singular value precision
672  TT[i - 1][j - 1] += M(j - 1, i - 1) / ww;
673  }
674  }
675  }
676  for(i = 1; i <= n; i++) {
677  for(j = 1; j <= n; j++) {
678  for(k = 1; k <= n; k++) {
679  II[i - 1][j - 1] += V(i - 1, k - 1) * TT[k - 1][j - 1];
680  }
681  }
682  }
683 }
684 
685 bool newton_fd(bool (*func)(fullVector<double> &, fullVector<double> &, void *),
686  fullVector<double> &x, void *data, double relax, double tolx)
687 {
688  const int MAXIT = 100;
689  const double EPS = 1.e-4;
690  const int N = x.size();
691 
692  fullMatrix<double> J(N, N);
693  fullVector<double> f(N), feps(N), dx(N);
694 
695  for(int iter = 0; iter < MAXIT; iter++) {
696  if(x.norm() > 1.e6) return false;
697  if(!func(x, f, data)) { return false; }
698 
699  bool isZero = false;
700  for(int k = 0; k < N; k++) {
701  if(f(k) == 0.)
702  isZero = true;
703  else
704  isZero = false;
705  if(isZero == false) break;
706  }
707  if(isZero) break;
708 
709  for(int j = 0; j < N; j++) {
710  double h = EPS * fabs(x(j));
711  if(h == 0.) h = EPS;
712  x(j) += h;
713  if(!func(x, feps, data)) { return false; }
714  for(int i = 0; i < N; i++) { J(i, j) = (feps(i) - f(i)) / h; }
715  x(j) -= h;
716  }
717 
718  if(N == 1)
719  dx(0) = f(0) / J(0, 0);
720  else if(!J.luSolve(f, dx)) {
721  return false;
722  }
723 
724  for(int i = 0; i < N; i++) x(i) -= relax * dx(i);
725 
726  if(dx.norm() < tolx) { return true; }
727  }
728  return false;
729 }
730 
731 /*
732  distance to triangle
733 
734  P(p) = p1 + t1 xi + t2 eta
735 
736  t1 = (p2-p1) ; t2 = (p3-p1) ;
737 
738  (P(p) - p) = d n
739 
740  (p1 + t1 xi + t2 eta - p) = d n
741  t1 xi + t2 eta - d n = p - p1
742 
743  | t1x t2x -nx | |xi | |px-p1x|
744  | t1y t2y -ny | |eta | = |py-p1y|
745  | t1z t2z -nz | |d | |pz-p1z|
746 
747  distance to segment
748 
749  P(p) = p1 + t (p2-p1)
750 
751  (p - P(p)) * (p2-p1) = 0
752  (p - p1 - t (p2-p1) ) * (p2-p1) = 0
753  - t ||p2-p1||^2 + (p-p1)(p2-p1) = 0
754 
755  t = (p-p1)*(p2-p1)/||p2-p1||^2
756 */
757 
758 void signedDistancePointTriangle(const SPoint3 &p1, const SPoint3 &p2,
759  const SPoint3 &p3, const SPoint3 &p, double &d,
760  SPoint3 &closePt)
761 {
762  SVector3 t1 = p2 - p1;
763  SVector3 t2 = p3 - p1;
764  SVector3 t3 = p3 - p2;
765  SVector3 n = crossprod(t1, t2);
766  n.normalize();
767  const double n2t1 = dot(t1, t1);
768  const double n2t2 = dot(t2, t2);
769  const double n2t3 = dot(t3, t3);
770  double mat[3][3] = {{t1.x(), t2.x(), -n.x()},
771  {t1.y(), t2.y(), -n.y()},
772  {t1.z(), t2.z(), -n.z()}};
773  double inv[3][3];
774  double det = inv3x3(mat, inv);
775  if(det == 0.0) return;
776 
777  double u, v;
778  SVector3 pp1 = p - p1;
779  u = (inv[0][0] * pp1.x() + inv[0][1] * pp1.y() + inv[0][2] * pp1.z());
780  v = (inv[1][0] * pp1.x() + inv[1][1] * pp1.y() + inv[1][2] * pp1.z());
781  d = (inv[2][0] * pp1.x() + inv[2][1] * pp1.y() + inv[2][2] * pp1.z());
782  double sign = (d > 0) ? 1. : -1.;
783  if(d == 0.) sign = 1.;
784  if(u >= 0. && v >= 0. && 1. - u - v >= 0.0) { // P(p) inside triangle
785  closePt = p1 + (p2 - p1) * u + (p3 - p1) * v;
786  }
787  else {
788  const double t12 = dot(pp1, t1) / n2t1;
789  const double t13 = dot(pp1, t2) / n2t2;
790  SVector3 pp2 = p - p2;
791  const double t23 = dot(pp2, t3) / n2t3;
792  d = 1.e10;
793  if(t12 >= 0 && t12 <= 1.) {
794  d = sign * std::min(fabs(d), p.distance(p1 + (p2 - p1) * t12));
795  closePt = p1 + (p2 - p1) * t12;
796  }
797  if(t13 >= 0 && t13 <= 1.) {
798  if(p.distance(p1 + (p3 - p1) * t13) < fabs(d))
799  closePt = p1 + (p3 - p1) * t13;
800  d = sign * std::min(fabs(d), p.distance(p1 + (p3 - p1) * t13));
801  }
802  if(t23 >= 0 && t23 <= 1.) {
803  if(p.distance(p2 + (p3 - p2) * t23) < fabs(d))
804  closePt = p2 + (p3 - p2) * t23;
805  d = sign * std::min(fabs(d), p.distance(p2 + (p3 - p2) * t23));
806  }
807  if(p.distance(p1) < fabs(d)) {
808  closePt = p1;
809  d = sign * std::min(fabs(d), p.distance(p1));
810  }
811  if(p.distance(p2) < fabs(d)) {
812  closePt = p2;
813  d = sign * std::min(fabs(d), p.distance(p2));
814  }
815  if(p.distance(p3) < fabs(d)) {
816  closePt = p3;
817  d = sign * std::min(fabs(d), p.distance(p3));
818  }
819  }
820 }
821 
822 void signedDistancesPointsTriangle(std::vector<double> &distances,
823  std::vector<SPoint3> &closePts,
824  const std::vector<SPoint3> &pts,
825  const SPoint3 &p1, const SPoint3 &p2,
826  const SPoint3 &p3)
827 {
828  const unsigned pts_size = pts.size();
829  distances.clear();
830  distances.resize(pts_size);
831  closePts.clear();
832  closePts.resize(pts_size);
833 
834  for(std::size_t i = 0; i < pts_size; ++i) distances[i] = 1.e22;
835 
836  SVector3 t1 = p2 - p1;
837  SVector3 t2 = p3 - p1;
838  SVector3 t3 = p3 - p2;
839  SVector3 n = crossprod(t1, t2);
840  n.normalize();
841  const double n2t1 = dot(t1, t1);
842  const double n2t2 = dot(t2, t2);
843  const double n2t3 = dot(t3, t3);
844  double mat[3][3] = {{t1.x(), t2.x(), -n.x()},
845  {t1.y(), t2.y(), -n.y()},
846  {t1.z(), t2.z(), -n.z()}};
847  double inv[3][3];
848  double det = inv3x3(mat, inv);
849  if(det == 0.0) return;
850 
851  for(std::size_t i = 0; i < pts.size(); i++) {
852  double d;
853  SPoint3 closePt;
854  const SPoint3 &p = pts[i];
855 
856  // signedDistancePointTrianglePrecomputed(p1, p2, p3, p, d, closePt);
857  double u, v;
858  SVector3 pp1 = p - p1;
859  u = (inv[0][0] * pp1.x() + inv[0][1] * pp1.y() + inv[0][2] * pp1.z());
860  v = (inv[1][0] * pp1.x() + inv[1][1] * pp1.y() + inv[1][2] * pp1.z());
861  d = (inv[2][0] * pp1.x() + inv[2][1] * pp1.y() + inv[2][2] * pp1.z());
862  double sign = (d > 0) ? 1. : -1.;
863  if(d == 0.) sign = 1.;
864  if(u >= 0 && v >= 0 && 1. - u - v >= 0.0) {
865  closePt = SPoint3(0., 0., 0.); // TO DO
866  }
867  else {
868  const double t12 = dot(pp1, t1) / n2t1;
869  const double t13 = dot(pp1, t2) / n2t2;
870  SVector3 pp2 = p - p2;
871  const double t23 = dot(pp2, t3) / n2t3;
872  d = 1.e10;
873  SPoint3 closePt;
874  if(t12 >= 0 && t12 <= 1.) {
875  d = sign * std::min(fabs(d), p.distance(p1 + (p2 - p1) * t12));
876  closePt = p1 + (p2 - p1) * t12;
877  }
878  if(t13 >= 0 && t13 <= 1.) {
879  if(p.distance(p1 + (p3 - p1) * t13) < fabs(d))
880  closePt = p1 + (p3 - p1) * t13;
881  d = sign * std::min(fabs(d), p.distance(p1 + (p3 - p1) * t13));
882  }
883  if(t23 >= 0 && t23 <= 1.) {
884  if(p.distance(p2 + (p3 - p2) * t23) < fabs(d))
885  closePt = p2 + (p3 - p2) * t23;
886  d = sign * std::min(fabs(d), p.distance(p2 + (p3 - p2) * t23));
887  }
888  if(p.distance(p1) < fabs(d)) {
889  closePt = p1;
890  d = sign * std::min(fabs(d), p.distance(p1));
891  }
892  if(p.distance(p2) < fabs(d)) {
893  closePt = p2;
894  d = sign * std::min(fabs(d), p.distance(p2));
895  }
896  if(p.distance(p3) < fabs(d)) {
897  closePt = p3;
898  d = sign * std::min(fabs(d), p.distance(p3));
899  }
900  }
901  // end signedDistance
902 
903  distances[i] = d;
904  closePts[i] = closePt;
905  }
906 }
907 
908 void signedDistancePointLine(const SPoint3 &p1, const SPoint3 &p2,
909  const SPoint3 &p, double &d, SPoint3 &closePt)
910 {
911  SVector3 v12 = p2 - p1;
912  SVector3 v1p = p - p1;
913  const double alpha = dot(v1p, v12) / dot(v12, v12);
914  if(alpha <= 0.)
915  closePt = p1;
916  else if(alpha >= 1.)
917  closePt = p2;
918  else
919  closePt = p1 + (p2 - p1) * alpha;
920  d = p.distance(closePt);
921 }
922 
923 void signedDistancesPointsLine(std::vector<double> &distances,
924  std::vector<SPoint3> &closePts,
925  const std::vector<SPoint3> &pts,
926  const SPoint3 &p1, const SPoint3 &p2)
927 {
928  distances.clear();
929  distances.resize(pts.size());
930  closePts.clear();
931  closePts.resize(pts.size());
932  for(std::size_t i = 0; i < pts.size(); i++) {
933  double d;
934  SPoint3 closePt;
935  const SPoint3 &p = pts[i];
936  signedDistancePointLine(p1, p2, p, d, closePt);
937  distances[i] = d;
938  closePts[i] = closePt;
939  }
940 }
941 
942 void changeReferential(const int direction, const SPoint3 &p,
943  const SPoint3 &closePt, const SPoint3 &p1,
944  const SPoint3 &p2, double *xp, double *yp,
945  double *otherp, double *x, double *y, double *other)
946 {
947  if(direction == 1) {
948  const SPoint3 &d1 = SPoint3(1.0, 0.0, 0.0);
949  const SPoint3 &d =
950  SPoint3(p2.x() - p1.x(), p2.y() - p1.y(), p2.z() - p1.z());
951  double norm = sqrt(d.x() * d.x() + d.y() * d.y() + d.z() * d.z());
952  const SPoint3 &dn = SPoint3(d.x() / norm, d.y() / norm, d.z() / norm);
953  const SPoint3 &d3 = SPoint3(d1.y() * dn.z() - d1.z() * dn.y(),
954  d1.z() * dn.x() - d1.x() * dn.z(),
955  d1.x() * dn.y() - d1.y() * dn.x());
956  norm = sqrt(d3.x() * d3.x() + d3.y() * d3.y() + d3.z() * d3.z());
957  const SPoint3 &d3n = SPoint3(d3.x() / norm, d3.y() / norm, d3.z() / norm);
958  const SPoint3 &d2 = SPoint3(d3n.y() * d1.z() - d3n.z() * d1.y(),
959  d3n.z() * d1.x() - d3n.x() * d1.z(),
960  d3n.x() * d1.y() - d3n.y() * d1.x());
961  norm = sqrt(d2.x() * d2.x() + d2.y() * d2.y() + d2.z() * d2.z());
962  const SPoint3 &d2n = SPoint3(d2.x() / norm, d2.y() / norm, d2.z() / norm);
963  *xp = p.x() * d1.x() + p.y() * d1.y() + p.z() * d1.z();
964  *yp = p.x() * d3n.x() + p.y() * d3n.y() + p.z() * d3n.z();
965  *otherp = p.x() * d2n.x() + p.y() * d2n.y() + p.z() * d2n.z();
966  *x = closePt.x() * d1.x() + closePt.y() * d1.y() + closePt.z() * d1.z();
967  *y = closePt.x() * d3n.x() + closePt.y() * d3n.y() + closePt.z() * d3n.z();
968  *other =
969  closePt.x() * d2n.x() + closePt.y() * d2n.y() + closePt.z() * d2n.z();
970  }
971  else {
972  const SPoint3 &d2 = SPoint3(0.0, 1.0, 0.0);
973  const SPoint3 &d =
974  SPoint3(p2.x() - p1.x(), p2.y() - p1.y(), p2.z() - p1.z());
975  double norm = sqrt(d.x() * d.x() + d.y() * d.y() + d.z() * d.z());
976  const SPoint3 &dn = SPoint3(d.x() / norm, d.y() / norm, d.z() / norm);
977  const SPoint3 &d3 = SPoint3(dn.y() * d2.z() - dn.z() * d2.y(),
978  dn.z() * d2.x() - dn.x() * d2.z(),
979  dn.x() * d2.y() - dn.y() * d2.x());
980  norm = sqrt(d3.x() * d3.x() + d3.y() * d3.y() + d3.z() * d3.z());
981  const SPoint3 &d3n = SPoint3(d3.x() / norm, d3.y() / norm, d3.z() / norm);
982  const SPoint3 &d1 = SPoint3(d2.y() * d3n.z() - d2.z() * d3n.y(),
983  d2.z() * d3n.x() - d2.x() * d3n.z(),
984  d2.x() * d3n.y() - d2.y() * d3n.x());
985  norm = sqrt(d1.x() * d1.x() + d1.y() * d1.y() + d1.z() * d1.z());
986  const SPoint3 &d1n = SPoint3(d1.x() / norm, d1.y() / norm, d1.z() / norm);
987  *xp = p.x() * d2.x() + p.y() * d2.y() + p.z() * d2.z();
988  *yp = p.x() * d3n.x() + p.y() * d3n.y() + p.z() * d3n.z();
989  *otherp = p.x() * d1n.x() + p.y() * d1n.y() + p.z() * d1n.z();
990  *x = closePt.x() * d2.x() + closePt.y() * d2.y() + closePt.z() * d2.z();
991  *y = closePt.x() * d3n.x() + closePt.y() * d3n.y() + closePt.z() * d3n.z();
992  *other =
993  closePt.x() * d1n.x() + closePt.y() * d1n.y() + closePt.z() * d1n.z();
994  }
995 }
996 
997 int computeDistanceRatio(const double &y, const double &yp, const double &x,
998  const double &xp, double *distance, const double &r1,
999  const double &r2)
1000 {
1001  double b;
1002  double a;
1003  if(y == yp) {
1004  b = -y;
1005  a = 0.0;
1006  }
1007  else {
1008  if(x == xp) {
1009  b = -x;
1010  a = 0.0;
1011  }
1012  else {
1013  b = (xp * y - x * yp) / (yp - y);
1014  if(yp == 0.0) { a = -(b + x) / y; }
1015  else {
1016  a = -(b + xp) / yp;
1017  }
1018  }
1019  }
1020  double ae;
1021  double be;
1022  double ce;
1023  double da = r1 * r1;
1024  double db = r2 * r2;
1025  if(y == yp) {
1026  ae = 1.0 / da;
1027  be = -(2 * x) / da;
1028  ce = (x * x / da) - 1.0;
1029  }
1030  else {
1031  if(x == xp) {
1032  ae = 1.0 / db;
1033  be = -(2.0 * y) / db;
1034  ce = (y * y / db) - 1.0;
1035  }
1036  else {
1037  if(fabs(a) < 0.00001) {
1038  ae = 1.0 / db;
1039  be = -(2.0 * y) / db;
1040  ce = (y * y / db) - 1.0;
1041  }
1042  else {
1043  double a2 = a * a;
1044  ae = (1.0 / da) + (1.0 / (db * a2));
1045  be = (2.0 * y) / (db * a) + (2.0 * b) / (a2 * db) - ((2.0 * x) / da);
1046  ce = (x * x) / da + (b * b) / (db * a2) + (2.0 * b * y) / (a * db) +
1047  (y * y / db) - 1.0;
1048  }
1049  }
1050  }
1051  double rho = be * be - 4 * ae * ce;
1052  double x1, x2, y1, y2, propdist;
1053  if(rho < 0) { return 1; }
1054  else {
1055  x1 = -(be + sqrt(rho)) / (2.0 * ae);
1056  x2 = (-be + sqrt(rho)) / (2.0 * ae);
1057  if(y == yp) {
1058  y1 = -b;
1059  y2 = -b;
1060  }
1061  else {
1062  if(x == xp) {
1063  y1 = x1;
1064  y2 = x2;
1065  x1 = -b;
1066  x2 = -b;
1067  }
1068  else {
1069  if(fabs(a) < 0.00001) {
1070  y1 = x1;
1071  y2 = x2;
1072  x1 = -b;
1073  x2 = -b;
1074  }
1075  else {
1076  y1 = -(b + x1) / a;
1077  y2 = -(b + x2) / a;
1078  }
1079  }
1080  }
1081  if(x1 == x2) {
1082  propdist = (y1 - y) / (yp - y);
1083  if(propdist < 0.0) { propdist = (y2 - y) / (yp - y); }
1084  }
1085  else {
1086  if(xp != x) {
1087  propdist = (x1 - x) / (xp - x);
1088  if(propdist < 0.0) { propdist = (x2 - x) / (xp - x); }
1089  }
1090  else {
1091  if(yp != y) {
1092  propdist = (y1 - y) / (yp - y);
1093  if(propdist < 0.0) { propdist = (y2 - y) / (yp - y); }
1094  }
1095  else {
1096  propdist = 0.01;
1097  }
1098  }
1099  }
1100  *distance = propdist;
1101  return 0;
1102  }
1103 }
1104 
1105 void signedDistancesPointsEllipsePoint(std::vector<double> &distances,
1106  std::vector<double> &distancesE,
1107  std::vector<int> &isInYarn,
1108  std::vector<SPoint3> &closePts,
1109  const std::vector<SPoint3> &pts,
1110  const SPoint3 &p1, const SPoint3 &p2,
1111  const double radius)
1112 {
1113  distances.clear();
1114  distances.resize(pts.size());
1115  distancesE.clear();
1116  distancesE.resize(pts.size());
1117  isInYarn.clear();
1118  isInYarn.resize(pts.size());
1119  closePts.clear();
1120  closePts.resize(pts.size());
1121  double d;
1122  for(std::size_t i = 0; i < pts.size(); i++) {
1123  SPoint3 closePt;
1124  const SPoint3 &p = pts[i];
1125  signedDistancePointLine(p1, p2, p, d, closePt);
1126  closePts[i] = closePt;
1127  distances[i] = d;
1128  if(d <= radius) {
1129  isInYarn[i] = 1;
1130  distancesE[i] = radius - d;
1131  }
1132  else {
1133  isInYarn[i] = 0;
1134  distancesE[i] = d - radius;
1135  }
1136  }
1137 }
1138 
1140  std::vector<double> &distances, std::vector<double> &distancesE,
1141  std::vector<int> &isInYarn, std::vector<SPoint3> &closePts,
1142  const std::vector<SPoint3> &pts, const SPoint3 &p1, const SPoint3 &p2,
1143  const double maxA, const double minA, const double maxB, const double minB,
1144  const int typeLevelSet)
1145 {
1146  distances.clear();
1147  distances.resize(pts.size());
1148  distancesE.clear();
1149  distancesE.resize(pts.size());
1150  isInYarn.clear();
1151  isInYarn.resize(pts.size());
1152  closePts.clear();
1153  closePts.resize(pts.size());
1154  double d;
1155  for(std::size_t i = 0; i < pts.size(); i++) {
1156  SPoint3 closePt;
1157  const SPoint3 &p = pts[i];
1158  signedDistancePointLine(p1, p2, p, d, closePt);
1159 
1160  distances[i] = d;
1161  closePts[i] = closePt;
1162  int direction = 0;
1163  if(!(p.x() == closePt.x() && p.y() == closePt.y() &&
1164  p.z() == closePt.z())) {
1165  double xp, yp, x, y, otherp, other, propdist;
1166  if(typeLevelSet == 2) {
1167  if(p1.x() == p2.x()) {
1168  direction = 1;
1169  if(fabs(closePt.x() - 0.0) < 0.00000001) isInYarn[i] = 1;
1170  if(fabs(closePt.x() - 2.2) < 0.00000001) isInYarn[i] = 1;
1171  }
1172  else {
1173  if(p1.y() == p2.y()) {
1174  direction = 2;
1175  if(fabs(closePt.y() - 0.0) < 0.00000001) isInYarn[i] = 6;
1176  if(fabs(closePt.y() - 2.2) < 0.00000001) isInYarn[i] = 6;
1177  }
1178  else {
1179  printf("Error %lf %lf\n", closePt.x(), closePt.y());
1180  }
1181  }
1182  }
1183  if(typeLevelSet == 1) {
1184  if(p1.x() == p2.x()) {
1185  direction = 1;
1186  // if (fabs(closePt.x() - 0.0) < 0.00000001) isInYarn[i] = 1;
1187  if(fabs(closePt.x() - 2.2) < 0.00000001) isInYarn[i] = 4;
1188  // if (fabs(closePt.x() - 4.4) < 0.00000001) isInYarn[i] = 2;
1189  // if (fabs(closePt.x() - 6.6) < 0.00000001) isInYarn[i] = 5;
1190  // if (fabs(closePt.x() - 8.8) < 0.00000001) isInYarn[i] = 3;
1191  // if (fabs(closePt.x() - 11.0) < 0.00000001) isInYarn[i] = 1;
1192  }
1193  else {
1194  if(p1.y() == p2.y()) {
1195  direction = 2;
1196  // if (fabs(closePt.y() - 0.0) < 0.00000001) isInYarn[i] = 6;
1197  if(fabs(closePt.y() - 2.2) < 0.00000001) isInYarn[i] = 7;
1198  // if (fabs(closePt.y() - 4.4) < 0.00000001) isInYarn[i] = 8;
1199  // if (fabs(closePt.y() - 6.6) < 0.00000001) isInYarn[i] = 9;
1200  // if (fabs(closePt.y() - 8.8) < 0.00000001) isInYarn[i] = 10;
1201  // if (fabs(closePt.y() - 11.0) < 0.00000001) isInYarn[i] = 6;
1202  }
1203  else {
1204  printf("Error %lf %lf\n", closePt.x(), closePt.y());
1205  }
1206  }
1207  }
1208  if(typeLevelSet == 4) {
1209  if(p1.x() == p2.x()) {
1210  direction = 1;
1211  if(fabs(closePt.x() - 0.0) < 0.00000001 && closePt.z() <= 0.35)
1212  isInYarn[i] = 1;
1213  if(fabs(closePt.x() - 2.2) < 0.00000001 && closePt.z() <= 0.35)
1214  isInYarn[i] = 4;
1215  if(fabs(closePt.x() - 4.4) < 0.00000001 && closePt.z() <= 0.35)
1216  isInYarn[i] = 2;
1217  if(fabs(closePt.x() - 6.6) < 0.00000001 && closePt.z() <= 0.35)
1218  isInYarn[i] = 5;
1219  if(fabs(closePt.x() - 8.8) < 0.00000001 && closePt.z() <= 0.35)
1220  isInYarn[i] = 3;
1221  if(fabs(closePt.x() - 11.0) < 0.00000001 && closePt.z() <= 0.35)
1222  isInYarn[i] = 1;
1223  if(fabs(closePt.x() - 0.0) < 0.00000001 && closePt.z() > 0.35)
1224  isInYarn[i] = 11;
1225  if(fabs(closePt.x() - 2.2) < 0.00000001 && closePt.z() > 0.35)
1226  isInYarn[i] = 14;
1227  if(fabs(closePt.x() - 4.4) < 0.00000001 && closePt.z() > 0.35)
1228  isInYarn[i] = 12;
1229  if(fabs(closePt.x() - 6.6) < 0.00000001 && closePt.z() > 0.35)
1230  isInYarn[i] = 15;
1231  if(fabs(closePt.x() - 8.8) < 0.00000001 && closePt.z() > 0.35)
1232  isInYarn[i] = 13;
1233  if(fabs(closePt.x() - 11.0) < 0.00000001 && closePt.z() > 0.35)
1234  isInYarn[i] = 11;
1235  }
1236  else {
1237  if(p1.y() == p2.y()) {
1238  direction = 2;
1239  if(fabs(closePt.y() - 0.0) < 0.00000001 && closePt.z() <= 0.35)
1240  isInYarn[i] = 6;
1241  if(fabs(closePt.y() - 2.2) < 0.00000001 && closePt.z() <= 0.35)
1242  isInYarn[i] = 7;
1243  if(fabs(closePt.y() - 4.4) < 0.00000001 && closePt.z() <= 0.35)
1244  isInYarn[i] = 8;
1245  if(fabs(closePt.y() - 6.6) < 0.00000001 && closePt.z() <= 0.35)
1246  isInYarn[i] = 9;
1247  if(fabs(closePt.y() - 8.8) < 0.00000001 && closePt.z() <= 0.35)
1248  isInYarn[i] = 10;
1249  if(fabs(closePt.y() - 11.0) < 0.00000001 && closePt.z() <= 0.35)
1250  isInYarn[i] = 6;
1251  if(fabs(closePt.y() - 0.0) < 0.00000001 && closePt.z() > 0.35)
1252  isInYarn[i] = 16;
1253  if(fabs(closePt.y() - 2.2) < 0.00000001 && closePt.z() > 0.35)
1254  isInYarn[i] = 17;
1255  if(fabs(closePt.y() - 4.4) < 0.00000001 && closePt.z() > 0.35)
1256  isInYarn[i] = 18;
1257  if(fabs(closePt.y() - 6.6) < 0.00000001 && closePt.z() > 0.35)
1258  isInYarn[i] = 19;
1259  if(fabs(closePt.y() - 8.8) < 0.00000001 && closePt.z() > 0.35)
1260  isInYarn[i] = 20;
1261  if(fabs(closePt.y() - 11.0) < 0.00000001 && closePt.z() > 0.35)
1262  isInYarn[i] = 16;
1263  }
1264  else {
1265  printf("Error %lf %lf\n", closePt.x(), closePt.y());
1266  }
1267  }
1268  }
1269  if(typeLevelSet == 3) {
1270  direction = 3;
1271  isInYarn[i] = 1;
1272  }
1273  if(typeLevelSet == 5) {
1274  if(p1.x() == p2.x()) {
1275  direction = 1;
1276  if(fabs(closePt.x() - 0.0) < 0.00000001) isInYarn[i] = 1;
1277  if(fabs(closePt.x() - 3.225) < 0.00000001) isInYarn[i] = 2;
1278  if(fabs(closePt.x() - 6.45) < 0.00000001) isInYarn[i] = 3;
1279  if(fabs(closePt.x() - 9.675) < 0.00000001) isInYarn[i] = 4;
1280  if(fabs(closePt.x() - 12.9) < 0.00000001) isInYarn[i] = 1;
1281  }
1282  else {
1283  if(p1.y() == p2.y()) {
1284  direction = 2;
1285  if(fabs(closePt.y() - 0.0) < 0.00000001) isInYarn[i] = 5;
1286  if(fabs(closePt.y() - 1.665) < 0.00000001) isInYarn[i] = 6;
1287  if(fabs(closePt.y() - 3.33) < 0.00000001) isInYarn[i] = 7;
1288  if(fabs(closePt.y() - 4.995) < 0.00000001) isInYarn[i] = 8;
1289  if(fabs(closePt.y() - 6.66) < 0.00000001) isInYarn[i] = 5;
1290  }
1291  else {
1292  printf("Error %lf %lf\n", closePt.x(), closePt.y());
1293  }
1294  }
1295  }
1296  changeReferential(direction, p, closePt, p1, p2, &xp, &yp, &otherp, &x,
1297  &y, &other);
1298  int result = 1;
1299  if(fabs(other - otherp) > 0.01) { result = 1; }
1300  else {
1301  if(direction == 1) {
1302  result = computeDistanceRatio(y, yp, x, xp, &propdist, maxA, minA);
1303  }
1304  else {
1305  if(direction == 2) {
1306  result = computeDistanceRatio(y, yp, x, xp, &propdist, maxB, minB);
1307  }
1308  }
1309  }
1310  if(result == 1) {
1311  distancesE[i] = 1.e10;
1312  isInYarn[i] = 0;
1313  }
1314  else {
1315  if(propdist < 1.0) {
1316  isInYarn[i] = 0;
1317  distancesE[i] = (1.0 / propdist) - 1.0;
1318  }
1319  else {
1320  distancesE[i] = (1.0 - (1.0 / propdist));
1321  }
1322  }
1323  }
1324  else {
1325  isInYarn[i] = 0;
1326  distancesE[i] = 1000000.0;
1327  }
1328  }
1329 }
1330 
1331 int intersection_segments(const SPoint2 &p1, const SPoint2 &p2,
1332  const SPoint2 &q1, const SPoint2 &q2, double x[2])
1333 {
1334  double xp_max = std::max(p1.x(), p2.x());
1335  double yp_max = std::max(p1.y(), p2.y());
1336  double xq_max = std::max(q1.x(), q2.x());
1337  double yq_max = std::max(q1.y(), q2.y());
1338 
1339  double xp_min = std::min(p1.x(), p2.x());
1340  double yp_min = std::min(p1.y(), p2.y());
1341  double xq_min = std::min(q1.x(), q2.x());
1342  double yq_min = std::min(q1.y(), q2.y());
1343  if(yq_min > yp_max || xq_min > xp_max || yq_max < yp_min || xq_max < xp_min) {
1344  return 0;
1345  }
1346  else {
1347  double A[2][2];
1348  A[0][0] = p2.x() - p1.x();
1349  A[0][1] = q1.x() - q2.x();
1350  A[1][0] = p2.y() - p1.y();
1351  A[1][1] = q1.y() - q2.y();
1352  double b[2] = {q1.x() - p1.x(), q1.y() - p1.y()};
1353  sys2x2(A, b, x);
1354  return (x[0] >= 0.0 && x[0] <= 1. && x[1] >= 0.0 && x[1] <= 1.);
1355  }
1356 }
1358 int intersection_segments(const SPoint3 &p1, const SPoint3 &p2,
1359  const SPoint3 &q1, const SPoint3 &q2, double x[2])
1360 {
1361  SVector3 v1(p2, p1), v2(q2, q1);
1362  double n1 = v1.norm();
1363  double n2 = v2.norm();
1364  double EPS = 1.e-10 * std::max(n1, n2);
1365  double A[2][2];
1366  A[0][0] = p2.x() - p1.x();
1367  A[0][1] = q1.x() - q2.x();
1368  A[1][0] = p2.y() - p1.y();
1369  A[1][1] = q1.y() - q2.y();
1370  double a[2] = {q1.x() - p1.x(), q1.y() - p1.y()};
1371  double B[2][2];
1372  B[0][0] = p2.z() - p1.z();
1373  B[0][1] = q1.z() - q2.z();
1374  B[1][0] = p2.y() - p1.y();
1375  B[1][1] = q1.y() - q2.y();
1376  double b[2] = {q1.z() - p1.z(), q1.y() - p1.y()};
1377  double C[2][2];
1378  C[0][0] = p2.z() - p1.z();
1379  C[0][1] = q1.z() - q2.z();
1380  C[1][0] = p2.x() - p1.x();
1381  C[1][1] = q1.x() - q2.x();
1382  double c[2] = {q1.z() - p1.z(), q1.x() - p1.x()};
1383  double detA = fabs(det2x2(A));
1384  double detB = fabs(det2x2(B));
1385  double detC = fabs(det2x2(C));
1386  // printf("%12.5E %12.5E %12.5E\n",detA,detB,detC);
1387  if(detA > detB && detA > detC)
1388  sys2x2(A, a, x);
1389  else if(detB > detA && detB > detC)
1390  sys2x2(B, b, x);
1391  else
1392  sys2x2(C, c, x);
1393  if(x[0] >= 0.0 && x[0] <= 1. && x[1] >= 0.0 && x[1] <= 1.) {
1394  SPoint3 x1(p1.x() * (1. - x[0]) + p2.x() * x[0],
1395  p1.y() * (1. - x[0]) + p2.y() * x[0],
1396  p1.z() * (1. - x[0]) + p2.z() * x[0]);
1397  SPoint3 x2(q1.x() * (1. - x[1]) + q2.x() * x[1],
1398  q1.y() * (1. - x[1]) + q2.y() * x[1],
1399  q1.z() * (1. - x[1]) + q2.z() * x[1]);
1400 
1401  SVector3 d(x2, x1);
1402  double nd = norm(d);
1403  if(nd > EPS) {
1404  x[0] = x[1] = 1.e22;
1405  return false;
1406  }
1407  return true;
1408  }
1409  return false;
1410 }
1411 
1412 void fillMeanPlane(double res[4], double t1[3], double t2[3],
1413  mean_plane &meanPlane)
1414 {
1415  for(int i = 0; i < 3; i++) meanPlane.plan[0][i] = t1[i];
1416  for(int i = 0; i < 3; i++) meanPlane.plan[1][i] = t2[i];
1417  for(int i = 0; i < 3; i++) meanPlane.plan[2][i] = res[i];
1418 
1419  meanPlane.a = res[0];
1420  meanPlane.b = res[1];
1421  meanPlane.c = res[2];
1422  meanPlane.d = res[3];
1423 
1424  meanPlane.x = meanPlane.y = meanPlane.z = 0.;
1425  if(fabs(meanPlane.a) >= fabs(meanPlane.b) &&
1426  fabs(meanPlane.a) >= fabs(meanPlane.c)) {
1427  meanPlane.x = meanPlane.d / meanPlane.a;
1428  }
1429  else if(fabs(meanPlane.b) >= fabs(meanPlane.a) &&
1430  fabs(meanPlane.b) >= fabs(meanPlane.c)) {
1431  meanPlane.y = meanPlane.d / meanPlane.b;
1432  }
1433  else {
1434  meanPlane.z = meanPlane.d / meanPlane.c;
1435  }
1436 }
1437 
1438 void computeMeanPlaneSimple(const std::vector<SPoint3> &points,
1439  mean_plane &meanPlane)
1440 {
1441  double xm = 0., ym = 0., zm = 0.;
1442  int ndata = points.size();
1443  int na = 3;
1444  for(int i = 0; i < ndata; i++) {
1445  xm += points[i].x();
1446  ym += points[i].y();
1447  zm += points[i].z();
1448  }
1449  xm /= (double)ndata;
1450  ym /= (double)ndata;
1451  zm /= (double)ndata;
1452 
1453  fullMatrix<double> U(ndata, na), V(na, na);
1454  fullVector<double> sigma(na);
1455  for(int i = 0; i < ndata; i++) {
1456  U(i, 0) = points[i].x() - xm;
1457  U(i, 1) = points[i].y() - ym;
1458  U(i, 2) = points[i].z() - zm;
1459  }
1460  U.svd(V, sigma);
1461  double res[4], svd[3];
1462  svd[0] = sigma(0);
1463  svd[1] = sigma(1);
1464  svd[2] = sigma(2);
1465  int min;
1466  if(fabs(svd[0]) < fabs(svd[1]) && fabs(svd[0]) < fabs(svd[2]))
1467  min = 0;
1468  else if(fabs(svd[1]) < fabs(svd[0]) && fabs(svd[1]) < fabs(svd[2]))
1469  min = 1;
1470  else
1471  min = 2;
1472  res[0] = V(0, min);
1473  res[1] = V(1, min);
1474  res[2] = V(2, min);
1475  norme(res);
1476 
1477  double ex[3], t1[3], t2[3];
1478 
1479  ex[0] = ex[1] = ex[2] = 0.0;
1480  if(res[0] == 0.)
1481  ex[0] = 1.0;
1482  else if(res[1] == 0.)
1483  ex[1] = 1.0;
1484  else
1485  ex[2] = 1.0;
1486 
1487  prodve(res, ex, t1);
1488  norme(t1);
1489  prodve(t1, res, t2);
1490  norme(t2);
1491 
1492  res[3] = (xm * res[0] + ym * res[1] + zm * res[2]);
1493 
1494  fillMeanPlane(res, t1, t2, meanPlane);
1495 }
1496 
1497 void projectPointToPlane(const SPoint3 &pt, SPoint3 &ptProj,
1498  const mean_plane &meanPlane)
1499 {
1500  double u = pt.x();
1501  double v = pt.y();
1502  double w = pt.z();
1503  double a = meanPlane.a;
1504  double b = meanPlane.b;
1505  double c = meanPlane.c;
1506  double d = meanPlane.d;
1507  double t0 = -(a * u + b * v + c * w + d) / (a * a + b * b + c * c);
1508 
1509  ptProj[0] = u + a * t0;
1510  ptProj[1] = v + b * t0;
1511  ptProj[2] = w + c * t0;
1512 }
1513 
1514 void projectPointsToPlane(const std::vector<SPoint3> &pts,
1515  std::vector<SPoint3> &ptsProj,
1516  const mean_plane &meanPlane)
1517 {
1518  ptsProj.resize(pts.size());
1519  for(std::size_t i = 0; i < pts.size(); i++) {
1520  projectPointToPlane(pts[i], ptsProj[i], meanPlane);
1521  }
1522 }
1523 
1524 void transformPointsIntoOrthoBasis(const std::vector<SPoint3> &ptsProj,
1525  std::vector<SPoint3> &pointsUV,
1526  const SPoint3 &ptCG,
1527  const mean_plane &meanPlane)
1528 {
1529  pointsUV.resize(ptsProj.size());
1530  SVector3 normal(meanPlane.a, meanPlane.b, meanPlane.c);
1531  SVector3 tangent, binormal;
1532  buildOrthoBasis(normal, tangent, binormal);
1533 
1534  for(std::size_t i = 0; i < ptsProj.size(); i++) {
1535  SVector3 pp(ptsProj[i][0] - ptCG[0], ptsProj[i][1] - ptCG[1],
1536  ptsProj[i][2] - ptCG[2]);
1537  pointsUV[i][0] = dot(pp, tangent);
1538  pointsUV[i][1] = dot(pp, binormal);
1539  pointsUV[i][2] = dot(pp, normal);
1540  }
1541 }
1542 
1544  void *data)
1545 {
1546  double *param = (double *)data;
1547  double x0 = param[0], x1 = param[1], y0 = param[2], y1 = param[3],
1548  ys = param[4];
1549  res(0) = (ys - 1 / x(0)) + 1 / x(0) * cosh(x(0) * (x0 - x(1))) - y0;
1550  res(1) = (ys - 1 / x(0)) + 1 / x(0) * cosh(x(0) * (x1 - x(1))) - y1;
1551  return true;
1552 }
1553 
1554 bool catenary(double x0, double x1, double y0, double y1, double ys, int N,
1555  double *yp)
1556 {
1557  // In the z=0 plane, catenary equation is y(x) = a + 1/b cosh(b(x-c))
1558  //
1559  // Three parameters a, b, c determined by imposing
1560  // - left point: y0 = y(x0) = a + 1/b cosh(b(x0-c))
1561  // - right point: y1 = y(x1) = a + 1/b cosh(b(x1-c))
1562  // - lowest point: ys = y(c) = a + 1/b , i.e. a = ys - 1/b
1563  //
1564  // Thus solve syst of 2 nl equations with 2 unknowns b and c:
1565  //
1566  // ys - 1/b + 1/b cosh(b(x0-c)) - y0 = 0
1567  // ys - 1/b + 1/b cosh(b(x1-c)) - y1 = 0
1568  double param[5] = {x0, x1, y0, y1, ys};
1569  fullVector<double> x(2);
1570  bool success = false, physical = true;
1571  double tolx = 1e-6 * fabs(x1 - x0);
1572  double toly = 1e-6 * fabs(std::max(y0, y1) - ys);
1573  if(x0 != x1) {
1574  x(0) = 1. / (x1 - x0);
1575  x(1) = (x0 + x1) / 2.;
1576  success = newton_fd(catenary_fct, x, param, 1, tolx);
1577  }
1578  if(success) {
1579  double a = ys - 1 / x(0);
1580  for(int i = 0; i < N; i++) {
1581  double r = x0 + (i + 1) * (x1 - x0) / (N + 1);
1582  yp[i] = a + 1 / x(0) * cosh(x(0) * (r - x(1)));
1583  if(yp[i] > std::max(y0, y1) + toly || yp[i] < ys - toly) {
1584  physical = false;
1585  break;
1586  }
1587  }
1588  }
1589  if(physical) return true;
1590 
1591  // could not solve: return linear interpolation
1592  for(int i = 0; i < N; i++) { yp[i] = y0 + (i + 1) * (y1 - y0) / (N + 1); }
1593  return false;
1594 }
crossprod
SVector3 crossprod(const SVector3 &a, const SVector3 &b)
Definition: SVector3.h:150
planarQuad_xyz2xy
void planarQuad_xyz2xy(double *x, double *y, double *z, double *xn, double *yn)
Definition: Numeric.cpp:376
inv3x3
double inv3x3(double mat[3][3], double inv[3][3])
Definition: Numeric.cpp:232
mean_plane::d
double d
Definition: Numeric.h:35
dot
double dot(const SVector3 &a, const SMetric3 &m, const SVector3 &b)
Definition: STensor3.h:71
mean_plane::c
double c
Definition: Numeric.h:35
circumCenterXY
void circumCenterXY(double *p1, double *p2, double *p3, double *res)
Definition: Numeric.cpp:317
distance
double distance(MVertex *v1, MVertex *v2)
Definition: MVertex.h:245
mean_plane::x
double x
Definition: Numeric.h:36
fullVector< double >
newton_fd
bool newton_fd(bool(*func)(fullVector< double > &, fullVector< double > &, void *), fullVector< double > &x, void *data, double relax, double tolx)
Definition: Numeric.cpp:685
SPoint2
Definition: SPoint2.h:12
float2char
char float2char(float f)
Definition: Numeric.cpp:452
Msg::Debug
static void Debug(const char *fmt,...)
Definition: GmshMessage.cpp:752
c
static double c(int i, int j, fullMatrix< double > &CA, const std::vector< SPoint3 > &P, const std::vector< SPoint3 > &Q)
Definition: discreteFrechetDistance.cpp:15
buildOrthoBasis
void buildOrthoBasis(SVector3 &normal, SVector3 &tangent, SVector3 &binormal)
Definition: SVector3.h:257
Msg::Error
static void Error(const char *fmt,...)
Definition: GmshMessage.cpp:482
fullVector::norm
scalar norm() const
det3x3
double det3x3(double mat[3][3])
Definition: Numeric.cpp:126
SPoint3
Definition: SPoint3.h:14
LegendrePolynomials::f
void f(int n, double u, double *val)
Definition: orthogonalBasis.cpp:77
signedDistancesPointsLine
void signedDistancesPointsLine(std::vector< double > &distances, std::vector< SPoint3 > &closePts, const std::vector< SPoint3 > &pts, const SPoint3 &p1, const SPoint3 &p2)
Definition: Numeric.cpp:923
fullVector::size
int size() const
Definition: fullMatrix.h:69
SVector3
Definition: SVector3.h:16
char2float
float char2float(char c)
Definition: Numeric.cpp:464
projectPointToPlane
void projectPointToPlane(const SPoint3 &pt, SPoint3 &ptProj, const mean_plane &meanPlane)
Definition: Numeric.cpp:1497
mean_plane::b
double b
Definition: Numeric.h:35
invert_singular_matrix3x3
void invert_singular_matrix3x3(double MM[3][3], double II[3][3])
Definition: Numeric.cpp:650
GmshMessage.h
mean_plane
Definition: Numeric.h:33
sys3x3
int sys3x3(double mat[3][3], double b[3], double res[3], double *det)
Definition: Numeric.cpp:147
signedDistancesPointsTriangle
void signedDistancesPointsTriangle(std::vector< double > &distances, std::vector< SPoint3 > &closePts, const std::vector< SPoint3 > &pts, const SPoint3 &p1, const SPoint3 &p2, const SPoint3 &p3)
Definition: Numeric.cpp:822
norm3
double norm3(double a[3])
Definition: Numeric.h:119
SPoint3::x
double x(void) const
Definition: SPoint3.h:125
fullMatrix< double >
signedDistancesPointsEllipseLine
void signedDistancesPointsEllipseLine(std::vector< double > &distances, std::vector< double > &distancesE, std::vector< int > &isInYarn, std::vector< SPoint3 > &closePts, const std::vector< SPoint3 > &pts, const SPoint3 &p1, const SPoint3 &p2, const double maxA, const double minA, const double maxB, const double minB, const int typeLevelSet)
Definition: Numeric.cpp:1139
myatan2
double myatan2(double a, double b)
Definition: Numeric.cpp:12
intersection_segments
int intersection_segments(const SPoint2 &p1, const SPoint2 &p2, const SPoint2 &q1, const SPoint2 &q2, double x[2])
Definition: Numeric.cpp:1331
det2x2
double det2x2(double mat[2][2])
Definition: Numeric.cpp:199
FindCubicRoots
void FindCubicRoots(const double coef[4], double real[3], double imag[3])
Definition: Numeric.cpp:574
signedDistancePointTriangle
void signedDistancePointTriangle(const SPoint3 &p1, const SPoint3 &p2, const SPoint3 &p3, const SPoint3 &p, double &d, SPoint3 &closePt)
Definition: Numeric.cpp:758
eigenvalue2x2
void eigenvalue2x2(double mat[2][2], double v[2])
Definition: Numeric.cpp:538
fillMeanPlane
void fillMeanPlane(double res[4], double t1[3], double t2[3], mean_plane &meanPlane)
Definition: Numeric.cpp:1412
angle_plan
double angle_plan(double v[3], double p1[3], double p2[3], double n[3])
Definition: Numeric.cpp:267
signedDistancesPointsEllipsePoint
void signedDistancesPointsEllipsePoint(std::vector< double > &distances, std::vector< double > &distancesE, std::vector< int > &isInYarn, std::vector< SPoint3 > &closePts, const std::vector< SPoint3 > &pts, const SPoint3 &p1, const SPoint3 &p2, const double radius)
Definition: Numeric.cpp:1105
matmat
void matmat(double mat1[3][3], double mat2[3][3], double res[3][3])
Definition: Numeric.cpp:54
a1
const double a1
Definition: GaussQuadratureHex.cpp:10
norm
void norm(const double *vec, double *norm)
Definition: gmshLevelset.cpp:202
Numeric.h
computeInnerRadiusForQuad
double computeInnerRadiusForQuad(double *x, double *y, int i)
Definition: Numeric.cpp:407
catenary_fct
static bool catenary_fct(fullVector< double > &x, fullVector< double > &res, void *data)
Definition: Numeric.cpp:1543
triangle_area2d
double triangle_area2d(double p0[2], double p1[2], double p2[2])
Definition: Numeric.cpp:309
mean_plane::a
double a
Definition: Numeric.h:35
mean_plane::z
double z
Definition: Numeric.h:36
SVector3::norm
double norm() const
Definition: SVector3.h:33
signedDistancePointLine
void signedDistancePointLine(const SPoint3 &p1, const SPoint3 &p2, const SPoint3 &p, double &d, SPoint3 &closePt)
Definition: Numeric.cpp:908
SPoint3::y
double y(void) const
Definition: SPoint3.h:127
prosca
double prosca(double const a[3], double const b[3])
Definition: Numeric.h:112
eigsort
void eigsort(double d[3])
Definition: Numeric.cpp:634
catenary
bool catenary(double x0, double x1, double y0, double y1, double ys, int N, double *yp)
Definition: Numeric.cpp:1554
computeDistanceRatio
int computeDistanceRatio(const double &y, const double &yp, const double &x, const double &xp, double *distance, const double &r1, const double &r2)
Definition: Numeric.cpp:997
fullMatrix::eig
bool eig(fullVector< double > &eigenValReal, fullVector< double > &eigenValImag, fullMatrix< scalar > &leftEigenVect, fullMatrix< scalar > &rightEigenVect, bool sortRealPart=false)
Definition: fullMatrix.h:727
prodve
void prodve(double a[3], double b[3], double c[3])
Definition: Numeric.h:105
SVector3::x
double x(void) const
Definition: SVector3.h:30
a2
const double a2
Definition: GaussQuadratureHex.cpp:12
norm2
double norm2(double a[3][3])
Definition: Numeric.cpp:38
triangle_area
double triangle_area(double p0[3], double p1[3], double p2[3])
Definition: Numeric.cpp:293
S
#define S
Definition: DefaultOptions.h:21
SPoint3::distance
double distance(const SPoint3 &p) const
Definition: SPoint3.h:176
trace2
double trace2(double mat[3][3])
Definition: Numeric.cpp:135
b1
const double b1
Definition: GaussQuadratureHex.cpp:14
mean_plane::y
double y
Definition: Numeric.h:36
inv2x2
double inv2x2(double mat[2][2], double inv[2][2])
Definition: Numeric.cpp:214
matvec
void matvec(double mat[3][3], double vec[3], double res[3])
Definition: Numeric.cpp:47
SPoint2::x
double x(void) const
Definition: SPoint2.h:86
projectPointsToPlane
void projectPointsToPlane(const std::vector< SPoint3 > &pts, std::vector< SPoint3 > &ptsProj, const mean_plane &meanPlane)
Definition: Numeric.cpp:1514
sys2x2
int sys2x2(double mat[2][2], double b[2], double res[2])
Definition: Numeric.cpp:101
mean_plane::plan
double plan[3][3]
Definition: Numeric.h:34
ComputeVonMises
double ComputeVonMises(double *V)
Definition: Numeric.cpp:496
SVector3::y
double y(void) const
Definition: SVector3.h:31
z
const double z
Definition: GaussQuadratureQuad.cpp:56
circumCenterXYZ
void circumCenterXYZ(double *p1, double *p2, double *p3, double *res, double *uv)
Definition: Numeric.cpp:342
det2x3
double det2x3(double mat[2][3])
Definition: Numeric.cpp:204
myacos
double myacos(double a)
Definition: Numeric.cpp:28
SVector3::z
double z(void) const
Definition: SVector3.h:32
direction
static void direction(Vertex *v1, Vertex *v2, double d[3])
Definition: Geo.cpp:218
normal3points
void normal3points(double x0, double y0, double z0, double x1, double y1, double z1, double x2, double y2, double z2, double n[3])
Definition: Numeric.cpp:76
eigenvalue
void eigenvalue(double mat[3][3], double v[3])
Definition: Numeric.cpp:550
fullMatrix::luSolve
bool luSolve(const fullVector< scalar > &rhs, fullVector< scalar > &result)
Definition: fullMatrix.h:603
SPoint3::z
double z(void) const
Definition: SPoint3.h:129
gradSimplex
void gradSimplex(double *x, double *y, double *z, double *v, double *grad)
Definition: Numeric.cpp:475
normal2points
void normal2points(double x0, double y0, double z0, double x1, double y1, double z1, double n[3])
Definition: Numeric.cpp:85
myasin
double myasin(double a)
Definition: Numeric.cpp:18
angle_02pi
double angle_02pi(double A3)
Definition: Numeric.cpp:255
norme
double norme(double a[3])
Definition: Numeric.h:123
changeReferential
void changeReferential(const int direction, const SPoint3 &p, const SPoint3 &closePt, const SPoint3 &p1, const SPoint3 &p2, double *xp, double *yp, double *otherp, double *x, double *y, double *other)
Definition: Numeric.cpp:942
transformPointsIntoOrthoBasis
void transformPointsIntoOrthoBasis(const std::vector< SPoint3 > &ptsProj, std::vector< SPoint3 > &pointsUV, const SPoint3 &ptCG, const mean_plane &meanPlane)
Definition: Numeric.cpp:1524
c1
const double c1
Definition: GaussQuadratureHex.cpp:16
trace3x3
double trace3x3(double mat[3][3])
Definition: Numeric.cpp:133
fullMatrix.h
sys3x3_with_tol
int sys3x3_with_tol(double mat[3][3], double b[3], double res[3], double *det)
Definition: Numeric.cpp:177
computeMeanPlaneSimple
void computeMeanPlaneSimple(const std::vector< SPoint3 > &points, mean_plane &meanPlane)
Definition: Numeric.cpp:1438
ComputeScalarRep
double ComputeScalarRep(int numComp, double *val, int tensorRep)
Definition: Numeric.cpp:506
SVector3::normalize
double normalize()
Definition: SVector3.h:38
SPoint2::y
double y(void) const
Definition: SPoint2.h:88