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

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