Milena (Olena)  User documentation 2.0a Id
 All Classes Namespaces Functions Variables Typedefs Enumerator Groups Pages
misc.hh
1 // Copyright (C) 2008, 2009 EPITA Research and Development Laboratory
2 // (LRDE)
3 //
4 // This file is part of Olena.
5 //
6 // Olena is free software: you can redistribute it and/or modify it under
7 // the terms of the GNU General Public License as published by the Free
8 // Software Foundation, version 2 of the License.
9 //
10 // Olena is distributed in the hope that it will be useful,
11 // but WITHOUT ANY WARRANTY; without even the implied warranty of
12 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13 // General Public License for more details.
14 //
15 // You should have received a copy of the GNU General Public License
16 // along with Olena. If not, see <http://www.gnu.org/licenses/>.
17 //
18 // As a special exception, you may use this file as part of a free
19 // software project without restriction. Specifically, if other files
20 // instantiate templates or use macros or inline functions from this
21 // file, or you compile this file and link it with other files to produce
22 // an executable, this file does not by itself cause the resulting
23 // executable to be covered by the GNU General Public License. This
24 // exception does not however invalidate any other reasons why the
25 // executable file might be covered by the GNU General Public License.
26 
27 // FIXME: Split this file.
28 // FIXME: Address license-related issues?
29 
30 #ifndef MILENA_APPS_MESH_SEGM_SKEL_MISC_HH
31 # define MILENA_APPS_MESH_SEGM_SKEL_MISC_HH
32 
33 # include <algorithm> // For std::swap.
34 # include <utility> // For std::pair.
35 
36 # include <mln/algebra/vec.hh>
37 # include <mln/algebra/mat.hh>
38 
39 # include <mln/norm/l2.hh>
40 
41 # include <mln/data/fill.hh>
42 
43 
50 // FIXME: We should turn `likely' and `unlikely' into functions.
51 #ifndef likely
52 # if !defined(__GNUC__) || (__GNUC__ == 2 && __GNUC_MINOR__ < 96)
53 # define likely(x) (x)
54 # define unlikely(x) (x)
55 # else
56 # define likely(x) (__builtin_expect((x), 1))
57 # define unlikely(x) (__builtin_expect((x), 0))
58 # endif
59 #endif
60 
61 namespace mln
62 {
63 
64  namespace algebra
65  {
66 
67  // From Trimesh.
68  // FIXME: Doc (input, return value, etc.).
69  // FIXME: Idea: Add reference to Numerical Recipes.
70 
76  template <unsigned N, typename T>
77  inline
78  bool
79  ldlt_decomp(mat<N, N, T>& A, vec<N, T>& rdiag)
80  {
81  vec<N - 1, T> v;
82  for (unsigned i = 0; i < N; ++i)
83  {
84  for (unsigned k = 0; k < i; ++k)
85  v[k] = A(i, k) * rdiag[k];
86  for (unsigned j = i; j < N; ++j)
87  {
88  T sum = A(i, j);
89  for (unsigned k = 0; k < i; k++)
90  sum -= v[k] * A(j, k);
91  if (i == j)
92  {
93  if (unlikely(sum <= T(0)))
94  return false;
95  rdiag[i] = T(1) / sum;
96  }
97  else
98  A(j, i) = sum;
99  }
100  }
101  return true;
102  }
103 
104  // From Trimesh.
105  // FIXME: Doc (input, return value, etc.).
106  // FIXME: Idea: Add reference to Numerical Recipes.
107 
109  template <unsigned N, typename T>
110  inline
111  void
112  ldlt_solve(const mat<N, N, T>& A, const vec<N, T>& rdiag,
113  const vec<N, T>& B, vec<N, T>& x)
114  {
115  for (unsigned i = 0; i < N; ++i)
116  {
117  T sum = B[i];
118  for (unsigned k = 0; k < i; ++k)
119  sum -= A(i, k) * x[k];
120  x[i] = sum * rdiag[i];
121  }
122  for (int i = N - 1; i >= 0; --i)
123  {
124  T sum = 0;
125  for (unsigned k = i + 1; k < N; ++k)
126  sum += A(k, i) * x[k];
127  x[i] -= sum * rdiag[i];
128  }
129  }
130 
131  } // end of namespace mln::algebra
132 
133 } // end of namespace mln
134 
135 
136 // FIXME: To be moved elsewhere (in mln/geom)?
137 
138 // Inspired from Trimesh.
139 
140 namespace mln
141 {
142 
143  namespace geom
144  {
145 
146  // FIXME: More doc (see Trimesh's).
147 
148  /* FIXME: We should restrict attached data to vertices in return
149  value. */
159  inline
162  {
163  // Shortcuts.
164  static const unsigned D = 2;
165  typedef space_2complex_geometry G;
166  typedef algebra::vec<3, float> vec3f;
167  typedef complex_image< D, G, vec3f > normal_t;
168 
169  normal_t normal(mesh);
170  data::fill(normal, literal::zero);
171 
173  // A neighborhood where neighbors are the set of 0-faces
174  // transitively adjacent to the reference point.
175  typedef mln::complex_m_face_neighborhood<D, G> adj_vertices_nbh_t;
176  adj_vertices_nbh_t adj_vertices_nbh;
177  mln_niter_(adj_vertices_nbh_t) adj_v(adj_vertices_nbh, f);
178  /* FIXME: We should be able to pas this value (m) either at the
179  construction of the neighborhood or at the construction of
180  the iterator. Alas, we can't inherit ctors (yet), so we
181  can't rely on a simple mixin approach. */
182  adj_v.iter().set_m(0);
183 
184  // Iterate on 2-faces (triangles).
185  for_all(f)
186  {
188  std::vector<mln_psite_(normal_t)> p;
189  p.reserve(3);
190  for_all(adj_v)
191  p.push_back(adj_v);
193  mln_assertion(p.size() == 3);
194 
195  /* FIXME: Is the explicit conversion to_site() really needed?
196  If not, we could provide a local shortcut like
197 
198  vec(p[0]) // Instead of p[0].to_site().front().to_vec()
199 
200  to shorten these lines. */
201  // FIXME: Factor as well?
202  vec3f a =
203  p[0].to_site().front().to_vec() - p[1].to_site().front().to_vec();
204  vec3f b =
205  p[1].to_site().front().to_vec() - p[2].to_site().front().to_vec();
206  vec3f c =
207  p[2].to_site().front().to_vec() - p[0].to_site().front().to_vec();
208 
209  // FIXME: Factor as well?
210  float l2a = norm::sqr_l2(a);
211  float l2b = norm::sqr_l2(b);
212  float l2c = norm::sqr_l2(c);
213  vec3f face_normal = algebra::vprod(a, b);
214 
215  normal(p[0]) += face_normal * (1.0f / (l2a * l2c));
216  normal(p[1]) += face_normal * (1.0f / (l2b * l2a));
217  normal(p[2]) += face_normal * (1.0f / (l2c * l2b));
218  }
219 
220  // Normalize normals. We don't need to iterate on all faces, just
221  // 0-faces.
223  for_all(v)
224  normal(v).normalize();
225 
226  return normal;
227  }
228 
229 
230  /* FIXME: We should restrict attached data to vertices in return
231  value. */
244  inline
245  std::pair<
248  >
250  {
251  // Shortcuts.
252  static const unsigned D = 2;
253  typedef space_2complex_geometry G;
254  typedef algebra::vec<3, float> vec3f;
255 
256  // Hold data for 2-faces only.
257  typedef complex_image< D, G, vec3f > corner_area_t;
258  // Hold data for 0-faces only.
259  typedef complex_image< D, G, float > point_area_t;
260  // Packed output.
261  typedef std::pair<corner_area_t, point_area_t> output_t;
262 
263  // Initialize output.
264  output_t output(mesh, mesh);
265  corner_area_t& corner_area = output.first;
266  point_area_t& point_area = output.second;
267  data::fill(corner_area, literal::zero);
268  data::fill(point_area, literal::zero);
269 
271  // A neighborhood where neighbors are the set of 0-faces
272  // transitively adjacent to the reference point.
273  typedef mln::complex_m_face_neighborhood<D, G> adj_vertices_nbh_t;
274  adj_vertices_nbh_t adj_vertices_nbh;
275  mln_niter_(adj_vertices_nbh_t) adj_v(adj_vertices_nbh, f);
276  /* FIXME: We should be able to pas this value (m) either at the
277  construction of the neighborhood or at the construction of
278  the iterator. Alas, we can't inherit ctors (yet), so we
279  can't rely on a simple mixin approach. */
280  adj_v.iter().set_m(0);
281 
282  for_all(f)
283  {
285  std::vector<mln_psite_(corner_area_t)> p;
286  p.reserve(3);
287  for_all(adj_v)
288  p.push_back(adj_v);
290  mln_assertion(p.size() == 3);
291 
292  // (Opposite) edge vectors.
293  algebra::vec<3, vec3f> e;
294  // FIXME: See above.
295  e.set
296  (p[2].to_site().front().to_vec() - p[1].to_site().front().to_vec(),
297  p[0].to_site().front().to_vec() - p[2].to_site().front().to_vec(),
298  p[1].to_site().front().to_vec() - p[0].to_site().front().to_vec());
299 
300  // Compute corner weights.
301  float area = 0.5f * norm::l2(algebra::vprod(e[0], e[1]));
302  vec3f sqr_norm;
303  sqr_norm.set(norm::sqr_l2(e[0]),
304  norm::sqr_l2(e[1]),
305  norm::sqr_l2(e[2]));
306  vec3f edge_weight;
307  edge_weight.set
308  (sqr_norm[0] * (sqr_norm[1] + sqr_norm[2] - sqr_norm[0]),
309  sqr_norm[1] * (sqr_norm[2] + sqr_norm[0] - sqr_norm[1]),
310  sqr_norm[2] * (sqr_norm[0] + sqr_norm[1] - sqr_norm[2]));
311 
312  if (edge_weight[0] <= 0.0f)
313  {
314  corner_area(f)[1] = -0.25f * sqr_norm[2] * area / (e[0] * e[2]);
315  corner_area(f)[2] = -0.25f * sqr_norm[1] * area / (e[0] * e[1]);
316  corner_area(f)[0] = area - corner_area(f)[1] - corner_area(f)[2];
317  }
318  else if (edge_weight[1] <= 0.0f)
319  {
320  corner_area(f)[2] = -0.25f * sqr_norm[0] * area / (e[1] * e[0]);
321  corner_area(f)[0] = -0.25f * sqr_norm[2] * area / (e[1] * e[2]);
322  corner_area(f)[1] = area - corner_area(f)[2] - corner_area(f)[0];
323  }
324  else if (edge_weight[2] <= 0.0f)
325  {
326  corner_area(f)[0] = -0.25f * sqr_norm[1] * area / (e[2] * e[1]);
327  corner_area(f)[1] = -0.25f * sqr_norm[0] * area / (e[2] * e[0]);
328  corner_area(f)[2] = area - corner_area(f)[0] - corner_area(f)[1];
329  }
330  else
331  {
332  float ewscale =
333  0.5f * area / (edge_weight[0] + edge_weight[1] + edge_weight[2]);
334  for (int i = 0; i < 3; ++i)
335  corner_area(f)[i] =
336  ewscale * (edge_weight[(i+1)%3] + edge_weight[(i+2)%3]);
337  }
338 
339  for (int i = 0; i < 3; ++i)
340  point_area(p[i]) += corner_area(f)[i];
341  }
342 
343  return output;
344  }
345 
346 
347  namespace internal
348  {
349 
365  static inline unsigned next(unsigned i) { return (i + 1) % 3; }
366  static inline unsigned prev(unsigned i) { return (i - 1) % 3; }
370 
371  inline
372  void
373  rot_coord_sys(const algebra::vec<3, float> &old_u,
374  const algebra::vec<3, float> &old_v,
375  const algebra::vec<3, float> &new_norm,
376  algebra::vec<3, float> &new_u,
377  algebra::vec<3, float> &new_v)
378  {
379  new_u = old_u;
380  new_v = old_v;
381  algebra::vec<3, float> old_norm = vprod(old_u, old_v);
382  float ndot = old_norm * new_norm;
383  if (unlikely(ndot <= -1.0f))
384  {
385  new_u = -new_u;
386  new_v = -new_v;
387  return;
388  }
389  algebra::vec<3, float> perp_old = new_norm - ndot * old_norm;
390  algebra::vec<3, float> dperp =
391  1.0f / (1 + ndot) * (old_norm + new_norm);
392  new_u -= dperp * (new_u * perp_old);
393  new_v -= dperp * (new_v * perp_old);
394  }
395 
396 
397  // FIXME: Add const to formals whenever we can.
398 
399  // Reproject a curvature tensor from the basis spanned by old_u and old_v
400  // (which are assumed to be unit-length and perpendicular) to the
401  // new_u, new_v basis.
402  inline
403  void
404  proj_curv(const algebra::vec<3, float>& old_u,
405  const algebra::vec<3, float>& old_v,
406  float old_ku, float old_kuv, float old_kv,
407  const algebra::vec<3, float>& new_u,
408  const algebra::vec<3, float>& new_v,
409  float& new_ku, float& new_kuv, float& new_kv)
410  {
411  algebra::vec<3, float> r_new_u, r_new_v;
412  rot_coord_sys(new_u, new_v, vprod(old_u, old_v), r_new_u, r_new_v);
413 
414  float u1 = r_new_u * old_u;
415  float v1 = r_new_u * old_v;
416  float u2 = r_new_v * old_u;
417  float v2 = r_new_v * old_v;
418  new_ku = old_ku * u1*u1 + old_kuv * (2.0f * u1*v1) + old_kv * v1*v1;
419  new_kuv = old_ku * u1*u2 + old_kuv * (u1*v2 + u2*v1) + old_kv * v1*v2;
420  new_kv = old_ku * u2*u2 + old_kuv * (2.0f * u2*v2) + old_kv * v2*v2;
421  }
422 
426  inline
427  void
428  diagonalize_curv(const algebra::vec<3, float>& old_u,
429  const algebra::vec<3, float>& old_v,
430  float ku, float kuv, float kv,
431  const algebra::vec<3, float>& new_norm,
432  algebra::vec<3, float>& pdir1,
433  algebra::vec<3, float>& pdir2,
434  float& k1, float& k2)
435  {
436  algebra::vec<3, float> r_old_u, r_old_v;
437  rot_coord_sys(old_u, old_v, new_norm, r_old_u, r_old_v);
438 
439  float c = 1, s = 0, tt = 0;
440  if (likely(kuv != 0.0f))
441  {
442  // Jacobi rotation to diagonalize.
443  float h = 0.5f * (kv - ku) / kuv;
444  tt = (h < 0.0f) ?
445  1.0f / (h - sqrt(1.0f + h*h)) :
446  1.0f / (h + sqrt(1.0f + h*h));
447  c = 1.0f / sqrt(1.0f + tt*tt);
448  s = tt * c;
449  }
450 
451  k1 = ku - tt * kuv;
452  k2 = kv + tt * kuv;
453 
454  if (fabs(k1) >= fabs(k2))
455  pdir1 = c*r_old_u - s*r_old_v;
456  else
457  {
458  std::swap(k1, k2);
459  pdir1 = s*r_old_u + c*r_old_v;
460  }
461  pdir2 = vprod(new_norm, pdir1);
462  }
463 
464  } // end of namespace mln::geom::internal
465 
466 
478  /* FIXME: We should restrict attached data to vertices in return
479  value. */
480  /* FIXME: Have a typedef for this return type? */
481  inline
482  std::pair<
483  complex_image< 2, mln::space_2complex_geometry, float >,
484  complex_image< 2, mln::space_2complex_geometry, float >
485  >
487  {
488  // Shortcuts.
489  static const unsigned D = 2;
490  typedef space_2complex_geometry G;
491  typedef algebra::vec<3, float> vec3f;
492  typedef algebra::mat<3, 3, float> mat3f;
493 
494  typedef complex_image< D, G, vec3f > corner_area_t;
495  typedef complex_image< D, G, float > point_area_t;
496 
497  // Normals at vertices.
499  normal_t normal = mesh_normal(mesh);
500 
501  // Areas ``belonging'' to a normal.
502  typedef complex_image< D, G, vec3f > corner_area_t;
503  typedef complex_image< D, G, float > point_area_t;
504  typedef std::pair<corner_area_t, point_area_t> corner_point_area_t;
505  corner_point_area_t corner_point_area = mesh_corner_point_area(mesh);
506  // Shortcuts.
507  corner_area_t& corner_area = corner_point_area.first;
508  point_area_t& point_area = corner_point_area.second;
509 
510  // Curvatures at each vertex.
512  typedef std::pair<curv_t, curv_t> output_t;
513  output_t output(mesh, mesh);
514  curv_t& curv1 = output.first;
515  curv_t& curv2 = output.second;
516  // ??
518  // Principal directions at each vertex.
521 
522  /*-------------------------------------------------.
523  | Set up an initial coordinate system per vertex. |
524  `-------------------------------------------------*/
525 
527  // A neighborhood where neighbors are the set of 0-faces
528  // transitively adjacent to the reference point.
529  typedef mln::complex_m_face_neighborhood<D, G> adj_vertices_nbh_t;
530  adj_vertices_nbh_t adj_vertices_nbh;
531  mln_niter_(adj_vertices_nbh_t) adj_v(adj_vertices_nbh, f);
532  /* FIXME: We should be able to pas this value (m) either at the
533  construction of the neighborhood or at the construction of
534  the iterator. Alas, we can't inherit ctors (yet), so we
535  can't rely on a simple mixin approach. */
536  adj_v.iter().set_m(0);
537 
538  for_all(f)
539  {
541  std::vector<mln_psite_(curv_t)> p;
542  p.reserve(3);
543  for_all(adj_v)
544  p.push_back(adj_v);
546  mln_assertion(p.size() == 3);
547 
548  // FIXME: See above.
549  pdir1(p[0]) =
550  p[1].to_site().front().to_vec() - p[0].to_site().front().to_vec();
551  pdir1(p[1]) =
552  p[2].to_site().front().to_vec() - p[1].to_site().front().to_vec();
553  pdir1(p[2]) =
554  p[0].to_site().front().to_vec() - p[2].to_site().front().to_vec();
555  }
556 
558  for_all(v)
559  {
560  pdir1(v) = algebra::vprod(pdir1(v), normal(v));
561  pdir1(v).normalize();
562  pdir2(v) = algebra::vprod(normal(v), pdir1(v));
563  }
564 
565  /*-----------------------------.
566  | Compute curvature per-face. |
567  `-----------------------------*/
568 
569  for_all(f)
570  {
571  std::vector<mln_psite_(curv_t)> p;
572  p.reserve(3);
573  for_all(adj_v)
574  p.push_back(adj_v);
575  mln_assertion(p.size() == 3);
576 
577  // (Opposite) edge vectors.
578  algebra::vec<3, vec3f> e;
579  // FIXME: See above.
580  e.set
581  (p[2].to_site().front().to_vec() - p[1].to_site().front().to_vec(),
582  p[0].to_site().front().to_vec() - p[2].to_site().front().to_vec(),
583  p[1].to_site().front().to_vec() - p[0].to_site().front().to_vec());
584 
585  // N-T-B coordinate system per face.
586  vec3f t = e[0];
587  t.normalize();
588  vec3f n = algebra::vprod(e[0], e[1]);
589  /* FIXME: Why don't we normalize N? Is it a normal vector by
590  construction? Or maybe we don't need it to be normal? */
591  vec3f b = algebra::vprod(n, t);
592  b.normalize();
593 
594  // Estimate curvature based on variation of normals along edges.
595  vec3f m = literal::zero;
596  mat3f w = literal::zero;
597  for (int j = 0; j < 3; ++j)
598  {
599  float u = e[j] * t;
600  float v = e[j] * b;
601  w(0, 0) += u * u;
602  w(0, 1) += u * v;
603  /* FIXME: Those two lines were commented in Trimesh's
604  original code.
605 
606  //w(1, 1) += v*v + u*u;
607  //w(1, 2) += u*v;
608 
609  */
610  w(2, 2) += v * v;
611  vec3f dn =
612  normal(p[internal::prev(j)]) - normal(p[internal::next(j)]);
613  float dnu = dn * t;
614  float dnv = dn * b;
615  m[0] += dnu*u;
616  m[1] += dnu*v + dnv*u;
617  m[2] += dnv*v;
618  }
619  w(1, 1) = w(0, 0) + w(2, 2);
620  w(1, 2) = w(0, 1);
621 
622  // Least squares solution.
623  vec3f diag;
624 #ifndef NDEBUG
625  bool ldlt_decomp_sucess_p = algebra::ldlt_decomp(w, diag);
626 #endif // ! NDEBUG
627  mln_assertion(ldlt_decomp_sucess_p);
628  algebra::ldlt_solve(w, diag, m, m);
629 
630  // Push it back out to the vertices
631  for (int j = 0; j < 3; ++j)
632  {
633  float c1, c12, c2;
634  internal::proj_curv(t, b, m[0], m[1], m[2],
635  pdir1(p[j]), pdir2(p[j]), c1, c12, c2);
636  float wt = corner_area(f)[j] / point_area(p[j]);
637  curv1(p[j]) += wt * c1;
638  curv12(p[j]) += wt * c12;
639  curv2(p[j]) += wt * c2;
640  }
641  }
642 
643  /*-------------------------------------------------------------.
644  | Compute principal directions and curvatures at each vertex. |
645  `-------------------------------------------------------------*/
646 
647  for_all(v)
648  internal::diagonalize_curv(pdir1(v), pdir2(v),
649  curv1(v), curv12(v), curv2(v),
650  normal(v), pdir1(v), pdir2(v),
651  curv1(v), curv2(v));
652 
653  return output;
654  }
655 
656  } // end of namespace mln::geom
657 
658 } // end of namespace mln
659 
660 #endif // ! MILENA_APPS_MESH_SEGM_SKEL_MISC_HH