H3C HEALPix library for PostgreSQL  (version 1.2)
h3c_poly.c
Go to the documentation of this file.
1 /*
2  Copyright (C) 2012 Gilles Landais (CDS)
3 
4  Author: Gilles Landais, Strasbourg astronomical Data Center (CDS)
5  Email: gilles.landais@unistra.fr
6 
7  This file is part of H3C.
8 
9  H3C is free software; you can redistribute it and/or modify
10  it under the terms of the GNU General Public License as published by
11  the Free Software Foundation; either version 2 of the License, or
12  (at your option) any later version.
13 
14  H3C is distributed in the hope that it will be useful,
15  but WITHOUT ANY WARRANTY; without even the implied warranty of
16  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17  GNU General Public License for more details.
18 
19  You should have received a copy of the GNU General Public License
20  along with H3C; if not, write to the Free Software
21  Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
22 */
23 
24 
55 #include "common.h"
56 
57 #if 0
58 void h3c_init_poly(h3c_poly *qp, int n)
59 {
60  qp->ra = malloc(n * sizeof(h3c_coord_t));
61  qp->dec = malloc(n * sizeof(h3c_coord_t));
62 
63  qp->ax = malloc(n * sizeof(h3c_coord_t));
64  qp->ay = malloc(n * sizeof(h3c_coord_t));
65 
66  qp->x = malloc(n * sizeof(h3c_coord_t));
67  qp->y = malloc(n * sizeof(h3c_coord_t));
68  qp->n = n;
69 }
70 #endif
71 static void h3c_prepare_poly(h3c_poly *qp)
72 {
73  int n = qp->n - 1 ;
74  int i;
75  h3c_coord_t *ax = qp->ax;
76  h3c_coord_t *ay = qp->ay;
77  h3c_coord_t *x = qp->x;
78  h3c_coord_t *y = qp->y;
79 
80  for(i = 0 ; i < n; i++)
81  {
82  ax[i] = x[i + 1] - x[i];
83  ay[i] = y[i + 1] - y[i];
84  }
85  ax[i] = x[0] - x[i];
86  ay[i] = y[0] - y[i];
87 }
88 
89 /* Cloned version of ang2ipix for outputting also the x,y on the cube face
90  * Coordinates on the cube face are x[-0.5,0.5] y[-0.5,0.5] */
91 static void ang2ipix_xy (/*struct q3c_prm *hprm,*/ h3c_coord_t ra, h3c_coord_t dec,
92  char *out_face_num, /*q3c_ipix_t *ipix, */h3c_coord_t *x_out,
93  h3c_coord_t *y_out)
94  /* ra in degrees, dec in degrees */
95  /* strictly 0<=ra<360 and -90<=dec<=90 */
96 {
97  h3c_coord_t x0 = 0,y0 = 0;
98  /*q3c_ipix_t nside = hprm->nside, *xbits = hprm->xbits,
99  *ybits = hprm->ybits, xi, yi, i1;*/
100  char face_num;
101  if (dec >= 90)
102  /* Poles */
103  {
104  face_num = 0;
105  x0 = H3C_HALF;
106  y0 = H3C_HALF;
107  *x_out = 0;
108  *y_out = 0;
109  goto END1;
110  }
111  else if (dec <= -90)
112  {
113  face_num = 5;
114  x0 = H3C_HALF;
115  y0 = H3C_HALF;
116  *x_out = 0;
117  *y_out = 0;
118  goto END1;
119  }
120 
121  face_num = fmod ((ra + 45) / 90, 4);
122  /*for equatorial pixels we'll have face_num from 1 to 4 */
123  x0 = tan (H3C_DEGRA * (ra - 90 * (h3c_coord_t)face_num));
124  y0 = tan (dec * H3C_DEGRA) /
125  cos (H3C_DEGRA * (ra - 90 * (h3c_coord_t)face_num));
126  face_num++;
127 
128  if (y0 > 1)
129  {
130  face_num = 0;
131  x0 = sin (H3C_DEGRA * ra) / tan (H3C_DEGRA * dec);
132  y0 = -1*cos (H3C_DEGRA * ra) / tan (H3C_DEGRA * dec);
133  }
134  else if (y0 < -1)
135  {
136  face_num = 5;
137  x0 = -1*sin (H3C_DEGRA * ra) / tan (H3C_DEGRA * dec);
138  y0 = -1*cos (H3C_DEGRA * ra) / tan (H3C_DEGRA * dec);
139  }
140 
141  *x_out = x0 / 2;
142  *y_out = y0 / 2;
143  x0 = (x0 + 1) / 2;
144  y0 = (y0 + 1) / 2;
145 
146  END1:
147 
148  /* Now I produce the final pixel value by converting x and y values to bitfields
149  * and combining them by interleaving, using the predefined arrays xbits and ybits
150  */
151 #if 0
152  xi = (q3c_ipix_t)(x0 * nside);
153  yi = (q3c_ipix_t)(y0 * nside);
154 
155  /* This two strings are written to handle the case of edges of base square */
156  if (xi == nside)
157  {
158  xi--;
159  }
160  if (yi == nside)
161  {
162  yi--;
163  }
164 
165  i1 = 1 << (Q3C_INTERLEAVED_NBITS);
166 
167 
168 #ifdef Q3C_INT4
169  *ipix = ((q3c_ipix_t)face_num) * nside * nside + xbits[xi % i1] +
170  ybits[yi % i1];
171  /*4byte computation*/
172 #endif /* Q3C_INT4 */
173 #ifdef Q3C_INT8
174  *ipix = ((q3c_ipix_t)face_num) * nside * nside + xbits[xi % i1] +
175  ybits[yi % i1] + (xbits[(xi >> Q3C_INTERLEAVED_NBITS) % i1] +
176  ybits[(yi >> Q3C_INTERLEAVED_NBITS) % i1]) * i1 * i1;
177  /*8byte computation*/
178 #endif /* Q3C_INT8 */
179 #endif
180  *out_face_num = face_num;
181 }
182 
183 static int h3c_check_point_in_poly(h3c_poly *qp, h3c_coord_t x0,
184  h3c_coord_t y0)
185 /* Implementation of the crossing algorithm */
186 {
187  int i, n = qp->n;
188  h3c_coord_t *y = qp->y;
189  h3c_coord_t *x = qp->x;
190  h3c_coord_t *ax = qp->ax;
191  h3c_coord_t *ay = qp->ay;
192  int result = !H3C_DISJUNCT;
193  for(i=0;i<n;i++)
194  {
195  if (((y0<=y[i])==(y0>y[(i+1)%n])) &&
196  ((x0-x[i])<(y0-y[i])*ax[i]/ay[i]))
197  {
198  result=!result;
199  }
200  }
201  return !result;
202 }
203 
204 static void h3c_get_minmax_poly(h3c_poly *qp, h3c_coord_t *xmin,
205  h3c_coord_t *xmax, h3c_coord_t *ymin,
206  h3c_coord_t *ymax)
207 {
208  int i;
209  const int n = qp->n;
210  h3c_coord_t *x = qp->x, *y = qp->y, t;
211  h3c_coord_t xmi, xma, ymi, yma;
212 
213  xmi = x[0];
214  xma = x[0];
215  ymi = y[0];
216  yma = y[0];
217 
218  for(i = 1; i < n; i++)
219  {
220  t = x[i];
221  if (t > xma)
222  {
223  xma = t;
224  }
225  else if (t < xmi)
226  {
227  xmi = t;
228  }
229 
230  t = y[i];
231 
232  if (t > yma)
233  {
234  yma = t;
235  }
236  else if (t < ymi)
237  {
238  ymi = t;
239  }
240  }
241 
242  *xmin = xmi;
243  *xmax = xma;
244  *ymin = ymi;
245  *ymax = yma;
246 }
247 
248 
249 static char h3c_get_facenum(h3c_coord_t ra, h3c_coord_t dec)
250  /* ra in degrees, dec in degrees */
251  /* strictly 0<=ra<360 and -90<=dec<=90 */
252 {
253  h3c_coord_t y0 = 0;
254  char face_num;
255 
256  if (dec >= 90)
257  /* Poles */
258  {
259  return 0;
260  }
261  else if (dec <= -90)
262  {
263  return 5;
264  }
265 
266  face_num = fmod ((ra + 45) / 90, 4);
267  /*for equatorial pixels we'll have face_num from 1 to 4 */
268 
269  y0 = tan(dec * H3C_DEGRA) /
270  cos(H3C_DEGRA * (ra - 90 * (h3c_coord_t)face_num));
271 
272  face_num++;
273 
274  if (y0 > 1)
275  {
276  return 0;
277  }
278  else if (y0 < -1)
279  {
280  return 5;
281  }
282  else
283  {
284  return face_num;
285  }
286 }
287 static char h3c_xy2facenum(h3c_coord_t x, h3c_coord_t y, char face_num0)
288 {
289  /* The input x, y should be >=-1 and <=1 */
290 
291 
292  h3c_coord_t ra = 0, dec = 0;
293  /* I do the initialization since gcc warn about probable not
294  * initialization of ra and dec
295  */
296 
297  /* This code have been cutted out from ipix2ang BEGIN */
298  if ((face_num0 >= 1) && (face_num0 <= 4))
299  {
300  ra = atan(x);
301  dec = H3C_RADEG * atan(y * cos(ra));
302  ra = ra * H3C_RADEG + ((h3c_coord_t)face_num0 - 1) * 90;
303  if (ra < 0)
304  {
305  ra += (h3c_coord_t)360;
306  }
307  }
308  else
309  {
310  if (face_num0 == 0)
311  {
312  ra = H3C_RADEG * atan2(x, -y);
313  dec = H3C_RADEG * atan(1 / sqrt(x * x + y * y));
314  if (ra < 0)
315  {
316  ra += (h3c_coord_t)360;
317  }
318  }
319  if (face_num0 == 5)
320  {
321  ra = H3C_RADEG * atan2(x, y);
322  dec = -H3C_RADEG * atan(1 / sqrt(x * x + y * y));
323  if (ra < 0)
324  {
325  ra += (h3c_coord_t)360;
326  }
327 
328  }
329  }
330  /* This code have been cutted out from ipix2ang END */
331 
332  return h3c_get_facenum(ra,dec);
333 }
334 
335 static char h3c_get_facenum_poly(h3c_poly *qp)
336 {
337  return h3c_get_facenum(qp->ra[0], qp->dec[0]);
338 }
339 
340 static void h3c_project_poly(h3c_poly *qp, char face_num)
341 {
342  h3c_coord_t ra1, dec1, tmp0;
343  h3c_coord_t *ra = qp->ra, *dec = qp->dec;
344  h3c_coord_t *x = qp->x, *y = qp->y, x0, y0;
345 
346  int i, n = qp->n;
347  if ((face_num > 0) && (face_num < 5))
348  {
349  face_num--; /* Just computation trick */
350  for (i = 0; i < n; i++)
351  {
352  ra1 = H3C_DEGRA * (ra[i] - 90 * (h3c_coord_t)face_num);
353  dec1 = H3C_DEGRA * dec[i];
354  x[i] = (tan(ra1)) / 2;
355  y[i] = (tan(dec1) / cos(ra1)) / 2;
356  }
357  /* Now x[i] and y[i] are coordinates on cube face [-0.5:0.5]x[-0.5:0.5] */
358  }
359  else if (face_num == 0)
360  {
361  for (i = 0; i < n; i++)
362  {
363  ra1 = H3C_DEGRA * ra[i];
364  dec1 = H3C_DEGRA * dec[i];
365 
366  tmp0 = 1 / tan(dec1);
367 #ifdef __USE_GNU
368  sincos(ra1, &x0, &y0);
369 #else
370  x0 = sin(ra1);
371  y0 = cos(ra1);
372 #endif
373  x0 *= tmp0;
374  y0 *= (-tmp0);
375  x[i] = x0 / 2;
376  y[i] = y0 / 2;
377  }
378  }
379  else
380  {
381  for (i = 0; i < n; i++)
382  {
383  ra1 = H3C_DEGRA * ra[i];
384  dec1 = H3C_DEGRA * dec[i];
385 
386  tmp0 = 1 / tan(dec1);
387 #ifdef __USE_GNU
388  sincos(ra1, &x0, &y0);
389 #else
390  x0 = sin(ra1);
391  y0 = cos(ra1);
392 #endif
393  x0 *= (-tmp0);
394  y0 *= (-tmp0);
395  x[i] = x0 / 2;
396  y[i] = y0 / 2;
397  }
398  }
399 }
400 #if 0
401 
402 
403 static char h3c_poly_intersection_check(h3c_poly *qp,
404  h3c_coord_t xl, h3c_coord_t xr,
405  h3c_coord_t yb, h3c_coord_t yt,
406  h3c_coord_t cur_size)
407 {
408  int i, n = qp->n;
409  h3c_coord_t *ax = qp->ax;
410  h3c_coord_t *ay = qp->ay;
411  h3c_coord_t *x = qp->x;
412  h3c_coord_t *y = qp->y;
413  h3c_coord_t txl, txr, tyb, tyt, axi, ayi, xi, yi, tmp, tmp1;
414  char ret = 0;
415  for( i = 0; i <n ; i++)
416  {
417  xi = x[i];
418  yi = y[i];
419  axi = ax[i];
420  ayi = ay[i];
421  txl = xl - xi;
422  txr = xr - xi;
423  tyb = yb - yi;
424  tyt = yt - yi;
425 
426  tmp = tyb / ayi;
427  tmp1 = axi * tmp - txl;
428  if ((tmp >= 0) && (tmp <= 1) && (tmp1 >= 0) && (tmp1 <= cur_size))
429  {
430  ret = 1;
431  break;
432  }
433 
434  tmp = tyt / ayi;
435  tmp1 = axi * tmp - txl;
436  if ((tmp >= 0) && (tmp <= 1) && (tmp1 >= 0) && (tmp1 <= cur_size))
437  {
438  ret = 1;
439  break;
440  }
441 
442  tmp = txl / axi;
443  tmp1 = ayi * tmp - tyb;
444  if ((tmp >= 0) && (tmp <= 1) && (tmp1 >= 0) && (tmp1 <= cur_size))
445  {
446  ret = 1;
447  break;
448  }
449 
450  tmp = txr / axi;
451  tmp1 = ayi * tmp - tyb;
452  if ((tmp >= 0) && (tmp <= 1) && (tmp1 >= 0) && (tmp1 <= cur_size))
453  {
454  ret = 1;
455  break;
456  }
457  }
458 
459  return ret;
460 }
461 
462 
463 
464 int h3c_poly_cover_check(h3c_poly *qp, h3c_coord_t xc_cur,
465  h3c_coord_t yc_cur, h3c_coord_t cur_size)
466 {
467  h3c_coord_t xl_cur, xr_cur, yb_cur, yt_cur;
468  int val;
469 
470 
471  /* Checking the intersection of ellipse and box
472  * The box parameters are set by variables xc_cur, yc_cur and cur_size
473  */
474  xl_cur = xc_cur - cur_size / 2; /* left */
475  xr_cur = xc_cur + cur_size / 2; /* right */
476  yb_cur = yc_cur - cur_size / 2; /* bottom */
477  yt_cur = yc_cur + cur_size / 2; /* top */
478 
479  /* Undef labels -- the labels when the current computed values dont allow
480  * to make the final decision about the intersection
481  */
482 
483  val = h3c_check_point_in_poly(qp, xl_cur, yb_cur);
484  if (val != H3C_DISJUNCT)
485  {
486  goto PARTUNDEF_CHECK01;
487  }
488 
489  val = h3c_check_point_in_poly(qp, xr_cur, yb_cur);
490  if (val != H3C_DISJUNCT)
491  {
492  return H3C_PARTIAL;
493  }
494 
495  val = h3c_check_point_in_poly(qp, xr_cur, yt_cur);
496  if (val != H3C_DISJUNCT)
497  {
498  return H3C_PARTIAL;
499  }
500 
501  val = h3c_check_point_in_poly(qp, xl_cur, yt_cur);
502 
503  if (val != H3C_DISJUNCT)
504  {
505  return H3C_PARTIAL;
506  }
507  else
508  {
509  if (h3c_poly_intersection_check(qp, xl_cur, xr_cur, yb_cur, yt_cur, cur_size)||
510  ((qp->x[0] > xl_cur) && (qp->x[0] < xr_cur) &&
511  (qp->y[0] > yb_cur) && (qp->y[0] < yt_cur)))
512  {
513  return H3C_PARTIAL;
514  }
515  else
516  {
517  return H3C_DISJUNCT;
518  }
519  }
520 
521 
522  PARTUNDEF_CHECK01:
523  val = h3c_check_point_in_poly(qp, xr_cur, yb_cur);
524  if (val == H3C_DISJUNCT)
525  {
526  return H3C_PARTIAL;
527  }
528 
529  //PARTUNDEF_CHECK11:
530  val = h3c_check_point_in_poly(qp, xr_cur, yt_cur);
531  if (val == H3C_DISJUNCT)
532  {
533  return H3C_PARTIAL;
534  }
535 
536 
537  //PARTUNDEF_CHECK10:
538  val = h3c_check_point_in_poly(qp, xl_cur, yt_cur);
539  if (val == H3C_DISJUNCT)
540  {
541  return H3C_PARTIAL;
542  }
543  else
544  {
545  return H3C_COVER;
546  }
547 }
548 #endif
549 
550 int h3c_check_sphere_point_in_poly(/*struct h3c_prm *hprm,*/ int n,
551  h3c_coord_t in_ra[], h3c_coord_t in_dec[],
552  h3c_coord_t ra0, h3c_coord_t dec0)
553 {
554  h3c_coord_t xmin,xmax,ymin, ymax;
555  static char faces[6], multi_flag;
556 /* h3c_ipix_t ipix;*/
557  h3c_coord_t points[4];
558  char face_num, face_num0, cur_face_num;
559  static h3c_coord_t x[3][H3C_MAX_N_POLY_VERTEX], y[3][H3C_MAX_N_POLY_VERTEX],
560  ax[3][H3C_MAX_N_POLY_VERTEX], ay[3][H3C_MAX_N_POLY_VERTEX], x0, y0;
561 
562  int face_count = -1, i;
563 
564  h3c_poly qp;
565 
566  ang2ipix_xy(/*hprm,*/ ra0, dec0, &cur_face_num,/* &ipix, */&x0, &y0);
567 
568  qp.ra = in_ra;
569  qp.dec = in_dec;
570  qp.n = n;
571 
572  face_num = h3c_get_facenum_poly(&qp);
573  faces[0] = face_num;
574 
575  qp.x = x[0];
576  qp.y = y[0];
577  qp.ax = ax[0];
578  qp.ay = ay[0];
579 
580  h3c_project_poly(&qp, face_num);
581  h3c_prepare_poly(&qp);
582 
583  h3c_get_minmax_poly(&qp, &xmin, &xmax, &ymin, &ymax);
584 
585  /* Now in a little bit ugly but fastest way I determine whether the ellipse
586  * intersect other faces or not, and if yes, I setup the array "points" to the
587  * multi_face loop.
588  */
589 
590  if (xmin < -H3C_HALF)
591  {
592  if (ymin < -H3C_HALF)
593  {
594  points[0] = xmax;
595  points[1] = ymin;
596  points[2] = xmin;
597  points[3] = ymax;
598  multi_flag = 2;
599  xmin = -H3C_HALF;
600  ymin = -H3C_HALF;
601  }
602  else
603  {
604  if (ymax > H3C_HALF)
605  {
606  points[0] = xmax;
607  points[1] = ymax;
608  points[2] = xmin;
609  points[3] = ymin;
610  multi_flag = 2;
611  xmin = -H3C_HALF;
612  ymax = H3C_HALF;
613  }
614  else
615  {
616  points[0] = xmin;
617  points[1] = (ymin + ymax) / 2;
618  multi_flag = 1;
619  xmin = -H3C_HALF;
620  }
621  }
622  }
623  else
624  {
625  if (xmax > H3C_HALF)
626  {
627  if (ymin < -H3C_HALF)
628  {
629  points[0] = xmin;
630  points[1] = ymin;
631  points[2] = xmax;
632  points[3] = ymax;
633  multi_flag = 2;
634  xmax = H3C_HALF;
635  ymin = -H3C_HALF;
636  }
637  else
638  {
639  if (ymax > H3C_HALF)
640  {
641  points[0] = xmin;
642  points[1] = ymax;
643  points[2] = xmax;
644  points[3] = ymax;
645  multi_flag = 2;
646  xmax = H3C_HALF;
647  ymax = H3C_HALF;
648  }
649  else
650  {
651  points[0] = xmax;
652  points[1] = (ymin + ymax) / 2;
653  multi_flag = 1;
654  xmax = H3C_HALF;
655  }
656  }
657  }
658  else
659  {
660  if (ymin < -H3C_HALF)
661  {
662  points[0] = (xmin + xmax) / 2;
663  points[1] = ymin;
664  multi_flag = 1;
665  ymin = -H3C_HALF;
666  }
667  else
668  {
669  if (ymax > H3C_HALF)
670  {
671  points[0] = (xmin + xmax) / 2;
672  points[1] = ymax;
673  multi_flag = 2;
674  ymax = H3C_HALF;
675  }
676  else
677  {
678  multi_flag = 0;
679  }
680  }
681  }
682  }
683 
684  face_num0 = face_num;
685 
686  for(face_count = 0; face_count <= multi_flag; face_count++)
687  {
688  /* This the beginning of the mega-loop over multiple faces */
689 
690  if (face_count > 0)
691  /* This "if" works when we pass through the secondary faces */
692  {
693  face_num = h3c_xy2facenum(2 * points[2 * face_count - 2],
694  2 * points[2 * face_count - 1], face_num0);
695 
696  faces[face_count] = face_num;
697 
698  qp.x = x[face_count];
699  qp.y = y[face_count];
700  qp.ax = ax[face_count];
701  qp.ay = ay[face_count];
702 
703  h3c_project_poly(&qp, faces[face_count]);
704  h3c_prepare_poly(&qp);
705  }
706  }
707 
708  for (i = 0; i <= multi_flag; i++)
709  {
710  if (faces[i] == cur_face_num)
711  {
712  face_count = i;
713  break;
714  }
715  }
716 
717  if (i == (multi_flag + 1))
718  {
719  return 0;
720  }
721 
722  qp.x = x[face_count];
723  qp.y = y[face_count];
724  qp.ax = ax[face_count];
725  qp.ay = ay[face_count];
726 
727  return h3c_check_point_in_poly(&qp, x0, y0);
728 }
729 
common.h