Main Page | Modules | Class Hierarchy | Alphabetical List | Class List | File List | Class Members | File Members | Related Pages

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.jpl.nasa.gov
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, 2005 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 "constants.h"
00056 #include "alm.h"
00057 #include "alm_map_tools.h"
00058 #include "alm_powspec_tools.h"
00059 
00060 using namespace std;
00061 
00062 planck_rng rng;
00063 
00064 void random_dir (pointing &ptg)
00065   {
00066   ptg.theta = acos(rng.rand_uni()*2-1);
00067   ptg.phi = rng.rand_uni()*twopi;
00068   }
00069 
00070 void check_ringnestring()
00071   {
00072   cout << "testing ring2nest(nest2ring(m))==m" << endl;
00073   for (int order=0; order<=13; ++order)
00074     {
00075     cout << "order = " << order << endl;
00076     Healpix_Base base (order,RING);
00077     for (int m=0; m<base.Npix(); ++m)
00078       {
00079       if (base.ring2nest(base.nest2ring(m))!=m)
00080         cout << "  PROBLEM: order = " << order << ", pixel = " << m << endl;
00081       }
00082     }
00083   }
00084 
00085 void check_pixangpix()
00086   {
00087   cout << "testing ang2pix(pix2ang(m))==m" << endl;
00088   for (int order=0; order<=13; ++order)
00089     {
00090     cout << "order = " << order << endl;
00091     Healpix_Base base1 (order,RING);
00092     Healpix_Base base2 (order,NEST);
00093     for (int m=0; m<base1.Npix(); ++m)
00094       {
00095       if (base1.ang2pix(base1.pix2ang(m))!=m)
00096         cout << "  PROBLEM: order = " << order << ", pixel = " << m << endl;
00097       if (base2.ang2pix(base2.pix2ang(m))!=m)
00098         cout << "  PROBLEM: order = " << order << ", pixel = " << m << endl;
00099       }
00100     }
00101   for (int nside=3; nside<1000; nside+=nside/2+1)
00102     {
00103     cout << "nside = " << nside << endl;
00104     Healpix_Base base (nside,RING,SET_NSIDE);
00105     for (int m=0; m<base.Npix(); ++m)
00106       {
00107       if (base.ang2pix(base.pix2ang(m))!=m)
00108         cout << "  PROBLEM: nside = " << nside << ", pixel = " << m << endl;
00109       }
00110     }
00111   }
00112 
00113 void check_angpixang()
00114   {
00115   cout << "testing pix2ang(ang2pix(ptg)) approx ptg" << endl;
00116   for (int nside=3; nside<1000; nside+=nside/2+1)
00117     {
00118     cout << "nside = " << nside << endl;
00119     Healpix_Base base (nside,RING,SET_NSIDE);
00120     double cosmaxang = cos (1.36*twopi/(8*nside));
00121     for (int m=0; m<100000; ++m)
00122       {
00123       pointing ptg;
00124       random_dir (ptg);
00125       if (dotprod(base.pix2ang(base.ang2pix(ptg)),ptg)<cosmaxang)
00126         cout << "  PROBLEM: nside = " << nside << ", ptg = " << ptg << endl;
00127       }
00128     }
00129   }
00130 
00131 void check_neighbors()
00132   {
00133   cout << "testing (pix2ang(neighbors(pix))) approx pix2ang(pix)" << endl;
00134   for (int order=0; order<=13; ++order)
00135     {
00136     cout << "order = " << order << endl;
00137     Healpix_Base base (order,NEST), base2(order,RING);
00138     double cosmaxang = cos (2.75*twopi/(8*base.Nside()));
00139     for (int m=0; m<100000; ++m)
00140       {
00141       int pix = rng.int_rand_uni()%base.Npix();
00142       fix_arr<int,8> nb;
00143       vec3 pixpt = base.pix2ang(pix);
00144       base.neighbors(pix,nb);
00145       int check=0;
00146       for (int n=0; n<8; ++n)
00147         {
00148         if (nb[n]<0) ++check;
00149         if ((nb[n]>=0) && (dotprod(base.pix2ang(nb[n]),pixpt)<cosmaxang))
00150           cout << "  PROBLEM: order = " << order << ", pix = " << pix << endl;
00151         }
00152       planck_assert((check<=1)||((order==0)&&(check<=2)),"too few neighbors");
00153       pixpt = base2.pix2ang(pix);
00154       base2.neighbors(pix,nb);
00155       for (int n=0; n<8; ++n)
00156         if ((nb[n]>=0) && (dotprod(base2.pix2ang(nb[n]),pixpt)<cosmaxang))
00157           cout << "  PROBLEM2: order = " << order << ", pix = " << pix << endl;
00158       }
00159     }
00160   for (int nside=3; nside<1000; nside+=nside/2+1)
00161     {
00162     cout << "nside = " << nside << endl;
00163     Healpix_Base base (nside,RING,SET_NSIDE);
00164     double cosmaxang = cos (2.75*twopi/(8*base.Nside()));
00165     for (int m=0; m<100000; ++m)
00166       {
00167       int pix = rng.int_rand_uni()%base.Npix();
00168       fix_arr<int,8> nb;
00169       vec3 pixpt = base.pix2ang(pix);
00170       base.neighbors(pix,nb);
00171       for (int n=0; n<8; ++n)
00172         if ((nb[n]>=0) && (dotprod(base.pix2ang(nb[n]),pixpt)<cosmaxang))
00173           cout << "  PROBLEM: nside = " << nside << ", pix = " << pix << endl;
00174       }
00175     }
00176   }
00177 
00178 void check_swap_scheme()
00179   {
00180   cout << "testing whether double swap_scheme() returns the original map"
00181        << endl << "(for orders 0 to 10)." << endl;
00182   for (int order=0; order<=10; ++order)
00183     {
00184     cout << "order = " << order << endl;
00185     Healpix_Map<unsigned char> map(order,NEST);
00186     for (int m=0; m<map.Npix(); ++m) map[m]=m&0xFF;
00187     map.swap_scheme();
00188     map.swap_scheme();
00189     for (int m=0; m<map.Npix(); ++m)
00190       if (map[m]!=(m&0xFF))
00191         cout << "  PROBLEM: order = " << order << ", pix = " << m << endl;
00192     }
00193   }
00194 
00195 void check_query_disc()
00196   {
00197   cout << "testing whether all pixels found by query_disc() really" << endl
00198        << "lie inside the disk (and vice versa)" << endl;
00199   for (int order=0; order<=5; ++order)
00200     {
00201     cout << "order = " << order << endl;
00202     Healpix_Map <bool> map (order,RING);
00203     map.fill(false);
00204     vector<int> list;
00205     for (int m=0; m<100000; ++m)
00206       {
00207       pointing ptg;
00208       random_dir (ptg);
00209       double rad = pi/1 * rng.rand_uni();
00210       map.query_disc(ptg,rad,list);
00211       vec3 vptg=ptg;
00212       double cosrad=cos(rad);
00213       for (unsigned int i=0; i<list.size(); ++i)
00214         map[list[i]] = true;
00215       for (int i=0; i<map.Npix(); ++i)
00216         {
00217         bool inside = dotprod(map.pix2ang(i),vptg)>cosrad;
00218         if (inside^map[i])
00219           cout << "  PROBLEM: order = " << order << ", ptg = " << ptg << endl;
00220         }
00221       for (unsigned int i=0; i<list.size(); ++i)
00222         map[list[i]] = false;
00223       }
00224     }
00225   }
00226 
00227 void helper_oop (int order)
00228   {
00229   Healpix_Map<double> map (order,RING), map2 (order,NEST), map3 (order,RING);
00230   for (int m=0; m<map.Npix(); ++m) map[m] = rng.rand_uni()+0.01;
00231   map2.Import(map);
00232   map3.Import(map2);
00233   for (int m=0; m<map.Npix(); ++m)
00234     if (!approx(map[m],map3[m],1e-12))
00235       cout << "PROBLEM: order = " << order << endl;
00236   }
00237 void helper_udgrade (int order, Healpix_Ordering_Scheme s1,
00238   Healpix_Ordering_Scheme s2)
00239   {
00240   Healpix_Map<double> map (order,s1), map2 (order+2,s2), map3 (order, s1);
00241   for (int m=0; m<map.Npix(); ++m) map[m] = rng.rand_uni()+0.01;
00242   map2.Import(map);
00243   map3.Import(map2);
00244   for (int m=0; m<map.Npix(); ++m)
00245     if (!approx(map[m],map3[m],1e-12))
00246       cout << "PROBLEM: order = " << order << endl;
00247   }
00248 void helper_udgrade2 (int nside)
00249   {
00250   Healpix_Map<double> map (nside,RING,SET_NSIDE), map2 (nside*3,RING,SET_NSIDE),
00251     map3 (nside, RING,SET_NSIDE);
00252   for (int m=0; m<map.Npix(); ++m) map[m] = rng.rand_uni()+0.01;
00253   map2.Import(map);
00254   map3.Import(map2);
00255   for (int m=0; m<map.Npix(); ++m)
00256     if (!approx(map[m],map3[m],1e-12))
00257       cout << "PROBLEM: nside = " << nside << endl;
00258   }
00259 
00260 void check_import()
00261   {
00262   cout << "testing out-of-place swapping" << endl;
00263   for (int order=0; order<=7; ++order)
00264     {
00265     cout << "order = " << order << endl;
00266     helper_oop(order);
00267     }
00268   cout << "testing downgrade(upgrade(map)) == map" << endl;
00269   for (int order=0; order<=7; ++order)
00270     {
00271     cout << "order = " << order << endl;
00272     cout << "RING, RING" << endl;
00273     helper_udgrade(order,RING,RING);
00274     cout << "RING, NEST" << endl;
00275     helper_udgrade(order,RING,NEST);
00276     cout << "NEST, NEST" << endl;
00277     helper_udgrade(order,NEST,NEST);
00278     cout << "NEST, RING" << endl;
00279     helper_udgrade(order,NEST,RING);
00280     }
00281   for (int nside=3; nside<500; nside+=nside/2+1)
00282     {
00283     cout << "nside = " << nside << endl;
00284     helper_udgrade2(nside);
00285     }
00286   }
00287 
00288 void check_average()
00289   {
00290   cout << "testing whether average(map) == average(downgraded map)" << endl;
00291   for (int order=1; order<=10; ++order)
00292     {
00293     cout << "order = " << order << endl;
00294     Healpix_Map<double> map (order,RING), map2(1,RING);
00295     double avg=0, avg2=0;
00296     for (int m=0; m<map.Npix(); ++m)
00297       {
00298       map[m] = rng.rand_uni()+0.01;
00299       avg+=map[m];
00300       }
00301     avg /= map.Npix();
00302     map2.Import(map);
00303     for (int m=0; m<map2.Npix(); ++m)
00304       avg2+=map2[m];
00305     avg2 /= map2.Npix();
00306 
00307     if (!approx(avg,avg2,1e-10))
00308       cout << "PROBLEM: order = " << order << " " << avg/avg2-1 << endl;
00309     }
00310   for (int nside=3; nside<1000; nside += nside/2+1)
00311     {
00312     cout << "nside = " << nside << endl;
00313     Healpix_Map<double> map (nside,RING,SET_NSIDE), map2(1,RING,SET_NSIDE);
00314     double avg=0, avg2=0;
00315     for (int m=0; m<map.Npix(); ++m)
00316       {
00317       map[m] = rng.rand_uni()+0.01;
00318       avg+=map[m];
00319       }
00320     avg /= map.Npix();
00321     map2.Import(map);
00322     for (int m=0; m<map2.Npix(); ++m)
00323       avg2+=map2[m];
00324     avg2 /= map2.Npix();
00325 
00326     if (!approx(avg,avg2,1e-10))
00327       cout << "PROBLEM: nside = " << nside << " " << avg/avg2-1 << endl;
00328     }
00329   }
00330 
00331 void random_alm (Alm<xcomplex<double> >&almT, Alm<xcomplex<double> >&almG,
00332   Alm<xcomplex<double> >&almC, int lmax, int mmax)
00333   {
00334   almT.Set(lmax,mmax); almG.Set(lmax,mmax); almC.Set(lmax,mmax);
00335 
00336   for (int l=0; l<=lmax; ++l)
00337     {
00338     almT(l,0).re=rng.rand_gauss(); almT(l,0).im=0.;
00339     almG(l,0).re=rng.rand_gauss(); almG(l,0).im=0.;
00340     almC(l,0).re=rng.rand_gauss(); almC(l,0).im=0.;
00341     }
00342 
00343   for (int m=1; m<=mmax; ++m)
00344     for (int l=m; l<=lmax; ++l)
00345       {
00346       almT(l,m).re=rng.rand_gauss(); almT(l,m).im=rng.rand_gauss();
00347       almG(l,m).re=rng.rand_gauss(); almG(l,m).im=rng.rand_gauss();
00348       almC(l,m).re=rng.rand_gauss(); almC(l,m).im=rng.rand_gauss();
00349       }
00350   almG(0,0)=almC(0,0)=almG(1,0)=almC(1,0)=almG(1,1)=almC(1,1)=0;
00351   }
00352 
00353 void random_alm (Alm<xcomplex<double> >&alm, int lmax, int mmax)
00354   {
00355   alm.Set(lmax,mmax);
00356 
00357   for (int l=0; l<=lmax; ++l)
00358     { alm(l,0).re=rng.rand_gauss(); alm(l,0).im=0.; }
00359 
00360   for (int m=1; m<=mmax; ++m)
00361     for (int l=m; l<=lmax; ++l)
00362       { alm(l,m).re=rng.rand_gauss(); alm(l,m).im=rng.rand_gauss(); }
00363   }
00364 
00365 void check_alm (const Alm<xcomplex<double> >&oalm,
00366   const Alm<xcomplex<double> >&alm, double epsilon)
00367   {
00368   for (int m=0; m<=alm.Mmax(); ++m)
00369     for (int l=m; l<=alm.Lmax(); ++l)
00370       {
00371       if (!abs_approx(oalm(l,m).re,alm(l,m).re,epsilon))
00372         cout << "Problemr " << l << " " << m << endl;
00373       if (!abs_approx(oalm(l,m).im,alm(l,m).im,epsilon))
00374         cout << "Problemi " << l << " " << m <<  endl;
00375       }
00376   }
00377 
00378 void check_alm2map2alm (int lmax, int mmax, int nside)
00379   {
00380   cout << "testing whether a_lm->map->a_lm returns original a_lm" << endl;
00381   cout << "lmax=" << lmax <<", mmax=" << mmax << ", nside=" << nside << endl;
00382   const double epsilon = 1e-8;
00383   Alm<xcomplex<double> > oalmT(lmax,mmax),almT(lmax,mmax),
00384     oalmG(lmax,mmax),almG(lmax,mmax),oalmC(lmax,mmax),almC(lmax,mmax);
00385   Healpix_Map<double> mapT(nside,RING,SET_NSIDE), mapQ(nside,RING,SET_NSIDE),
00386     mapU(nside,RING,SET_NSIDE);
00387 
00388   random_alm(oalmT,oalmG,oalmC,lmax,mmax);
00389   alm2map(oalmT,mapT);
00390   map2alm_iter2(mapT,almT,1e-12,1e-12);
00391   check_alm (oalmT, almT, epsilon);
00392 
00393   alm2map_pol(oalmT,oalmG,oalmC,mapT,mapQ,mapU);
00394   map2alm_pol_iter2(mapT,mapQ,mapU,almT,almG,almC,1e-12,1e-12);
00395   check_alm (oalmT, almT, epsilon);
00396   check_alm (oalmG, almG, epsilon);
00397   check_alm (oalmC, almC, epsilon);
00398   }
00399 
00400 void check_smooth_alm ()
00401   {
00402   cout << "testing whether unsmooth(smooth(a_lm)) returns a_lm" << endl;
00403   const double epsilon = 1e-14;
00404   const double fwhm_arcmin = 100;
00405   const int lmax=300, mmax=300;
00406   Alm<xcomplex<double> > oalmT(lmax,mmax),almT(lmax,mmax),
00407     oalmG(lmax,mmax),almG(lmax,mmax),oalmC(lmax,mmax),almC(lmax,mmax);
00408 
00409   random_alm(oalmT,oalmG,oalmC,lmax,mmax);
00410 
00411   almT=oalmT; almG=oalmG; almC=oalmC;
00412   smooth_with_Gauss (almT, fwhm_arcmin);
00413   smooth_with_Gauss (almT, -fwhm_arcmin);
00414   check_alm (oalmT, almT, epsilon);
00415   almT=oalmT;
00416   smooth_with_Gauss (almT, almG, almC, fwhm_arcmin);
00417   smooth_with_Gauss (almT, almG, almC, -fwhm_arcmin);
00418   check_alm (oalmT, almT, epsilon);
00419   check_alm (oalmG, almG, epsilon);
00420   check_alm (oalmC, almC, epsilon);
00421   }
00422 
00423 void check_rot_alm ()
00424   {
00425   cout << "testing whether rot^-1(rot(a_lm)) returns a_lm" << endl;
00426   const double epsilon = 2e-13;
00427   const int lmax=300;
00428   Alm<xcomplex<double> > oalm(lmax,lmax),alm(lmax,lmax);
00429 
00430   random_alm(oalm,lmax,lmax);
00431 
00432   alm=oalm;
00433   rotate_alm (alm,3,4,5);
00434   rotate_alm (alm,-5,-4,-3);
00435   check_alm (oalm, alm, epsilon);
00436   }
00437 
00438 int main()
00439   {
00440   check_smooth_alm();
00441   check_rot_alm();
00442   check_alm2map2alm(620,620,256);
00443   check_alm2map2alm(620,2,256);
00444   check_average();
00445   check_import();
00446   check_neighbors();
00447   check_pixangpix();
00448   check_angpixang();
00449   check_ringnestring();
00450   check_query_disc();
00451   check_swap_scheme();
00452   }

Generated on Fri Jul 8 09:37:15 2005 for Healpix C++