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

Generated on Thu Oct 8 14:48:52 2015 for Healpix C++