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

alm_map_tools.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) 2003, 2004, 2005 Max-Planck-Society
00029  *  Author: Martin Reinecke
00030  */
00031 
00032 #include "alm_map_tools.h"
00033 #include "alm.h"
00034 #include "healpix_map.h"
00035 #include "fftpack_support.h"
00036 #include "ylmgen.h"
00037 #include "xcomplex.h"
00038 
00039 using namespace std;
00040 
00041 namespace {
00042 
00043 void init_lam_fact_1d (int m, arr<double> &lam_fact)
00044   {
00045   for (int l=m; l<lam_fact.size(); ++l)
00046     lam_fact[l] = (l<2) ? 0. : 2*sqrt((2*l+1.)/(2*l-1.) * (l*l-m*m));
00047   }
00048 
00049 void init_lam_fact_deriv_1d (int m, arr<double> &lam_fact)
00050   {
00051   lam_fact[m]=0;
00052   for (int l=m+1; l<lam_fact.size(); ++l)
00053     lam_fact[l] = sqrt((2*l+1.)/(2*l-1.) * (l*l-m*m));
00054   }
00055 
00056 void init_normal_l (arr<double> &normal_l)
00057   {
00058   for (int l=0; l<normal_l.size(); ++l)
00059     normal_l[l] = (l<2) ? 0. : sqrt(1./((l+2.)*(l+1.)*l*(l-1.)));
00060   }
00061 
00062 void get_chunk_info (int nrings, int &nchunks, int &chunksize)
00063   {
00064   nchunks = nrings/max(100,nrings/10) + 1;
00065   chunksize = (nrings+nchunks-1)/nchunks;
00066   }
00067 
00068 void fill_work (const xcomplex<double> *datain, int nph, int mmax,
00069   bool shifted, const arr<xcomplex<double> > &shiftarr,
00070   arr<xcomplex<double> > &work)
00071   {
00072   for (int m=1; m<nph; ++m) work[m]=0;
00073   work[0]=datain[0];
00074 
00075   int cnt1=0, cnt2=nph;
00076   for (int m=1; m<=mmax; ++m)
00077     {
00078     if (++cnt1==nph) cnt1=0;
00079     if (--cnt2==-1) cnt2=nph-1;
00080     xcomplex<double> tmp = shifted ? (datain[m]*shiftarr[m]) : datain[m];
00081     work[cnt1] += tmp;
00082     work[cnt2] += conj(tmp);
00083     }
00084   }
00085 
00086 void read_work (const arr<xcomplex<double> >& work, int nph, int mmax,
00087   bool shifted, const arr<xcomplex<double> > &shiftarr,
00088   xcomplex<double> *dataout)
00089   {
00090   int cnt2=0;
00091   for (int m=0; m<=mmax; ++m)
00092     {
00093     dataout[m] = work[cnt2];
00094     if (++cnt2==nph) cnt2=0;
00095     }
00096   if (shifted)
00097     for (int m=0; m<=mmax; ++m) dataout[m] *= shiftarr[m];
00098   }
00099 
00100 void recalc_map2alm (int nph, int mmax, rfft &plan,
00101   arr<xcomplex<double> > &shiftarr)
00102   {
00103   if (plan.size() == nph) return;
00104   plan.Set (nph);
00105   double f1 = pi/nph;
00106   for (int m=0; m<=mmax; ++m)
00107     {
00108     if (m<nph)
00109       shiftarr[m].Set (cos(m*f1),-sin(m*f1));
00110     else
00111       shiftarr[m]=-shiftarr[m-nph];
00112     }
00113   }
00114 
00115 template<typename T> void fft_map2alm (int nph, int mmax, bool shifted,
00116   double weight, rfft &plan, T *mapN, T *mapS,
00117   xcomplex<double> *phas_n, xcomplex<double> *phas_s,
00118   const arr<xcomplex<double> > &shiftarr, arr<xcomplex<double> > &work)
00119   {
00120   for (int m=0; m<nph; ++m) work[m] = mapN[m]*weight;
00121   plan.forward_c(work);
00122   read_work (work, nph, mmax, shifted, shiftarr, phas_n);
00123   if (mapN!=mapS)
00124     {
00125     for (int m=0; m<nph; ++m) work[m] = mapS[m]*weight;
00126     plan.forward_c(work);
00127     read_work (work, nph, mmax, shifted, shiftarr, phas_s);
00128     }
00129   else
00130     for (int m=0; m<=mmax; ++m) phas_s[m]=0;
00131   }
00132 
00133 void recalc_alm2map (int nph, int mmax, rfft &plan,
00134   arr<xcomplex<double> > &shiftarr)
00135   {
00136   if (plan.size() == nph) return;
00137   plan.Set (nph);
00138   double f1 = pi/nph;
00139   for (int m=0; m<=mmax; ++m)
00140     {
00141     if (m<nph)
00142       shiftarr[m].Set (cos(m*f1),sin(m*f1));
00143     else
00144       shiftarr[m]=-shiftarr[m-nph];
00145     }
00146   }
00147 
00148 template<typename T> void fft_alm2map (int nph, int mmax, bool shifted,
00149   rfft &plan, T *mapN, T *mapS, xcomplex<double> *b_north,
00150   xcomplex<double> *b_south, const arr<xcomplex<double> > &shiftarr,
00151   arr<xcomplex<double> > &work)
00152   {
00153   fill_work (b_north, nph, mmax, shifted, shiftarr, work);
00154   plan.backward_c(work);
00155   for (int m=0; m<nph; ++m) mapN[m] = work[m].re;
00156   if (mapN==mapS) return;
00157   fill_work (b_south, nph, mmax, shifted, shiftarr, work);
00158   plan.backward_c(work);
00159   for (int m=0; m<nph; ++m) mapS[m] = work[m].re;
00160   }
00161 
00162 } // namespace
00163 
00164 template<typename T> void map2alm (const Healpix_Map<T> &map,
00165   Alm<xcomplex<T> > &alm, const arr<double> &weight, bool add_alm)
00166   {
00167   planck_assert (map.Scheme()==RING, "map2alm: map must be in RING scheme");
00168   planck_assert (weight.size()>=2*map.Nside(),
00169     "map2alm: weight array has too few entries");
00170 
00171   int lmax = alm.Lmax(), mmax = alm.Mmax(), nside = map.Nside();
00172 
00173   int nchunks, chunksize;
00174   get_chunk_info(2*nside,nchunks,chunksize);
00175 
00176   arr2<xcomplex<double> > phas_n(chunksize,mmax+1), phas_s(chunksize,mmax+1);
00177   arr<double> cth(chunksize), sth(chunksize);
00178   double normfact = pi/(3*nside*nside);
00179 
00180   if (!add_alm) alm.SetToZero();
00181 
00182   for (int chunk=0; chunk<nchunks; ++chunk)
00183     {
00184     int llim=chunk*chunksize, ulim=min(llim+chunksize,2*nside);
00185 
00186 #pragma omp parallel
00187 {
00188     arr<xcomplex<double> > shiftarr(mmax+1), work(4*nside);
00189     rfft plan;
00190 
00191     int ith;
00192 #pragma omp for schedule(dynamic,1)
00193     for (ith=llim; ith<ulim; ++ith)
00194       {
00195       int istart_north, istart_south, nph;
00196       bool shifted;
00197       map.get_ring_info (ith+1,istart_north,nph,cth[ith-llim],sth[ith-llim],
00198                          shifted);
00199       istart_south = 12*nside*nside - istart_north - nph;
00200 
00201       recalc_map2alm (nph, mmax, plan, shiftarr);
00202       fft_map2alm (nph, mmax, shifted, weight[ith]*normfact, plan,
00203         &map[istart_north], &map[istart_south], phas_n[ith-llim],
00204         phas_s[ith-llim], shiftarr, work);
00205       }
00206 } // end of parallel region
00207 
00208 #pragma omp parallel
00209 {
00210     Ylmgen generator(lmax,mmax,1e-30);
00211     arr<double> Ylm;
00212     arr<xcomplex<double> > alm_tmp(lmax+1);
00213     int m;
00214 #pragma omp for schedule(dynamic,1)
00215     for (m=0; m<=mmax; ++m)
00216       {
00217       for (int l=m; l<=lmax; ++l) alm_tmp[l].Set(0.,0.);
00218       for (int ith=0; ith<ulim-llim; ++ith)
00219         {
00220         int l;
00221         generator.get_Ylm(cth[ith],sth[ith],m,Ylm,l);
00222         if (l<=lmax)
00223           {
00224           xcomplex<double> p1 = phas_n[ith][m]+phas_s[ith][m],
00225                            p2 = phas_n[ith][m]-phas_s[ith][m];
00226 
00227           if ((l-m)&1) goto middle;
00228 start:    alm_tmp[l].re += p1.re*Ylm[l]; alm_tmp[l].im += p1.im*Ylm[l];
00229           if (++l>lmax) goto end;
00230 middle:   alm_tmp[l].re += p2.re*Ylm[l]; alm_tmp[l].im += p2.im*Ylm[l];
00231           if (++l<=lmax) goto start;
00232 end:      ;
00233           }
00234         }
00235       xcomplex<T> *palm = alm.mstart(m);
00236       for (int l=m; l<=lmax; ++l)
00237         { palm[l].re += alm_tmp[l].re; palm[l].im += alm_tmp[l].im; }
00238       }
00239 } // end of parallel region
00240     }
00241   }
00242 
00243 template void map2alm (const Healpix_Map<float> &map,
00244   Alm<xcomplex<float> > &alm, const arr<double> &weight,
00245   bool add_alm);
00246 template void map2alm (const Healpix_Map<double> &map,
00247   Alm<xcomplex<double> > &alm, const arr<double> &weight,
00248   bool add_alm);
00249 
00250 template<typename T> void map2alm_iter (const Healpix_Map<T> &map,
00251   Alm<xcomplex<T> > &alm, int num_iter, const arr<double> &weight)
00252   {
00253   map2alm(map,alm,weight);
00254   for (int iter=1; iter<=num_iter; ++iter)
00255     {
00256     Healpix_Map<T> map2(map.Nside(),map.Scheme(),SET_NSIDE);
00257     alm2map(alm,map2);
00258     for (int m=0; m<map.Npix(); ++m)
00259       map2[m] = map[m]-map2[m];
00260     map2alm(map2,alm,weight,true);
00261     }
00262   }
00263 
00264 template void map2alm_iter (const Healpix_Map<float> &map,
00265   Alm<xcomplex<float> > &alm, int num_iter,
00266   const arr<double> &weight);
00267 template void map2alm_iter (const Healpix_Map<double> &map,
00268   Alm<xcomplex<double> > &alm, int num_iter,
00269   const arr<double> &weight);
00270 
00271 template<typename T> void map2alm_iter2 (const Healpix_Map<T> &map,
00272   Alm<xcomplex<T> > &alm, double err_abs, double err_rel)
00273   {
00274   double x_err_abs=1./err_abs, x_err_rel=1./err_rel;
00275   arr<double> wgt(2*map.Nside());
00276   wgt.fill(1);
00277   Healpix_Map<T> map2(map);
00278   alm.SetToZero();
00279   while(true)
00280     {
00281     map2alm(map2,alm,wgt,true);
00282     alm2map(alm,map2);
00283     double errmeasure=0;
00284     for (int m=0; m<map.Npix(); ++m)
00285       {
00286       double err = abs(map[m]-map2[m]);
00287       double rel = (map[m]!=0) ? abs(err/map[m]) : 1e300;
00288       errmeasure = max(errmeasure,min(err*x_err_abs,rel*x_err_rel));
00289       map2[m] = map[m]-map2[m];
00290       }
00291 cout << "map error measure: " << errmeasure << endl;
00292     if (errmeasure<1) break;
00293     }
00294   }
00295 
00296 template void map2alm_iter2 (const Healpix_Map<double> &map,
00297   Alm<xcomplex<double> > &alm, double err_abs, double err_rel);
00298 
00299 #define SETUP_LAMBDA \
00300   const double t1  = lam_lm1m*lam_fact[l]; \
00301   const double a_w = (l-m2)*two_on_s2 + l*(l-1); \
00302   const double a_x = twocth*(l-1)*lam_lm; \
00303   const double lambda_w = a_w*lam_lm - t1*c_on_s2; \
00304   const double lambda_x = m_on_s2 * (a_x-t1);
00305 
00306 #define MAP2ALM_POL_MACRO(T1,Q1,Q2,U1,U2) \
00307   { \
00308   double lam_lm1m=lam_lm; \
00309   lam_lm=Ylm[l]; \
00310   alm_tmp[l][0].re += T1.re*lam_lm; alm_tmp[l][0].im += T1.im*lam_lm; \
00311   SETUP_LAMBDA \
00312   alm_tmp[l][1].re += Q1.re*lambda_w - U2.im*lambda_x; \
00313   alm_tmp[l][1].im += Q1.im*lambda_w + U2.re*lambda_x; \
00314   alm_tmp[l][2].re += U1.re*lambda_w + Q2.im*lambda_x; \
00315   alm_tmp[l][2].im += U1.im*lambda_w - Q2.re*lambda_x; \
00316   }
00317 
00318 template<typename T> void map2alm_pol
00319   (const Healpix_Map<T> &mapT,
00320    const Healpix_Map<T> &mapQ,
00321    const Healpix_Map<T> &mapU,
00322    Alm<xcomplex<T> > &almT,
00323    Alm<xcomplex<T> > &almG,
00324    Alm<xcomplex<T> > &almC,
00325    const arr<double> &weightT,
00326    const arr<double> &weightQ,
00327    const arr<double> &weightU,
00328    bool add_alm)
00329   {
00330   planck_assert (mapT.Scheme()==RING,
00331     "map2alm_pol: maps must be in RING scheme");
00332   planck_assert (mapT.conformable(mapQ) && mapT.conformable(mapU),
00333     "map2alm_pol: maps are not conformable");
00334   planck_assert (almT.conformable(almG) && almT.conformable(almC),
00335     "map2alm_pol: a_lms are not conformable");
00336   planck_assert ((weightT.size()>=2*mapT.Nside()) &&
00337     (weightQ.size()>=2*mapT.Nside()) && (weightU.size()>=2*mapT.Nside()),
00338     "map2alm_pol: at least one weight array has too few entries");
00339 
00340   int lmax = almT.Lmax(), mmax = almT.Mmax(), nside = mapT.Nside();
00341 
00342   arr<double> normal_l (lmax+1);
00343   init_normal_l (normal_l);
00344 
00345   int nchunks, chunksize;
00346   get_chunk_info(2*nside,nchunks,chunksize);
00347 
00348   arr2<xcomplex<double> > phas_nT(chunksize,mmax+1),phas_sT(chunksize,mmax+1),
00349                           phas_nQ(chunksize,mmax+1),phas_sQ(chunksize,mmax+1),
00350                           phas_nU(chunksize,mmax+1),phas_sU(chunksize,mmax+1);
00351 
00352   arr<double> cth(chunksize), sth(chunksize);
00353   double normfact = pi/(3*nside*nside);
00354 
00355   if (!add_alm)
00356     { almT.SetToZero(); almG.SetToZero(); almC.SetToZero(); }
00357 
00358   for (int chunk=0; chunk<nchunks; ++chunk)
00359     {
00360     int llim=chunk*chunksize, ulim=min(llim+chunksize,2*nside);
00361 
00362 #pragma omp parallel
00363 {
00364     arr<xcomplex<double> > shiftarr(mmax+1), work(4*nside);
00365     rfft plan;
00366 
00367     int ith;
00368 #pragma omp for schedule(dynamic,1)
00369     for (ith=llim; ith<ulim; ++ith)
00370       {
00371       int istart_north, istart_south, nph;
00372       bool shifted;
00373       mapT.get_ring_info (ith+1,istart_north,nph,cth[ith-llim],sth[ith-llim],
00374                           shifted);
00375       istart_south = 12*nside*nside - istart_north - nph;
00376 
00377       recalc_map2alm (nph, mmax, plan, shiftarr);
00378       fft_map2alm (nph, mmax, shifted, weightT[ith]*normfact, plan,
00379         &mapT[istart_north], &mapT[istart_south], phas_nT[ith-llim],
00380         phas_sT[ith-llim], shiftarr, work);
00381       fft_map2alm (nph, mmax, shifted, weightQ[ith]*normfact, plan,
00382         &mapQ[istart_north], &mapQ[istart_south], phas_nQ[ith-llim],
00383         phas_sQ[ith-llim], shiftarr, work);
00384       fft_map2alm (nph, mmax, shifted, weightU[ith]*normfact, plan,
00385         &mapU[istart_north], &mapU[istart_south], phas_nU[ith-llim],
00386         phas_sU[ith-llim], shiftarr, work);
00387       }
00388 } // end of parallel region
00389 
00390 #pragma omp parallel
00391 {
00392     Ylmgen generator(lmax,mmax,1e-30);
00393     arr<double> Ylm;
00394     arr<double> lam_fact(lmax+1);
00395     arr<xcomplex<double>[3] > alm_tmp(lmax+1);
00396     int m;
00397 #pragma omp for schedule(dynamic,1)
00398     for (m=0; m<=mmax; ++m)
00399       {
00400       init_lam_fact_1d (m,lam_fact);
00401 
00402       for (int l=m; l<alm_tmp.size(); ++l)
00403         alm_tmp[l][0]=alm_tmp[l][1]=alm_tmp[l][2] = 0;
00404 
00405       for (int ith=0; ith<ulim-llim; ++ith)
00406         {
00407         int l;
00408         generator.get_Ylm(cth[ith],sth[ith],m,Ylm,l);
00409         if (l<=lmax)
00410           {
00411           double one_on_s2 = 1/(sth[ith]*sth[ith]);
00412           double c_on_s2 = cth[ith] * one_on_s2;
00413           double two_on_s2 = 2*one_on_s2;
00414           double twocth = 2*cth[ith];
00415           int m2 = m*m;
00416           double m_on_s2 = m*one_on_s2;
00417 
00418           xcomplex<double> T1 = phas_nT[ith][m]+phas_sT[ith][m],
00419                            T2 = phas_nT[ith][m]-phas_sT[ith][m],
00420                            Q1 = phas_nQ[ith][m]+phas_sQ[ith][m],
00421                            Q2 = phas_nQ[ith][m]-phas_sQ[ith][m],
00422                            U1 = phas_nU[ith][m]+phas_sU[ith][m],
00423                            U2 = phas_nU[ith][m]-phas_sU[ith][m];
00424 
00425           double lam_lm = 0;
00426           if ((l-m)&1) goto middle;
00427 start:    MAP2ALM_POL_MACRO(T1,Q1,Q2,U1,U2)
00428           if (++l>lmax) goto end;
00429 middle:   MAP2ALM_POL_MACRO(T2,Q2,Q1,U2,U1)
00430           if (++l<=lmax) goto start;
00431 end:      ;
00432           }
00433         }
00434       xcomplex<T> *palmT=almT.mstart(m), *palmG=almG.mstart(m),
00435                   *palmC=almC.mstart(m);
00436       for (int l=m;l<=lmax;++l)
00437         {
00438         palmT[l].re += alm_tmp[l][0].re;
00439         palmT[l].im += alm_tmp[l][0].im;
00440         palmG[l].re += alm_tmp[l][1].re*normal_l[l];
00441         palmG[l].im += alm_tmp[l][1].im*normal_l[l];
00442         palmC[l].re += alm_tmp[l][2].re*normal_l[l];
00443         palmC[l].im += alm_tmp[l][2].im*normal_l[l];
00444         }
00445       }
00446 } // end of parallel region
00447     }
00448   }
00449 
00450 template void map2alm_pol
00451   (const Healpix_Map<float> &mapT,
00452    const Healpix_Map<float> &mapQ,
00453    const Healpix_Map<float> &mapU,
00454    Alm<xcomplex<float> > &almT,
00455    Alm<xcomplex<float> > &almG,
00456    Alm<xcomplex<float> > &almC,
00457    const arr<double> &weightT,
00458    const arr<double> &weightQ,
00459    const arr<double> &weightU,
00460    bool add_alm);
00461 template void map2alm_pol
00462   (const Healpix_Map<double> &mapT,
00463    const Healpix_Map<double> &mapQ,
00464    const Healpix_Map<double> &mapU,
00465    Alm<xcomplex<double> > &almT,
00466    Alm<xcomplex<double> > &almG,
00467    Alm<xcomplex<double> > &almC,
00468    const arr<double> &weightT,
00469    const arr<double> &weightQ,
00470    const arr<double> &weightU,
00471    bool add_alm);
00472 
00473 template<typename T> void map2alm_pol_iter
00474   (const Healpix_Map<T> &mapT,
00475    const Healpix_Map<T> &mapQ,
00476    const Healpix_Map<T> &mapU,
00477    Alm<xcomplex<T> > &almT,
00478    Alm<xcomplex<T> > &almG,
00479    Alm<xcomplex<T> > &almC,
00480    int num_iter,
00481    const arr<double> &weightT,
00482    const arr<double> &weightQ,
00483    const arr<double> &weightU)
00484   {
00485   map2alm_pol(mapT,mapQ,mapU,almT,almG,almC,weightT,weightQ,weightU);
00486   for (int iter=1; iter<=num_iter; ++iter)
00487     {
00488     Healpix_Map<T> mapT2(mapT.Nside(),mapT.Scheme(),SET_NSIDE),
00489                    mapQ2(mapT.Nside(),mapT.Scheme(),SET_NSIDE),
00490                    mapU2(mapT.Nside(),mapT.Scheme(),SET_NSIDE);
00491 
00492     alm2map_pol(almT,almG,almC,mapT2,mapQ2,mapU2);
00493     for (int m=0; m<mapT.Npix(); ++m)
00494       {
00495       mapT2[m] = mapT[m]-mapT2[m];
00496       mapQ2[m] = mapQ[m]-mapQ2[m];
00497       mapU2[m] = mapU[m]-mapU2[m];
00498       }
00499     map2alm_pol(mapT2,mapQ2,mapU2,almT,almG,almC,weightT,weightQ,weightU,true);
00500     }
00501   }
00502 
00503 template void map2alm_pol_iter
00504   (const Healpix_Map<float> &mapT,
00505    const Healpix_Map<float> &mapQ,
00506    const Healpix_Map<float> &mapU,
00507    Alm<xcomplex<float> > &almT,
00508    Alm<xcomplex<float> > &almG,
00509    Alm<xcomplex<float> > &almC,
00510    int num_iter,
00511    const arr<double> &weightT,
00512    const arr<double> &weightQ,
00513    const arr<double> &weightU);
00514 template void map2alm_pol_iter
00515   (const Healpix_Map<double> &mapT,
00516    const Healpix_Map<double> &mapQ,
00517    const Healpix_Map<double> &mapU,
00518    Alm<xcomplex<double> > &almT,
00519    Alm<xcomplex<double> > &almG,
00520    Alm<xcomplex<double> > &almC,
00521    int num_iter,
00522    const arr<double> &weightT,
00523    const arr<double> &weightQ,
00524    const arr<double> &weightU);
00525 
00526 template<typename T> void map2alm_pol_iter2
00527   (const Healpix_Map<T> &mapT,
00528    const Healpix_Map<T> &mapQ,
00529    const Healpix_Map<T> &mapU,
00530    Alm<xcomplex<T> > &almT,
00531    Alm<xcomplex<T> > &almG,
00532    Alm<xcomplex<T> > &almC,
00533    double err_abs, double err_rel)
00534   {
00535   arr<double> wgt(2*mapT.Nside());
00536   wgt.fill(1);
00537   Healpix_Map<T> mapT2(mapT), mapQ2(mapQ), mapU2(mapU);
00538   almT.SetToZero(); almG.SetToZero(); almC.SetToZero();
00539   while(true)
00540     {
00541     map2alm_pol(mapT2,mapQ2,mapU2,almT,almG,almC,wgt,wgt,wgt,true);
00542     alm2map_pol(almT,almG,almC,mapT2,mapQ2,mapU2);
00543     double errmeasure=0;
00544     for (int m=0; m<mapT.Npix(); ++m)
00545       {
00546       double err = abs(mapT[m]-mapT2[m]);
00547       double rel = (mapT[m]!=0) ? abs(err/mapT[m]) : 1e300;
00548       errmeasure = max(errmeasure,min(err/err_abs,rel/err_rel));
00549       mapT2[m] = mapT[m]-mapT2[m];
00550       err = abs(mapQ[m]-mapQ2[m]);
00551       rel = (mapQ[m]!=0) ? abs(err/mapQ[m]) : 1e300;
00552       errmeasure = max(errmeasure,min(err/err_abs,rel/err_rel));
00553       mapQ2[m] = mapQ[m]-mapQ2[m];
00554       err = abs(mapU[m]-mapU2[m]);
00555       rel = (mapU[m]!=0) ? abs(err/mapU[m]) : 1e300;
00556       errmeasure = max(errmeasure,min(err/err_abs,rel/err_rel));
00557       mapU2[m] = mapU[m]-mapU2[m];
00558       }
00559 cout << "map error measure: " << errmeasure << endl;
00560     if (errmeasure<1) break;
00561     }
00562   }
00563 
00564 template void map2alm_pol_iter2
00565   (const Healpix_Map<double> &mapT,
00566    const Healpix_Map<double> &mapQ,
00567    const Healpix_Map<double> &mapU,
00568    Alm<xcomplex<double> > &almT,
00569    Alm<xcomplex<double> > &almG,
00570    Alm<xcomplex<double> > &almC,
00571    double err_abs, double err_rel);
00572 
00573 template<typename T> void alm2map (const Alm<xcomplex<T> > &alm,
00574   Healpix_Map<T> &map)
00575   {
00576   int lmax = alm.Lmax(), mmax = alm.Mmax(), nside = map.Nside();
00577 
00578   int nchunks, chunksize;
00579   get_chunk_info(2*nside,nchunks,chunksize);
00580 
00581   arr2<xcomplex<double> > b_north(chunksize,mmax+1), b_south(chunksize,mmax+1);
00582   arr<double> cth(chunksize),sth(chunksize);
00583 
00584   for (int chunk=0; chunk<nchunks; ++chunk)
00585     {
00586     int llim=chunk*chunksize, ulim=min(llim+chunksize,2*nside);
00587     for (int ith=llim; ith<ulim; ++ith)
00588       {
00589       int nph, istart_north;
00590       bool shifted;
00591       map.get_ring_info (ith+1,istart_north,nph, cth[ith-llim],sth[ith-llim],
00592                          shifted);
00593       }
00594 
00595 #pragma omp parallel
00596 {
00597     Ylmgen generator(lmax,mmax,1e-30);
00598     arr<double> Ylm;
00599     arr<xcomplex<double> > alm_tmp(lmax+1);
00600     int m;
00601 #pragma omp for schedule(dynamic,1)
00602     for (m=0; m<=mmax; ++m)
00603       {
00604       for (int l=m; l<=lmax; ++l)
00605         alm_tmp[l]=alm(l,m);
00606 
00607       for (int ith=0; ith<ulim-llim; ++ith)
00608         {
00609         int l;
00610         generator.get_Ylm(cth[ith],sth[ith],m,Ylm,l);
00611         if (l<=lmax)
00612           {
00613           xcomplex<double> p1=0, p2=0;
00614 
00615           if ((l-m)&1) goto middle;
00616 start:    p1.re += alm_tmp[l].re*Ylm[l]; p1.im += alm_tmp[l].im*Ylm[l];
00617           if (++l>lmax) goto end;
00618 middle:   p2.re += alm_tmp[l].re*Ylm[l]; p2.im += alm_tmp[l].im*Ylm[l];
00619           if (++l<=lmax) goto start;
00620 end:      b_north[ith][m] = p1+p2; b_south[ith][m] = p1-p2;
00621           }
00622         else
00623           {
00624           b_north[ith][m] = b_south[ith][m] = 0;
00625           }
00626         }
00627       }
00628 } // end of parallel region
00629 
00630 #pragma omp parallel
00631 {
00632     arr<xcomplex<double> > shiftarr(mmax+1), work(4*nside);
00633     rfft plan;
00634     int ith;
00635 #pragma omp for schedule(dynamic,1)
00636     for (ith=llim; ith<ulim; ++ith)
00637       {
00638       int istart_north, istart_south, nph;
00639       double dum1, dum2;
00640       bool shifted;
00641       map.get_ring_info (ith+1,istart_north,nph,dum1,dum2,shifted);
00642       istart_south = map.Npix()-istart_north-nph;
00643       recalc_alm2map (nph, mmax, plan, shiftarr);
00644       fft_alm2map (nph, mmax, shifted, plan, &map[istart_north],
00645         &map[istart_south], b_north[ith-llim], b_south[ith-llim],
00646         shiftarr, work);
00647       }
00648 } // end of parallel region
00649     }
00650   }
00651 
00652 template void alm2map (const Alm<xcomplex<double> > &alm,
00653   Healpix_Map<double> &map);
00654 template void alm2map (const Alm<xcomplex<float> > &alm,
00655   Healpix_Map<float> &map);
00656 
00657 #define ALM2MAP_POL_MACRO(bT1,bQ1,bQ2,bU1,bU2) \
00658   { \
00659   double lam_lm1m = lam_lm; \
00660   lam_lm = Ylm[l]; \
00661   SETUP_LAMBDA \
00662   bT1.re+=alm_tmp[l][0].re*lam_lm;   bT1.im+=alm_tmp[l][0].im*lam_lm; \
00663   bQ1.re+=alm_tmp[l][1].re*lambda_w; bQ1.im+=alm_tmp[l][1].im*lambda_w; \
00664   bQ2.re-=alm_tmp[l][2].im*lambda_x; bQ2.im+=alm_tmp[l][2].re*lambda_x; \
00665   bU1.re-=alm_tmp[l][2].re*lambda_w; bU1.im-=alm_tmp[l][2].im*lambda_w; \
00666   bU2.re-=alm_tmp[l][1].im*lambda_x; bU2.im+=alm_tmp[l][1].re*lambda_x; \
00667   }
00668 
00669 template<typename T> void alm2map_pol
00670   (const Alm<xcomplex<T> > &almT,
00671    const Alm<xcomplex<T> > &almG,
00672    const Alm<xcomplex<T> > &almC,
00673    Healpix_Map<T> &mapT,
00674    Healpix_Map<T> &mapQ,
00675    Healpix_Map<T> &mapU)
00676   {
00677   planck_assert (mapT.Scheme()==RING,
00678     "alm2map_pol: maps must be in RING scheme");
00679   planck_assert (mapT.conformable(mapQ) && mapT.conformable(mapU),
00680     "alm2map_pol: maps are not conformable");
00681   planck_assert (almT.conformable(almG) && almT.conformable(almC),
00682     "alm2map_pol: a_lms are not conformable");
00683 
00684   int lmax = almT.Lmax(), mmax = almT.Mmax(), nside = mapT.Nside();
00685 
00686   arr<double> normal_l (lmax+1);
00687   init_normal_l (normal_l);
00688 
00689   int nchunks, chunksize;
00690   get_chunk_info(2*nside,nchunks,chunksize);
00691 
00692   arr2<xcomplex<double> >
00693     b_north_T(chunksize,mmax+1), b_south_T(chunksize,mmax+1),
00694     b_north_Q(chunksize,mmax+1), b_south_Q(chunksize,mmax+1),
00695     b_north_U(chunksize,mmax+1), b_south_U(chunksize,mmax+1);
00696 
00697   arr<double> cth(chunksize),sth(chunksize);
00698 
00699   for (int chunk=0; chunk<nchunks; ++chunk)
00700     {
00701     int llim=chunk*chunksize, ulim=min(llim+chunksize,2*nside);
00702     for (int ith=llim; ith<ulim; ++ith)
00703       {
00704       int nph, istart_north;
00705       bool shifted;
00706       mapT.get_ring_info (ith+1,istart_north,nph, cth[ith-llim],sth[ith-llim],
00707                           shifted);
00708       }
00709 
00710 #pragma omp parallel
00711 {
00712     Ylmgen generator(lmax,mmax,1e-30);
00713     arr<double> Ylm;
00714     arr<double> lam_fact (lmax+1);
00715     arr<xcomplex<double>[3]> alm_tmp(lmax+1);
00716     int m;
00717 #pragma omp for schedule(dynamic,1)
00718     for (m=0; m<=mmax; ++m)
00719       {
00720       int m2 = m*m;
00721       init_lam_fact_1d (m,lam_fact);
00722       for (int l=m; l<=lmax; ++l)
00723         {
00724         alm_tmp[l][0] = almT(l,m);
00725         alm_tmp[l][1] = almG(l,m)*(-normal_l[l]);
00726         alm_tmp[l][2] = almC(l,m)*(-normal_l[l]);
00727         }
00728       for (int ith=0; ith<ulim-llim; ++ith)
00729         {
00730         int l;
00731         generator.get_Ylm(cth[ith],sth[ith],m,Ylm,l);
00732         if (l<=lmax)
00733           {
00734           double one_on_s2 = 1/(sth[ith]*sth[ith]);
00735           double c_on_s2 = cth[ith] * one_on_s2;
00736           double two_on_s2 = 2*one_on_s2;
00737           double m_on_s2 = m*one_on_s2;
00738           double twocth = 2*cth[ith];
00739 
00740           xcomplex<double> bT1=0, bT2=0, bQ1=0, bQ2=0, bU1=0, bU2=0;
00741 
00742           double lam_lm = 0;
00743           if ((l-m)&1) goto middle;
00744 start:    ALM2MAP_POL_MACRO(bT1,bQ1,bQ2,bU1,bU2)
00745           if (++l>lmax) goto end;
00746 middle:   ALM2MAP_POL_MACRO(bT2,bQ2,bQ1,bU2,bU1)
00747           if (++l<=lmax) goto start;
00748 end:      b_north_T[ith][m] = bT1+bT2; b_south_T[ith][m] = bT1-bT2;
00749           b_north_Q[ith][m] =-bQ1-bQ2; b_south_Q[ith][m] =-bQ1+bQ2;
00750           b_north_U[ith][m] = bU1+bU2; b_south_U[ith][m] = bU1-bU2;
00751           }
00752         else
00753           {
00754           b_north_T[ith][m] = b_south_T[ith][m] = 0;
00755           b_north_Q[ith][m] = b_south_Q[ith][m] = 0;
00756           b_north_U[ith][m] = b_south_U[ith][m] = 0;
00757           }
00758         }
00759       }
00760 } // end of parallel region
00761 
00762 #pragma omp parallel
00763 {
00764     arr<xcomplex<double> > shiftarr(mmax+1), work(4*nside);
00765     rfft plan;
00766     int ith;
00767 #pragma omp for schedule(dynamic,1)
00768     for (ith=llim; ith<ulim; ++ith)
00769       {
00770       int nph, istart_north, istart_south;
00771       double dum1, dum2;
00772       bool shifted;
00773       mapT.get_ring_info (ith+1,istart_north,nph,dum1,dum2,shifted);
00774       istart_south = mapT.Npix()-istart_north-nph;
00775       recalc_alm2map (nph, mmax, plan, shiftarr);
00776       fft_alm2map (nph, mmax, shifted, plan, &mapT[istart_north],
00777         &mapT[istart_south], b_north_T[ith-llim], b_south_T[ith-llim],
00778         shiftarr, work);
00779       fft_alm2map (nph, mmax, shifted, plan, &mapQ[istart_north],
00780         &mapQ[istart_south], b_north_Q[ith-llim], b_south_Q[ith-llim],
00781         shiftarr, work);
00782       fft_alm2map (nph, mmax, shifted, plan, &mapU[istart_north],
00783         &mapU[istart_south], b_north_U[ith-llim], b_south_U[ith-llim],
00784         shiftarr, work);
00785       }
00786 } // end of parallel region
00787     }
00788   }
00789 
00790 template void alm2map_pol (const Alm<xcomplex<double> > &almT,
00791                            const Alm<xcomplex<double> > &almG,
00792                            const Alm<xcomplex<double> > &almC,
00793                            Healpix_Map<double> &mapT,
00794                            Healpix_Map<double> &mapQ,
00795                            Healpix_Map<double> &mapU);
00796 
00797 template void alm2map_pol (const Alm<xcomplex<float> > &almT,
00798                            const Alm<xcomplex<float> > &almG,
00799                            const Alm<xcomplex<float> > &almC,
00800                            Healpix_Map<float> &mapT,
00801                            Healpix_Map<float> &mapQ,
00802                            Healpix_Map<float> &mapU);
00803 
00804 #define ALM2MAP_DER1_MACRO(b1,bdth1,bdph1) \
00805   { \
00806   double lam_lm1m = lam_lm; \
00807   lam_lm = Ylm[l]; \
00808   const double t1 = alm_tmp[l].re*lam_lm; \
00809   const double t2 = alm_tmp[l].im*lam_lm; \
00810   const double t3 = l*cotanth; \
00811   const double t4 = one_on_s*lam_lm1m*lam_fact[l]; \
00812   b1.re+=t1; b1.im+=t2; \
00813   bdth1.re+=t3*t1-t4*alm_tmp[l].re; bdth1.im+=t3*t2-t4*alm_tmp[l].im; \
00814   bdph1.re-=m*t2; bdph1.im+=m*t1; \
00815   }
00816 
00817 template<typename T> void alm2map_der1
00818   (const Alm<xcomplex<T> > &alm,
00819    Healpix_Map<T> &map,
00820    Healpix_Map<T> &mapdth,
00821    Healpix_Map<T> &mapdph)
00822   {
00823   planck_assert (map.Scheme()==RING,
00824     "alm2map_der1: maps must be in RING scheme");
00825   planck_assert (map.conformable(mapdth) && map.conformable(mapdph),
00826     "alm2map_der1: maps are not conformable");
00827 
00828   int lmax = alm.Lmax(), mmax = alm.Mmax(), nside = map.Nside();
00829 
00830   int nchunks, chunksize;
00831   get_chunk_info(2*nside,nchunks,chunksize);
00832 
00833   arr2<xcomplex<double> >
00834     b_north(chunksize,mmax+1), b_south(chunksize,mmax+1),
00835     b_north_dth(chunksize,mmax+1), b_south_dth(chunksize,mmax+1),
00836     b_north_dph(chunksize,mmax+1), b_south_dph(chunksize,mmax+1);
00837 
00838   arr<double> cth(chunksize),sth(chunksize);
00839 
00840   for (int chunk=0; chunk<nchunks; ++chunk)
00841     {
00842     int llim=chunk*chunksize, ulim=min(llim+chunksize,2*nside);
00843     for (int ith=llim; ith<ulim; ++ith)
00844       {
00845       int nph, istart_north;
00846       bool shifted;
00847       map.get_ring_info (ith+1,istart_north,nph, cth[ith-llim],sth[ith-llim],
00848                          shifted);
00849       }
00850 
00851 #pragma omp parallel
00852 {
00853     Ylmgen generator(lmax,mmax,1e-30);
00854     arr<double> Ylm;
00855     arr<double> lam_fact (lmax+1);
00856     int m;
00857 #pragma omp for schedule(dynamic,1)
00858     for (m=0; m<=mmax; ++m)
00859       {
00860       const xcomplex<T> *alm_tmp=alm.mstart(m);
00861       init_lam_fact_deriv_1d (m,lam_fact);
00862       for (int ith=0; ith<ulim-llim; ++ith)
00863         {
00864         int l;
00865         generator.get_Ylm(cth[ith],sth[ith],m,Ylm,l);
00866         if (l<=lmax)
00867           {
00868           double one_on_s = 1/sth[ith];
00869           double cotanth = cth[ith]*one_on_s;
00870 
00871           xcomplex<double> b1=0, b2=0, bdth1=0, bdth2=0, bdph1=0, bdph2=0;
00872 
00873           double lam_lm = 0;
00874           if ((l-m)&1) goto middle;
00875 start:    ALM2MAP_DER1_MACRO(b1,bdth1,bdph1)
00876           if (++l>lmax) goto end;
00877 middle:   ALM2MAP_DER1_MACRO(b2,bdth2,bdph2)
00878           if (++l<=lmax) goto start;
00879 end:      b_north[ith][m] = b1+b2; b_south[ith][m] = b1-b2;
00880           b_north_dth[ith][m] = bdth1+bdth2; b_south_dth[ith][m] = bdth2-bdth1;
00881           b_north_dph[ith][m] = (bdph1+bdph2)*one_on_s;
00882           b_south_dph[ith][m] = (bdph1-bdph2)*one_on_s;
00883           }
00884         else
00885           {
00886           b_north[ith][m] = b_south[ith][m] = 0;
00887           b_north_dth[ith][m] = b_south_dth[ith][m] = 0;
00888           b_north_dph[ith][m] = b_south_dph[ith][m] = 0;
00889           }
00890         }
00891       }
00892 } // end of parallel region
00893 
00894 #pragma omp parallel
00895 {
00896     arr<xcomplex<double> > shiftarr(mmax+1), work(4*nside);
00897     rfft plan;
00898     int ith;
00899 #pragma omp for schedule(dynamic,1)
00900     for (ith=llim; ith<ulim; ++ith)
00901       {
00902       int nph, istart_north, istart_south;
00903       double dum1, dum2;
00904       bool shifted;
00905       map.get_ring_info (ith+1,istart_north,nph,dum1,dum2,shifted);
00906       istart_south = map.Npix()-istart_north-nph;
00907       recalc_alm2map (nph, mmax, plan, shiftarr);
00908       fft_alm2map (nph, mmax, shifted, plan, &map[istart_north],
00909         &map[istart_south], b_north[ith-llim], b_south[ith-llim],
00910         shiftarr, work);
00911       fft_alm2map (nph, mmax, shifted, plan, &mapdth[istart_north],
00912         &mapdth[istart_south], b_north_dth[ith-llim], b_south_dth[ith-llim],
00913         shiftarr, work);
00914       fft_alm2map (nph, mmax, shifted, plan, &mapdph[istart_north],
00915         &mapdph[istart_south], b_north_dph[ith-llim], b_south_dph[ith-llim],
00916         shiftarr, work);
00917       }
00918 } // end of parallel region
00919     }
00920   }
00921 
00922 template void alm2map_der1 (const Alm<xcomplex<double> > &alm,
00923                             Healpix_Map<double> &map,
00924                             Healpix_Map<double> &map_dth,
00925                             Healpix_Map<double> &map_dph);
00926 
00927 template void alm2map_der1 (const Alm<xcomplex<float> > &alm,
00928                             Healpix_Map<float> &map,
00929                             Healpix_Map<float> &map_dth,
00930                             Healpix_Map<float> &map_dph);

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