healpix_base.h

Go to the documentation of this file.
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 /*! \file healpix_base.h
00028  *  Copyright (C) 2003-2012 Max-Planck-Society
00029  *  \author Martin Reinecke
00030  */
00031 
00032 #ifndef HEALPIX_BASE_H
00033 #define HEALPIX_BASE_H
00034 
00035 #include <vector>
00036 #include "healpix_tables.h"
00037 #include "pointing.h"
00038 #include "arr.h"
00039 #include "rangeset.h"
00040 
00041 /*! Functionality related to the HEALPix pixelisation. */
00042 template<typename I> class T_Healpix_Base: public Healpix_Tables
00043   {
00044   protected:
00045     /*! The order of the map; -1 for nonhierarchical map. */
00046     int order_;
00047     /*! The N_side parameter of the map; 0 if not allocated. */
00048     I nside_;
00049     I npface_, ncap_, npix_;
00050     double fact1_, fact2_;
00051     /*! The map's ordering scheme. */
00052     Healpix_Ordering_Scheme scheme_;
00053 
00054     /*! Returns the number of the next ring to the north of \a z=cos(theta).
00055         It may return 0; in this case \a z lies north of all rings. */
00056     inline I ring_above (double z) const;
00057     void in_ring (I iz, double phi0, double dphi, rangeset<I> &pixset) const;
00058 
00059     template<typename I2> void query_multidisc (const arr<vec3> &norm,
00060       const arr<double> &rad, int fact, rangeset<I2> &pixset) const;
00061 
00062     void query_multidisc_general (const arr<vec3> &norm, const arr<double> &rad,
00063       bool inclusive, const std::vector<int> &cmds, rangeset<I> &pixset) const;
00064 
00065     void query_strip_internal (double theta1, double theta2, bool inclusive,
00066       rangeset<I> &pixset) const;
00067 
00068     inline I spread_bits (int v) const;
00069     inline int compress_bits (I v) const;
00070 
00071     I xyf2nest(int ix, int iy, int face_num) const;
00072     void nest2xyf(I pix, int &ix, int &iy, int &face_num) const;
00073     I xyf2ring(int ix, int iy, int face_num) const;
00074     void ring2xyf(I pix, int &ix, int &iy, int &face_num) const;
00075 
00076     I loc2pix (double z, double phi, double sth, bool have_sth) const;
00077     void pix2loc (I pix, double &z, double &phi, double &sth, bool &have_sth)
00078       const;
00079 
00080     void xyf2loc(double x, double y, int face, double &z, double &ph,
00081       double &sth, bool &have_sth) const;
00082 
00083     I nest_peano_helper (I pix, int dir) const;
00084 
00085     typedef I (T_Healpix_Base::*swapfunc)(I pix) const;
00086 
00087   public:
00088     static const int order_max;
00089 
00090     /*! Calculates the map order from its \a N_side parameter.
00091         Returns -1 if \a nside is not a power of 2.
00092         \param nside the \a N_side parameter */
00093     static int nside2order (I nside);
00094     /*! Calculates the \a N_side parameter from the number of pixels.
00095         \param npix the number of pixels */
00096     static I npix2nside (I npix);
00097     /*! Constructs an unallocated object. */
00098     T_Healpix_Base ();
00099     /*! Constructs an object with a given \a order and the ordering
00100         scheme \a scheme. */
00101     T_Healpix_Base (int order, Healpix_Ordering_Scheme scheme)
00102       { Set (order, scheme); }
00103     /*! Constructs an object with a given \a nside and the ordering
00104         scheme \a scheme. The \a nside_dummy parameter must be set to
00105         SET_NSIDE. */
00106     T_Healpix_Base (I nside, Healpix_Ordering_Scheme scheme, const nside_dummy)
00107       { SetNside (nside, scheme); }
00108 
00109     /*! Adjusts the object to \a order and \a scheme. */
00110     void Set (int order, Healpix_Ordering_Scheme scheme);
00111     /*! Adjusts the object to \a nside and \a scheme. */
00112     void SetNside (I nside, Healpix_Ordering_Scheme scheme);
00113 
00114     /*! Returns the z-coordinate of the ring \a ring. This also works
00115         for the (not really existing) rings 0 and 4*nside. */
00116     double ring2z (I ring) const;
00117     /*! Returns the number of the ring in which \a pix lies. */
00118     I pix2ring (I pix) const;
00119 
00120     I xyf2pix(int ix, int iy, int face_num) const
00121       {
00122       return (scheme_==RING) ?
00123         xyf2ring(ix,iy,face_num) : xyf2nest(ix,iy,face_num);
00124       }
00125     void pix2xyf(I pix, int &ix, int &iy, int &face_num) const
00126       {
00127       (scheme_==RING) ?
00128         ring2xyf(pix,ix,iy,face_num) : nest2xyf(pix,ix,iy,face_num);
00129       }
00130 
00131     /*! Translates a pixel number from NEST to RING. */
00132     I nest2ring (I pix) const;
00133     /*! Translates a pixel number from RING to NEST. */
00134     I ring2nest (I pix) const;
00135     /*! Translates a pixel number from NEST to its Peano index. */
00136     I nest2peano (I pix) const;
00137     /*! Translates a pixel number from its Peano index to NEST. */
00138     I peano2nest (I pix) const;
00139 
00140     /*! Returns the number of the pixel which contains the angular coordinates
00141         (\a z:=cos(theta), \a phi).
00142         \note This method is inaccurate near the poles at high resolutions. */
00143     I zphi2pix (double z, double phi) const
00144       { return loc2pix(z,phi,0.,false); }
00145 
00146     /*! Returns the number of the pixel which contains the angular coordinates
00147         \a ang. */
00148     I ang2pix (const pointing &ang) const
00149       {
00150       const double pi_=3.141592653589793238462643383279502884197;
00151       planck_assert((ang.theta>=0)&&(ang.theta<=pi_),"invalid theta value");
00152       return ((ang.theta<0.01) || (ang.theta > 3.14159-0.01)) ?
00153         loc2pix(cos(ang.theta),ang.phi,sin(ang.theta),true) :
00154         loc2pix(cos(ang.theta),ang.phi,0.,false);
00155       }
00156     /*! Returns the number of the pixel which contains the vector \a vec
00157         (\a vec is normalized if necessary). */
00158     I vec2pix (const vec3 &vec) const
00159       {
00160       double xl = 1./vec.Length();
00161       double phi = safe_atan2(vec.y,vec.x);
00162       double nz = vec.z*xl;
00163       if (std::abs(nz)>0.99)
00164         return loc2pix (nz,phi,sqrt(vec.x*vec.x+vec.y*vec.y)*xl,true);
00165       else
00166         return loc2pix (nz,phi,0,false);
00167       }
00168 
00169     /*! Returns the angular coordinates (\a z:=cos(theta), \a phi) of the center
00170         of the pixel with number \a pix.
00171         \note This method is inaccurate near the poles at high resolutions. */
00172     void pix2zphi (I pix, double &z, double &phi) const
00173       {
00174       bool dum_b;
00175       double dum_d;
00176       pix2loc(pix,z,phi,dum_d,dum_b);
00177       }
00178 
00179     /*! Returns the angular coordinates of the center of the pixel with
00180         number \a pix. */
00181     pointing pix2ang (I pix) const
00182       {
00183       double z, phi, sth;
00184       bool have_sth;
00185       pix2loc (pix,z,phi,sth,have_sth);
00186       return have_sth ? pointing(atan2(sth,z),phi) : pointing(acos(z),phi);
00187       }
00188     /*! Returns the vector to the center of the pixel with number \a pix. */
00189     vec3 pix2vec (I pix) const
00190       {
00191       double z, phi, sth;
00192       bool have_sth;
00193       pix2loc (pix,z,phi,sth,have_sth);
00194       if (have_sth)
00195         return vec3(sth*cos(phi),sth*sin(phi),z);
00196       else
00197         {
00198         vec3 res;
00199         res.set_z_phi (z, phi);
00200         return res;
00201         }
00202       }
00203 
00204     template<typename I2> void query_disc_internal (pointing ptg, double radius,
00205       int fact, rangeset<I2> &pixset) const;
00206 
00207     /*! Returns the range set of all pixels whose centers lie within the disk
00208         defined by \a dir and \a radius.
00209         \param dir the angular coordinates of the disk center
00210         \param radius the radius (in radians) of the disk
00211         \param pixset a \a rangeset object containing the indices of all pixels
00212            whose centers lie inside the disk
00213         \note This method is more efficient in the RING scheme. */
00214     void query_disc (pointing ptg, double radius, rangeset<I> &pixset) const;
00215     /*! Returns the range set of all pixels which overlap with the disk
00216         defined by \a dir and \a radius.
00217         \param dir the angular coordinates of the disk center
00218         \param radius the radius (in radians) of the disk
00219         \param pixset a \a rangeset object containing the indices of all pixels
00220            overlapping with the disk.
00221         \param fact The overlapping test will be done at the resolution
00222            \a fact*nside. For NESTED ordering, \a fact must be a power of 2,
00223            else it can be any positive integer. A typical choice would be 4.
00224         \note This method may return some pixels which don't overlap with
00225            the disk at all. The higher \a fact is chosen, the fewer false
00226            positives are returned, at the cost of increased run time.
00227         \note This method is more efficient in the RING scheme. */
00228     void query_disc_inclusive (pointing ptg, double radius, rangeset<I> &pixset,
00229       int fact=1) const;
00230 
00231     /*! \deprecated Please use the version based on \a rangeset */
00232     void query_disc (const pointing &dir, double radius,
00233       std::vector<I> &listpix) const
00234       {
00235       rangeset<I> pixset;
00236       query_disc(dir,radius,pixset);
00237       pixset.toVector(listpix);
00238       }
00239     /*! \deprecated Please use the version based on \a rangeset */
00240     void query_disc_inclusive (const pointing &dir, double radius,
00241       std::vector<I> &listpix, int fact=1) const
00242       {
00243       rangeset<I> pixset;
00244       query_disc_inclusive(dir,radius,pixset,fact);
00245       pixset.toVector(listpix);
00246       }
00247 
00248     template<typename I2> void query_polygon_internal
00249       (const std::vector<pointing> &vertex, int fact,
00250       rangeset<I2> &pixset) const;
00251 
00252     /*! Returns a range set of pixels whose centers lie within the convex
00253         polygon defined by the \a vertex array.
00254         \param vertex array containing the vertices of the polygon.
00255         \param pixset a \a rangeset object containing the indices of all pixels
00256            whose centers lie inside the polygon
00257         \note This method is more efficient in the RING scheme. */
00258     void query_polygon (const std::vector<pointing> &vertex,
00259       rangeset<I> &pixset) const;
00260 
00261     /*! Returns a range set of pixels which overlap with the convex
00262         polygon defined by the \a vertex array.
00263         \param vertex array containing the vertices of the polygon.
00264         \param pixset a \a rangeset object containing the indices of all pixels
00265            overlapping with the polygon.
00266         \param fact The overlapping test will be done at the resolution
00267            \a fact*nside. For NESTED ordering, \a fact must be a power of 2,
00268            else it can be any positive integer. A typical choice would be 4.
00269         \note This method may return some pixels which don't overlap with
00270            the polygon at all. The higher \a fact is chosen, the fewer false
00271            positives are returned, at the cost of increased run time.
00272         \note This method is more efficient in the RING scheme. */
00273     void query_polygon_inclusive (const std::vector<pointing> &vertex,
00274       rangeset<I> &pixset, int fact=1) const;
00275 
00276     /*! Returns a range set of pixels whose centers lie within the colatitude
00277         range defined by \a theta1 and \a theta2 (if \a inclusive==false), or
00278         which overlap with this region (if \a inclusive==true). If
00279         \a theta1<theta2, the region between both angles is considered,
00280         otherwise the regions \a 0<theta<theta2 and \a theta1<theta<pi.
00281         \param theta1 first colatitude
00282         \param theta2 second colatitude
00283         \param inclusive if \a false, return the exact set of pixels whose
00284            pixels centers lie within the region; if \a true, return all pixels
00285            that overlap with the region. */
00286     void query_strip (double theta1, double theta2, bool inclusive,
00287       rangeset<I> &pixset) const;
00288 
00289     /*! Returns useful information about a given ring of the map.
00290         \param ring the ring number (the number of the first ring is 1)
00291         \param startpix the number of the first pixel in the ring
00292                (NOTE: this is always given in the RING numbering scheme!)
00293         \param ringpix the number of pixels in the ring
00294         \param costheta the cosine of the colatitude of the ring
00295         \param sintheta the sine of the colatitude of the ring
00296         \param shifted if \a true, the center of the first pixel is not at
00297                \a phi=0 */
00298     void get_ring_info (I ring, I &startpix, I &ringpix,
00299       double &costheta, double &sintheta, bool &shifted) const;
00300     /*! Returns useful information about a given ring of the map.
00301         \param ring the ring number (the number of the first ring is 1)
00302         \param startpix the number of the first pixel in the ring
00303                (NOTE: this is always given in the RING numbering scheme!)
00304         \param ringpix the number of pixels in the ring
00305         \param theta the colatitude (in radians) of the ring
00306         \param shifted if \a true, the center of the first pixel is not at
00307                \a phi=0 */
00308     void get_ring_info2 (I ring, I &startpix, I &ringpix,
00309       double &theta, bool &shifted) const;
00310     /*! Returns useful information about a given ring of the map.
00311         \param ring the ring number (the number of the first ring is 1)
00312         \param startpix the number of the first pixel in the ring
00313                (NOTE: this is always given in the RING numbering scheme!)
00314         \param ringpix the number of pixels in the ring
00315         \param shifted if \a true, the center of the first pixel is not at
00316                \a phi=0 */
00317     void get_ring_info_small (I ring, I &startpix, I &ringpix,
00318         bool &shifted) const;
00319     /*! Returns the neighboring pixels of \a pix in \a result.
00320         On exit, \a result contains (in this order)
00321         the pixel numbers of the SW, W, NW, N, NE, E, SE and S neighbor
00322         of \a pix. If a neighbor does not exist (this can only be the case
00323         for the W, N, E and S neighbors), its entry is set to -1.
00324 
00325         \note This method works in both RING and NEST schemes, but is
00326           considerably faster in the NEST scheme. */
00327     void neighbors (I pix, fix_arr<I,8> &result) const;
00328     /*! Returns interpolation information for the direction \a ptg.
00329         The surrounding pixels are returned in \a pix, their corresponding
00330         weights in \a wgt.
00331         \note This method works in both RING and NEST schemes, but is
00332           considerably faster in the RING scheme. */
00333     void get_interpol (const pointing &ptg, fix_arr<I,4> &pix,
00334                        fix_arr<double,4> &wgt) const;
00335 
00336     /*! Returns the order parameter of the object. */
00337     int Order() const { return order_; }
00338     /*! Returns the \a N_side parameter of the object. */
00339     I Nside() const { return nside_; }
00340     /*! Returns the number of pixels of the object. */
00341     I Npix() const { return npix_; }
00342     /*! Returns the ordering scheme of the object. */
00343     Healpix_Ordering_Scheme Scheme() const { return scheme_; }
00344 
00345     /*! Returns \a true, if both objects have the same nside and scheme,
00346         else \a false. */
00347     bool conformable (const T_Healpix_Base &other) const
00348       { return ((nside_==other.nside_) && (scheme_==other.scheme_)); }
00349 
00350     /*! Swaps the contents of two Healpix_Base objects. */
00351     void swap (T_Healpix_Base &other);
00352 
00353     /*! Returns the maximum angular distance (in radian) between any pixel
00354         center and its corners. */
00355     double max_pixrad() const;
00356 
00357     /*! Returns the maximum angular distance (in radian) between any pixel
00358         center and its corners in a given ring. */
00359     double max_pixrad(I ring) const;
00360 
00361     /*! Returns a set of points along the boundary of the given pixel.
00362         \a step=1 gives 4 points on the corners. The first point corresponds
00363         to the northernmost corner, the subsequent points follow the pixel
00364         boundary through west, south and east corners.
00365         \param pix pixel index number
00366         \param step the number of returned points is 4*step. */
00367     void boundaries (I pix, tsize step, std::vector<vec3> &out) const;
00368 
00369     arr<int> swap_cycles() const;
00370   };
00371 
00372 /*! T_Healpix_Base for Nside up to 2^13. */
00373 typedef T_Healpix_Base<int> Healpix_Base;
00374 /*! T_Healpix_Base for Nside up to 2^29. */
00375 typedef T_Healpix_Base<int64> Healpix_Base2;
00376 
00377 #endif

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