GLAMERDOC++
Gravitational Lensing Code Library
Loading...
Searching...
No Matches
concave_hull.h
1//
2// concave_hull.h
3// GLAMER
4//
5// Created by bmetcalf on 23/02/16.
6//
7//
8
9#ifndef concave_hull_h
10#define concave_hull_h
11
12#include <set>
13#include "point.h"
14#include "geometry.h"
15#include "image_processing.h"
16#include "simpleTreeVec.h"
17
18namespace Utilities{
19
20template<typename T,typename P>
21Point_2d subtract(T& p1,P& p2){
22 return Point_2d(p1[0] - p2[0], p1[1] - p2[1]);
23}
24template<typename P>
25double crossD(P &O,P &A,P &B){
26 return (A[0] - O[0]) * (B[1] - O[1]) - (A[1] - O[1]) * (B[0] - O[0]);
27}
28
29
31template <typename T>
32size_t RemoveIntersections(std::vector<T> &curve){
33
34 if(curve.size() <=3) return 0;
35
36 size_t N = curve.size(),count = 0;
37 T tmp;
38
39 curve.push_back(curve[0]);
40
41 for(size_t i=0;i<N-2;++i){
42 for(size_t j=i+2;j<N;++j){
43 if(Utilities::Geometry::intersect(curve[i].x,curve[i+1].x,curve[j].x,curve[j+1].x)){
44
45 size_t k=i+1,l=j;
46 while(k < l){
47 std::swap(curve[k] , curve[l]);
48 ++k;
49 --l;
50 }
51 ++count;
52 }
53 }
54 }
55
56 //assert(curve[0]==curve.back());
57 curve.pop_back();
58
59 return count;
60}
61
66std::vector<Point_2d> TighterHull(const std::vector<Point_2d> &v);
67
69std::vector<Point_2d> TightestHull(const std::vector<Point_2d> &v);
70
71//template <typename T>
72//std::vector<T> TightHull(const std::vector<T> &curve){
73//
74// if(curve.size() <=3) return curve;
75//
76// size_t N = curve.size(),count = 0;
77// T tmp;
78//
79// // find left most point
80// long init = 0;
81// double pmin = curve[0][0];
82// for(int i=1 ; i<N ; ++i){
83// if(curve[i][0] < pmin){
84// init = i;
85// pmin = curve[i][0];
86// }
87// }
88//
89// // orientation == 1 clockwise
90// Geometry::CYCLIC cyc(N);
91// int orientation = sign( (curve[cyc[init-1]] - curve[init])^(curve[cyc[init+1]] - curve[init]) );
92// if(orientation == 0){
93// orientation = sign( curve[cyc[init+1]][1] - curve[init][1]);
94// }
95//
96// orientation *= -1; // orientation == 1 is counter clockwise
97//
98// // make a copy where first one is the leftmost
99// std::vector<T> hull(N);
100// for(size_t i=0;i<N;++i){
101// hull[i] = curve[ cyc[ orientation * i + init] ];
102// }
103//
104// for(size_t i=0;i<N-2;++i){
105//
106// for(size_t j=i+2;j<N;++j){
107// if(Utilities::Geometry::intersect(hull[i].x,hull[ (i+1) % N ].x,hull[j].x,hull[ (j+1) % N ].x)){
108//
109//
110// long ii = cyc[i-1];
111// //if(i==0) ii = N-1;
112// //else ii = i-1;
113//
114// Point_2d x = hull[i] - hull[ii];
115// x.unitize();
116//
117// if(
118// //( x^(hull[ (j+1) % N ] - hull[i]) ) / (hull[ (j+1) % N ] - hull[i]).length()
119// atan2( x^(hull[ (j+1) % N ] - hull[i]) , x*(hull[ (j+1) % N ] - hull[i]) )
120// <
121// atan2( x^(hull[j] - hull[i]) , x*(hull[j] - hull[i]) )
122// //( x^(hull[j] - hull[i]) ) / (hull[j] - hull[i]).length()
123// ){
124//
125// // inner loop - remove all points between i and j+1
126//
127// int n = j - i;
128// int k = j+1;
129// while(k < N){
130// //std::swap(hull[k],hull[k-n]);
131// hull[k-n] = hull[k];
132// ++k;
133// }
134// N -= n;
135// }else{
136//
137// // outer loop - put them in reverse order
138// size_t k=i+1,l=j;
139// while(k < l){
140// std::swap(hull[k] , hull[l]);
141// ++k;
142// --l;
143// }
144// }
145// ++count;
146// }
147// }
148// }
149//
150// //assert(hull[0]==hull.back());
151// hull.resize(N);
152//
153// assert(N>2);
154// return hull;
155//}
156
159 const Point_2d &x2,
160 const Point_2d &x3,
162
164Point_2d RandomInConvexPoly(const std::vector<Point_2d> &pp,
166std::vector<Point_2d> RandomInConvexPoly(const std::vector<Point_2d> &pp,
167 int N,
169
171Point_2d RandomInPoly(std::vector<Point_2d> &pp,
173std::vector<Point_2d> RandomInPoly(std::vector<Point_2d> &pp,
174 int N,
176
178Point_2d RandomNearPoly(std::vector<Point_2d> &pp
179 ,double R
181std::vector<Point_2d> RandomNearPoly(std::vector<Point_2d> &pp
182 ,int N
183 ,double R
185
186
187
200template <typename P>
201void find_boundaries(std::vector<bool> &bitmap // = true inside
202 ,long nx // number of pixels in x direction
203 ,std::vector<std::vector<P> > &points
204 ,std::vector<bool> &hits_edge
205 ,bool add_to_vector=false
206 ,bool outer_only=false
207 ){
208
209 size_t n = bitmap.size();
210 long ny = n/nx;
211
212 if(n != nx*ny){
213 std::cerr << "Wrong sizes in Utilities::find_boundaries." << std::endl;
214 throw std::invalid_argument("invalid size");
215 }
216
217 std::vector<bool> not_used(n,true);
218
219 // pad edge of field with bitmap=false
220 for(size_t i=0 ; i<nx ; ++i) bitmap[i]=false;
221 size_t j = nx*(ny-1);
222 for(size_t i=0 ; i<nx ; ++i) bitmap[i + j]=false;
223 for(size_t i=0 ; i<ny ; ++i) bitmap[i*nx]=false;
224 j = nx-1;
225 for(size_t i=0 ; i<ny ; ++i) bitmap[j + i*nx]=false;
226
227 std::list<std::list<Point_2d>> contours;
228
229 if(!add_to_vector){
230 hits_edge.resize(0);
231 }
232
233 bool done = false;
234 long kfirst_in_bound = -1;
235 while(!done){
236 // find first cell in edge
237 size_t k=0;
238 int type;
239 for( k = kfirst_in_bound + 1 ; k < n - nx ; ++k){
240 if(k % nx != nx-1){ // one less cells than points
241 type = 0;
242 if(bitmap[k] ) type +=1;
243 if(bitmap[k+1]) type += 10;
244 if(bitmap[k + nx]) type += 100;
245 if(bitmap[k + nx + 1]) type += 1000;
246
247 if(type > 0
248 && type != 1111
249 && not_used[k]
250 ) break;
251 }
252 }
253
254 kfirst_in_bound = k;
255
256 if(k == n-nx){
257 done=true;
258 }else{ // found an edge
259
260 contours.resize(contours.size() + 1);
261 std::list<Point_2d> &contour = contours.back();
262 hits_edge.push_back(false);
263
264
265 int type;
266 int face_in=0;
267 size_t n_edge = 0;
268
269 // follow edge until we return to the first point
270 while(k != kfirst_in_bound || n_edge==0){
271
272 if(n_edge >= n){ // infinite loop, output debugging data
273 std::cerr << "Too many points in Utilities::find_boundaries()." << std::endl;
274 std::cerr << "kfirst_in_bound " << kfirst_in_bound << std::endl;
275 std::cerr << " contour is output to boundary_error_file.csv and bitmap_error_file.csv" << std::endl;
276 {
277 std::ofstream file("bitmap_error_file.csv");
278 file << "in,x,y" << std::endl;
279 for(size_t i=0 ; i<n ; ++i){
280 file << bitmap[i] << "," << i%nx << "," << i/nx << std::endl;
281 }
282 }
283
284 {
285 std::ofstream file("boundary_error_file.csv");
286 file << "contour,x,y" << std::endl;
287 int i = 0;
288 for(auto &v : contours){
289 for(Point_2d &p : v){
290 file << i << "," << p[0] << "," << p[1] << std::endl;
291 }
292 ++i;
293 }
294 }
295 throw std::runtime_error("caught in loop.");
296 }
297
298 if(k%nx == 0 || k%nx == nx-2) hits_edge.back() = true;
299 if(k/nx == 0 || k/nx == ny-2) hits_edge.back() = true;
300
301 not_used[k] = false;
302
303 ++n_edge;
304 type = 0;
305 // find type of cell
306 if(bitmap[k]) type +=1;
307 if(bitmap[k+1]) type += 10;
308 if(bitmap[k + nx]) type += 100;
309 if(bitmap[k + nx + 1]) type += 1000;
310
311 if(type == 0 || type == 1111){ // all in or all out
312 throw std::runtime_error("off edge!!");
313 }else if(type == 1 || type == 1110){ // lower left only
314
315 if(face_in==0){
316 contour.push_back( Point_2d( k%nx + 0.5
317 , k/nx ) );
318 //contour.push_back( (i_points[k] + i_points[k+1]) / 2 );
319 face_in=1;
320 k -= nx;
321 }else{
322 contour.push_back( Point_2d( k%nx
323 , k/nx + 0.5 ) );
324 //contour.push_back( (i_points[k] + i_points[k+nx]) / 2 );
325 face_in=2;
326 k -= 1;
327 }
328
329 }else if(type == 10 || type == 1101){ // lower right only
330
331 if(face_in==2){
332 contour.push_back( Point_2d( k%nx + 0.5
333 , k/nx ) );
334 //contour.push_back( (i_points[k] + i_points[k+1]) / 2 );
335 face_in=1;
336 k -= nx;
337 }else{
338 contour.push_back( Point_2d( k%nx + 1
339 , k/nx + 0.5) );
340 //contour.push_back( (i_points[k+nx+1] + i_points[k+1]) / 2 );
341 face_in=0;
342 k += 1;
343 }
344
345 }else if(type == 100 || type == 1011){ // upper left only
346
347 if(face_in==0){
348 contour.push_back( Point_2d( k%nx + 0.5
349 , k/nx + 1 ) );
350 //contour.push_back( (i_points[k+nx] + i_points[k+nx+1]) / 2 );
351 face_in=3;
352 k += nx;
353 }else{
354 contour.push_back( Point_2d( k%nx
355 , k/nx + 0.5 ) );
356 //contour.push_back( (i_points[k] + i_points[k+nx]) / 2 );
357 face_in=2;
358 k -= 1;
359 }
360
361 }else if(type == 1000 || type == 111){ // upper right only
362
363 if(face_in==1){
364 contour.push_back( Point_2d( k%nx + 1
365 , k/nx + 0.5 ) );
366 //contour.push_back( (i_points[k+1] + i_points[k+nx+1]) / 2 );
367 face_in=0;
368 k += 1;
369 }else{
370 contour.push_back( Point_2d( k%nx + 0.5
371 , k/nx + 1 ) );
372 //contour.push_back( (i_points[k+nx] + i_points[k+nx+1]) / 2 );
373 face_in=3;
374 k += nx;
375 }
376
377 }else if(type == 11 || type == 1100){ // lower two
378
379 if(face_in==0){
380 contour.push_back( Point_2d( k%nx + 1
381 , k/nx + 0.5 ) );
382 //contour.push_back( (i_points[k+1] + i_points[k+nx+1]) / 2 );
383 k += 1;
384 }else{
385 contour.push_back( Point_2d( k%nx
386 , k/nx + 0.5 ) );
387 //contour.push_back( (i_points[k] + i_points[k+nx]) / 2 );
388 face_in = 2;
389 k -= 1;
390 }
391
392 }else if(type == 1010 || type == 101){ // right two
393
394 if(face_in==1){
395 contour.push_back( Point_2d( k%nx + 0.5
396 , k/nx ) );
397 //contour.push_back( (i_points[k] + i_points[k+1]) / 2 );
398 k -= nx;
399 }else{
400 contour.push_back( Point_2d( k%nx + 0.5
401 , k/nx + 1 ) );
402 //contour.push_back( (i_points[k+nx] + i_points[k+nx+1]) / 2 );
403 face_in = 3;
404 k += nx;
405 }
406
407 }else if(type == 1001){ // lower left upper right
408
409 if(face_in==0){
410 contour.push_back( Point_2d( k%nx + 0.5
411 , k/nx + 1 ) );
412 //contour.push_back( (i_points[k+nx] + i_points[k+nx+1]) / 2 );
413 face_in=3;
414 k += nx;
415 }else if(face_in==1){
416 contour.push_back( Point_2d( k%nx
417 , k/nx + 0.5 ) );
418 //contour.push_back( (i_points[k] + i_points[k+nx]) / 2 );
419 face_in=2;
420 k -= 1;
421 }else if(face_in==2){
422 contour.push_back( Point_2d( k%nx + 0.5
423 , k/nx ) );
424 //contour.push_back( (i_points[k] + i_points[k+1]) / 2 );
425 face_in=1;
426 k -= nx;
427 }else{
428 contour.push_back( Point_2d( k%nx + 1
429 , k/nx + 0.5 ) );
430 //contour.push_back( (i_points[k+nx + 1] + i_points[k+1]) / 2 );
431 face_in=0;
432 k += 1;
433 }
434
435 }else if(type == 110){ // upper left lower right
436
437 if(face_in==0){
438 contour.push_back( Point_2d( k%nx + 0.5
439 , k/nx ) );
440 //contour.push_back( (i_points[k] + i_points[k+1]) / 2 );
441 face_in=1;
442 k -= nx;
443 }else if(face_in==1){
444 contour.push_back( Point_2d( k%nx + 1
445 , k/nx + 0.5 ) );
446 //contour.push_back( (i_points[k+1] + i_points[k+nx+1]) / 2 );
447 face_in=0;
448 k += 1;
449 }else if(face_in==2){
450 contour.push_back( Point_2d( k%nx + 0.5
451 , k/nx + 1 ) );
452 //contour.push_back( (i_points[k + nx] + i_points[k+nx+1]) / 2 );
453 face_in=3;
454 k += nx;
455 }else{
456 contour.push_back( Point_2d( k%nx
457 , k/nx + 0.5 ) );
458 //contour.push_back( (i_points[k] + i_points[k+nx]) / 2 );
459 face_in=2;
460 k -= 1;
461 }
462 }
463 }
464
465 // special case of the diamand with a hole
466 //if(n_edge == 12 && bitmap[k + nx + 1] ){
467 // n_edge=0;
468 // face_in = 2;
469 //}
470
471 }
472 if(outer_only) break;
473 }
474
475 long offset = 0;
476 if(!add_to_vector){
477 points.resize(contours.size());
478 }else{
479 offset = points.size();
480 points.resize(points.size() + contours.size());
481 }
482
483 // copy lists of points into vectors
484 int i=0;
485 for(auto &c: contours){
486 points[offset+i].resize(c.size());
487 size_t j=0;
488 for(auto &p: c){
489 points[offset+i][j] = p;
490 ++j;
491 }
492 ++i;
493 }
494}
495
497void find_islands(std::vector<bool> &bitmap // = true inside
498 ,long nx // number of pixels in x direction
499 ,std::vector<std::vector<long> > &indexes
500 ,std::vector<bool> &hits_edge
501 ,bool add_to_vector=false
502 );
503
505double interior_mass(const std::vector<Point_2d> &alpha
506 ,const std::vector<Point_2d> &x);
507
509double interior_mass(const std::vector<RAY> &rays);
510
512template<typename T>
513std::vector<T> convex_hull(const std::vector<T> &PP)
514{
515
516 if(PP.size() <= 3){
517 return PP;
518 }
519
520 std::vector<T> P=PP;
521
522 size_t n = P.size();
523 size_t k = 0;
524 std::vector<T> hull(2*n);
525
526 // Sort points lexicographically
527 std::sort(P.begin(), P.end(),
528 [](T p1,T p2){
529 if(p1[0]==p2[0]) return p1[1] < p2[1];
530 return p1[0] < p2[0];});
531
532
533 // Build lower hull
534 for (size_t i = 0; i < n; i++) {
535 while (k >= 2 && crossD(hull[k-2], hull[k-1], P[i]) <= 0){
536 k--;
537 }
538 hull[k++] = P[i];
539 }
540
541 // Build upper hull
542 for (long i = n-2, t = k+1; i >= 0; i--) {
543 while (k >= t && crossD(hull[k-2], hull[k-1], P[i]) <= 0){
544 k--;
545 assert(k > 1);
546 }
547 hull[k++] = P[i];
548 }
549
550
551 hull.resize(k);
552 hull.pop_back();
553
554 return hull;
555}
556
559template<typename T>
560void convex_hull(const std::vector<T> &P,std::vector<size_t> &hull)
561{
562
563 size_t n = P.size();
564
565 if(n <= 3){
566 hull.resize(n);
567 size_t i = 0;
568 for(size_t &d : hull) d = i++;
569 return;
570 }
571
572 //std::vector<T> hull(n + 1);
573 hull.resize(n + 1);
574 size_t k = 0;
575
576 // Sort points lexicographically
577 std::vector<size_t> sort_index(P.size());
578 {
579 int i=0;
580 for(size_t &a : sort_index) a = i++;
581 }
582 std::sort(sort_index.begin(), sort_index.end(),
583 [&P](size_t i1,size_t i2){
584 if(P[i1][0]==P[i2][0]) return P[i1][1] < P[i2][1];
585 return P[i1][0] < P[i2][0];});
586
587
588 // Build lower hull
589 for (size_t i = 0; i < n; i++) {
590 while (k >= 2 && crossD(P[hull[k-2]], P[hull[k-1]], P[sort_index[i]]) <= 0){
591 k--;
592 }
593 hull[k++] = sort_index[i];
594 }
595
596 // Build upper hull
597 for (long i = n-2, t = k+1; i >= 0; i--) {
598 while (k >= t && crossD(P[hull[k-2]], P[hull[k-1]],P[sort_index[i]]) <= 0){
599 k--;
600 assert(k > 1);
601 }
602 hull[k++] = sort_index[i];
603 }
604
605 hull.resize(k);
606 hull.pop_back();
607
608 return;
609}
610
611
612struct Edge{
613 Edge(size_t i,double l):length(l),index(i){};
614 double length = 0;
615 size_t index = 0;
616};
632template<typename T>
633void concave(std::vector<T> &init_points
634 ,std::vector<T> &hull_out,double scale)
635{
636
637 //typedef typename InputIt::value_type point;
638
639 bool TEST = false;
640
641 // find the convex hull
642 std::vector<T> hull = convex_hull(init_points);
643
644 if(init_points.size() == hull.size()) return;
645
646 std::list<Edge> edges;
647 std::list<T> leftovers;
648 hull.push_back(hull[0]);
649
650 double tmp;
651
652 // find pairs of hull points that are further appart than scale
653 for(int i=0;i<hull.size()-1;++i){
654 tmp = (subtract(hull[i],hull[i+1])).length();
655 if(tmp > scale) edges.emplace_back(i,tmp);// .push_back(std::pair<size_t,double>(i,tmp));
656 }
657
658 if(edges.size() == 0) return;
659
660 if(TEST){
661 PosType cent[] ={0.5,0.5};
662 PixelMap<float> map(cent,256, 1.0/256);
663
664 map.drawgrid(10,0.5);
665
666 map.drawPoints(init_points,0.01,1.0);
667 map.drawCurve(hull,2);
668
669 map.printFITS("!test_concave.fits");
670 }
671
672 // sort edges by length
673 edges.sort([](const Edge &p1,const Edge &p2){return p1.length > p2.length;});
674 assert(edges.front().length >= edges.back().length);
675
676
677 { // make a list of points that are not already in the convex hull
678
679 hull.pop_back();
680 std::vector<size_t> index(hull.size());
681 size_t i = 0;
682 for(size_t &ind: index) ind = i++;
683 std::sort(index.begin(),index.end(),
684 [&hull](size_t p1,size_t p2){
685 if(hull[p1][0] == hull[p2][0]) return hull[p1][1] < hull[p2][1];
686 return hull[p1][0] < hull[p2][0];});
687
688 size_t j = 0;
689
690 for(auto pit = init_points.begin(); pit != init_points.end() ; ++pit){
691
692 if((hull[index[j]][0] == (*pit)[0])*(hull[index[j]][1] == (*pit)[1])){
693 ++j;
694 }else{
695 leftovers.push_back(*pit);
696 }
697
698 }
699 assert(hull.size() + leftovers.size() == init_points.size());
700 hull.push_back(hull[0]);
701 }
702
703 typename std::list<T>::iterator p,nextpoint;
704 //auto p = leftovers.begin();
705
706 double minarea,area,co1;//,co2;
707 size_t i;
708 Point_2d vo,v1,v2;
709
710 //int count = 0;
711 while(edges.front().length > scale ){
712
713 if(leftovers.size() == 0) break;
714
715 // find the point which if added to this edge would change the area least
716 minarea = HUGE_VAL;
717 i = edges.front().index;
718
719 for(p = leftovers.begin() ; p != leftovers.end() ; ++p){
720
721 v1 = subtract(*p,hull[i]);
722
723 if( edges.front().length > v1.length()){
724
725 vo = subtract(hull[i+1], hull[i]);
726 //v2 = subtract(*p,hull[i+1]);
727 area = vo^v1;
728 co1 = v1*vo;
729 //co2 = v2*vo;
730
731 if( co1 > 0 && (area > 0)*(area < minarea) ){
732 // if(co1 > 0 && area > 0 && index_length.front().second > v1.length()){
733 minarea = area;
734 nextpoint = p;
735 }
736 }
737 }
738
739 if(minarea == HUGE_VAL){
740 // if there is no acceptable point for this edge continue with the second longest edge
741 edges.pop_front();
742 continue;
743 }
744
745 // insert new point into hull
746 T tmp = *nextpoint;
747 leftovers.erase(nextpoint);
748 hull.insert(hull.begin() + i + 1,tmp);
749
750 // update index_length list
751 Edge new1 = edges.front();
752 edges.pop_front();
753 new1.length = (subtract(hull[i],hull[i+1])).length();
754 Edge new2(i+1,(subtract(hull[i+1],hull[i+2])).length());
755
756 for(auto &p : edges) if(p.index > i) ++(p.index);
757
758 auto p = edges.begin();
759 if(new1.length > scale){
760 for(p = edges.begin() ; p != edges.end() ; ++p){
761 if((*p).length < new1.length){
762 edges.insert(p,new1);
763 break;
764 }
765 }
766 if(p==edges.end()) edges.push_back(new1);
767 }
768 if(new2.length > scale){
769 for(p = edges.begin(); p != edges.end() ; ++p){
770 if((*p).length < new2.length){
771 edges.insert(p,new2);
772 break;
773 }
774 }
775 if(p==edges.end()) edges.push_back(new2);
776 }
777
778
779 if(TEST){
780 PosType cent[] ={0.5,0.5};
781 PixelMap<float> map(cent,256, 1.0/256);
782
783 map.drawgrid(10,0.5);
784
785 map.drawPoints(init_points,0.03,1.0);
786 map.drawCurve(hull,2);
787
788 map.printFITS("!test_concave.fits");
789 }
790 }
791
792 hull.pop_back();
793
794 //Utilities::RemoveIntersections(hull);
795
796 std::swap(hull,hull_out);
797
798 return;
799}
800
801template<typename T>
802std::vector<T> concave2(std::vector<T> &init_points,double scale)
803{
804
805 //typedef typename InputIt::value_type point;
806
807 bool TEST = false;
808
809 // find the convex hull
810 std::vector<T> hull = convex_hull(init_points);
811
812 if(init_points.size() == hull.size()) return hull;
813
814 std::list<Edge> edges;
815 std::list<T> leftovers;
816 hull.push_back(hull[0]);
817
818 double tmp;
819
820 // find pairs of hull points that are further appart than scale
821 for(int i=0;i<hull.size()-1;++i){
822 tmp = (subtract(hull[i],hull[i+1])).length();
823 if(tmp > scale) edges.emplace_back(i,tmp);// .push_back(std::pair<size_t,double>(i,tmp));
824 }
825
826 if(edges.size() == 0) return hull;
827
828 if(TEST){
829 PosType cent[] ={0.5,0.5};
830 PixelMap<float> map(cent,256, 1.0/256);
831
832 map.drawgrid(10,0.5);
833
834 map.drawPoints(init_points,0.01,1.0);
835 map.drawCurve(hull,2);
836
837 map.printFITS("!test_concave.fits");
838 }
839
840 // sort edges by length
841 edges.sort([](const Edge &p1,const Edge &p2){return p1.length > p2.length;});
842 assert(edges.front().length >= edges.back().length);
843
844
845 { // make a list of points that are not already in the convex hull
846
847 hull.pop_back();
848 std::vector<size_t> index(hull.size());
849 size_t i = 0;
850 for(size_t &ind: index) ind = i++;
851 std::sort(index.begin(),index.end(),
852 [&hull](size_t p1,size_t p2){
853 if(hull[p1][0] == hull[p2][0]) return hull[p1][1] < hull[p2][1];
854 return hull[p1][0] < hull[p2][0];});
855
856 size_t j = 0;
857
858 for(auto pit = init_points.begin(); pit != init_points.end() ; ++pit){
859
860 if((hull[index[j]][0] == (*pit)[0])*(hull[index[j]][1] == (*pit)[1])){
861 ++j;
862 }else{
863 leftovers.push_back(*pit);
864 }
865
866 }
867 assert(hull.size() + leftovers.size() == init_points.size());
868 hull.push_back(hull[0]);
869 }
870
871 typename std::list<T>::iterator p,nextpoint;
872 //auto p = leftovers.begin();
873
874 double minarea,area,co1;
875 size_t i;
876 Point_2d vo,v1,v2;
877
878 while(edges.front().length > scale ){
879
880 if(leftovers.size() == 0) break;
881
882 // find the point which if added to this edge would change the area least
883 minarea = HUGE_VAL;
884 i = edges.front().index;
885
886 for(p = leftovers.begin() ; p != leftovers.end() ; ++p){
887
888 v1 = subtract(*p,hull[i]);
889
890 if( edges.front().length > v1.length()){
891
892 vo = subtract(hull[i+1], hull[i]);
893 //v2 = subtract(*p,hull[i+1]);
894 area = vo^v1;
895 co1 = v1*vo;
896 //co2 = v2*vo;
897
898 if( co1 > 0 && (area > 0)*(area < minarea) ){
899 // if(co1 > 0 && area > 0 && index_length.front().second > v1.length()){
900 minarea = area;
901 nextpoint = p;
902 }
903 }
904 }
905
906 if(minarea == HUGE_VAL){
907 // if there is no acceptable point for this edge continue with the second longest edge
908 edges.pop_front();
909 continue;
910 }
911
912 // insert new point into hull
913 T tmp = *nextpoint;
914 leftovers.erase(nextpoint);
915 hull.insert(hull.begin() + i + 1,tmp);
916
917 // update index_length list
918 Edge new1 = edges.front();
919 edges.pop_front();
920 new1.length = (subtract(hull[i],hull[i+1])).length();
921 Edge new2(i+1,(subtract(hull[i+1],hull[i+2])).length());
922
923 for(auto &p : edges) if(p.index > i) ++(p.index);
924
925 auto p = edges.begin();
926 if(new1.length > scale){
927 for(p = edges.begin() ; p != edges.end() ; ++p){
928 if((*p).length < new1.length){
929 edges.insert(p,new1);
930 break;
931 }
932 }
933 if(p==edges.end()) edges.push_back(new1);
934 }
935 if(new2.length > scale){
936 for(p = edges.begin(); p != edges.end() ; ++p){
937 if((*p).length < new2.length){
938 edges.insert(p,new2);
939 break;
940 }
941 }
942 if(p==edges.end()) edges.push_back(new2);
943 }
944
945
946 if(TEST){
947 PosType cent[] ={0.5,0.5};
948 PixelMap<float> map(cent,256, 1.0/256);
949
950 map.drawgrid(10,0.5);
951
952 map.drawPoints(init_points,0.03,1.0);
953 map.drawCurve(hull,2);
954
955 map.printFITS("!test_concave.fits");
956 }
957 }
958
959 hull.pop_back();
960
961 //Utilities::RemoveIntersections(hull);
962
963 return hull;
964}
965
966// shortest distance between P and the segment (S1,S2)
967template <typename Ptype>
968double distance_to_segment(const Ptype &P
969 ,const Ptype &S1
970 ,const Ptype &S2
971 ){
972
973 Ptype D = S2-S1;
974 double s = (P-S1)*D / D.length_sqr();
975 if(isnan(s) || s<=0){
976 return (P-S1).length();
977 }else if(s>=1){
978 return (P-S2).length();
979 }else{
980 return (S1 + D*s - P).length();
981 }
982}
983template <typename Ptype>
984double distance_to_segment(const Ptype &P
985 ,const Ptype &S1
986 ,const Ptype &S2
987 ,Ptype &closest_point
988 ){
989
990 Ptype D = S2-S1;
991 double s = (P-S1)*D / D.length_sqr();
992 if(isnan(s) || s<=0){
993 closest_point = S1;
994 }else if(s>=1){
995 closest_point = S2;
996 }else{
997 closest_point = S1 + D*s;
998 }
999 return (closest_point - P).length();
1000}
1001
1002template <typename Ptype>
1003bool segments_cross(const Ptype &a1,const Ptype &a2
1004 ,const Ptype &b1,const Ptype &b2){
1005 Ptype db= b2 - b1;
1006 Ptype da= a2 - a1;
1007 Ptype d1= b1 - a1;
1008
1009 double tmp = (db^da);
1010 if(tmp==0) return false; // parallel case
1011
1012 double B = (da^d1) / tmp;
1013 if( (B>1) || (B<0)) return false;
1014 B = (db^d1) / tmp;
1015 if( (B>1) || (B<0)) return false;
1016 return true;
1017}
1018
1020template <typename Ptype>
1021bool inCurve(const Ptype &x,const std::vector<Ptype> &H){
1022 if(H.size() < 3) return false;
1023 Point_2d p1,p2;
1024 long w=0;
1025 size_t n=H.size();
1026 for(size_t i=0 ; i<n ; ++i){
1027 size_t j = (i+1)%n;
1028 if(H[j][1] != H[i][1]){ // ignore horizontal segments
1029 p1[1] = H[i][1] - x[1];
1030 p2[1] = H[j][1] - x[1];
1031 if( p1[1]*p2[1] < 0 ){
1032 p1[0] = H[i][0] - x[0];
1033 p2[0] = H[j][0] - x[0];
1034
1035 if(p1[0]>=0 && p2[0]>=0){
1036 w += 2*( p1[1] >= 0 ) - 1;
1037 }else if(p1[0] > 0 || p2[0] > 0){
1038 double s = p1[0] - (p2[0]-p1[0])*p1[1]/(p2[1]-p1[1]);
1039 if(s >= 0){
1040 w += 2*( p1[1] >= 0 ) - 1;
1041 }
1042 }
1043 }else if(p2[1]==0){ // second point on lines
1044 p2[0] = H[j][0] - x[0];
1045 // this needs to be handled separately because of the
1046 // possibility that the curve touches the line a one vertex
1047 // but does not cross it.
1048 if(p2[0] >= 0){
1049
1050 size_t k = (j+1)%n;
1051 while( H[k][1] - x[1] == 0 && k!=i) k = (k+1)%n;
1052
1053 if(k!=i && p1[1]*( H[k][1] - x[1] ) <= 0){
1054 w += 2*( p1[1] >= 0 ) - 1;
1055 }
1056 }
1057 }
1058 }
1059 }
1060 return w != 0;
1061}
1062
1063template <typename Ptype>
1064bool inhull(PosType x[],const std::vector<Ptype> &H){
1065
1066 size_t n = H.size();
1067 if(n <=2) return false;
1068 long w=0;
1069 Ptype dH,dD;
1070 double B,A;
1071 for(size_t i=0 ; i<n-1 ; ++i){
1072 dH = H[i+1]-H[i];
1073 if(dH[1] == 0){ // // horizontal segment
1074 if(x[1] == H[i][1]){
1075 B = (x[0]-H[i][0])/dH[0];
1076 if( B>=0 && B<=1 ) return true; // point on boundary
1077 }
1078 }else{
1079 dD[0] = x[0]-H[i][0];
1080 dD[1] = x[1]-H[i][1];
1081 A = dD[1]/dH[1];
1082 if( A > 1 || A <= 0) continue;
1083 B = (dH^dD)/dH[1];
1084 if(B==0){ // on the boundary
1085 return true; // boundaries are always in
1086 }else if(B>0){
1087 if(dH[1] > 0) ++w;
1088 else --w;
1089 }
1090 }
1091 }
1092
1093 return abs(w);
1094}
1095
1096template <>
1097inline bool inhull<Point *>(PosType x[],const std::vector<Point *> &H){
1098
1099 size_t n = H.size();
1100 if(n <=2) return false;
1101 long w=0;
1102 Point_2d dH,dD;
1103 double B,A;
1104 for(size_t i=0 ; i<n-1 ; ++i){
1105 dH = *H[i+1]-*H[i];
1106 if(dH[1] == 0){ // // horizontal segment
1107 if(x[1] == H[i]->x[1]){
1108 B = (x[0]-H[i]->x[0])/dH[0];
1109 if( B>=0 && B<=1 ) return true; // point on boundary
1110 }
1111 }else{
1112 dD[0] = x[0]-H[i]->x[0];
1113 dD[1] = x[1]-H[i]->x[1];
1114 A = dD[1]/dH[1];
1115 if( A > 1 || A <= 0) continue;
1116 B = (dH^dD)/dH[1];
1117 if(B==0){ // on the boundary
1118 return true; // boundaries are always in
1119 }else if(B>0){
1120 if(dH[1] > 0) ++w;
1121 else --w;
1122 }
1123 }
1124 }
1125
1126 return w;
1127}
1128
1130template <>
1131inline bool inhull<RAY>(PosType x[],const std::vector<RAY> &H){
1132
1133 size_t n = H.size();
1134 if(n <=2) return false;
1135 long w=0;
1136 Point_2d dH,dD;
1137 double B,A;
1138 for(size_t i=0 ; i<n-1 ; ++i){
1139 dH = H[i+1].x-H[i].x;
1140 if(dH[1] == 0){ // // horizontal segment
1141 if(x[1] == H[i].x[1]){
1142 B = (x[0]-H[i].x[0])/dH[0];
1143 if( B>=0 && B<=1 ) return true; // point on boundary
1144 }
1145 }else{
1146 dD[0] = x[0]-H[i].x[0];
1147 dD[1] = x[1]-H[i].x[1];
1148 A = dD[1]/dH[1];
1149 if( A > 1 || A <= 0) continue;
1150 B = (dH^dD)/dH[1];
1151 if(B==0){ // on the boundary
1152 return true; // boundaries are always in
1153 }else if(B>0){
1154 if(dH[1] > 0) ++w;
1155 else --w;
1156 }
1157 }
1158 }
1159
1160 return w;
1161}
1162
1163template <>
1164inline bool inhull<PosType *>(PosType x[],const std::vector<PosType *> &H){
1165
1166 size_t n = H.size();
1167 if(n <=2) return false;
1168 long w=0;
1169 Point_2d dH,dD;
1170 double B,A;
1171 for(size_t i=0 ; i<n-1 ; ++i){
1172 dH[0] = H[i+1][0]-H[i][0];
1173 dH[1] = H[i+1][1]-H[i][1];
1174 if(dH[1] == 0){ // // horizontal segment
1175 if(x[1] == H[i][1]){
1176 B = (x[0]-H[i][0])/dH[0];
1177 if( B>=0 && B<=1 ) return true; // point on boundary
1178 }
1179 }else{
1180 dD[0] = x[0]-H[i][0];
1181 dD[1] = x[1]-H[i][1];
1182 A = dD[1]/dH[1];
1183 if( A > 1 || A <= 0) continue;
1184 B = (dH^dD)/dH[1];
1185 if(B==0){ // on the boundary
1186 return true; // boundaries are always in
1187 }else if(B>0){
1188 if(dH[1] > 0) ++w;
1189 else --w;
1190 }
1191 }
1192 }
1193
1194 return w;
1195}
1196
1197
1198/*** \brief Calculate the k nearest neighbors concave hull.
1199
1200 This algorithem is guarenteed to find a curve that serounds an island of points, but not all islands unless check==true.
1201
1202 If it at first fails with the k input, k will increase. The final k relaces the input k.
1203
1204 check determines whether a final check is done to determine whether all the raminaing points are inside the hull. Otherwise it is possible to get multiple disconnected clusters and the hull will surround only one.
1205 */
1206
1207template <typename Ptype>
1208std::vector<Ptype> concaveK(std::vector<Ptype> &points,int &k,bool check=true)
1209{
1210 //std::cout << "finding hull .... ";
1211 if(points.size() <= 3){
1212 return points;
1213 }
1214
1215 if(k < 3) k =3;
1216
1217 size_t npoints = points.size();
1218
1219 // {
1220 // tree.pop(7);
1221 // RandomNumbers_NR ran(123);
1222 // std::vector<double> radii;
1223 // std::vector<size_t> neighbors;
1224 // Point_2d point;
1225 // for(int i = 0 ; i<100; ++i){
1226 // point[0] = (1-2*ran());
1227 // point[1] = (1-2*ran());
1228 // tree.NearestNeighbors(point.x,k,radii,neighbors);
1229 // }
1230 // for(auto p : points){
1231 // tree.NearestNeighbors(p.x,k,radii,neighbors);
1232 // }
1233 // }
1234 std::vector<Ptype> hull;
1235
1236 double xmin = points[0][0];
1237 size_t first_point_index = 0;
1238 for(size_t i=0 ; i<npoints ; ++i){
1239 if(points[i][0] < xmin){
1240 xmin = points[i][0];
1241 first_point_index = i;
1242 }
1243 }
1244
1245 Ptype firstpoint = points[first_point_index];
1246
1247 std::vector<size_t> remaining_index;
1248 bool segmented = true;
1249 while(segmented){
1250 bool found=false;
1251 while(!found){ // if hull leaves points out repeat with larger k
1252 hull.resize(0);
1253 TreeSimpleVec<Ptype> tree(points.data(),npoints,2);
1254
1255 if(k>=points.size()){
1256 hull = convex_hull(points);
1257 return hull;
1258 }
1259
1260 hull.push_back(firstpoint);
1261 // point is not popper frum the free so that it can be found at the end
1262
1263 std::vector<double> radii;
1264 std::vector<size_t> neighbors;
1265 Ptype hull_end = hull.back();
1266
1267 Ptype v1;
1268 Ptype last_point = Ptype(hull[0][0],hull[0][1]-1);
1269 Ptype v2;
1270
1271 size_t new_index;
1272 double theta_max=0;
1273 Ptype new_point;
1274 std::vector<double> thetas;
1275 std::vector<int> sorted_index(k);
1276
1277 while((hull[0] != hull.back() && found) || hull.size()==1 ){
1278
1279 tree.NearestNeighbors(hull_end.x,k,radii,neighbors);
1280
1281 long Nneighbors = neighbors.size(); // incase there are not k left
1282 thetas.resize(Nneighbors);
1283 sorted_index.resize(Nneighbors);
1284
1285 v1 = hull_end - last_point;
1286 theta_max = 0;
1287 for(int i=0; i<Nneighbors ; ++i){
1288 v2 = points[neighbors[i]] - hull_end;
1289 double cross = v1^v2,dot = (v1*v2);
1290 if(cross == 0 && dot <= 0) thetas[i] = -PI; // prevents backtracking at beginning
1291 else thetas[i] = atan2( cross , dot );
1292 sorted_index[i] = i;
1293 }
1294
1295 std::sort(sorted_index.begin(),sorted_index.end()
1296 ,[&thetas](int i,int j){return thetas[i] > thetas[j];} );
1297
1298 int trial=0;
1299 bool intersect = true;
1300 while(intersect && trial < Nneighbors){
1301
1302 new_index = neighbors[ sorted_index[trial] ];
1303 new_point = points[ new_index ];
1304
1305 if(hull.size() <= 3) intersect = false;
1306
1307 // check that new edge doesn't cross earlier edge
1308 for(long i = hull.size() - 2 ; i > 1 ; --i){
1309 intersect = segments_cross(hull[i],hull[i-1]
1310 ,hull.back(),new_point);
1311
1312 if(intersect){
1313 break;
1314 }
1315 }
1316
1317 if(new_index == first_point_index && !intersect){
1318 // check that neighbors are inside hull that would be closed
1319 // this is to prevent pre-mature closing of the loop
1320 hull.push_back(new_point);
1321 for(size_t i : neighbors){
1322 if(i != new_index && !inCurve(points[i],hull)){
1323 intersect = true;
1324 break;
1325 }
1326 }
1327 hull.pop_back();
1328 }
1329
1330 ++trial;
1331 }
1332
1333 if(intersect){
1334 // std::cout << "Intersection ..." << k << std::endl;
1335 //
1336 // std::ofstream logfile("testpoint.csv");
1337 // for(auto &p : points){
1338 // logfile << p[0] <<","<< p[1] << std::endl;
1339 // }
1340 // logfile.close();
1341 //
1342 // logfile.open("testhull.csv");
1343 // for(auto &p : hull){
1344 // logfile << p[0] <<","<< p[1] << std::endl;
1345 // }
1346 // logfile << new_point[0] <<","<< new_point[1] << std::endl;
1347 // logfile.close();
1348
1349 k *= 2;
1350 found = false;
1351 }else{
1352 hull.push_back(new_point);
1353 last_point = hull_end;
1354 hull_end = hull.back();
1355 //tree.print();
1356 tree.pop(new_index);
1357 //tree.print();
1358 found = true;
1359
1360 // {
1361 // Ptype center;
1362 // PixelMap map(center.x, 512, 3. / 512.);
1363 //
1364 // map.drawCurve(hull,1);
1365 // map.drawPoints(points,0,2);
1366 // map.printFITS("!concavek_test"+ std::to_string(f) +".fits");
1367 // std::cout << "Test plot : concavek_test" << f++ << ".fits" << std::endl;
1368 // }
1369 }
1370 }
1371 remaining_index = tree.get_index();
1372 }
1373 // test if all remaining points are in the hull
1374
1375
1376 segmented = false;
1377 if(check){
1378 for(auto i : remaining_index){
1379 if(!inCurve(points[i],hull)){
1380 segmented = true;
1381 k *= 2;
1382 // std::cout << "point outside ... ";
1383 // std::ofstream logfile("testpoint.csv");
1384 // for(auto &p : points){
1385 // logfile << p[0] <<","<< p[1] << std::endl;
1386 // }
1387 // logfile.close();
1388 //
1389 // logfile.open("testhull.csv");
1390 // for(auto &p : hull){
1391 // logfile << p[0] <<","<< p[1] << std::endl;
1392 // }
1393 // logfile.close();
1394 //
1395
1396 break;
1397 }
1398 }
1399 }
1400 }
1401
1402 //std::cout << "found hull" << std::endl;
1403 //hull.pop_back();
1404
1405 return hull;
1406}
1407
1408
1409template <typename Ptype>
1410void testconcaveK(){
1411 std::vector<Ptype> points(200);
1412 RandomNumbers_NR ran(2312);
1413
1414 for(Ptype &p : points){
1415 p[0] = (1-ran());
1416 p[1] = (1-2*ran());
1417 }
1418
1419 double x=-1;
1420 for(int i=0 ; i< points.size()/4 ; ++i){
1421 points[i][0] = 1;
1422 points[i][1] = x;
1423 x += 1/25.;
1424 }
1425 x=-1;
1426 for(int i=points.size()/4 ; i< 2*points.size()/4 ; ++i){
1427 points[i][0] = -1;
1428 points[i][1] = x;
1429 x += 1/25.;
1430 }
1431 x=-1;
1432 for(int i=2*points.size()/4 ; i< 3*points.size()/4 ; ++i){
1433 points[i][0] = x;
1434 points[i][1] = 1;
1435 x += 1/25.;
1436 }
1437 x=-1;
1438 for(int i=3*points.size()/4 ; i< points.size() ; ++i){
1439 points[i][0] = x;
1440 points[i][1] = -1;
1441 x += 1/25.;
1442 }
1443
1444 // for(int i=0 ; i< points.size()/2 ; ++i){
1445 // points[i][0] = (1-ran()) - 2;
1446 // points[i][1] = (1-2*ran());
1447 // }
1448 //
1449 //
1450 // points[0][0] = 0; points[0][1] = -1;
1451 // points[1][0] = -1; points[1][1] = 0;
1452 // points[2][0] = 0; points[2][1] = 1;
1453 // points[3][0] = 1; points[3][1] = 0;
1454 //
1455 // points[4][0] = -1; points[4][1] = 1;
1456 // points[5][0] = 1; points[5][1] = 1;
1457
1458 int k = 5;
1459 std::vector<Ptype> hull = concaveK(points,k);
1460 Ptype center;
1461 PixelMap<float> map(center.x, 512, 3. / 512.);
1462
1463 map.drawCurve(hull,1);
1464 map.drawPoints(points,0,2);
1465 map.printFITS("!concavek_test.fits");
1466 std::cout << "Test plot : concavek_test.fits" << std::endl;
1467}
1468
1470Point_2d line_intersection(const Point_2d &v1,const Point_2d &v2,
1471 const Point_2d &w1,const Point_2d &w2);
1472
1473
1475bool circleIntersetsCurve(const Point_2d &x,double r,const std::vector<Point_2d> &v);
1476
1479bool circleOverlapsCurve(const Point_2d &x,double r,const std::vector<Point_2d> &v);
1480
1489std::vector<Point_2d> envelope(const std::vector<Point_2d> &v
1490 ,const std::vector<Point_2d> &w);
1491
1497std::vector<Point_2d> envelope2(const std::vector<Point_2d> &v
1498 ,const std::vector<Point_2d> &w);
1499
1507std::vector<std::vector<Point_2d> > thicken_poly(
1508 const std::vector<Point_2d> &v
1509 ,double R);
1510
1511}
1512#endif /* concave_hull_h */
Image structure that can be manipulated and exported to/from fits files.
Definition pixelmap.h:42
void printFITS(std::string filename, bool Xflip=false, bool ctype=false, bool verbose=false)
Output the pixel map as a fits file.
Definition pixelmap.cpp:1318
void drawgrid(int N, PosType value)
draw a grid on the image that divides the each demension into N cells
Definition pixelmap.cpp:1717
A tree for doing quick searches in multidimensional space. A pointer to an array of objects type T is...
Definition simpleTreeVec.h:22
void pop(IndexType index_of_point)
pop this point out of the tree
Definition simpleTreeVec.h:61
std::vector< IndexType > get_index() const
returns the index numbers of the remaining points which are less than the original if pop() was used.
Definition simpleTreeVec.h:105
void NearestNeighbors(D *ray, int Nneighbors, std::vector< D > &radii, std::vector< IndexType > &neighbors)
Finds the nearest N neighbors and puts their index numbers in an array, also returns the distance to ...
Definition simpleTreeVec.h:464
This is a class for generating random numbers. It simplifies and fool proofs initialization and allow...
Definition utilities_slsim.h:1059
Definition utilities.h:39
size_t RemoveIntersections(std::vector< T > &curve)
removes the intersections of the curve
Definition concave_hull.h:32
Point_2d RandomInPoly(std::vector< Point_2d > &pp, Utilities::RandomNumbers_NR &ran)
return a point within a polygon that doesn't need to be convex
Definition curve_routines.cpp:3667
void find_boundaries(std::vector< bool > &bitmap, long nx, std::vector< std::vector< P > > &points, std::vector< bool > &hits_edge, bool add_to_vector=false, bool outer_only=false)
finds ordered boundaries to regions where bitmap == true
Definition concave_hull.h:201
double interior_mass(const std::vector< Point_2d > &alpha, const std::vector< Point_2d > &x)
this returns area within the curve x average kappa iwithin the curve
Definition curve_routines.cpp:2078
std::vector< Point_2d > envelope(const std::vector< Point_2d > &v, const std::vector< Point_2d > &w)
Find a curve that is made up of segments from v and w that surrounds them and does not self intersect...
Definition curve_routines.cpp:3126
void concave(std::vector< T > &init_points, std::vector< T > &hull_out, double scale)
Creates the concave hull of a group of 2 dimensional points by the shrink-wrap algorithm.
Definition concave_hull.h:633
std::vector< Point_2d > TighterHull(const std::vector< Point_2d > &v)
Definition curve_routines.cpp:3275
bool inhull< RAY >(PosType x[], const std::vector< RAY > &H)
finds in x is within the curve discribed by the H[].x points ie image points
Definition concave_hull.h:1131
std::vector< Point_2d > envelope2(const std::vector< Point_2d > &v, const std::vector< Point_2d > &w)
Find a curve that is made up of segments from v and w that surrounds them and does not self intersect...
Definition curve_routines.cpp:3067
std::vector< std::vector< Point_2d > > thicken_poly(const std::vector< Point_2d > &v, double R)
return the boundaries of the region that is within R of the curve v
Definition curve_routines.cpp:2966
Point_2d RandomInTriangle(const Point_2d &x1, const Point_2d &x2, const Point_2d &x3, Utilities::RandomNumbers_NR &ran)
returns random point within a trinagle
Definition curve_routines.cpp:3600
std::vector< Ptype > concaveK(std::vector< Ptype > &points, int &k, bool check=true)
Definition concave_hull.h:1208
std::vector< T > convex_hull(const std::vector< T > &PP)
Returns a vector of points on the convex hull in counter-clockwise order.
Definition concave_hull.h:513
std::vector< Point_2d > TightestHull(const std::vector< Point_2d > &v)
Finds a concave envelope for an arbitrary closed curve. This is done by gridding and then finding poi...
Definition curve_routines.cpp:3519
Point_2d RandomInConvexPoly(const std::vector< Point_2d > &pp, Utilities::RandomNumbers_NR &ran)
return a point within a convex polygon
Definition curve_routines.cpp:3614
void find_islands(std::vector< bool > &bitmap, long nx, std::vector< std::vector< long > > &indexes, std::vector< bool > &hits_edge, bool add_to_vector=false)
find the indexes of areas with bitmap[]=true broken up into diconnected islands
Definition curve_routines.cpp:2670
bool circleOverlapsCurve(const Point_2d &x, double r, const std::vector< Point_2d > &v)
Returns true if there is any overlap between a curve and a circle. This includea cases where one comp...
Definition curve_routines.cpp:2955
Point_2d line_intersection(const Point_2d &v1, const Point_2d &v2, const Point_2d &w1, const Point_2d &w2)
return the interesetion point of line defined by v1,v2 and w1,w2
Definition curve_routines.cpp:2928
Point_2d RandomNearPoly(std::vector< Point_2d > &pp, double R, Utilities::RandomNumbers_NR &ran)
return a point within distance R of a polygon, i.e. the center of a cicle that intersects the interio...
Definition curve_routines.cpp:3737
bool circleIntersetsCurve(const Point_2d &x, double r, const std::vector< Point_2d > &v)
returns true if a circle of radius r around the point x intersects with the curve v....
Definition curve_routines.cpp:2936
bool inCurve(const Ptype &x, const std::vector< Ptype > &H)
returns true if x is within the polygon H
Definition concave_hull.h:1021
Class for representing points or vectors in 2 dimensions. Not that the dereferencing operator is over...
Definition point.h:48
PosType length() const
length
Definition point.h:134