source: Sophya/trunk/SophyaLib/SkyMap/HEALPixUtils.cc@ 4042

Last change on this file since 4042 was 2985, checked in by ansari, 19 years ago

1/ Suite codage HEALPix-NEST , test avec transforme Ylm OK
2/ Petites modifs et corrections diverses (ajout SizeIndex2Resol() par exemple)

Reza 21 Juin 2006

File size: 22.3 KB
Line 
1// utilitaires de pixelisation HEALPix
2#include "sopnamsp.h"
3#include "HEALPixUtils.h"
4#include <iostream>
5#include <math.h>
6//#include <complex>
7#include "tvector.h"
8#include "smathconst.h"
9extern "C"
10{
11#include <stdio.h>
12#include <stdlib.h>
13#include <unistd.h>
14}
15
16using namespace SOPHYA;
17
18//////////////////////////////////////////////////////////////////////////
19//
20// ------------- Classe PIXELS_XY -----------------------
21//
22class PIXELS_XY
23{
24
25public :
26
27static PIXELS_XY& instance();
28
29NDataBlock<int_4> pix2x_;
30NDataBlock<int_4> pix2y_;
31NDataBlock<int_4> x2pix_;
32NDataBlock<int_4> y2pix_;
33
34private :
35
36PIXELS_XY();
37void mk_pix2xy();
38void mk_xy2pix();
39};
40
41
42
43//*******************************************************************
44// Class PIXELS_XY
45// Construction des tableaux necessaires a la traduction des indices RING en
46// indices NESTED (ou l'inverse)
47//*******************************************************************
48
49PIXELS_XY::PIXELS_XY()
50{
51 pix2x_.ReSize(1024);
52 pix2x_.Reset();
53 pix2y_.ReSize(1024);
54 pix2y_.Reset();
55 x2pix_.ReSize(128);
56 x2pix_.Reset();
57 y2pix_.ReSize(128);
58 y2pix_.Reset();
59 mk_pix2xy();
60 mk_xy2pix();
61}
62
63// Instance unique de la classe PIXELS_XY
64static PIXELS_XY * _singleton = NULL;
65
66PIXELS_XY& PIXELS_XY::instance()
67{
68 if (_singleton == NULL) _singleton = new PIXELS_XY ;
69 return (*_singleton);
70}
71
72void PIXELS_XY::mk_pix2xy()
73{
74 /*
75 ==================================================
76 subroutine mk_pix2xy
77 ==================================================
78 c constructs the array giving x and y in the face from pixel number
79 c for the nested (quad-cube like) ordering of pixels
80 c
81 c the bits corresponding to x and y are interleaved in the pixel number
82 c one breaks up the pixel number by even and odd bits
83 ==================================================
84 */
85 // tranlated from FORTRAN (Gorski) to C, by B. Revenu, revised Guy Le Meur
86 // (16/12/98)
87
88 int kpix, jpix, IX, IY, IP, ID;
89
90 for(kpix = 0; kpix < 1024; kpix++)
91 {
92 jpix = kpix;
93 IX = 0;
94 IY = 0;
95 IP = 1 ;// ! bit position (in x and y)
96 while( jpix!=0 )
97 { // ! go through all the bits
98 ID=jpix%2;// ! bit value (in kpix), goes in ix
99 jpix = jpix/2;
100 IX = ID*IP+IX;
101
102 ID=jpix%2;// ! bit value (in kpix), goes in iy
103 jpix = jpix/2;
104 IY = ID*IP+IY;
105
106 IP = 2*IP;// ! next bit (in x and y)
107 }
108 pix2x_(kpix) = IX;// ! in 0,31
109 pix2y_(kpix) = IY;// ! in 0,31
110 }
111}
112
113void PIXELS_XY::mk_xy2pix()
114{
115 /*
116 =================================================
117 subroutine mk_xy2pix
118 =================================================
119 c sets the array giving the number of the pixel lying in (x,y)
120 c x and y are in {1,128}
121 c the pixel number is in {0,128**2-1}
122 c
123 c if i-1 = sum_p=0 b_p * 2^p
124 c then ix = sum_p=0 b_p * 4^p
125 c iy = 2*ix
126 c ix + iy in {0, 128**2 -1}
127 =================================================
128 */
129 // tranlated from FORTRAN (Gorski) to C, by B. Revenu, revised Guy Le Meur
130 // (16/12/98)
131
132 int K,IP,I,J,ID;
133 for(I = 1; I <= 128; I++)
134 {
135 J = I-1;// !pixel numbers
136 K = 0;//
137 IP = 1;//
138 truc : if( J==0 )
139 {
140 x2pix_(I-1) = K;
141 y2pix_(I-1) = 2*K;
142 }
143 else
144 {
145 ID = J%2; //RzModFloor (int)fmod(J,2);
146 J = J/2;
147 K = IP*ID+K;
148 IP = IP*4;
149 goto truc;
150 }
151 }
152}
153
154
155
156int_4 HEALPix::nest2ring(int_4 nside, int_4 ipnest)
157{
158 /*
159 ====================================================
160 subroutine nest2ring(nside, ipnest, ipring)
161 ====================================================
162 c conversion from NESTED to RING pixel number
163 ====================================================
164 */
165 // tranlated from FORTRAN (Gorski) to C, by B. Revenu, revised Guy Le Meur
166 // (16/12/98)
167
168 const PIXELS_XY& PXY= PIXELS_XY::instance();
169
170 int npix, npface, face_num, ncap, n_before;
171 int ipf, ip_low, ip_trunc, ip_med, ip_hi;
172 int ix, iy, jrt, jr, nr, jpt, jp, kshift, nl4;
173 int ns_max=8192;
174 int jrll[12]={2, 2, 2, 2, 3, 3, 3, 3, 4, 4, 4, 4};
175 int jpll[12]={1, 3, 5, 7, 0, 2, 4, 6, 1, 3, 5, 7};
176
177 if( nside<1 || nside>ns_max ) {
178 char buff[64];
179 sprintf(buff,"HEALPix::nest2ring nside(=%d) out of range ", nside);
180 throw RangeCheckError(PExcLongMessage(buff));
181 }
182 npix = 12 * nside* nside;
183 if( ipnest<0 || ipnest>npix-1 ) {
184 char buff[64];
185 sprintf(buff,"HEALPix::nest2ring ipnest(=%d) out of range ", ipnest);
186 throw RangeCheckError(PExcLongMessage(buff));
187 }
188
189 ncap = 2* nside*( nside-1);// ! number of points in the North Polar cap
190 nl4 = 4* nside;
191
192 //c finds the face, and the number in the face
193 npface = nside* nside;
194 //cccccc ip = ipnest - 1 ! in {0,npix-1}
195
196 face_num = ipnest/npface;// ! face number in {0,11}
197 ipf =ipnest%npface;// ! pixel number in the face {0,npface-1}
198 //c finds the x,y on the face (starting from the lowest corner)
199 //c from the pixel number
200 ip_low=ipf%1024; // ! content of the last 10 bits
201 ip_trunc = ipf/1024; // ! truncation of the last 10 bits
202 ip_med=ip_trunc%1024; // ! content of the next 10 bits
203 ip_hi = ip_trunc/1024;// ! content of the high weight 10 bits
204
205 ix = 1024*PXY.pix2x_(ip_hi)+32*PXY.pix2x_(ip_med)+PXY.pix2x_(ip_low);
206 iy = 1024*PXY.pix2y_(ip_hi)+32*PXY.pix2y_(ip_med)+PXY.pix2y_(ip_low);
207
208 //c transforms this in (horizontal, vertical) coordinates
209 jrt = ix + iy;// ! 'vertical' in {0,2*(nside-1)}
210 jpt = ix - iy;// ! 'horizontal' in {-nside+1,nside-1}
211
212 //c computes the z coordinate on the sphere
213 // jr = jrll[face_num+1]*nside - jrt - 1;// ! ring number in {1,4*nside-1}
214 jr = jrll[face_num]*nside - jrt - 1;
215 nr = nside;// ! equatorial region (the most frequent)
216 n_before = ncap + nl4 * (jr - nside);
217 kshift=(jr - nside)%2;
218 if( jr<nside ) {//then ! north pole region
219 nr = jr;
220 n_before = 2 * nr * (nr - 1);
221 kshift = 0;
222 }
223 else if( jr>3*nside ) {//then ! south pole region
224 nr = nl4 - jr;
225 n_before = npix - 2 * (nr + 1) * nr;
226 kshift = 0;
227 }
228
229 //c computes the phi coordinate on the sphere, in [0,2Pi]
230 jp = (jpll[face_num]*nr + jpt + 1 + kshift)/2;// ! 'phi' number in the ring in {1,4*nr}
231
232 if( jp>nl4 ) jp = jp - nl4;
233 if( jp<1 ) jp = jp + nl4;
234
235 //RzDel int aux=n_before + jp - 1;
236 return (n_before + jp - 1);// ! in {0, npix-1}
237}
238
239
240int_4 HEALPix::ring2nest(int_4 nside, int_4 ipring)
241{
242 /*
243 ==================================================
244 subroutine ring2nest(nside, ipring, ipnest)
245 ==================================================
246 c conversion from RING to NESTED pixel number
247 ==================================================
248 */
249 // tranlated from FORTRAN (Gorski) to C, by B. Revenu, revised Guy Le Meur
250 // (16/12/98)
251
252 const PIXELS_XY& PXY= PIXELS_XY::instance();
253
254 double fihip, hip;
255 int npix, nl2, nl4, ncap, ip, iphi, ipt, ipring1;
256 int kshift, face_num, nr;
257 int irn, ire, irm, irs, irt, ifm , ifp;
258 int ix, iy, ix_low, ix_hi, iy_low, iy_hi, ipf;
259 int ns_max(8192);
260
261 // coordinate of the lowest corner of each face
262 int jrll[12]={2, 2, 2, 2, 3, 3, 3, 3, 4, 4, 4, 4};// ! in unit of nside
263 int jpll[12]={1, 3, 5, 7, 0, 2, 4, 6, 1, 3, 5, 7};//! in unit of nside/2
264
265 if( nside<1 || nside>ns_max ) {
266 char buff[64];
267 sprintf(buff,"HEALPix::ring2nest nside(=%d) out of range ", nside);
268 throw RangeCheckError(PExcLongMessage(buff));
269 }
270 npix = 12 * nside*nside;
271 if( ipring<0 || ipring>npix-1 ) {
272 char buff[64];
273 sprintf(buff,"HEALPix::ring2nest ipring(=%d) out of range ", ipring);
274 throw RangeCheckError(PExcLongMessage(buff));
275 }
276
277 nl2 = 2*nside;
278 nl4 = 4*nside;
279 npix = 12*nside*nside;// ! total number of points
280 ncap = 2*nside*(nside-1);// ! points in each polar cap, =0 for nside =1
281 ipring1 = ipring + 1;
282
283 //c finds the ring number, the position of the ring and the face number
284 if( ipring1<=ncap ) {//then
285
286 hip = ipring1/2.;
287 fihip = floor ( hip );
288 irn = (int)floor( sqrt( hip - sqrt(fihip) ) ) + 1;// ! counted from North pole
289 iphi = ipring1 - 2*irn*(irn - 1);
290
291 kshift = 0;
292 nr = irn ;// ! 1/4 of the number of points on the current ring
293 face_num = (iphi-1) / irn;// ! in {0,3}
294 }
295 else if( ipring1<=nl2*(5*nside+1) ) {//then
296
297 ip = ipring1 - ncap - 1;
298 irn = ip / nl4 + nside; //RzModFloor (int)floor( ip / nl4 ) + nside; ! counted from North pole
299 iphi = ip%nl4 + 1; //RzModFloor (int)fmod(ip,nl4) + 1;
300
301 kshift = (irn+nside)%2; //RzModFloor (int)fmod(irn+nside,2); ! 1 if irn+nside is odd, 0 otherwise
302 nr = nside;
303 ire = irn - nside + 1;// ! in {1, 2*nside +1}
304 irm = nl2 + 2 - ire;
305 ifm = (iphi - ire/2 + nside -1) / nside;// ! face boundary
306 ifp = (iphi - irm/2 + nside -1) / nside;
307 if( ifp==ifm ) {//then ! faces 4 to 7
308 face_num = ifp%4+4; // RzModFloor (int)fmod(ifp,4) + 4;
309 }
310 else if( ifp + 1==ifm ) {//then ! (half-)faces 0 to 3
311 face_num = ifp;
312 }
313 else if( ifp - 1==ifm ) {//then ! (half-)faces 8 to 11
314 face_num = ifp + 7;
315 }
316 }
317 else {
318
319 ip = npix - ipring1 + 1;
320 hip = ip/2.;
321 fihip = floor ( hip );
322 irs = (int)floor( sqrt( hip - sqrt(fihip) ) ) + 1;// ! counted from South pole
323 iphi = 4*irs + 1 - (ip - 2*irs*(irs-1));
324
325 kshift = 0;
326 nr = irs;
327 irn = nl4 - irs;
328 face_num = (iphi-1) / irs + 8;// ! in {8,11}
329 }
330
331 //c finds the (x,y) on the face
332 irt = irn - jrll[face_num]*nside + 1;// ! in {-nside+1,0}
333 ipt = 2*iphi - jpll[face_num]*nr - kshift - 1;// ! in {-nside+1,nside-1}
334
335
336 if( ipt>=nl2 ) ipt = ipt - 8*nside;// ! for the face #4
337
338 ix = (ipt - irt ) / 2;
339 iy = -(ipt + irt ) / 2;
340
341 ix_low = ix%128; //RzModFloor (int)fmod(ix,128);
342 ix_hi = ix/128;
343 iy_low = iy%128; //RzModFloor (int)fmod(iy,128);
344 iy_hi = iy/128;
345 ipf=(PXY.x2pix_(ix_hi)+PXY.y2pix_(iy_hi))*(128*128)+(PXY.x2pix_(ix_low)+PXY.y2pix_(iy_low));
346
347 return (ipf + face_num* nside *nside);// ! in {0, 12*nside**2 - 1}
348}
349
350int_4 HEALPix::ang2pix_ring(int_4 nside, double theta, double phi)
351{
352 /*
353 ==================================================
354 c gives the pixel number ipix (RING)
355 c corresponding to angles theta and phi
356 c==================================================
357 */
358 // tranlated from FORTRAN (Gorski) to C, by B. Revenu, revised Guy Le Meur
359 // (16/12/98)
360
361 int nl2, nl4, ncap, npix, jp, jm, ipix1;
362 double z, za, tt, tp, tmp;
363 int ir, ip, kshift;
364
365 double piover2(Pi/2.);
366 double twopi(2.*Pi);
367 double z0(2./3.);
368 int ns_max(8192);
369
370 if( nside<1 || nside>ns_max ) {
371 char buff[64];
372 sprintf(buff,"HEALPix::ang2pix_ring nside(=%d) out of range ", nside);
373 throw RangeCheckError(PExcLongMessage(buff));
374 }
375
376 if( theta<0. || theta>Pi) {
377 char buff[64];
378 sprintf(buff,"HEALPix::ang2pix_ring theta(=%g) out of range ", theta);
379 throw RangeCheckError(PExcLongMessage(buff));
380 }
381
382 z = cos(theta);
383 za = fabs(z);
384 if( phi >= twopi) phi = phi - twopi;
385 if (phi < 0.) phi = phi + twopi;
386 tt = phi / piover2;// ! in [0,4)
387
388 nl2 = 2*nside;
389 nl4 = 4*nside;
390 ncap = nl2*(nside-1);// ! number of pixels in the north polar cap
391 npix = 12*nside*nside;
392
393 if( za <= z0 ) {
394
395 jp = (int)floor(nside*(0.5 + tt - z*0.75));// ! index of ascending edge line
396 jm = (int)floor(nside*(0.5 + tt + z*0.75));// ! index of descending edge line
397
398 ir = nside + 1 + jp - jm;// ! in {1,2n+1} (ring number counted from z=2/3)
399 kshift = 0;
400 //RzModFloor if (fmod(ir,2)==0.) kshift = 1; ! kshift=1 if ir even, 0 otherwise
401 if ((ir%2)==0) kshift = 1;// ! kshift=1 if ir even, 0 otherwise
402
403 //RzModFloor ip = (int)floor( ( jp+jm - nside + kshift + 1 ) / 2 ) + 1; ! in {1,4n}
404 ip = ( jp+jm - nside + kshift + 1 )/2 + 1;
405 if( ip>nl4 ) ip = ip - nl4;
406
407 ipix1 = ncap + nl4*(ir-1) + ip ;
408 }
409 else {
410
411 tp = tt - floor(tt);// !MOD(tt,1.d0)
412 tmp = sqrt( 3.*(1. - za) );
413
414 jp = (int)floor( nside * tp * tmp );// ! increasing edge line index
415 jm = (int)floor( nside * (1. - tp) * tmp );// ! decreasing edge line index
416
417 ir = jp + jm + 1;// ! ring number counted from the closest pole
418 ip = (int)floor( tt * ir ) + 1;// ! in {1,4*ir}
419 if( ip>4*ir ) ip = ip - 4*ir;
420
421 ipix1 = 2*ir*(ir-1) + ip;
422 if( z<=0. ) {
423 ipix1 = npix - 2*ir*(ir+1) + ip;
424 }
425 }
426 return (ipix1 - 1);// ! in {0, npix-1}
427}
428
429int_4 HEALPix::ang2pix_nest(int_4 nside, double theta, double phi)
430{
431 /*
432 ==================================================
433 subroutine ang2pix_nest(nside, theta, phi, ipix)
434 ==================================================
435 c gives the pixel number ipix (NESTED)
436 c corresponding to angles theta and phi
437 c
438 c the computation is made to the highest resolution available (nside=8192)
439 c and then degraded to that required (by integer division)
440 c this doesn't cost more, and it makes sure
441 c that the treatement of round-off will be consistent
442 c for every resolution
443 ==================================================
444 */
445 // tranlated from FORTRAN (Gorski) to C, by B. Revenu, revised Guy Le Meur
446 // (16/12/98)
447
448 const PIXELS_XY& PXY= PIXELS_XY::instance();
449
450 double z, za, z0, tt, tp, tmp;
451 int face_num,jp,jm;
452 int ifp, ifm;
453 int ix, iy, ix_low, ix_hi, iy_low, iy_hi, ipf, ntt;
454 double piover2(Pi/2.), twopi(2.*Pi);
455 int ns_max(8192);
456
457 if( nside<1 || nside>ns_max ) {
458 char buff[64];
459 sprintf(buff,"HEALPix:::ang2pix_nest nside(=%d) out of range ", nside);
460 throw RangeCheckError(PExcLongMessage(buff));
461 }
462 if( theta<0 || theta>Pi ) {
463 char buff[64];
464 sprintf(buff,"HEALPix:::ang2pix_nest theta(=%g) out of range ", theta);
465 throw RangeCheckError(PExcLongMessage(buff));
466 }
467 z = cos(theta);
468 za = fabs(z);
469 z0 = 2./3.;
470 if( phi>=twopi ) phi = phi - twopi;
471 if( phi<0. ) phi = phi + twopi;
472 tt = phi / piover2;// ! in [0,4[
473 if( za<=z0 ) { // then ! equatorial region
474
475 //(the index of edge lines increase when the longitude=phi goes up)
476 jp = (int)floor(ns_max*(0.5 + tt - z*0.75));// ! ascending edge line index
477 jm = (int)floor(ns_max*(0.5 + tt + z*0.75));// ! descending edge line index
478
479 //c finds the face
480 ifp = jp / ns_max;// ! in {0,4}
481 ifm = jm / ns_max;
482 if( ifp==ifm ) face_num = (ifp%4)+4; //RzModFloor (int)fmod(ifp,4) + 4; then ! faces 4 to 7
483 else if( ifp<ifm ) face_num = ifp%4; //RzModFloor (int)fmod(ifp,4); (half-)faces 0 to 3
484 else face_num = (ifm%4) + 8; //RzModFloor (int)fmod(ifm,4) + 8;! (half-)faces 8 to 11
485
486 ix = jm%ns_max; //RzModFloor (int)fmod(jm, ns_max);
487 iy = ns_max - (jp%ns_max) - 1;//RzModFloor ns_max - (int)fmod(jp, ns_max) - 1;
488 }
489 else { //! polar region, za > 2/3
490
491 ntt = (int)floor(tt);
492 if( ntt>=4 ) ntt = 3;
493 tp = tt - ntt;
494 tmp = sqrt( 3.*(1. - za) );// ! in ]0,1]
495
496 //(the index of edge lines increase when distance from the closest pole goes up)
497 jp = (int)floor(ns_max*tp*tmp); // ! line going toward the pole as phi increases
498 jm = (int)floor(ns_max*(1.-tp)*tmp); // ! that one goes away of the closest pole
499 jp = (int)min(ns_max-1, jp);// ! for points too close to the boundary
500 jm = (int)min(ns_max-1, jm);
501
502 // finds the face and pixel's (x,y)
503 if( z>=0 ) {
504 face_num = ntt;// ! in {0,3}
505 ix = ns_max - jm - 1;
506 iy = ns_max - jp - 1;
507 }
508 else {
509 face_num = ntt + 8;// ! in {8,11}
510 ix = jp;
511 iy = jm;
512 }
513 }
514
515 ix_low = (ix%128); //RzModFloor (int)fmod(ix,128);
516 ix_hi = ix/128;
517 iy_low = (iy%128); //RzModFloor ((int)fmod(iy,128);
518 iy_hi = iy/128;
519 ipf= (PXY.x2pix_(ix_hi)+PXY.y2pix_(iy_hi))*(128*128)+(PXY.x2pix_(ix_low)+PXY.y2pix_(iy_low));
520 // ipf = ipf / pow(ns_max/nside,2.);// ! in {0, nside**2 - 1}
521 // return ( ipf + face_num*pow(nside,2));// ! in {0, 12*nside**2 - 1}
522 // $CHECK$ Reza 25/10/99 , pow remplace par *
523 ipf = ipf / ((ns_max/nside)*(ns_max/nside));
524 return (ipf + face_num*nside*nside);
525}
526
527void HEALPix::pix2ang_ring(int_4 nside,int_4 ipix,double& theta,double& phi)
528{
529 /*
530 ===================================================
531 c gives theta and phi corresponding to pixel ipix (RING)
532 c for a parameter nside
533 ===================================================
534 */
535 // tranlated from FORTRAN (Gorski) to C, by B. Revenu, revised Guy Le Meur
536 // (16/12/98)
537
538 int nl2, nl4, npix, ncap, iring, iphi, ip, ipix1;
539 double fact1, fact2, fodd, hip, fihip;
540
541 int ns_max(8192);
542
543 if( nside<1 || nside>ns_max ) {
544 char buff[64];
545 sprintf(buff,"HEALPix::pix2ang_ring nside(=%d) out of range ", nside);
546 throw RangeCheckError(PExcLongMessage(buff));
547 }
548 npix = 12*nside*nside; // ! total number of points
549 if( ipix<0 || ipix>npix-1 ) {
550 char buff[64];
551 sprintf(buff,"HEALPix::pix2ang_ring ipix(=%d) out of range ", ipix);
552 throw RangeCheckError(PExcLongMessage(buff));
553 }
554
555 ipix1 = ipix + 1; // in {1, npix}
556 nl2 = 2*nside;
557 nl4 = 4*nside;
558 ncap = 2*nside*(nside-1);// ! points in each polar cap, =0 for nside =1
559 fact1 = 1.5*nside;
560 fact2 = 3.0*nside*nside;
561
562 if( ipix1 <= ncap ) { //! North Polar cap -------------
563
564 hip = ipix1/2.;
565 fihip = floor(hip);
566 iring = (int)floor( sqrt( hip - sqrt(fihip) ) ) + 1;// ! counted from North pole
567 iphi = ipix1 - 2*iring*(iring - 1);
568
569 theta = acos( 1. - iring*iring / fact2 );
570 phi = ((double)iphi - 0.5) * Pi/(2.*iring);
571 // cout << theta << " " << phi << endl;
572 }
573 else if( ipix1 <= nl2*(5*nside+1) ) {//then ! Equatorial region ------
574
575 ip = ipix1 - ncap - 1;
576 iring = ip/nl4 + nside; //RzModFloor (int)floor( ip / nl4 ) + nside; ! counted from North pole
577 iphi = ip%nl4 + 1;
578
579 fodd = 0.5 * (1 + (iring+nside)%2 );// ! 1 if iring+nside is odd, 1/2 otherwise
580 theta = acos( (nl2 - iring) / fact1 );
581 phi = ((double)iphi - fodd) * Pi /(2.*nside);
582 }
583 else {//! South Polar cap -----------------------------------
584
585 ip = npix - ipix1 + 1;
586 hip = ip/2.;
587 fihip = floor(hip);
588 iring = (int)floor( sqrt( hip - sqrt(fihip) ) ) + 1;// ! counted from South pole
589 iphi = (int)(4.*iring + 1 - (ip - 2.*iring*(iring-1)));
590
591 theta = acos( -1. + iring*iring / fact2 );
592 phi = ((double)iphi - 0.5) * Pi/(2.*iring);
593 // cout << theta << " " << phi << endl;
594 }
595}
596
597void HEALPix::pix2ang_nest(int_4 nside,int_4 ipix,double& theta,double& phi)
598{
599 /*
600 ==================================================
601 subroutine pix2ang_nest(nside, ipix, theta, phi)
602 ==================================================
603 c gives theta and phi corresponding to pixel ipix (NESTED)
604 c for a parameter nside
605 ==================================================
606 */
607 // tranlated from FORTRAN (Gorski) to C, by B. Revenu, revised Guy Le Meur
608 // (16/12/98)
609
610 const PIXELS_XY& PXY= PIXELS_XY::instance();
611
612 int npix, npface, face_num;
613 int ipf, ip_low, ip_trunc, ip_med, ip_hi;
614 int ix, iy, jrt, jr, nr, jpt, jp, kshift, nl4;
615 double z, fn, fact1, fact2;
616 double piover2(Pi/2.);
617 int ns_max(8192);
618
619 // ! coordinate of the lowest corner of each face
620 int jrll[12]={2, 2, 2, 2, 3, 3, 3, 3, 4, 4, 4, 4};//! in unit of nside
621 int jpll[12]={1, 3, 5, 7, 0, 2, 4, 6, 1, 3, 5, 7};// ! in unit of nside/2
622
623 if( nside<1 || nside>ns_max ) {
624 char buff[64];
625 sprintf(buff,"HEALPix::pix2ang_nest nside(=%d) out of range ", nside);
626 throw RangeCheckError(PExcLongMessage(buff));
627 }
628 npix = 12 * nside*nside;
629 if( ipix<0 || ipix>npix-1 ) {
630 char buff[64];
631 sprintf(buff,"HEALPix::pix2ang_nest ipix(=%d) out of range ", ipix);
632 throw RangeCheckError(PExcLongMessage(buff));
633 }
634
635 fn = 1.*nside;
636 fact1 = 1./(3.*fn*fn);
637 fact2 = 2./(3.*fn);
638 nl4 = 4*nside;
639
640 //c finds the face, and the number in the face
641 npface = nside*nside;
642
643 face_num = ipix/npface;// ! face number in {0,11}
644 ipf = ipix%npface; //RzModFloor (int)fmod(ipix,npface); ! pixel number in the face {0,npface-1}
645
646 //c finds the x,y on the face (starting from the lowest corner)
647 //c from the pixel number
648 ip_low = ipf%1024; //RzModFloor (int)fmod(ipf,1024); ! content of the last 10 bits
649 ip_trunc = ipf/1024 ;// ! truncation of the last 10 bits
650 ip_med = ip_trunc%1024; //RzModFloor (int)fmod(ip_trunc,1024); ! content of the next 10 bits
651 ip_hi = ip_trunc/1024 ;//! content of the high weight 10 bits
652
653 ix = 1024*PXY.pix2x_(ip_hi)+32*PXY.pix2x_(ip_med)+PXY.pix2x_(ip_low);
654 iy = 1024*PXY.pix2y_(ip_hi)+32*PXY.pix2y_(ip_med)+PXY.pix2y_(ip_low);
655
656 //c transforms this in (horizontal, vertical) coordinates
657 jrt = ix + iy;// ! 'vertical' in {0,2*(nside-1)}
658 jpt = ix - iy;// ! 'horizontal' in {-nside+1,nside-1}
659
660 //c computes the z coordinate on the sphere
661 // jr = jrll[face_num+1]*nside - jrt - 1;// ! ring number in {1,4*nside-1}
662 jr = jrll[face_num]*nside - jrt - 1;
663 nr = nside;// ! equatorial region (the most frequent)
664 z = (2*nside-jr)*fact2;
665 kshift = (jr-nside) % 2; //RzModFloor (int)fmod(jr - nside, 2);
666 if( jr<nside ) { //then ! north pole region
667 nr = jr;
668 z = 1. - nr*nr*fact1;
669 kshift = 0;
670 }
671 else {
672 if( jr>3*nside ) {// then ! south pole region
673 nr = nl4 - jr;
674 z = - 1. + nr*nr*fact1;
675 kshift = 0;
676 }
677 }
678 theta = acos(z);
679
680 //c computes the phi coordinate on the sphere, in [0,2Pi]
681 // jp = (jpll[face_num+1]*nr + jpt + 1 + kshift)/2;// ! 'phi' number in the ring in {1,4*nr}
682 jp = (jpll[face_num]*nr + jpt + 1 + kshift)/2;
683 if( jp>nl4 ) jp = jp - nl4;
684 if( jp<1 ) jp = jp + nl4;
685 phi = (jp - (kshift+1)*0.5) * (piover2 / nr);
686}
687
688
689/*!
690 \brief return the size index value corresponding to resolution \b res
691 Computes the size index m providing a resolution better or equal the
692 specified value \b res (in radian)
693
694 - m : size-index
695 - number of pixels for 4Pi : 12*m*m
696
697 \sa SOPHYA::Angle
698*/
699int_4 HEALPix::ResolToSizeIndex(double res)
700{
701 int_4 m = 1;
702 double cres = sqrt(M_PI/3.); // 4PI/12
703 while (cres > 1.001*res) {
704 m *= 2;
705 cres = sqrt(4.*M_PI/(12*m*m));
706 }
707 return m;
708}
709
710//! return the pixel resolution for a given size index (nside)
711double HEALPix::SizeIndexToResol(int_4 m)
712{
713 double res = 4.*M_PI/(double)(12*m*m);
714 return sqrt(res);
715}
Note: See TracBrowser for help on using the repository browser.