hpxtest.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) 2004-2012 Max-Planck-Society
00029  *  Author: Martin Reinecke
00030  */
00031 
00032 /*
00033 
00034 Candidates for testing the validity of the Healpix routines:
00035 
00036 - done: ang2pix(pix2ang(i)) = i for all Healpix_Bases
00037 - done: pix2ang(ang2pix(ptg)) dot ptg > 1-delta for all Healpix_Bases
00038 - done: ring2nest(nest2ring(i)) = i for all hierarchical Healpix_Bases
00039 - done: downgrade(upgrade(map)) = map for all maps
00040 - done: map and downgraded map should have same average
00041 - done: alm2map(map2alm(map)) approx map (same for pol)
00042 - partly done: neighbor tests
00043 - powspec -> alm -> powspec (should produce similar powspecs, also for pol)
00044 - done: two swap_schemes() should produce original map
00045 - done: query_disc tests (dot products etc.)
00046 - a_lms: test Set(), Scale(), Add(), alm(l,m) = alm.mstart(m)[l], etc.
00047 
00048 */
00049 
00050 #include <iostream>
00051 #include "healpix_base.h"
00052 #include "healpix_map.h"
00053 #include "arr.h"
00054 #include "planck_rng.h"
00055 #include "lsconstants.h"
00056 #include "alm.h"
00057 #include "alm_healpix_tools.h"
00058 #include "alm_powspec_tools.h"
00059 #include "geom_utils.h"
00060 #include "walltimer.h"
00061 #include "announce.h"
00062 
00063 using namespace std;
00064 
00065 const int nsamples = 1000000;
00066 
00067 planck_rng rng;
00068 
00069 namespace {
00070 
00071 void random_dir (pointing &ptg)
00072   {
00073   ptg.theta = acos(rng.rand_uni()*2-1);
00074   ptg.phi = rng.rand_uni()*twopi;
00075   }
00076 void random_zphi (double &z, double &phi)
00077   {
00078   z = rng.rand_uni()*2-1;
00079   phi = rng.rand_uni()*twopi;
00080   }
00081 
00082 template<typename I> string bname()
00083   { return string("(basetype: ")+type2typename<I>()+")"; }
00084 
00085 template<typename I> void check_ringnestring()
00086   {
00087   cout << "testing ring2nest(nest2ring(m))==m " << bname<I>() << endl;
00088   for (int order=0; order<=T_Healpix_Base<I>::order_max; ++order)
00089     {
00090     T_Healpix_Base<I> base (order,RING);
00091     for (int m=0; m<nsamples; ++m)
00092       {
00093       I pix = I(rng.rand_uni()*base.Npix());
00094       if (base.ring2nest(base.nest2ring(pix))!=pix)
00095         cout << "  PROBLEM: order = " << order << ", pixel = " << pix << endl;
00096       }
00097     }
00098   }
00099 
00100 template<typename I> void check_nestpeanonest()
00101   {
00102   cout << "testing peano2nest(nest2peano(m))==m " << bname<I>() << endl;
00103   for (int order=0; order<=T_Healpix_Base<I>::order_max; ++order)
00104     {
00105     T_Healpix_Base<I> base (order,NEST);
00106     for (int m=0; m<nsamples; ++m)
00107       {
00108       I pix = I(rng.rand_uni()*base.Npix());
00109       if (base.peano2nest(base.nest2peano(pix))!=pix)
00110         cout << "  PROBLEM: order = " << order << ", pixel = " << pix << endl;
00111       }
00112     }
00113   }
00114 
00115 template<typename I> void check_pixzphipix()
00116   {
00117   cout << "testing zphi2pix(pix2zphi(m))==m " << bname<I>() << endl;
00118   int omax=T_Healpix_Base<I>::order_max;
00119   for (int order=0; order<=omax; ++order)
00120     {
00121     T_Healpix_Base<I> base1 (order,RING), base2 (order,NEST);
00122     for (int m=0; m<nsamples; ++m)
00123       {
00124       double z,phi;
00125       I pix = I(rng.rand_uni()*base1.Npix());
00126       base1.pix2zphi(pix,z,phi);
00127       if (base1.zphi2pix(z,phi)!=pix)
00128         cout << "  PROBLEM: order = " << order << ", pixel = " << pix << endl;
00129       base2.pix2zphi(pix,z,phi);
00130       if (base2.zphi2pix(z,phi)!=pix)
00131         cout << "  PROBLEM: order = " << order << ", pixel = " << pix << endl;
00132       }
00133     }
00134   for (I nside=3; nside<(I(1)<<omax); nside+=nside/2+1)
00135     {
00136     T_Healpix_Base<I> base (nside,RING,SET_NSIDE);
00137     for (int m=0; m<nsamples; ++m)
00138       {
00139       double z,phi;
00140       I pix = I(rng.rand_uni()*base.Npix());
00141       base.pix2zphi(pix,z,phi);
00142       if (base.zphi2pix(z,phi)!=pix)
00143         cout << "  PROBLEM: nside = " << nside << ", pixel = " << pix << endl;
00144       }
00145     }
00146   }
00147 
00148 template<typename I> void check_zphipixzphi()
00149   {
00150   cout << "testing pix2zphi(zphi2pix(ptg)) approx zphi " << bname<I>() << endl;
00151   int omax=T_Healpix_Base<I>::order_max;
00152   for (int order=0; order<=omax; ++order)
00153     {
00154     T_Healpix_Base<I> base1 (order,NEST), base2 (order,RING);
00155     double mincos = min (cos(base1.max_pixrad()),0.999999999999999);
00156     for (int m=0; m<nsamples; ++m)
00157       {
00158       double z,phi,z2,phi2;
00159       random_zphi (z,phi);
00160       base1.pix2zphi(base1.zphi2pix(z,phi),z2,phi2);
00161       if (cosdist_zphi(z,phi,z2,phi2)<mincos)
00162         cout << "  PROBLEM: order = " << order
00163              << ", zphi = " << z << ", " << phi << endl;
00164       base2.pix2zphi(base2.zphi2pix(z,phi),z2,phi2);
00165       if (cosdist_zphi(z,phi,z2,phi2)<mincos)
00166         cout << "  PROBLEM: order = " << order
00167              << ", zphi = " << z << ", " << phi << endl;
00168       }
00169     }
00170   for (int nside=3; nside<(I(1)<<omax); nside+=nside/2+1)
00171     {
00172     T_Healpix_Base<I> base (nside,RING,SET_NSIDE);
00173     double mincos = min (cos(base.max_pixrad()),0.999999999999999);
00174     for (int m=0; m<nsamples; ++m)
00175       {
00176       double z,phi,z2,phi2;
00177       random_zphi (z,phi);
00178       base.pix2zphi(base.zphi2pix(z,phi),z2,phi2);
00179       if (cosdist_zphi(z,phi,z2,phi2)<mincos)
00180         cout << "  PROBLEM: nside = " << nside
00181              << ", zphi = " << z << ", " << phi << endl;
00182       }
00183     }
00184   }
00185 
00186 template<typename I> void check_pixangpix()
00187   {
00188   cout << "testing ang2pix(pix2ang(m))==m " << bname<I>() << endl;
00189   int omax=T_Healpix_Base<I>::order_max;
00190   for (int order=0; order<=omax; ++order)
00191     {
00192     T_Healpix_Base<I> base1 (order,RING), base2 (order,NEST);
00193     for (int m=0; m<nsamples; ++m)
00194       {
00195       I pix = I(rng.rand_uni()*base1.Npix());
00196       if (base1.ang2pix(base1.pix2ang(pix))!=pix)
00197         cout << "  PROBLEM: order = " << order << ", pixel = " << pix << endl;
00198       if (base2.ang2pix(base2.pix2ang(pix))!=pix)
00199         cout << "  PROBLEM: order = " << order << ", pixel = " << pix << endl;
00200       }
00201     }
00202   for (I nside=3; nside<(I(1)<<omax); nside+=nside/2+1)
00203     {
00204     T_Healpix_Base<I> base (nside,RING,SET_NSIDE);
00205     for (int m=0; m<nsamples; ++m)
00206       {
00207       I pix = I(rng.rand_uni()*base.Npix());
00208       if (base.ang2pix(base.pix2ang(pix))!=pix)
00209         cout << "  PROBLEM: nside = " << nside << ", pixel = " << pix << endl;
00210       }
00211     }
00212   }
00213 
00214 template<typename I> void check_neighbors()
00215   {
00216   cout << "testing neighbor function " << bname<I>() << endl;
00217   int omax=T_Healpix_Base<I>::order_max;
00218   for (int order=0; order<=omax; ++order)
00219     {
00220     T_Healpix_Base<I> base (order,NEST), base2(order,RING);
00221     double maxang = 2.01*base.max_pixrad();
00222     for (int m=0; m<nsamples/10; ++m)
00223       {
00224       I pix = I(rng.rand_uni()*base.Npix());
00225       fix_arr<I,8> nb,nb2;
00226       vec3 pixpt = base.pix2vec(pix);
00227       base.neighbors(pix,nb);
00228       base2.neighbors(base.nest2ring(pix),nb2);
00229       for (int n=0; n<8; ++n)
00230         if (nb[n]<0)
00231           planck_assert(nb2[n]<0,"neighbor inconsistency");
00232         else
00233           planck_assert(base.nest2ring(nb[n])==nb2[n],"neighbor inconsistency");
00234       sort(&nb[0],&nb[0]+8);
00235       int check=0;
00236       for (int n=0; n<8; ++n)
00237         {
00238         if (nb[n]<0)
00239           ++check;
00240         else
00241           {
00242           if (v_angle(base.pix2vec(nb[n]),pixpt)>maxang)
00243             cout << " PROBLEM: order = " << order << ", pix = " << pix << endl;
00244           if ((n>0) && (nb[n]==nb[n-1]))
00245             cout << " PROBLEM: order = " << order << ", pix = " << pix << endl;
00246           }
00247         }
00248       planck_assert((check<=1)||((order==0)&&(check<=2)),"too few neighbors");
00249       }
00250     }
00251   for (I nside=3; nside<(I(1)<<omax); nside+=nside/2+1)
00252     {
00253     T_Healpix_Base<I> base (nside,RING,SET_NSIDE);
00254     double maxang = 2.01*base.max_pixrad();
00255     for (int m=0; m<nsamples/10; ++m)
00256       {
00257       I pix = I(rng.rand_uni()*base.Npix());
00258       fix_arr<I,8> nb;
00259       vec3 pixpt = base.pix2vec(pix);
00260       base.neighbors(pix,nb);
00261       for (int n=0; n<8; ++n)
00262         if ((nb[n]>=0) && (v_angle(base.pix2vec(nb[n]),pixpt)>maxang))
00263           cout << "  PROBLEM: nside = " << nside << ", pix = " << pix << endl;
00264       }
00265     }
00266   }
00267 
00268 void check_swap_scheme()
00269   {
00270   cout << "testing whether double swap_scheme() returns the original map"
00271        << endl << "(for orders 0 to 10)." << endl;
00272   for (int order=0; order<=10; ++order)
00273     {
00274     Healpix_Map<uint8> map(order,NEST);
00275     for (int m=0; m<map.Npix(); ++m) map[m]=uint8(m&0xFF);
00276     map.swap_scheme();
00277     map.swap_scheme();
00278     for (int m=0; m<map.Npix(); ++m)
00279       if (map[m]!=(m&0xFF))
00280         cout << "  PROBLEM: order = " << order << ", pix = " << m << endl;
00281     }
00282   }
00283 
00284 void check_query_disc_strict (Healpix_Ordering_Scheme scheme)
00285   {
00286   cout << "testing whether all pixels found by query_disc() really" << endl
00287        << "lie inside the disk (and vice versa)" << endl;
00288   cout << "Ordering scheme: " << (scheme==RING ? "RING" : "NEST") << endl;
00289   for (int order=0; order<=5; ++order)
00290     {
00291     Healpix_Map<bool> map (order,scheme);
00292     map.fill(false);
00293     Healpix_Map<vec3> vmap(order,scheme);
00294     for (int m=0; m<vmap.Npix(); ++m)
00295       vmap[m]=vmap.pix2vec(m);
00296     rangeset<int> pixset;
00297     for (int m=0; m<100000; ++m)
00298       {
00299       pointing ptg;
00300       random_dir (ptg);
00301       double rad = pi/1 * rng.rand_uni();
00302       map.query_disc(ptg,rad,pixset);
00303       vec3 vptg=ptg;
00304       double cosrad=cos(rad);
00305       for (tsize j=0; j<pixset.size(); ++j)
00306         for (int i=pixset.ivbegin(j); i<pixset.ivend(j); ++i)
00307           map[i] = true;
00308       for (int i=0; i<map.Npix(); ++i)
00309         {
00310         bool inside = dotprod(vmap[i],vptg)>cosrad;
00311         if (inside^map[i])
00312           cout << "  PROBLEM: order = " << order << ", ptg = " << ptg << endl;
00313         }
00314       for (tsize j=0; j<pixset.size(); ++j)
00315         for (int i=pixset.ivbegin(j); i<pixset.ivend(j); ++i)
00316           map[i] = false;
00317       }
00318     }
00319   }
00320 
00321 template<typename I>void check_query_disc()
00322   {
00323   cout << "checking query_disc() " << bname<I>() << endl;
00324   int omax=min(20,T_Healpix_Base<I>::order_max);
00325   for (int order=0; order<=omax; ++order)
00326     {
00327     T_Healpix_Base<I> rbase (order,RING), nbase (order,NEST);
00328     rangeset<I> pr,pn;
00329     int niter=max(1,min(1000,100000>>order));
00330     for (int m=0; m<niter; ++m)
00331       {
00332       pointing ptg;
00333       random_dir (ptg);
00334       double rad = pi/1 * rng.rand_uni();
00335       rbase.query_disc(ptg,rad,pr);
00336       nbase.query_disc(ptg,rad,pn);
00337       if (pn.nval()!=pr.nval())
00338         cout << "PROBLEM: set size mismatch" << endl;
00339       if (order<5) // takes too long for high resolutions
00340         for (tsize i=0; i<pr.size(); ++i)
00341           for (tsize j=pr.ivbegin(i);j<pr.ivend(i); ++j)
00342             if (!pn.contains(rbase.ring2nest(j)))
00343               cout << "PROBLEM: set mismatch" << endl;
00344 
00345       for (tsize fct=8; fct>0; fct>>=1)
00346         {
00347         rangeset<I> pri,pni;
00348         rbase.query_disc_inclusive(ptg,rad,pri,fct);
00349         nbase.query_disc_inclusive(ptg,rad,pni,fct);
00350         if (pni.nval()>pri.nval())
00351           cout << "PROBLEM: RING inclusive < NEST inclusive" << endl;
00352         if (pni.nval()<pri.nval())
00353           cout << "Warning: RING inclusive > NEST inclusive "
00354                << order << " " << fct << endl;
00355         if (order<5) // takes too long for high resolutions
00356           {
00357           for (tsize i=0; i<pni.size(); ++i)
00358             for (tsize j=pni.ivbegin(i);j<pni.ivend(i); ++j)
00359               if (!pri.contains(rbase.nest2ring(j)))
00360                 cout << "PROBLEM: NEST inclusive not in RING inclusive" << endl;
00361           for (tsize i=0; i<pn.size(); ++i)
00362             for (tsize j=pn.ivbegin(i);j<pn.ivend(i); ++j)
00363               if (!pni.contains(j))
00364                 cout << "PROBLEM: Inclusive not superset of normal" << endl;
00365           }
00366         }
00367       }
00368     }
00369   }
00370 template<typename I>void check_query_polygon()
00371   {
00372   cout << "checking query_polygon() " << bname<I>() << endl;
00373   int omax=min(20,T_Healpix_Base<I>::order_max);
00374   for (int order=0; order<=omax; ++order)
00375     {
00376     T_Healpix_Base<I> rbase (order,RING), nbase (order,NEST);
00377     rangeset<I> pixset;
00378     int niter=max(1,min(1000,100000>>order));
00379     for (int m=0; m<niter; ++m)
00380       {
00381       vector<pointing> corner(3);
00382       random_dir(corner[0]); random_dir(corner[1]); random_dir(corner[2]);
00383       rbase.query_polygon(corner,pixset);
00384       I nval = pixset.nval();
00385       nbase.query_polygon(corner,pixset);
00386       if (nval!=pixset.nval())
00387         cout << "  PROBLEM: number of pixels different: "
00388              << nval << " vs. " << pixset.nval() << endl;
00389       rbase.query_polygon_inclusive(corner,pixset,4);
00390       I nv1=pixset.nval();
00391       nbase.query_polygon_inclusive(corner,pixset,4);
00392       I nv2=pixset.nval();
00393       if (nv1<nv2)
00394         cout << "  PROBLEM: inclusive(RING)<inclusive(NEST): "
00395              << nv1 << " vs. " << nv2 << endl;
00396       if (nv2<nval)
00397         cout << "  PROBLEM: inclusive(NEST)<non-inclusive: "
00398              << nv2 << " vs. " << nval << endl;
00399       }
00400     }
00401   }
00402 
00403 void helper_oop (int order)
00404   {
00405   Healpix_Map<double> map (order,RING), map2 (order,NEST), map3 (order,RING);
00406   for (int m=0; m<map.Npix(); ++m) map[m] = rng.rand_uni()+0.01;
00407   map2.Import(map);
00408   map3.Import(map2);
00409   for (int m=0; m<map.Npix(); ++m)
00410     if (!approx(map[m],map3[m],1e-12))
00411       cout << "PROBLEM: order = " << order << endl;
00412   }
00413 void helper_udgrade (int order, Healpix_Ordering_Scheme s1,
00414   Healpix_Ordering_Scheme s2)
00415   {
00416   Healpix_Map<double> map (order,s1), map2 (order+2,s2), map3 (order, s1);
00417   for (int m=0; m<map.Npix(); ++m) map[m] = rng.rand_uni()+0.01;
00418   map2.Import(map);
00419   map3.Import(map2);
00420   for (int m=0; m<map.Npix(); ++m)
00421     if (!approx(map[m],map3[m],1e-12))
00422       cout << "PROBLEM: order = " << order << endl;
00423   }
00424 void helper_udgrade2 (int nside)
00425   {
00426   Healpix_Map<double> map (nside,RING,SET_NSIDE), map2 (nside*3,RING,SET_NSIDE),
00427     map3 (nside, RING,SET_NSIDE);
00428   for (int m=0; m<map.Npix(); ++m) map[m] = rng.rand_uni()+0.01;
00429   map2.Import(map);
00430   map3.Import(map2);
00431   for (int m=0; m<map.Npix(); ++m)
00432     if (!approx(map[m],map3[m],1e-12))
00433       cout << "PROBLEM: nside = " << nside << endl;
00434   }
00435 
00436 void check_import()
00437   {
00438   cout << "testing out-of-place swapping" << endl;
00439   for (int order=0; order<=7; ++order)
00440     helper_oop(order);
00441   cout << "testing downgrade(upgrade(map)) == map" << endl;
00442   for (int order=0; order<=7; ++order)
00443     {
00444     helper_udgrade(order,RING,RING);
00445     helper_udgrade(order,RING,NEST);
00446     helper_udgrade(order,NEST,NEST);
00447     helper_udgrade(order,NEST,RING);
00448     }
00449   for (int nside=3; nside<500; nside+=nside/2+1)
00450     helper_udgrade2(nside);
00451   }
00452 
00453 void check_average()
00454   {
00455   cout << "testing whether average(map) == average(downgraded map)" << endl;
00456   for (int order=1; order<=10; ++order)
00457     {
00458     Healpix_Map<double> map (order,RING), map2(1,RING);
00459     for (int m=0; m<map.Npix(); ++m)
00460       map[m] = rng.rand_uni()+2;
00461     map2.Import(map);
00462     double avg=map.average(), avg2=map2.average();
00463     if (!approx(avg,avg2,1e-12))
00464       cout << "PROBLEM: order = " << order << " " << avg/avg2-1 << endl;
00465     }
00466   for (int nside=3; nside<1000; nside += nside/2+1)
00467     {
00468     Healpix_Map<double> map (nside,RING,SET_NSIDE), map2(1,RING,SET_NSIDE);
00469     for (int m=0; m<map.Npix(); ++m)
00470       map[m] = rng.rand_uni()+2;
00471     map2.Import(map);
00472     double avg=map.average(), avg2=map2.average();
00473     if (!approx(avg,avg2,1e-12))
00474       cout << "PROBLEM: nside = " << nside << " " << avg/avg2-1 << endl;
00475     }
00476   }
00477 
00478 void random_alm (Alm<xcomplex<double> >&almT, Alm<xcomplex<double> >&almG,
00479   Alm<xcomplex<double> >&almC, int lmax, int mmax)
00480   {
00481   almT.Set(lmax,mmax); almG.Set(lmax,mmax); almC.Set(lmax,mmax);
00482 
00483   for (int l=0; l<=lmax; ++l)
00484     {
00485     almT(l,0).re=rng.rand_gauss(); almT(l,0).im=0.;
00486     almG(l,0).re=rng.rand_gauss(); almG(l,0).im=0.;
00487     almC(l,0).re=rng.rand_gauss(); almC(l,0).im=0.;
00488     }
00489 
00490   for (int m=1; m<=mmax; ++m)
00491     for (int l=m; l<=lmax; ++l)
00492       {
00493       almT(l,m).re=rng.rand_gauss(); almT(l,m).im=rng.rand_gauss();
00494       almG(l,m).re=rng.rand_gauss(); almG(l,m).im=rng.rand_gauss();
00495       almC(l,m).re=rng.rand_gauss(); almC(l,m).im=rng.rand_gauss();
00496       }
00497   almG(0,0)=almC(0,0)=almG(1,0)=almC(1,0)=almG(1,1)=almC(1,1)=0;
00498   }
00499 
00500 void random_alm (Alm<xcomplex<double> >&alm, int lmax, int mmax)
00501   {
00502   alm.Set(lmax,mmax);
00503 
00504   for (int l=0; l<=lmax; ++l)
00505     { alm(l,0).re=rng.rand_gauss(); alm(l,0).im=0.; }
00506 
00507   for (int m=1; m<=mmax; ++m)
00508     for (int l=m; l<=lmax; ++l)
00509       { alm(l,m).re=rng.rand_gauss(); alm(l,m).im=rng.rand_gauss(); }
00510   }
00511 
00512 void check_alm (const Alm<xcomplex<double> >&oalm,
00513   const Alm<xcomplex<double> >&alm, double epsilon)
00514   {
00515   for (int m=0; m<=alm.Mmax(); ++m)
00516     for (int l=m; l<=alm.Lmax(); ++l)
00517       {
00518       if (!abs_approx(oalm(l,m).re,alm(l,m).re,epsilon))
00519         cout << "Problemr " << l << " " << m << endl;
00520       if (!abs_approx(oalm(l,m).im,alm(l,m).im,epsilon))
00521         cout << "Problemi " << l << " " << m <<  endl;
00522       }
00523   }
00524 
00525 void check_alm2map2alm (int lmax, int mmax, int nside)
00526   {
00527   cout << "testing whether a_lm->map->a_lm returns original a_lm" << endl;
00528   cout << "lmax=" << lmax <<", mmax=" << mmax << ", nside=" << nside << endl;
00529   const double epsilon = 1e-8;
00530   Alm<xcomplex<double> > oalmT(lmax,mmax),almT(lmax,mmax),
00531     oalmG(lmax,mmax),almG(lmax,mmax),oalmC(lmax,mmax),almC(lmax,mmax);
00532   Healpix_Map<double> mapT(nside,RING,SET_NSIDE), mapQ(nside,RING,SET_NSIDE),
00533     mapU(nside,RING,SET_NSIDE);
00534 
00535   random_alm(oalmT,oalmG,oalmC,lmax,mmax);
00536   alm2map(oalmT,mapT);
00537   map2alm_iter2(mapT,almT,1e-12,1e-12);
00538   check_alm (oalmT, almT, epsilon);
00539 
00540   alm2map_spin(oalmG,oalmC,mapQ,mapU,1);
00541   map2alm_spin_iter2(mapQ,mapU,almG,almC,1,1e-12,1e-12);
00542   check_alm (oalmG, almG, epsilon);
00543   check_alm (oalmC, almC, epsilon);
00544 
00545   alm2map_pol(oalmT,oalmG,oalmC,mapT,mapQ,mapU);
00546   map2alm_pol_iter2(mapT,mapQ,mapU,almT,almG,almC,1e-12,1e-12);
00547   check_alm (oalmT, almT, epsilon);
00548   check_alm (oalmG, almG, epsilon);
00549   check_alm (oalmC, almC, epsilon);
00550   }
00551 
00552 void check_smooth_alm ()
00553   {
00554   cout << "testing whether unsmooth(smooth(a_lm)) returns a_lm" << endl;
00555   const double epsilon = 1e-14;
00556   const double fwhm = 100.*arcmin2rad;
00557   const int lmax=300, mmax=300;
00558   Alm<xcomplex<double> > oalmT(lmax,mmax),almT(lmax,mmax),
00559     oalmG(lmax,mmax),almG(lmax,mmax),oalmC(lmax,mmax),almC(lmax,mmax);
00560 
00561   random_alm(oalmT,oalmG,oalmC,lmax,mmax);
00562 
00563   almT=oalmT; almG=oalmG; almC=oalmC;
00564   smoothWithGauss (almT, fwhm);
00565   smoothWithGauss (almT, -fwhm);
00566   check_alm (oalmT, almT, epsilon);
00567   almT=oalmT;
00568   smoothWithGauss (almT, almG, almC, fwhm);
00569   smoothWithGauss (almT, almG, almC, -fwhm);
00570   check_alm (oalmT, almT, epsilon);
00571   check_alm (oalmG, almG, epsilon);
00572   check_alm (oalmC, almC, epsilon);
00573   }
00574 
00575 void check_rot_alm ()
00576   {
00577   cout << "testing whether rot^-1(rot(a_lm)) returns a_lm" << endl;
00578   const double epsilon = 2e-13;
00579   const int lmax=300;
00580   Alm<xcomplex<double> > oalm(lmax,lmax),alm(lmax,lmax);
00581 
00582   random_alm(oalm,lmax,lmax);
00583 
00584   alm=oalm;
00585   rotate_alm (alm,3,4,5);
00586   rotate_alm (alm,-5,-4,-3);
00587   check_alm (oalm, alm, epsilon);
00588   }
00589 
00590 void check_isqrt()
00591   {
00592   cout << "testing whether isqrt() works reliably" << endl;
00593   uint64 val=uint64(0xF234)<<16, valsq=val*val;
00594   if (isqrt(valsq)!=val) cout << "PROBLEM1" << endl;
00595   if (isqrt(valsq-1)!=val-1) cout << "PROBLEM2" << endl;
00596   }
00597 
00598 void check_pix2ang_acc()
00599   {
00600   cout << "testing accuracy of pix2ang at the poles" << endl;
00601   for (int m=0; m<=29;++m)
00602     {
00603     Healpix_Base2 base(m,RING);
00604     if (base.pix2ang(1).theta==0.)
00605       cout << "PROBLEM: order " << m << endl;
00606     }
00607   }
00608 
00609 const int nsteps=1000000;
00610 
00611 template<typename I>void perf_neighbors(const string &name,
00612   Healpix_Ordering_Scheme scheme)
00613   {
00614   tsize cnt=0;
00615   wallTimers.start(name);
00616   int omax=T_Healpix_Base<I>::order_max;
00617   for (int order=0; order<=omax; ++order)
00618     {
00619     T_Healpix_Base<I> base (order,scheme);
00620     I dpix=max(base.Npix()/nsteps,I(1));
00621     fix_arr<I,8> nres;
00622     for (I pix=0; pix<base.Npix(); pix+=dpix)
00623       {
00624       base.neighbors(pix,nres);
00625       ++cnt;
00626       }
00627     }
00628   wallTimers.stop(name);
00629   cout << name << ": " << cnt/wallTimers.acc(name)*1e-6 << "MOps/s" << endl;
00630   }
00631 template<typename I>void perf_pix2ang(const string &name,
00632   Healpix_Ordering_Scheme scheme, double &dummy)
00633   {
00634   tsize cnt=0;
00635   wallTimers.start(name);
00636   int omax=T_Healpix_Base<I>::order_max;
00637   for (int order=0; order<=omax; ++order)
00638     {
00639     T_Healpix_Base<I> base (order,scheme);
00640     I dpix=max(base.Npix()/nsteps,I(1));
00641     for (I pix=0; pix<base.Npix(); pix+=dpix)
00642       {
00643       pointing p(base.pix2ang(pix));
00644       dummy+=p.theta+p.phi;
00645       ++cnt;
00646       }
00647     }
00648   wallTimers.stop(name);
00649   cout << name << ": " << cnt/wallTimers.acc(name)*1e-6 << "MOps/s" << endl;
00650   }
00651 template<typename I>void perf_pix2vec(const string &name,
00652   Healpix_Ordering_Scheme scheme, double &dummy)
00653   {
00654   tsize cnt=0;
00655   wallTimers.start(name);
00656   int omax=T_Healpix_Base<I>::order_max;
00657   for (int order=0; order<=omax; ++order)
00658     {
00659     T_Healpix_Base<I> base (order,scheme);
00660     I dpix=max(base.Npix()/nsteps,I(1));
00661     for (I pix=0; pix<base.Npix(); pix+=dpix)
00662       {
00663       vec3 v(base.pix2vec(pix));
00664       dummy+=v.x+v.y+v.z;
00665       ++cnt;
00666       }
00667     }
00668   wallTimers.stop(name);
00669   cout << name << ": " << cnt/wallTimers.acc(name)*1e-6 << "MOps/s" << endl;
00670   }
00671 template<typename I>void perf_pix2zphi(const string &name,
00672   Healpix_Ordering_Scheme scheme, double &dummy)
00673   {
00674   tsize cnt=0;
00675   wallTimers.start(name);
00676   int omax=T_Healpix_Base<I>::order_max;
00677   for (int order=0; order<=omax; ++order)
00678     {
00679     T_Healpix_Base<I> base (order,scheme);
00680     I dpix=max(base.Npix()/nsteps,I(1));
00681     double z,phi;
00682     for (I pix=0; pix<base.Npix(); pix+=dpix)
00683       {
00684       base.pix2zphi(pix,z,phi);
00685       dummy+=z+phi;
00686       ++cnt;
00687       }
00688     }
00689   wallTimers.stop(name);
00690   cout << name << ": " << cnt/wallTimers.acc(name)*1e-6 << "MOps/s" << endl;
00691   }
00692 template<typename I>void perf_zphi2pix(const string &name,
00693   Healpix_Ordering_Scheme scheme, double &dummy)
00694   {
00695   tsize cnt=0;
00696   wallTimers.start(name);
00697   int omax=T_Healpix_Base<I>::order_max;
00698   double dz=2./sqrt(nsteps);
00699   double dph=twopi/sqrt(nsteps);
00700   for (int order=0; order<=omax; ++order)
00701     {
00702     T_Healpix_Base<I> base (order,scheme);
00703     for (double z=-1; z<1; z+=dz)
00704       for (double phi=0; phi<twopi; phi+=dph)
00705         {
00706         dummy+=base.zphi2pix(z,phi);
00707         ++cnt;
00708         }
00709     }
00710   wallTimers.stop(name);
00711   cout << name << ": " << cnt/wallTimers.acc(name)*1e-6 << "MOps/s" << endl;
00712   }
00713 template<typename I>void perf_ang2pix(const string &name,
00714   Healpix_Ordering_Scheme scheme, double &dummy)
00715   {
00716   tsize cnt=0;
00717   wallTimers.start(name);
00718   int omax=T_Healpix_Base<I>::order_max;
00719   double dth=pi/sqrt(nsteps);
00720   double dph=twopi/sqrt(nsteps);
00721   for (int order=0; order<=omax; ++order)
00722     {
00723     T_Healpix_Base<I> base (order,scheme);
00724     for (double theta=0; theta<pi; theta+=dth)
00725       for (double phi=0; phi<twopi; phi+=dph)
00726         {
00727         dummy+=base.ang2pix(pointing(theta+1e-15*phi,phi));
00728         ++cnt;
00729         }
00730     }
00731   wallTimers.stop(name);
00732   cout << name << ": " << cnt/wallTimers.acc(name)*1e-6 << "MOps/s" << endl;
00733   }
00734 template<typename I>void perf_ring2nest(const string &name,double &dummy)
00735   {
00736   tsize cnt=0;
00737   wallTimers.start(name);
00738   int omax=T_Healpix_Base<I>::order_max;
00739   for (int order=0; order<=omax; ++order)
00740     {
00741     T_Healpix_Base<I> base (order,RING);
00742     I dpix=max(base.Npix()/nsteps,I(1));
00743     for (I pix=0; pix<base.Npix(); pix+=dpix)
00744       {
00745       dummy+=base.ring2nest(pix);
00746       ++cnt;
00747       }
00748     }
00749   wallTimers.stop(name);
00750   cout << name << ": " << cnt/wallTimers.acc(name)*1e-6 << "MOps/s" << endl;
00751   }
00752 template<typename I>void perf_nest2ring(const string &name,double &dummy)
00753   {
00754   tsize cnt=0;
00755   wallTimers.start(name);
00756   int omax=T_Healpix_Base<I>::order_max;
00757   for (int order=0; order<=omax; ++order)
00758     {
00759     T_Healpix_Base<I> base (order,RING);
00760     I dpix=max(base.Npix()/nsteps,I(1));
00761     for (I pix=0; pix<base.Npix(); pix+=dpix)
00762       {
00763       dummy+=base.nest2ring(pix);
00764       ++cnt;
00765       }
00766     }
00767   wallTimers.stop(name);
00768   cout << name << ": " << cnt/wallTimers.acc(name)*1e-6 << "MOps/s" << endl;
00769   }
00770 template<typename I>void perf_peano2nest(const string &name,double &dummy)
00771   {
00772   tsize cnt=0;
00773   wallTimers.start(name);
00774   int omax=T_Healpix_Base<I>::order_max;
00775   for (int order=0; order<=omax; ++order)
00776     {
00777     T_Healpix_Base<I> base (order,NEST);
00778     I dpix=max(base.Npix()/nsteps,I(1));
00779     for (I pix=0; pix<base.Npix(); pix+=dpix)
00780       {
00781       dummy+=base.peano2nest(pix);
00782       ++cnt;
00783       }
00784     }
00785   wallTimers.stop(name);
00786   cout << name << ": " << cnt/wallTimers.acc(name)*1e-6 << "MOps/s" << endl;
00787   }
00788 template<typename I>void perf_nest2peano(const string &name,double &dummy)
00789   {
00790   tsize cnt=0;
00791   wallTimers.start(name);
00792   int omax=T_Healpix_Base<I>::order_max;
00793   for (int order=0; order<=omax; ++order)
00794     {
00795     T_Healpix_Base<I> base (order,NEST);
00796     I dpix=max(base.Npix()/nsteps,I(1));
00797     for (I pix=0; pix<base.Npix(); pix+=dpix)
00798       {
00799       dummy+=base.nest2peano(pix);
00800       ++cnt;
00801       }
00802     }
00803   wallTimers.stop(name);
00804   cout << name << ": " << cnt/wallTimers.acc(name)*1e-6 << "MOps/s" << endl;
00805   }
00806 template<typename I>void perf_query_disc(const string &name,
00807   Healpix_Ordering_Scheme scheme, double &dummy)
00808   {
00809   tsize cnt=0;
00810   T_Healpix_Base<I> base(1024,scheme,SET_NSIDE);
00811   wallTimers.start(name);
00812   for (int m=0; m<1000; ++m)
00813     {
00814     rangeset<I> pix;
00815     base.query_disc(vec3(1,0,0),halfpi/9,pix);
00816     dummy+=pix.size();
00817     ++cnt;
00818     }
00819   wallTimers.stop(name);
00820   cout << name << ": " << cnt/wallTimers.acc(name)*1e-3 << "kOps/s" << endl;
00821   }
00822 template<typename I>void perf_query_triangle(const string &name,
00823   Healpix_Ordering_Scheme scheme, double &dummy)
00824   {
00825   tsize cnt=0;
00826   T_Healpix_Base<I> base(1024,scheme,SET_NSIDE);
00827   vector<pointing> corner;
00828   corner.push_back(vec3(1,0.01,0.01));
00829   corner.push_back(vec3(0.01,1,0.01));
00830   corner.push_back(vec3(0.01,0.01,1));
00831   wallTimers.start(name);
00832   for (int m=0; m<1000; ++m)
00833     {
00834     rangeset<I> pix;
00835     base.query_polygon(corner,pix);
00836     dummy+=pix.size();
00837     ++cnt;
00838     }
00839   wallTimers.stop(name);
00840   cout << name << ": " << cnt/wallTimers.acc(name)*1e-3 << "kOps/s" << endl;
00841   }
00842 template<typename I>void perf_query_polygon(const string &name,
00843   Healpix_Ordering_Scheme scheme, double &dummy)
00844   {
00845   tsize cnt=0;
00846   T_Healpix_Base<I> base(1024,scheme,SET_NSIDE);
00847   vector<pointing> corner;
00848   corner.push_back(vec3(1,0.01,0.01));
00849   corner.push_back(vec3(1,1,-0.3));
00850   corner.push_back(vec3(0.01,1,0.01));
00851   corner.push_back(vec3(0.01,0.01,1));
00852   wallTimers.start(name);
00853   for (int m=0; m<1000; ++m)
00854     {
00855     rangeset<I> pix;
00856     base.query_polygon(corner,pix);
00857     dummy+=pix.size();
00858     ++cnt;
00859     }
00860   wallTimers.stop(name);
00861   cout << name << ": " << cnt/wallTimers.acc(name)*1e-3 << "kOps/s" << endl;
00862   }
00863 
00864 void perftest()
00865   {
00866   double dummy=0;
00867   cout << "Measuring performance of Healpix_Base methods." << endl;
00868   perf_pix2zphi<int>   ("pix2zphi (RING):int  ",RING,dummy);
00869   perf_pix2zphi<int>   ("pix2zphi (NEST):int  ",NEST,dummy);
00870   perf_pix2zphi<int64> ("pix2zphi (RING):int64",RING,dummy);
00871   perf_pix2zphi<int64> ("pix2zphi (NEST):int64",NEST,dummy);
00872   perf_zphi2pix<int>   ("zphi2pix (RING):int  ",RING,dummy);
00873   perf_zphi2pix<int>   ("zphi2pix (NEST):int  ",NEST,dummy);
00874   perf_zphi2pix<int64> ("zphi2pix (RING):int64",RING,dummy);
00875   perf_zphi2pix<int64> ("zphi2pix (NEST):int64",NEST,dummy);
00876   perf_pix2ang<int>    ("pix2ang  (RING):int  ",RING,dummy);
00877   perf_pix2ang<int>    ("pix2ang  (NEST):int  ",NEST,dummy);
00878   perf_pix2ang<int64>  ("pix2ang  (RING):int64",RING,dummy);
00879   perf_pix2ang<int64>  ("pix2ang  (NEST):int64",NEST,dummy);
00880   perf_ang2pix<int>    ("ang2pix  (RING):int  ",RING,dummy);
00881   perf_ang2pix<int>    ("ang2pix  (NEST):int  ",NEST,dummy);
00882   perf_ang2pix<int64>  ("ang2pix  (RING):int64",RING,dummy);
00883   perf_ang2pix<int64>  ("ang2pix  (NEST):int64",NEST,dummy);
00884   perf_pix2vec<int>    ("pix2vec  (RING):int  ",RING,dummy);
00885   perf_pix2vec<int>    ("pix2vec  (NEST):int  ",NEST,dummy);
00886   perf_pix2vec<int64>  ("pix2vec  (RING):int64",RING,dummy);
00887   perf_pix2vec<int64>  ("pix2vec  (NEST):int64",NEST,dummy);
00888   perf_neighbors<int>  ("neighbors(NEST):int  ",NEST);
00889   perf_neighbors<int>  ("neighbors(RING):int  ",RING);
00890   perf_neighbors<int64>("neighbors(NEST):int64",NEST);
00891   perf_neighbors<int64>("neighbors(RING):int64",RING);
00892   perf_ring2nest<int>  ("ring2nest      :int  ",dummy);
00893   perf_ring2nest<int64>("ring2nest      :int64",dummy);
00894   perf_nest2ring<int>  ("nest2ring      :int  ",dummy);
00895   perf_nest2ring<int64>("nest2ring      :int64",dummy);
00896   perf_peano2nest<int>  ("peano2nest     :int  ",dummy);
00897   perf_peano2nest<int64>("peano2nest     :int64",dummy);
00898   perf_nest2peano<int>  ("nest2peano     :int  ",dummy);
00899   perf_nest2peano<int64>("nest2peano     :int64",dummy);
00900   perf_query_disc<int>      ("query_disc    (RING):int  ",RING,dummy);
00901   perf_query_disc<int>      ("query_disc    (NEST):int  ",NEST,dummy);
00902   perf_query_disc<int64>    ("query_disc    (RING):int64",RING,dummy);
00903   perf_query_disc<int64>    ("query_disc    (NEST):int64",NEST,dummy);
00904   perf_query_triangle<int>  ("query_triangle(RING):int  ",RING,dummy);
00905   perf_query_triangle<int>  ("query_triangle(NEST):int  ",NEST,dummy);
00906   perf_query_triangle<int64>("query_triangle(RING):int64",RING,dummy);
00907   perf_query_triangle<int64>("query_triangle(NEST):int64",NEST,dummy);
00908   perf_query_polygon<int>   ("query_polygon (RING):int  ",RING,dummy);
00909   perf_query_polygon<int>   ("query_polygon (NEST):int  ",NEST,dummy);
00910   perf_query_polygon<int64> ("query_polygon (RING):int64",RING,dummy);
00911   perf_query_polygon<int64> ("query_polygon (NEST):int64",NEST,dummy);
00912 
00913   if (dummy<0) cout << dummy << endl;
00914   }
00915 } // unnamed namespace
00916 
00917 int main(int argc, const char **argv)
00918   {
00919   module_startup ("hpxtest",argc,argv,1,"");
00920   perftest();
00921   check_isqrt();
00922   check_pix2ang_acc();
00923   check_smooth_alm();
00924   check_rot_alm();
00925   check_alm2map2alm(620,620,256);
00926   check_alm2map2alm(620,2,256);
00927   check_average();
00928   check_import();
00929   check_ringnestring<int>();
00930   check_ringnestring<int64>();
00931   check_nestpeanonest<int>();
00932   check_nestpeanonest<int64>();
00933   check_pixzphipix<int>();
00934   check_pixzphipix<int64>();
00935   check_zphipixzphi<int>();
00936   check_zphipixzphi<int64>();
00937   check_pixangpix<int>();
00938   check_pixangpix<int64>();
00939   check_neighbors<int>();
00940   check_neighbors<int64>();
00941   check_swap_scheme();
00942   check_query_disc_strict(RING);
00943   check_query_disc_strict(NEST);
00944   check_query_disc<int>();
00945   check_query_disc<int64>();
00946   check_query_polygon<int>();
00947   check_query_polygon<int64>();
00948   }

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