/* * This file is part of libsharp. * * libsharp is free software; you can redistribute it and/or modify * it under the terms of the GNU General Public License as published by * the Free Software Foundation; either version 2 of the License, or * (at your option) any later version. * * libsharp is distributed in the hope that it will be useful, * but WITHOUT ANY WARRANTY; without even the implied warranty of * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the * GNU General Public License for more details. * * You should have received a copy of the GNU General Public License * along with libsharp; if not, write to the Free Software * Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA */ /* * libsharp is being developed at the Max-Planck-Institut fuer Astrophysik * and financially supported by the Deutsches Zentrum fuer Luft- und Raumfahrt * (DLR). */ /* \file sharp_testsuite.c * * Copyright (C) 2012-2013 Max-Planck-Society * \author Martin Reinecke */ #include #include #ifdef USE_MPI #include "mpi.h" #include "sharp_mpi.h" #endif #ifdef _OPENMP #include #endif #include "sharp.h" #include "sharp_internal.h" #include "sharp_geomhelpers.h" #include "sharp_almhelpers.h" #include "c_utils.h" #include "sharp_announce.h" #include "memusage.h" #include "sharp_vecsupport.h" typedef complex double dcmplx; int ntasks, mytask; static double drand (double min, double max, int *state) { *state = (((*state) * 1103515245) + 12345) & 0x7fffffff; return min + (max-min)*(*state)/(0x7fffffff+1.0); } static void random_alm (dcmplx *alm, sharp_alm_info *helper, int spin, int cnt) { #pragma omp parallel { int mi; #pragma omp for schedule (dynamic,100) for (mi=0;minm; ++mi) { int m=helper->mval[mi]; int state=1234567*cnt+8912*m; // random seed for (int l=m;l<=helper->lmax; ++l) { if ((lnm; i+=ntasks,++nmnew) { ainfo->mval[nmnew]=ainfo->mval[i]; ainfo->mvstart[nmnew]=ofs-ainfo->mval[nmnew]; ofs+=ainfo->lmax-ainfo->mval[nmnew]+1; } ainfo->nm=nmnew; } static void reduce_geom_info(sharp_geom_info *ginfo) { int npairsnew=0; ptrdiff_t ofs = 0; for (int i=mytask; inpairs; i+=ntasks,++npairsnew) { ginfo->pair[npairsnew]=ginfo->pair[i]; ginfo->pair[npairsnew].r1.ofs=ofs; ofs+=ginfo->pair[npairsnew].r1.nph; ginfo->pair[npairsnew].r2.ofs=ofs; if (ginfo->pair[npairsnew].r2.nph>0) ofs+=ginfo->pair[npairsnew].r2.nph; } ginfo->npairs=npairsnew; } #endif static ptrdiff_t get_nalms(const sharp_alm_info *ainfo) { ptrdiff_t res=0; for (int i=0; inm; ++i) res += ainfo->lmax-ainfo->mval[i]+1; return res; } static ptrdiff_t get_npix(const sharp_geom_info *ginfo) { ptrdiff_t res=0; for (int i=0; inpairs; ++i) { res += ginfo->pair[i].r1.nph; if (ginfo->pair[i].r2.nph>0) res += ginfo->pair[i].r2.nph; } return res; } static double *get_sqsum_and_invert (dcmplx **alm, ptrdiff_t nalms, int ncomp) { double *sqsum=RALLOC(double,ncomp); for (int i=0; imaxdiff) maxdiff=sqr; } maxdiff=sqrt(maxdiff); #ifdef USE_MPI MPI_Allreduce(&sum,&sumtot,1,MPI_DOUBLE,MPI_SUM,MPI_COMM_WORLD); MPI_Allreduce(&sqsum[i],&sqsumtot,1,MPI_DOUBLE,MPI_SUM,MPI_COMM_WORLD); MPI_Allreduce(&maxdiff,&maxdifftot,1,MPI_DOUBLE,MPI_MAX,MPI_COMM_WORLD); #else sumtot=sum; sqsumtot=sqsum[i]; maxdifftot=maxdiff; #endif sumtot=sqrt(sumtot/nalms_tot); sqsumtot=sqrt(sqsumtot/nalms_tot); (*err_abs)[i]=maxdifftot; (*err_rel)[i]=sumtot/sqsumtot; } } static int good_fft_size(int n) { if (n<=6) return n; int bestfac=2*n; for (int f2=1; f2=n) bestfac=f235; return bestfac; } static void get_infos (const char *gname, int lmax, int *mmax, int *gpar1, int *gpar2, sharp_geom_info **ginfo, sharp_alm_info **ainfo) { UTIL_ASSERT(lmax>=0,"lmax must not be negative"); if (*mmax<0) *mmax=lmax; UTIL_ASSERT(*mmax<=lmax,"mmax larger than lmax"); if (mytask==0) printf ("lmax: %d, mmax: %d\n",lmax,*mmax); sharp_make_triangular_alm_info(lmax,*mmax,1,ainfo); #ifdef USE_MPI reduce_alm_info(*ainfo); #endif if (strcmp(gname,"healpix")==0) { if (*gpar1<1) *gpar1=lmax/2; if (*gpar1==0) ++(*gpar1); sharp_make_healpix_geom_info (*gpar1, 1, ginfo); if (mytask==0) printf ("HEALPix grid, nside=%d\n",*gpar1); } else if (strcmp(gname,"gauss")==0) { if (*gpar1<1) *gpar1=lmax+1; if (*gpar2<1) *gpar2=2*(*mmax)+1; sharp_make_gauss_geom_info (*gpar1, *gpar2, 0., 1, *gpar2, ginfo); if (mytask==0) printf ("Gauss-Legendre grid, nlat=%d, nlon=%d\n",*gpar1,*gpar2); } else if (strcmp(gname,"fejer1")==0) { if (*gpar1<1) *gpar1=2*lmax+1; if (*gpar2<1) *gpar2=2*(*mmax)+1; sharp_make_fejer1_geom_info (*gpar1, *gpar2, 0., 1, *gpar2, ginfo); if (mytask==0) printf ("Fejer1 grid, nlat=%d, nlon=%d\n",*gpar1,*gpar2); } else if (strcmp(gname,"fejer2")==0) { if (*gpar1<1) *gpar1=2*lmax+1; if (*gpar2<1) *gpar2=2*(*mmax)+1; sharp_make_fejer2_geom_info (*gpar1, *gpar2, 0., 1, *gpar2, ginfo); if (mytask==0) printf ("Fejer2 grid, nlat=%d, nlon=%d\n",*gpar1,*gpar2); } else if (strcmp(gname,"cc")==0) { if (*gpar1<1) *gpar1=2*lmax+1; if (*gpar2<1) *gpar2=2*(*mmax)+1; sharp_make_cc_geom_info (*gpar1, *gpar2, 0., 1, *gpar2, ginfo); if (mytask==0) printf("Clenshaw-Curtis grid, nlat=%d, nlon=%d\n",*gpar1,*gpar2); } else if (strcmp(gname,"smallgauss")==0) { int nlat=*gpar1, nlon=*gpar2; if (nlat<1) nlat=lmax+1; if (nlon<1) nlon=2*(*mmax)+1; *gpar1=nlat; *gpar2=nlon; sharp_make_gauss_geom_info (nlat, nlon, 0., 1, nlon, ginfo); ptrdiff_t npix_o=get_npix(*ginfo); size_t ofs=0; for (int i=0; i<(*ginfo)->npairs; ++i) { sharp_ringpair *pair=&((*ginfo)->pair[i]); int pring=1+2*sharp_get_mlim(lmax,0,pair->r1.sth,pair->r1.cth); if (pring>nlon) pring=nlon; pring=good_fft_size(pring); pair->r1.nph=pring; pair->r1.weight*=nlon*1./pring; pair->r1.ofs=ofs; ofs+=pring; if (pair->r2.nph>0) { pair->r2.nph=pring; pair->r2.weight*=nlon*1./pring; pair->r2.ofs=ofs; ofs+=pring; } } if (mytask==0) { ptrdiff_t npix=get_npix(*ginfo); printf("Small Gauss grid, nlat=%d, npix=%ld, savings=%.2f%%\n", nlat,(long)npix,(npix_o-npix)*100./npix_o); } } else UTIL_FAIL("unknown grid geometry"); #ifdef USE_MPI reduce_geom_info(*ginfo); #endif } static void check_sign_scale(void) { int lmax=50; int mmax=lmax; sharp_geom_info *tinfo; int nrings=lmax+1; int ppring=2*lmax+2; ptrdiff_t npix=(ptrdiff_t)nrings*ppring; sharp_make_gauss_geom_info (nrings, ppring, 0., 1, ppring, &tinfo); /* flip theta to emulate the "old" Gaussian grid geometry */ for (int i=0; inpairs; ++i) { const double pi=3.141592653589793238462643383279502884197; tinfo->pair[i].r1.cth=-tinfo->pair[i].r1.cth; tinfo->pair[i].r2.cth=-tinfo->pair[i].r2.cth; tinfo->pair[i].r1.theta=pi-tinfo->pair[i].r1.theta; tinfo->pair[i].r2.theta=pi-tinfo->pair[i].r2.theta; } sharp_alm_info *alms; sharp_make_triangular_alm_info(lmax,mmax,1,&alms); ptrdiff_t nalms = ((mmax+1)*(mmax+2))/2 + (mmax+1)*(lmax-mmax); for (int ntrans=1; ntrans<10; ++ntrans) { double **map; ALLOC2D(map,double,2*ntrans,npix); dcmplx **alm; ALLOC2D(alm,dcmplx,2*ntrans,nalms); for (int i=0; i<2*ntrans; ++i) for (int j=0; j=9,"usage: grid lmax mmax geom1 geom2 spin ntrans"); int lmax=atoi(argv[3]); int mmax=atoi(argv[4]); int gpar1=atoi(argv[5]); int gpar2=atoi(argv[6]); int spin=atoi(argv[7]); int ntrans=atoi(argv[8]); if (mytask==0) printf("Testing map analysis accuracy.\n"); if (mytask==0) printf("spin=%d, ntrans=%d\n", spin, ntrans); sharp_geom_info *ginfo; sharp_alm_info *ainfo; get_infos (argv[2], lmax, &mmax, &gpar1, &gpar2, &ginfo, &ainfo); int ncomp = ntrans*((spin==0) ? 1 : 2); double t_a2m=1e30, t_m2a=1e30; unsigned long long op_a2m, op_m2a; double *err_abs,*err_rel; double t_acc=0; int nrpt=0; while(1) { ++nrpt; double ta2m2, tm2a2; do_sht (ginfo, ainfo, spin, ntrans, 0, &err_abs, &err_rel, &ta2m2, &tm2a2, &op_a2m, &op_m2a); if (ta2m22.) { if (mytask==0) printf("Best of %d runs\n",nrpt); break; } DEALLOC(err_abs); DEALLOC(err_rel); } if (mytask==0) printf("wall time for alm2map: %fs\n",t_a2m); if (mytask==0) printf("Performance: %fGFLOPs/s\n",1e-9*op_a2m/t_a2m); if (mytask==0) printf("wall time for map2alm: %fs\n",t_m2a); if (mytask==0) printf("Performance: %fGFLOPs/s\n",1e-9*op_m2a/t_m2a); if (mytask==0) for (int i=0; i=9,"usage: grid lmax mmax geom1 geom2 spin ntrans"); int lmax=atoi(argv[3]); int mmax=atoi(argv[4]); int gpar1=atoi(argv[5]); int gpar2=atoi(argv[6]); int spin=atoi(argv[7]); int ntrans=atoi(argv[8]); if (mytask==0) printf("Testing map analysis accuracy.\n"); if (mytask==0) printf("spin=%d, ntrans=%d\n", spin, ntrans); sharp_geom_info *ginfo; sharp_alm_info *ainfo; get_infos (argv[2], lmax, &mmax, &gpar1, &gpar2, &ginfo, &ainfo); double ta2m_auto=1e30, tm2a_auto=1e30, ta2m_min=1e30, tm2a_min=1e30; unsigned long long opa2m_min=0, opm2a_min=0; int nvmin_a2m=-1, nvmin_m2a=-1; for (int nv=0; nv<=6; ++nv) { int ntries=0; double tacc=0; do { double t_a2m, t_m2a; unsigned long long op_a2m, op_m2a; double *err_abs,*err_rel; do_sht (ginfo, ainfo, spin, ntrans, nv, &err_abs, &err_rel, &t_a2m, &t_m2a, &op_a2m, &op_m2a); DEALLOC(err_abs); DEALLOC(err_rel); tacc+=t_a2m+t_m2a; ++ntries; if (nv==0) { if (t_a2m=2,"need at least one command line argument"); if (strcmp(argv[1],"acctest")==0) sharp_acctest(); else if (strcmp(argv[1],"test")==0) sharp_test(argc,argv); else if (strcmp(argv[1],"bench")==0) sharp_bench(argc,argv); else UTIL_FAIL("unknown command"); #ifdef USE_MPI MPI_Finalize(); #endif return 0; }