healpix_base.cc

00001 /*
00002  *  This file is part of Healpix_cxx.
00003  *
00004  *  Healpix_cxx is free software; you can redistribute it and/or modify
00005  *  it under the terms of the GNU General Public License as published by
00006  *  the Free Software Foundation; either version 2 of the License, or
00007  *  (at your option) any later version.
00008  *
00009  *  Healpix_cxx is distributed in the hope that it will be useful,
00010  *  but WITHOUT ANY WARRANTY; without even the implied warranty of
00011  *  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00012  *  GNU General Public License for more details.
00013  *
00014  *  You should have received a copy of the GNU General Public License
00015  *  along with Healpix_cxx; if not, write to the Free Software
00016  *  Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA  02110-1301  USA
00017  *
00018  *  For more information about HEALPix, see http://healpix.sourceforge.net
00019  */
00020 
00021 /*
00022  *  Healpix_cxx is being developed at the Max-Planck-Institut fuer Astrophysik
00023  *  and financially supported by the Deutsches Zentrum fuer Luft- und Raumfahrt
00024  *  (DLR).
00025  */
00026 
00027 /*
00028  *  Copyright (C) 2003-2012 Max-Planck-Society
00029  *  Author: Martin Reinecke
00030  */
00031 
00032 #include "healpix_base.h"
00033 #include "geom_utils.h"
00034 #include "lsconstants.h"
00035 
00036 using namespace std;
00037 
00038 template<> const int T_Healpix_Base<int  >::order_max=13;
00039 template<> const int T_Healpix_Base<int64>::order_max=29;
00040 
00041 template<typename I> int T_Healpix_Base<I>::nside2order (I nside)
00042   {
00043   planck_assert (nside>I(0), "invalid value for Nside");
00044   return ((nside)&(nside-1)) ? -1 : ilog2(nside);
00045   }
00046 template<typename I> I T_Healpix_Base<I>::npix2nside (I npix)
00047   {
00048   I res=isqrt(npix/I(12));
00049   planck_assert (npix==res*res*I(12), "invalid value for npix");
00050   return res;
00051   }
00052 
00053 template<typename I> I T_Healpix_Base<I>::ring_above (double z) const
00054   {
00055   double az=abs(z);
00056   if (az<=twothird) // equatorial region
00057     return I(nside_*(2-1.5*z));
00058   I iring = I(nside_*sqrt(3*(1-az)));
00059   return (z>0) ? iring : 4*nside_-iring-1;
00060   }
00061 
00062 namespace {
00063 
00064 /* Short note on the "zone":
00065    zone = 0: pixel lies completely outside the queried shape
00066           1: pixel may overlap with the shape, pixel center is outside
00067           2: pixel center is inside the shape, but maybe not the complete pixel
00068           3: pixel lies completely inside the shape */
00069 
00070 template<typename I, typename I2> inline void check_pixel (int o, int order_,
00071   int omax, int zone, rangeset<I2> &pixset, I pix, vector<pair<I,int> > &stk,
00072   bool inclusive, int &stacktop)
00073   {
00074   if (zone==0) return;
00075 
00076   if (o<order_)
00077     {
00078     if (zone>=3)
00079       {
00080       int sdist=2*(order_-o); // the "bit-shift distance" between map orders
00081       pixset.append(pix<<sdist,(pix+1)<<sdist); // output all subpixels
00082       }
00083     else // (zone>=1)
00084       for (int i=0; i<4; ++i)
00085         stk.push_back(make_pair(4*pix+3-i,o+1)); // add children
00086     }
00087   else if (o>order_) // this implies that inclusive==true
00088     {
00089     if (zone>=2) // pixel center in shape
00090       {
00091       pixset.append(pix>>(2*(o-order_))); // output the parent pixel at order_
00092       stk.resize(stacktop); // unwind the stack
00093       }
00094     else // (zone>=1): pixel center in safety range
00095       {
00096       if (o<omax) // check sublevels
00097         for (int i=0; i<4; ++i) // add children in reverse order
00098           stk.push_back(make_pair(4*pix+3-i,o+1));
00099       else // at resolution limit
00100         {
00101         pixset.append(pix>>(2*(o-order_))); // output the parent pixel at order_
00102         stk.resize(stacktop); // unwind the stack
00103         }
00104       }
00105     }
00106   else // o==order_
00107     {
00108     if (zone>=2)
00109       pixset.append(pix);
00110     else if (inclusive) // and (zone>=1)
00111       {
00112       if (order_<omax) // check sublevels
00113         {
00114         stacktop=stk.size(); // remember current stack position
00115         for (int i=0; i<4; ++i) // add children in reverse order
00116           stk.push_back(make_pair(4*pix+3-i,o+1));
00117         }
00118       else // at resolution limit
00119         pixset.append(pix); // output the pixel
00120       }
00121     }
00122   }
00123 
00124 template<typename I> bool check_pixel_ring (const T_Healpix_Base<I> &b1,
00125   const T_Healpix_Base<I> &b2, I pix, I nr, I ipix1, int fct,
00126   double cz, double cphi, double cosrp2, I cpix)
00127   {
00128   if (pix>=nr) pix-=nr;
00129   if (pix<0) pix+=nr;
00130   pix+=ipix1;
00131   if (pix==cpix) return false; // disk center in pixel => overlap
00132   int px,py,pf;
00133   b1.pix2xyf(pix,px,py,pf);
00134   for (int i=0; i<fct-1; ++i) // go along the 4 edges
00135     {
00136     I ox=fct*px, oy=fct*py;
00137     double pz,pphi;
00138     b2.pix2zphi(b2.xyf2pix(ox+i,oy,pf),pz,pphi);
00139     if (cosdist_zphi(pz,pphi,cz,cphi)>cosrp2) // overlap
00140       return false;
00141     b2.pix2zphi(b2.xyf2pix(ox+fct-1,oy+i,pf),pz,pphi);
00142     if (cosdist_zphi(pz,pphi,cz,cphi)>cosrp2) // overlap
00143       return false;
00144     b2.pix2zphi(b2.xyf2pix(ox+fct-1-i,oy+fct-1,pf),pz,pphi);
00145     if (cosdist_zphi(pz,pphi,cz,cphi)>cosrp2) // overlap
00146       return false;
00147     b2.pix2zphi(b2.xyf2pix(ox,oy+fct-1-i,pf),pz,pphi);
00148     if (cosdist_zphi(pz,pphi,cz,cphi)>cosrp2) // overlap
00149       return false;
00150     }
00151   return true;
00152   }
00153 
00154 } // unnamed namespace
00155 
00156 template<typename I> template<typename I2>
00157   void T_Healpix_Base<I>::query_disc_internal
00158   (pointing ptg, double radius, int fact, rangeset<I2> &pixset) const
00159   {
00160   bool inclusive = (fact!=0);
00161   pixset.clear();
00162   ptg.normalize();
00163 
00164   if (scheme_==RING)
00165     {
00166     I fct=1;
00167     if (inclusive)
00168       {
00169       planck_assert (((I(1)<<order_max)/nside_)>=fact,
00170         "invalid oversampling factor");
00171       fct = fact;
00172       }
00173     T_Healpix_Base b2;
00174     double rsmall, rbig;
00175     if (fct>1)
00176       {
00177       b2.SetNside(fct*nside_,RING);
00178       rsmall = radius+b2.max_pixrad();
00179       rbig = radius+max_pixrad();
00180       }
00181     else
00182       rsmall = rbig = inclusive ? radius+max_pixrad() : radius;
00183 
00184     if (rsmall>=pi)
00185       { pixset.append(0,npix_); return; }
00186 
00187     rbig = min(pi,rbig);
00188 
00189     double cosrsmall = cos(rsmall);
00190     double cosrbig = cos(rbig);
00191 
00192     double z0 = cos(ptg.theta);
00193     double xa = 1./sqrt((1-z0)*(1+z0));
00194 
00195     I cpix=zphi2pix(z0,ptg.phi);
00196 
00197     double rlat1 = ptg.theta - rsmall;
00198     double zmax = cos(rlat1);
00199     I irmin = ring_above (zmax)+1;
00200 
00201     if ((rlat1<=0) && (irmin>1)) // north pole in the disk
00202       {
00203       I sp,rp; bool dummy;
00204       get_ring_info_small(irmin-1,sp,rp,dummy);
00205       pixset.append(0,sp+rp);
00206       }
00207 
00208     if ((fct>1) && (rlat1>0)) irmin=max(I(1),irmin-1);
00209 
00210     double rlat2 = ptg.theta + rsmall;
00211     double zmin = cos(rlat2);
00212     I irmax = ring_above (zmin);
00213 
00214     if ((fct>1) && (rlat2<pi)) irmax=min(4*nside_-1,irmax+1);
00215 
00216     for (I iz=irmin; iz<=irmax; ++iz)
00217       {
00218       double z=ring2z(iz);
00219       double x = (cosrbig-z*z0)*xa;
00220       double ysq = 1-z*z-x*x;
00221       double dphi = (ysq<=0) ? pi-1e-15 : atan2(sqrt(ysq),x);
00222       I nr, ipix1;
00223       bool shifted;
00224       get_ring_info_small(iz,ipix1,nr,shifted);
00225       double shift = shifted ? 0.5 : 0.;
00226 
00227       I ipix2 = ipix1 + nr - 1; // highest pixel number in the ring
00228 
00229       I ip_lo = ifloor<I>(nr*inv_twopi*(ptg.phi-dphi) - shift)+1;
00230       I ip_hi = ifloor<I>(nr*inv_twopi*(ptg.phi+dphi) - shift);
00231 
00232       if (fct>1)
00233         {
00234         while ((ip_lo<=ip_hi) && check_pixel_ring
00235                (*this,b2,ip_lo,nr,ipix1,fct,z0,ptg.phi,cosrsmall,cpix))
00236           ++ip_lo;
00237         while ((ip_hi>ip_lo) && check_pixel_ring
00238                (*this,b2,ip_hi,nr,ipix1,fct,z0,ptg.phi,cosrsmall,cpix))
00239           --ip_hi;
00240         }
00241 
00242       if (ip_lo<=ip_hi)
00243         {
00244         if (ip_hi>=nr)
00245           { ip_lo-=nr; ip_hi-=nr; }
00246         if (ip_lo<0)
00247           {
00248           pixset.append(ipix1,ipix1+ip_hi+1);
00249           pixset.append(ipix1+ip_lo+nr,ipix2+1);
00250           }
00251         else
00252           pixset.append(ipix1+ip_lo,ipix1+ip_hi+1);
00253         }
00254       }
00255     if ((rlat2>=pi) && (irmax+1<4*nside_)) // south pole in the disk
00256       {
00257       I sp,rp; bool dummy;
00258       get_ring_info_small(irmax+1,sp,rp,dummy);
00259       pixset.append(sp,npix_);
00260       }
00261     }
00262   else // scheme_==NEST
00263     {
00264     if (radius>=pi) // disk covers the whole sphere
00265       { pixset.append(0,npix_); return; }
00266 
00267     int oplus = 0;
00268     if (inclusive)
00269       {
00270       planck_assert ((I(1)<<(order_max-order_))>=fact,
00271         "invalid oversampling factor");
00272       planck_assert ((fact&(fact-1))==0,
00273         "oversampling factor must be a power of 2");
00274       oplus=ilog2(fact);
00275       }
00276     int omax=order_+oplus; // the order up to which we test
00277 
00278     vec3 vptg(ptg);
00279     arr<T_Healpix_Base<I> > base(omax+1);
00280     arr<double> crpdr(omax+1), crmdr(omax+1);
00281     double cosrad=cos(radius);
00282     for (int o=0; o<=omax; ++o) // prepare data at the required orders
00283       {
00284       base[o].Set(o,NEST);
00285       double dr=base[o].max_pixrad(); // safety distance
00286       crpdr[o] = (radius+dr>pi) ? -1. : cos(radius+dr);
00287       crmdr[o] = (radius-dr<0.) ?  1. : cos(radius-dr);
00288       }
00289     vector<pair<I,int> > stk; // stack for pixel numbers and their orders
00290     stk.reserve(12+3*omax); // reserve maximum size to avoid reallocation
00291     for (int i=0; i<12; ++i) // insert the 12 base pixels in reverse order
00292       stk.push_back(make_pair(I(11-i),0));
00293 
00294     int stacktop=0; // a place to save a stack position
00295 
00296     while (!stk.empty()) // as long as there are pixels on the stack
00297       {
00298       // pop current pixel number and order from the stack
00299       I pix=stk.back().first;
00300       int o=stk.back().second;
00301       stk.pop_back();
00302 
00303       double z,phi;
00304       base[o].pix2zphi(pix,z,phi);
00305       // cosine of angular distance between pixel center and disk center
00306       double cangdist=cosdist_zphi(vptg.z,ptg.phi,z,phi);
00307 
00308       if (cangdist>crpdr[o])
00309         {
00310         int zone = (cangdist<cosrad) ? 1 : ((cangdist<=crmdr[o]) ? 2 : 3);
00311 
00312         check_pixel (o, order_, omax, zone, pixset, pix, stk, inclusive,
00313           stacktop);
00314         }
00315       }
00316     }
00317   }
00318 
00319 template<typename I> void T_Healpix_Base<I>::query_disc
00320   (pointing ptg, double radius, rangeset<I> &pixset) const
00321   {
00322   query_disc_internal (ptg, radius, 0, pixset);
00323   }
00324 
00325 template<typename I> void T_Healpix_Base<I>::query_disc_inclusive
00326   (pointing ptg, double radius, rangeset<I> &pixset, int fact) const
00327   {
00328   planck_assert(fact>0,"fact must be a positive integer");
00329   if ((sizeof(I)<8) && (((I(1)<<order_max)/nside_)<fact))
00330     {
00331     T_Healpix_Base<int64> base2(nside_,scheme_,SET_NSIDE);
00332     base2.query_disc_internal(ptg,radius,fact,pixset);
00333     return;
00334     }
00335   query_disc_internal (ptg, radius, fact, pixset);
00336   }
00337 
00338 template<typename I> template<typename I2>
00339   void T_Healpix_Base<I>::query_multidisc (const arr<vec3> &norm,
00340   const arr<double> &rad, int fact, rangeset<I2> &pixset) const
00341   {
00342   bool inclusive = (fact!=0);
00343   tsize nv=norm.size();
00344   planck_assert(nv==rad.size(),"inconsistent input arrays");
00345   pixset.clear();
00346 
00347   if (scheme_==RING)
00348     {
00349     I fct=1;
00350     if (inclusive)
00351       {
00352       planck_assert (((I(1)<<order_max)/nside_)>=fact,
00353         "invalid oversampling factor");
00354       fct = fact;
00355       }
00356     T_Healpix_Base b2;
00357     double rpsmall, rpbig;
00358     if (fct>1)
00359       {
00360       b2.SetNside(fct*nside_,RING);
00361       rpsmall = b2.max_pixrad();
00362       rpbig = max_pixrad();
00363       }
00364     else
00365       rpsmall = rpbig = inclusive ? max_pixrad() : 0;
00366 
00367     I irmin=1, irmax=4*nside_-1;
00368     vector<double> z0,xa,cosrsmall,cosrbig;
00369     vector<pointing> ptg;
00370     vector<I> cpix;
00371     for (tsize i=0; i<nv; ++i)
00372       {
00373       double rsmall=rad[i]+rpsmall;
00374       if (rsmall<pi)
00375         {
00376         double rbig=min(pi,rad[i]+rpbig);
00377         pointing pnt=pointing(norm[i]);
00378         cosrsmall.push_back(cos(rsmall));
00379         cosrbig.push_back(cos(rbig));
00380         double cth=cos(pnt.theta);
00381         z0.push_back(cth);
00382         if (fct>1) cpix.push_back(zphi2pix(cth,pnt.phi));
00383         xa.push_back(1./sqrt((1-cth)*(1+cth)));
00384         ptg.push_back(pnt);
00385 
00386         double rlat1 = pnt.theta - rsmall;
00387         double zmax = cos(rlat1);
00388         I irmin_t = (rlat1<=0) ? 1 : ring_above (zmax)+1;
00389 
00390         if ((fct>1) && (rlat1>0)) irmin_t=max(I(1),irmin_t-1);
00391 
00392         double rlat2 = pnt.theta + rsmall;
00393         double zmin = cos(rlat2);
00394         I irmax_t = (rlat2>=pi) ? 4*nside_-1 : ring_above (zmin);
00395 
00396         if ((fct>1) && (rlat2<pi)) irmax_t=min(4*nside_-1,irmax_t+1);
00397 
00398         if (irmax_t < irmax) irmax=irmax_t;
00399         if (irmin_t > irmin) irmin=irmin_t;
00400         }
00401       }
00402 
00403     for (I iz=irmin; iz<=irmax; ++iz)
00404       {
00405       double z=ring2z(iz);
00406       I ipix1,nr;
00407       bool shifted;
00408       get_ring_info_small(iz,ipix1,nr,shifted);
00409       double shift = shifted ? 0.5 : 0.;
00410       rangeset<I2> tr;
00411       tr.append(ipix1,ipix1+nr);
00412       for (tsize j=0; j<z0.size(); ++j)
00413         {
00414         double x = (cosrbig[j]-z*z0[j])*xa[j];
00415         double ysq = 1.-z*z-x*x;
00416         double dphi = (ysq<=0) ? pi-1e-15 : atan2(sqrt(ysq),x);
00417         I ip_lo = ifloor<I>(nr*inv_twopi*(ptg[j].phi-dphi) - shift)+1;
00418         I ip_hi = ifloor<I>(nr*inv_twopi*(ptg[j].phi+dphi) - shift);
00419         if (fct>1)
00420           {
00421           while ((ip_lo<=ip_hi) && check_pixel_ring
00422             (*this,b2,ip_lo,nr,ipix1,fct,z0[j],ptg[j].phi,cosrsmall[j],cpix[j]))
00423             ++ip_lo;
00424           while ((ip_hi>ip_lo) && check_pixel_ring
00425             (*this,b2,ip_hi,nr,ipix1,fct,z0[j],ptg[j].phi,cosrsmall[j],cpix[j]))
00426             --ip_hi;
00427           }
00428         if (ip_hi>=nr)
00429           { ip_lo-=nr; ip_hi-=nr;}
00430         if (ip_lo<0)
00431           tr.remove(ipix1+ip_hi+1,ipix1+ip_lo+nr);
00432         else
00433           tr.intersect(ipix1+ip_lo,ipix1+ip_hi+1);
00434         }
00435       pixset.append(tr);
00436       }
00437     }
00438   else // scheme_ == NEST
00439     {
00440     int oplus = 0;
00441     if (inclusive)
00442       {
00443       planck_assert ((I(1)<<(order_max-order_))>=fact,
00444         "invalid oversampling factor");
00445       planck_assert ((fact&(fact-1))==0,
00446         "oversampling factor must be a power of 2");
00447       oplus=ilog2(fact);
00448       }
00449     int omax=order_+oplus; // the order up to which we test
00450 
00451     // TODO: ignore all disks with radius>=pi
00452 
00453     arr<T_Healpix_Base<I> > base(omax+1);
00454     arr3<double> crlimit(omax+1,nv,3);
00455     for (int o=0; o<=omax; ++o) // prepare data at the required orders
00456       {
00457       base[o].Set(o,NEST);
00458       double dr=base[o].max_pixrad(); // safety distance
00459       for (tsize i=0; i<nv; ++i)
00460         {
00461         crlimit(o,i,0) = (rad[i]+dr>pi) ? -1. : cos(rad[i]+dr);
00462         crlimit(o,i,1) = (o==0) ? cos(rad[i]) : crlimit(0,i,1);
00463         crlimit(o,i,2) = (rad[i]-dr<0.) ?  1. : cos(rad[i]-dr);
00464         }
00465       }
00466 
00467     vector<pair<I,int> > stk; // stack for pixel numbers and their orders
00468     stk.reserve(12+3*omax); // reserve maximum size to avoid reallocation
00469     for (int i=0; i<12; ++i) // insert the 12 base pixels in reverse order
00470       stk.push_back(make_pair(I(11-i),0));
00471 
00472     int stacktop=0; // a place to save a stack position
00473 
00474     while (!stk.empty()) // as long as there are pixels on the stack
00475       {
00476       // pop current pixel number and order from the stack
00477       I pix=stk.back().first;
00478       int o=stk.back().second;
00479       stk.pop_back();
00480 
00481       vec3 pv(base[o].pix2vec(pix));
00482 
00483       tsize zone=3;
00484       for (tsize i=0; i<nv; ++i)
00485         {
00486         double crad=dotprod(pv,norm[i]);
00487         for (tsize iz=0; iz<zone; ++iz)
00488           if (crad<crlimit(o,i,iz))
00489             if ((zone=iz)==0) goto bailout;
00490         }
00491 
00492       check_pixel (o, order_, omax, zone, pixset, pix, stk, inclusive,
00493         stacktop);
00494       bailout:;
00495       }
00496     }
00497   }
00498 
00499 template<typename I> void T_Healpix_Base<I>::query_multidisc_general
00500   (const arr<vec3> &norm, const arr<double> &rad, bool inclusive,
00501   const vector<int> &cmds, rangeset<I> &pixset) const
00502   {
00503   tsize nv=norm.size();
00504   planck_assert(nv==rad.size(),"inconsistent input arrays");
00505   pixset.clear();
00506 
00507   if (scheme_==RING)
00508     {
00509     planck_fail ("not yet implemented");
00510     }
00511   else // scheme_ == NEST
00512     {
00513     int oplus=inclusive ? 2 : 0;
00514     int omax=min(order_max,order_+oplus); // the order up to which we test
00515 
00516     // TODO: ignore all disks with radius>=pi
00517 
00518     arr<T_Healpix_Base<I> > base(omax+1);
00519     arr3<double> crlimit(omax+1,nv,3);
00520     for (int o=0; o<=omax; ++o) // prepare data at the required orders
00521       {
00522       base[o].Set(o,NEST);
00523       double dr=base[o].max_pixrad(); // safety distance
00524       for (tsize i=0; i<nv; ++i)
00525         {
00526         crlimit(o,i,0) = (rad[i]+dr>pi) ? -1. : cos(rad[i]+dr);
00527         crlimit(o,i,1) = (o==0) ? cos(rad[i]) : crlimit(0,i,1);
00528         crlimit(o,i,2) = (rad[i]-dr<0.) ?  1. : cos(rad[i]-dr);
00529         }
00530       }
00531 
00532     vector<pair<I,int> > stk; // stack for pixel numbers and their orders
00533     stk.reserve(12+3*omax); // reserve maximum size to avoid reallocation
00534     for (int i=0; i<12; ++i) // insert the 12 base pixels in reverse order
00535       stk.push_back(make_pair(I(11-i),0));
00536 
00537     int stacktop=0; // a place to save a stack position
00538     arr<tsize> zone(nv);
00539 
00540     vector<tsize> zstk; zstk.reserve(cmds.size());
00541 
00542     while (!stk.empty()) // as long as there are pixels on the stack
00543       {
00544       // pop current pixel number and order from the stack
00545       I pix=stk.back().first;
00546       int o=stk.back().second;
00547       stk.pop_back();
00548 
00549       vec3 pv(base[o].pix2vec(pix));
00550 
00551       for (tsize i=0; i<nv; ++i)
00552         {
00553         zone[i]=3;
00554         double crad=dotprod(pv,norm[i]);
00555         for (tsize iz=0; iz<zone[i]; ++iz)
00556           if (crad<crlimit(o,i,iz))
00557             zone[i]=iz;
00558         }
00559 
00560       for (tsize i=0; i<cmds.size(); ++i)
00561         {
00562         tsize tmp;
00563         switch (cmds[i])
00564           {
00565           case -1: // union
00566             tmp=zstk.back(); zstk.pop_back();
00567             zstk.back() = max(zstk.back(),tmp);
00568             break;
00569           case -2: // intersection
00570             tmp=zstk.back(); zstk.pop_back();
00571             zstk.back() = min(zstk.back(),tmp);
00572             break;
00573           default: // add value
00574             zstk.push_back(zone[cmds[i]]);
00575           }
00576         }
00577       planck_assert(zstk.size()==1,"inconsistent commands");
00578       tsize zn=zstk[0]; zstk.pop_back();
00579 
00580       check_pixel (o, order_, omax, zn, pixset, pix, stk, inclusive,
00581         stacktop);
00582       }
00583     }
00584   }
00585 
00586 template<> inline int T_Healpix_Base<int>::spread_bits (int v) const
00587   { return utab[v&0xff] | (utab[(v>>8)&0xff]<<16); }
00588 template<> inline int64 T_Healpix_Base<int64>::spread_bits (int v) const
00589   {
00590   return  int64(utab[ v     &0xff])      | (int64(utab[(v>> 8)&0xff])<<16)
00591        | (int64(utab[(v>>16)&0xff])<<32) | (int64(utab[(v>>24)&0xff])<<48);
00592   }
00593 
00594 template<> inline int T_Healpix_Base<int>::compress_bits (int v) const
00595   {
00596   int raw = (v&0x5555) | ((v&0x55550000)>>15);
00597   return ctab[raw&0xff] | (ctab[raw>>8]<<4);
00598   }
00599 template<> inline int T_Healpix_Base<int64>::compress_bits (int64 v) const
00600   {
00601   int64 raw = v&0x5555555555555555ull;
00602   raw|=raw>>15;
00603   return ctab[ raw     &0xff]      | (ctab[(raw>> 8)&0xff]<< 4)
00604       | (ctab[(raw>>32)&0xff]<<16) | (ctab[(raw>>40)&0xff]<<20);
00605   }
00606 
00607 template<typename I> inline void T_Healpix_Base<I>::nest2xyf (I pix, int &ix,
00608   int &iy, int &face_num) const
00609   {
00610   face_num = pix>>(2*order_);
00611   pix &= (npface_-1);
00612   ix = compress_bits(pix);
00613   iy = compress_bits(pix>>1);
00614   }
00615 
00616 template<typename I> inline I T_Healpix_Base<I>::xyf2nest (int ix, int iy,
00617   int face_num) const
00618   { return (I(face_num)<<(2*order_)) + spread_bits(ix) + (spread_bits(iy)<<1); }
00619 
00620 template<typename I> void T_Healpix_Base<I>::ring2xyf (I pix, int &ix, int &iy,
00621   int &face_num) const
00622   {
00623   I iring, iphi, kshift, nr;
00624   I nl2 = 2*nside_;
00625 
00626   if (pix<ncap_) // North Polar cap
00627     {
00628     iring = (1+isqrt(1+2*pix))>>1; //counted from North pole
00629     iphi  = (pix+1) - 2*iring*(iring-1);
00630     kshift = 0;
00631     nr = iring;
00632     face_num=(iphi-1)/nr;
00633     }
00634   else if (pix<(npix_-ncap_)) // Equatorial region
00635     {
00636     I ip = pix - ncap_;
00637     I tmp = (order_>=0) ? ip>>(order_+2) : ip/(4*nside_);
00638     iring = tmp+nside_;
00639     iphi = ip-tmp*4*nside_ + 1;
00640     kshift = (iring+nside_)&1;
00641     nr = nside_;
00642     I ire = iring-nside_+1,
00643       irm = nl2+2-ire;
00644     I ifm = iphi - ire/2 + nside_ -1,
00645       ifp = iphi - irm/2 + nside_ -1;
00646     if (order_>=0)
00647       { ifm >>= order_; ifp >>= order_; }
00648     else
00649       { ifm /= nside_; ifp /= nside_; }
00650     face_num = (ifp==ifm) ? (ifp|4) : ((ifp<ifm) ? ifp : (ifm+8));
00651     }
00652   else // South Polar cap
00653     {
00654     I ip = npix_ - pix;
00655     iring = (1+isqrt(2*ip-1))>>1; //counted from South pole
00656     iphi  = 4*iring + 1 - (ip - 2*iring*(iring-1));
00657     kshift = 0;
00658     nr = iring;
00659     iring = 2*nl2-iring;
00660     face_num = 8 + (iphi-1)/nr;
00661     }
00662 
00663   I irt = iring - (jrll[face_num]*nside_) + 1;
00664   I ipt = 2*iphi- jpll[face_num]*nr - kshift -1;
00665   if (ipt>=nl2) ipt-=8*nside_;
00666 
00667   ix =  (ipt-irt) >>1;
00668   iy = (-ipt-irt) >>1;
00669   }
00670 
00671 template<typename I> I T_Healpix_Base<I>::xyf2ring (int ix, int iy,
00672   int face_num) const
00673   {
00674   I nl4 = 4*nside_;
00675   I jr = (jrll[face_num]*nside_) - ix - iy  - 1;
00676 
00677   I nr, kshift, n_before;
00678 
00679   bool shifted;
00680   get_ring_info_small(jr,n_before,nr,shifted);
00681   nr>>=2;
00682   kshift=1-shifted;
00683   I jp = (jpll[face_num]*nr + ix - iy + 1 + kshift) / 2;
00684   planck_assert(jp<=4*nr,"must not happen");
00685   if (jp<1) jp+=nl4; // assumption: if this triggers, then nl4==4*nr
00686 
00687   return n_before + jp - 1;
00688   }
00689 
00690 template<typename I> T_Healpix_Base<I>::T_Healpix_Base ()
00691   : order_(-1), nside_(0), npface_(0), ncap_(0), npix_(0),
00692     fact1_(0), fact2_(0), scheme_(RING) {}
00693 
00694 template<typename I> void T_Healpix_Base<I>::Set (int order,
00695   Healpix_Ordering_Scheme scheme)
00696   {
00697   planck_assert ((order>=0)&&(order<=order_max), "bad order");
00698   order_  = order;
00699   nside_  = I(1)<<order;
00700   npface_ = nside_<<order_;
00701   ncap_   = (npface_-nside_)<<1;
00702   npix_   = 12*npface_;
00703   fact2_  = 4./npix_;
00704   fact1_  = (nside_<<1)*fact2_;
00705   scheme_ = scheme;
00706   }
00707 template<typename I> void T_Healpix_Base<I>::SetNside (I nside,
00708   Healpix_Ordering_Scheme scheme)
00709   {
00710   order_  = nside2order(nside);
00711   planck_assert ((scheme!=NEST) || (order_>=0),
00712     "SetNside: nside must be power of 2 for nested maps");
00713   nside_  = nside;
00714   npface_ = nside_*nside_;
00715   ncap_   = (npface_-nside_)<<1;
00716   npix_   = 12*npface_;
00717   fact2_  = 4./npix_;
00718   fact1_  = (nside_<<1)*fact2_;
00719   scheme_ = scheme;
00720   }
00721 
00722 template<typename I> double T_Healpix_Base<I>::ring2z (I ring) const
00723   {
00724   if (ring<nside_)
00725     return 1 - ring*ring*fact2_;
00726   if (ring <=3*nside_)
00727     return (2*nside_-ring)*fact1_;
00728   ring=4*nside_ - ring;
00729   return ring*ring*fact2_ - 1;
00730   }
00731 
00732 template<typename I> I T_Healpix_Base<I>::pix2ring (I pix) const
00733   {
00734   if (scheme_==RING)
00735     {
00736     if (pix<ncap_) // North Polar cap
00737       return (1+I(isqrt(1+2*pix)))>>1; // counted from North pole
00738     else if (pix<(npix_-ncap_)) // Equatorial region
00739       return (pix-ncap_)/(4*nside_) + nside_; // counted from North pole
00740     else // South Polar cap
00741       return 4*nside_-((1+I(isqrt(2*(npix_-pix)-1)))>>1);
00742     }
00743   else
00744     {
00745     int face_num, ix, iy;
00746     nest2xyf(pix,ix,iy,face_num);
00747     return (I(jrll[face_num])<<order_) - ix - iy - 1;
00748     }
00749   }
00750 
00751 template<typename I> I T_Healpix_Base<I>::nest2ring (I pix) const
00752   {
00753   planck_assert(order_>=0, "hierarchical map required");
00754   int ix, iy, face_num;
00755   nest2xyf (pix, ix, iy, face_num);
00756   return xyf2ring (ix, iy, face_num);
00757   }
00758 
00759 template<typename I> I T_Healpix_Base<I>::ring2nest (I pix) const
00760   {
00761   planck_assert(order_>=0, "hierarchical map required");
00762   int ix, iy, face_num;
00763   ring2xyf (pix, ix, iy, face_num);
00764   return xyf2nest (ix, iy, face_num);
00765   }
00766 
00767 template<typename I> inline I T_Healpix_Base<I>::nest_peano_helper
00768   (I pix, int dir) const
00769   {
00770   int face = pix>>(2*order_);
00771   I result = 0;
00772   uint8 path = peano_face2path[dir][face];
00773 
00774   for (int shift=2*order_-2; shift>=0; shift-=2)
00775     {
00776     uint8 spix = uint8((pix>>shift) & 0x3);
00777     result <<= 2;
00778     result |= peano_subpix[dir][path][spix];
00779     path=peano_subpath[dir][path][spix];
00780     }
00781 
00782   return result + (I(peano_face2face[dir][face])<<(2*order_));
00783   }
00784 
00785 template<typename I> I T_Healpix_Base<I>::nest2peano (I pix) const
00786   { return nest_peano_helper(pix,0); }
00787 
00788 template<typename I> I T_Healpix_Base<I>::peano2nest (I pix) const
00789   { return nest_peano_helper(pix,1); }
00790 
00791 template<typename I> I T_Healpix_Base<I>::loc2pix (double z, double phi,
00792   double sth, bool have_sth) const
00793   {
00794   double za = abs(z);
00795   double tt = fmodulo(phi*inv_halfpi,4.0); // in [0,4)
00796 
00797   if (scheme_==RING)
00798     {
00799     if (za<=twothird) // Equatorial region
00800       {
00801       I nl4 = 4*nside_;
00802       double temp1 = nside_*(0.5+tt);
00803       double temp2 = nside_*z*0.75;
00804       I jp = I(temp1-temp2); // index of  ascending edge line
00805       I jm = I(temp1+temp2); // index of descending edge line
00806 
00807       // ring number counted from z=2/3
00808       I ir = nside_ + 1 + jp - jm; // in {1,2n+1}
00809       I kshift = 1-(ir&1); // kshift=1 if ir even, 0 otherwise
00810 
00811       I t1 = jp+jm-nside_+kshift+1+nl4+nl4;
00812       I ip = (order_>0) ?
00813         (t1>>1)&(nl4-1) : ((t1>>1)%nl4); // in {0,4n-1}
00814 
00815       return ncap_ + (ir-1)*nl4 + ip;
00816       }
00817     else  // North & South polar caps
00818       {
00819       double tp = tt-I(tt);
00820       double tmp = ((za<0.99)||(!have_sth)) ?
00821                    nside_*sqrt(3*(1-za)) :
00822                    nside_*sth/sqrt((1.+za)/3.);
00823 
00824       I jp = I(tp*tmp); // increasing edge line index
00825       I jm = I((1.0-tp)*tmp); // decreasing edge line index
00826 
00827       I ir = jp+jm+1; // ring number counted from the closest pole
00828       I ip = I(tt*ir); // in {0,4*ir-1}
00829       planck_assert((ip>=0)&&(ip<4*ir),"must not happen");
00830       //ip %= 4*ir;
00831 
00832       return (z>0)  ?  2*ir*(ir-1) + ip  :  npix_ - 2*ir*(ir+1) + ip;
00833       }
00834     }
00835   else // scheme_ == NEST
00836     {
00837     if (za<=twothird) // Equatorial region
00838       {
00839       double temp1 = nside_*(0.5+tt);
00840       double temp2 = nside_*(z*0.75);
00841       I jp = I(temp1-temp2); // index of  ascending edge line
00842       I jm = I(temp1+temp2); // index of descending edge line
00843       I ifp = jp >> order_;  // in {0,4}
00844       I ifm = jm >> order_;
00845       int face_num = (ifp==ifm) ? (ifp|4) : ((ifp<ifm) ? ifp : (ifm+8));
00846 
00847       int ix = jm & (nside_-1),
00848           iy = nside_ - (jp & (nside_-1)) - 1;
00849       return xyf2nest(ix,iy,face_num);
00850       }
00851     else // polar region, za > 2/3
00852       {
00853       int ntt = min(3,int(tt));
00854       double tp = tt-ntt;
00855       double tmp = ((za<0.99)||(!have_sth)) ?
00856                    nside_*sqrt(3*(1-za)) :
00857                    nside_*sth/sqrt((1.+za)/3.);
00858 
00859       I jp = I(tp*tmp); // increasing edge line index
00860       I jm = I((1.0-tp)*tmp); // decreasing edge line index
00861       jp=min(jp,nside_-1); // for points too close to the boundary
00862       jm=min(jm,nside_-1);
00863       return (z>=0) ?
00864         xyf2nest(nside_-jm -1,nside_-jp-1,ntt) : xyf2nest(jp,jm,ntt+8);
00865       }
00866     }
00867   }
00868 
00869 template<typename I> void T_Healpix_Base<I>::pix2loc (I pix, double &z,
00870   double &phi, double &sth, bool &have_sth) const
00871   {
00872   have_sth=false;
00873   if (scheme_==RING)
00874     {
00875     if (pix<ncap_) // North Polar cap
00876       {
00877       I iring = (1+I(isqrt(1+2*pix)))>>1; // counted from North pole
00878       I iphi  = (pix+1) - 2*iring*(iring-1);
00879 
00880       double tmp=(iring*iring)*fact2_;
00881       z = 1.0 - tmp;
00882       if (z>0.99) { sth=sqrt(tmp*(2.-tmp)); have_sth=true; }
00883       phi = (iphi-0.5) * halfpi/iring;
00884       }
00885     else if (pix<(npix_-ncap_)) // Equatorial region
00886       {
00887       I nl4 = 4*nside_;
00888       I ip  = pix - ncap_;
00889       I tmp = (order_>=0) ? ip>>(order_+2) : ip/nl4;
00890       I iring = tmp + nside_,
00891         iphi = ip-nl4*tmp+1;
00892       // 1 if iring+nside is odd, 1/2 otherwise
00893       double fodd = ((iring+nside_)&1) ? 1 : 0.5;
00894 
00895       z = (2*nside_-iring)*fact1_;
00896       phi = (iphi-fodd) * pi*0.75*fact1_;
00897       }
00898     else // South Polar cap
00899       {
00900       I ip = npix_ - pix;
00901       I iring = (1+I(isqrt(2*ip-1)))>>1; // counted from South pole
00902       I iphi  = 4*iring + 1 - (ip - 2*iring*(iring-1));
00903 
00904       double tmp=(iring*iring)*fact2_;
00905       z = tmp - 1.0;
00906       if (z<-0.99) { sth=sqrt(tmp*(2.-tmp)); have_sth=true; }
00907       phi = (iphi-0.5) * halfpi/iring;
00908       }
00909     }
00910   else
00911     {
00912     int face_num, ix, iy;
00913     nest2xyf(pix,ix,iy,face_num);
00914 
00915     I jr = (I(jrll[face_num])<<order_) - ix - iy - 1;
00916 
00917     I nr;
00918     if (jr<nside_)
00919       {
00920       nr = jr;
00921       double tmp=(nr*nr)*fact2_;
00922       z = 1 - tmp;
00923       if (z>0.99) { sth=sqrt(tmp*(2.-tmp)); have_sth=true; }
00924       }
00925     else if (jr > 3*nside_)
00926       {
00927       nr = nside_*4-jr;
00928       double tmp=(nr*nr)*fact2_;
00929       z = tmp - 1;
00930       if (z<-0.99) { sth=sqrt(tmp*(2.-tmp)); have_sth=true; }
00931       }
00932     else
00933       {
00934       nr = nside_;
00935       z = (2*nside_-jr)*fact1_;
00936       }
00937 
00938     I tmp=I(jpll[face_num])*nr+ix-iy;
00939     planck_assert(tmp<8*nr,"must not happen");
00940     if (tmp<0) tmp+=8*nr;
00941     phi = (nr==nside_) ? 0.75*halfpi*tmp*fact1_ :
00942                          (0.5*halfpi*tmp)/nr;
00943     }
00944   }
00945 
00946 template<typename I> template<typename I2>
00947   void T_Healpix_Base<I>::query_polygon_internal
00948   (const vector<pointing> &vertex, int fact, rangeset<I2> &pixset) const
00949   {
00950   bool inclusive = (fact!=0);
00951   tsize nv=vertex.size();
00952   tsize ncirc = inclusive ? nv+1 : nv;
00953   planck_assert(nv>=3,"not enough vertices in polygon");
00954   arr<vec3> vv(nv);
00955   for (tsize i=0; i<nv; ++i)
00956     vv[i]=vertex[i].to_vec3();
00957   arr<vec3> normal(ncirc);
00958   int flip=0;
00959   for (tsize i=0; i<nv; ++i)
00960     {
00961     normal[i]=crossprod(vv[i],vv[(i+1)%nv]);
00962     double hnd=dotprod(normal[i],vv[(i+2)%nv]);
00963     planck_assert(abs(hnd)>1e-10,"degenerate corner");
00964     if (i==0)
00965       flip = (hnd<0.) ? -1 : 1;
00966     else
00967       planck_assert(flip*hnd>0,"polygon is not convex");
00968     normal[i]*=flip/normal[i].Length();
00969     }
00970   arr<double> rad(ncirc,halfpi);
00971   if (inclusive)
00972     {
00973     double cosrad;
00974     find_enclosing_circle (vv, normal[nv], cosrad);
00975     rad[nv]=acos(cosrad);
00976     }
00977   query_multidisc(normal,rad,fact,pixset);
00978   }
00979 
00980 template<typename I> void T_Healpix_Base<I>::query_polygon
00981   (const vector<pointing> &vertex, rangeset<I> &pixset) const
00982   {
00983   query_polygon_internal(vertex, 0, pixset);
00984   }
00985 
00986 template<typename I> void T_Healpix_Base<I>::query_polygon_inclusive
00987   (const vector<pointing> &vertex, rangeset<I> &pixset, int fact) const
00988   {
00989   planck_assert(fact>0,"fact must be a positive integer");
00990   if ((sizeof(I)<8) && (((I(1)<<order_max)/nside_)<fact))
00991     {
00992     T_Healpix_Base<int64> base2(nside_,scheme_,SET_NSIDE);
00993     base2.query_polygon_internal(vertex,fact,pixset);
00994     return;
00995     }
00996   query_polygon_internal(vertex, fact, pixset);
00997   }
00998 
00999 template<typename I> void T_Healpix_Base<I>::query_strip_internal
01000   (double theta1, double theta2, bool inclusive, rangeset<I> &pixset) const
01001   {
01002   if (scheme_==RING)
01003     {
01004     I ring1 = max(I(1),1+ring_above(cos(theta1))),
01005       ring2 = min(4*nside_-1,ring_above(cos(theta2)));
01006     if (inclusive)
01007       {
01008       ring1 = max(I(1),ring1-1);
01009       ring2 = min(4*nside_-1,ring2+1);
01010       }
01011 
01012     I sp1,rp1,sp2,rp2;
01013     bool dummy;
01014     get_ring_info_small(ring1,sp1,rp1,dummy);
01015     get_ring_info_small(ring2,sp2,rp2,dummy);
01016     I pix1 = sp1,
01017       pix2 = sp2+rp2;
01018     if (pix1<=pix2) pixset.append(pix1,pix2);
01019     }
01020   else
01021     planck_fail("query_strip not yet implemented for NESTED");
01022   }
01023 
01024 template<typename I> void T_Healpix_Base<I>::query_strip (double theta1,
01025   double theta2, bool inclusive, rangeset<I> &pixset) const
01026   {
01027   pixset.clear();
01028 
01029   if (theta1<theta2)
01030     query_strip_internal(theta1,theta2,inclusive,pixset);
01031   else
01032     {
01033     query_strip_internal(0.,theta2,inclusive,pixset);
01034     rangeset<I> ps2;
01035     query_strip_internal(theta1,pi,inclusive,ps2);
01036     pixset.append(ps2);
01037     }
01038   }
01039 
01040 template<typename I> inline void T_Healpix_Base<I>::get_ring_info_small
01041   (I ring, I &startpix, I &ringpix, bool &shifted) const
01042   {
01043   if (ring < nside_)
01044     {
01045     shifted = true;
01046     ringpix = 4*ring;
01047     startpix = 2*ring*(ring-1);
01048     }
01049   else if (ring < 3*nside_)
01050     {
01051     shifted = ((ring-nside_) & 1) == 0;
01052     ringpix = 4*nside_;
01053     startpix = ncap_ + (ring-nside_)*ringpix;
01054     }
01055   else
01056     {
01057     shifted = true;
01058     I nr= 4*nside_-ring;
01059     ringpix = 4*nr;
01060     startpix = npix_-2*nr*(nr+1);
01061     }
01062   }
01063 
01064 template<typename I> void T_Healpix_Base<I>::get_ring_info (I ring, I &startpix,
01065   I &ringpix, double &costheta, double &sintheta, bool &shifted) const
01066   {
01067   I northring = (ring>2*nside_) ? 4*nside_-ring : ring;
01068   if (northring < nside_)
01069     {
01070     double tmp = northring*northring*fact2_;
01071     costheta = 1 - tmp;
01072     sintheta = sqrt(tmp*(2-tmp));
01073     ringpix = 4*northring;
01074     shifted = true;
01075     startpix = 2*northring*(northring-1);
01076     }
01077   else
01078     {
01079     costheta = (2*nside_-northring)*fact1_;
01080     sintheta = sqrt((1+costheta)*(1-costheta));
01081     ringpix = 4*nside_;
01082     shifted = ((northring-nside_) & 1) == 0;
01083     startpix = ncap_ + (northring-nside_)*ringpix;
01084     }
01085   if (northring != ring) // southern hemisphere
01086     {
01087     costheta = -costheta;
01088     startpix = npix_ - startpix - ringpix;
01089     }
01090   }
01091 template<typename I> void T_Healpix_Base<I>::get_ring_info2 (I ring,
01092   I &startpix, I &ringpix, double &theta, bool &shifted) const
01093   {
01094   I northring = (ring>2*nside_) ? 4*nside_-ring : ring;
01095   if (northring < nside_)
01096     {
01097     double tmp = northring*northring*fact2_;
01098     double costheta = 1 - tmp;
01099     double sintheta = sqrt(tmp*(2-tmp));
01100     theta = atan2(sintheta,costheta);
01101     ringpix = 4*northring;
01102     shifted = true;
01103     startpix = 2*northring*(northring-1);
01104     }
01105   else
01106     {
01107     theta = acos((2*nside_-northring)*fact1_);
01108     ringpix = 4*nside_;
01109     shifted = ((northring-nside_) & 1) == 0;
01110     startpix = ncap_ + (northring-nside_)*ringpix;
01111     }
01112   if (northring != ring) // southern hemisphere
01113     {
01114     theta = pi-theta;
01115     startpix = npix_ - startpix - ringpix;
01116     }
01117   }
01118 
01119 template<typename I> void T_Healpix_Base<I>::neighbors (I pix,
01120   fix_arr<I,8> &result) const
01121   {
01122   int ix, iy, face_num;
01123   (scheme_==RING) ?
01124     ring2xyf(pix,ix,iy,face_num) : nest2xyf(pix,ix,iy,face_num);
01125 
01126   const I nsm1 = nside_-1;
01127   if ((ix>0)&&(ix<nsm1)&&(iy>0)&&(iy<nsm1))
01128     {
01129     if (scheme_==RING)
01130       for (int m=0; m<8; ++m)
01131         result[m] = xyf2ring(ix+nb_xoffset[m],iy+nb_yoffset[m],face_num);
01132     else
01133       {
01134       I fpix = I(face_num)<<(2*order_),
01135         px0=spread_bits(ix  ), py0=spread_bits(iy  )<<1,
01136         pxp=spread_bits(ix+1), pyp=spread_bits(iy+1)<<1,
01137         pxm=spread_bits(ix-1), pym=spread_bits(iy-1)<<1;
01138 
01139       result[0] = fpix+pxm+py0; result[1] = fpix+pxm+pyp;
01140       result[2] = fpix+px0+pyp; result[3] = fpix+pxp+pyp;
01141       result[4] = fpix+pxp+py0; result[5] = fpix+pxp+pym;
01142       result[6] = fpix+px0+pym; result[7] = fpix+pxm+pym;
01143       }
01144     }
01145   else
01146     {
01147     for (int i=0; i<8; ++i)
01148       {
01149       int x=ix+nb_xoffset[i], y=iy+nb_yoffset[i];
01150       int nbnum=4;
01151       if (x<0)
01152         { x+=nside_; nbnum-=1; }
01153       else if (x>=nside_)
01154         { x-=nside_; nbnum+=1; }
01155       if (y<0)
01156         { y+=nside_; nbnum-=3; }
01157       else if (y>=nside_)
01158         { y-=nside_; nbnum+=3; }
01159 
01160       int f = nb_facearray[nbnum][face_num];
01161       if (f>=0)
01162         {
01163         int bits = nb_swaparray[nbnum][face_num>>2];
01164         if (bits&1) x=nside_-x-1;
01165         if (bits&2) y=nside_-y-1;
01166         if (bits&4) std::swap(x,y);
01167         result[i] = (scheme_==RING) ? xyf2ring(x,y,f) : xyf2nest(x,y,f);
01168         }
01169       else
01170         result[i] = -1;
01171       }
01172     }
01173   }
01174 
01175 template<typename I> void T_Healpix_Base<I>::get_interpol (const pointing &ptg,
01176   fix_arr<I,4> &pix, fix_arr<double,4> &wgt) const
01177   {
01178   planck_assert((ptg.theta>=0)&&(ptg.theta<=pi),"invalid theta value");
01179   double z = cos (ptg.theta);
01180   I ir1 = ring_above(z);
01181   I ir2 = ir1+1;
01182   double theta1, theta2, w1, tmp, dphi;
01183   I sp,nr;
01184   bool shift;
01185   I i1,i2;
01186   if (ir1>0)
01187     {
01188     get_ring_info2 (ir1, sp, nr, theta1, shift);
01189     dphi = twopi/nr;
01190     tmp = (ptg.phi/dphi - .5*shift);
01191     i1 = (tmp<0) ? I(tmp)-1 : I(tmp);
01192     w1 = (ptg.phi-(i1+.5*shift)*dphi)/dphi;
01193     i2 = i1+1;
01194     if (i1<0) i1 +=nr;
01195     if (i2>=nr) i2 -=nr;
01196     pix[0] = sp+i1; pix[1] = sp+i2;
01197     wgt[0] = 1-w1; wgt[1] = w1;
01198     }
01199   if (ir2<(4*nside_))
01200     {
01201     get_ring_info2 (ir2, sp, nr, theta2, shift);
01202     dphi = twopi/nr;
01203     tmp = (ptg.phi/dphi - .5*shift);
01204     i1 = (tmp<0) ? I(tmp)-1 : I(tmp);
01205     w1 = (ptg.phi-(i1+.5*shift)*dphi)/dphi;
01206     i2 = i1+1;
01207     if (i1<0) i1 +=nr;
01208     if (i2>=nr) i2 -=nr;
01209     pix[2] = sp+i1; pix[3] = sp+i2;
01210     wgt[2] = 1-w1; wgt[3] = w1;
01211     }
01212 
01213   if (ir1==0)
01214     {
01215     double wtheta = ptg.theta/theta2;
01216     wgt[2] *= wtheta; wgt[3] *= wtheta;
01217     double fac = (1-wtheta)*0.25;
01218     wgt[0] = fac; wgt[1] = fac; wgt[2] += fac; wgt[3] +=fac;
01219     pix[0] = (pix[2]+2)&3;
01220     pix[1] = (pix[3]+2)&3;
01221     }
01222   else if (ir2==4*nside_)
01223     {
01224     double wtheta = (ptg.theta-theta1)/(pi-theta1);
01225     wgt[0] *= (1-wtheta); wgt[1] *= (1-wtheta);
01226     double fac = wtheta*0.25;
01227     wgt[0] += fac; wgt[1] += fac; wgt[2] = fac; wgt[3] =fac;
01228     pix[2] = ((pix[0]+2)&3)+npix_-4;
01229     pix[3] = ((pix[1]+2)&3)+npix_-4;
01230     }
01231   else
01232     {
01233     double wtheta = (ptg.theta-theta1)/(theta2-theta1);
01234     wgt[0] *= (1-wtheta); wgt[1] *= (1-wtheta);
01235     wgt[2] *= wtheta; wgt[3] *= wtheta;
01236     }
01237 
01238   if (scheme_==NEST)
01239     for (tsize m=0; m<pix.size(); ++m)
01240       pix[m] = ring2nest(pix[m]);
01241   }
01242 
01243 template<typename I> void T_Healpix_Base<I>::swap (T_Healpix_Base &other)
01244   {
01245   std::swap(order_,other.order_);
01246   std::swap(nside_,other.nside_);
01247   std::swap(npface_,other.npface_);
01248   std::swap(ncap_,other.ncap_);
01249   std::swap(npix_,other.npix_);
01250   std::swap(fact1_,other.fact1_);
01251   std::swap(fact2_,other.fact2_);
01252   std::swap(scheme_,other.scheme_);
01253   }
01254 
01255 template<typename I> double T_Healpix_Base<I>::max_pixrad() const
01256   {
01257   vec3 va,vb;
01258   va.set_z_phi (2./3., pi/(4*nside_));
01259   double t1 = 1.-1./nside_;
01260   t1*=t1;
01261   vb.set_z_phi (1-t1/3, 0);
01262   return v_angle(va,vb);
01263   }
01264 
01265 template<typename I> double T_Healpix_Base<I>::max_pixrad(I ring) const
01266   {
01267   if (ring>=2*nside_) ring=4*nside_-ring;
01268   double z=ring2z(ring), z_up=(ring>1) ? ring2z(ring-1) : 1.;
01269   vec3 mypos, uppos;
01270   uppos.set_z_phi(z_up,0);
01271   if (ring<=nside_)
01272     {
01273     mypos.set_z_phi(z,pi/(4*ring));
01274     return v_angle(mypos,uppos);
01275     }
01276   mypos.set_z_phi(z,0);
01277   double vdist=v_angle(mypos,uppos);
01278   double hdist=sqrt(1.-z*z)*pi/(4*nside_);
01279   return max(hdist,vdist);
01280   }
01281 
01282 template<typename I> void T_Healpix_Base<I>::xyf2loc (double x, double y,
01283   int face, double &z, double &phi, double &sth, bool &have_sth) const
01284   {
01285   have_sth = false;
01286   double jr = jrll[face] - x - y;
01287   double nr;
01288   if (jr<1)
01289     {
01290     nr = jr;
01291     double tmp = nr*nr/3.;
01292     z = 1 - tmp;
01293     if (z > 0.99)
01294       {
01295       sth = std::sqrt(tmp*(2.0-tmp));
01296       have_sth = true;
01297       }
01298     }
01299   else if (jr>3)
01300     {
01301     nr = 4-jr;
01302     double tmp = nr*nr/3.;
01303     z = tmp - 1;
01304     if (z<-0.99)
01305       {
01306       sth = std::sqrt(tmp*(2.-tmp));
01307       have_sth = true;
01308       }
01309     }
01310   else
01311     {
01312     nr = 1;
01313     z = (2-jr)*2./3.;
01314     }
01315 
01316   double tmp=jpll[face]*nr+x-y;
01317   if (tmp<0) tmp+=8;
01318   if (tmp>=8) tmp-=8;
01319   phi = (nr<1e-15) ? 0 : (0.5*halfpi*tmp)/nr;
01320   }
01321 
01322 namespace {
01323 
01324 vec3 locToVec3 (double z, double phi, double sth, bool have_sth)
01325   {
01326   if (have_sth)
01327     return vec3(sth*cos(phi),sth*sin(phi),z);
01328   else
01329     {
01330     vec3 res;
01331     res.set_z_phi (z, phi);
01332     return res;
01333     }
01334   }
01335 
01336 } // unnamed namespace
01337 
01338 template<typename I> void T_Healpix_Base<I>::boundaries(I pix, tsize step,
01339   vector<vec3> &out) const
01340   {
01341   out.resize(4*step);
01342   int ix, iy, face;
01343   pix2xyf(pix, ix, iy, face);
01344   double dc = 0.5 / nside_;
01345   double xc = (ix + 0.5)/nside_, yc = (iy + 0.5)/nside_;
01346   double d = 1.0/(step*nside_);
01347   for (tsize i=0; i<step; ++i)
01348     {
01349     double z, phi, sth;
01350     bool have_sth;
01351     xyf2loc(xc+dc-i*d, yc+dc, face, z, phi, sth, have_sth);
01352     out[i] = locToVec3(z, phi, sth, have_sth);
01353     xyf2loc(xc-dc, yc+dc-i*d, face, z, phi, sth, have_sth);
01354     out[i+step] = locToVec3(z, phi, sth, have_sth);
01355     xyf2loc(xc-dc+i*d, yc-dc, face, z, phi, sth, have_sth);
01356     out[i+2*step] = locToVec3(z, phi, sth, have_sth);
01357     xyf2loc(xc+dc, yc-dc+i*d, face, z, phi, sth, have_sth);
01358     out[i+3*step] = locToVec3(z, phi, sth, have_sth);
01359     }
01360   }
01361 
01362 template<typename I> arr<int> T_Healpix_Base<I>::swap_cycles() const
01363   {
01364   planck_assert(order_>=0, "need hierarchical map");
01365   planck_assert(order_<=13, "map too large");
01366   arr<int> result(swap_clen[order_]);
01367   tsize ofs=0;
01368   for (int m=0; m<order_;++m) ofs+=swap_clen[m];
01369   for (tsize m=0; m<result.size();++m) result[m]=swap_cycle[m+ofs];
01370   return result;
01371   }
01372 
01373 template class T_Healpix_Base<int>;
01374 template class T_Healpix_Base<int64>;

Generated on Wed Apr 24 11:31:18 2013 for Healpix C++