udgrade_harmonic_cxx_module.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.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) 2012 Max-Planck-Society
00029  *  Author: Martin Reinecke
00030  */
00031 
00032 #include "paramfile.h"
00033 #include "alm.h"
00034 #include "xcomplex.h"
00035 #include "healpix_data_io.h"
00036 #include "healpix_map.h"
00037 #include "healpix_map_fitsio.h"
00038 #include "alm_healpix_tools.h"
00039 #include "levels_facilities.h"
00040 #include "announce.h"
00041 #include "lsconstants.h"
00042 
00043 using namespace std;
00044 
00045 namespace {
00046 
00047 template<typename T> void udgrade_harmonic_cxx (paramfile &params)
00048   {
00049   string infile = params.template find<string>("infile");
00050   string outfile = params.template find<string>("outfile");
00051   int nlmax = params.template find<int>("nlmax");
00052   int nside = params.template find<int>("nside");
00053   int nside_pixwin_in = params.template find<int>("nside_pixwin_in",0);
00054   planck_assert (nside_pixwin_in>=0,"nside_pixwin_in must be >= 0");
00055   int nside_pixwin_out = params.template find<int>("nside_pixwin_out",0);
00056   planck_assert (nside_pixwin_out>=0,"nside_pixwin_out must be >= 0");
00057   int num_iter = params.template find<int>("iter_order",0);
00058   bool polarisation = params.template find<bool>("polarisation",false);
00059 
00060   string datadir;
00061   if ((nside_pixwin_in>0) || (nside_pixwin_out>0))
00062     datadir = params.template find<string>("healpix_data");
00063 
00064   if (!polarisation)
00065     {
00066     Healpix_Map<T> map;
00067     read_Healpix_map_from_fits(infile,map,1,2);
00068 
00069     arr<double> weight;
00070     get_ring_weights (params,map.Nside(),weight);
00071 
00072     Alm<xcomplex<T> > alm(nlmax,nlmax);
00073     double avg=map.average();
00074     map.Add(T(-avg));
00075     if (map.Scheme()==NEST) map.swap_scheme();
00076     map2alm_iter(map,alm,num_iter,weight);
00077 
00078     arr<double> temp(nlmax+1);
00079     if (nside_pixwin_in>0)
00080       {
00081       read_pixwin(datadir,nside_pixwin_in,temp);
00082       for (int l=0; l<=nlmax; ++l)
00083         temp[l] = 1/temp[l];
00084       alm.ScaleL (temp);
00085       }
00086     if (nside_pixwin_out>0)
00087       {
00088       read_pixwin(datadir,nside_pixwin_out,temp);
00089       alm.ScaleL (temp);
00090       }
00091 
00092     alm(0,0) += T(avg*sqrt(fourpi));
00093     map.SetNside(nside,RING);
00094     alm2map (alm,map);
00095 
00096     write_Healpix_map_to_fits (outfile,map,planckType<T>());
00097     }
00098   else
00099     {
00100     Healpix_Map<T> mapT, mapQ, mapU;
00101     read_Healpix_map_from_fits(infile,mapT,mapQ,mapU,2);
00102 
00103     arr<double> weight;
00104     get_ring_weights (params,mapT.Nside(),weight);
00105 
00106     Alm<xcomplex<T> > almT(nlmax,nlmax), almG(nlmax,nlmax), almC(nlmax,nlmax);
00107     double avg=mapT.average();
00108     mapT.Add(T(-avg));
00109     if (mapT.Scheme()==NEST) mapT.swap_scheme();
00110     if (mapQ.Scheme()==NEST) mapQ.swap_scheme();
00111     if (mapU.Scheme()==NEST) mapU.swap_scheme();
00112     map2alm_pol_iter
00113       (mapT,mapQ,mapU,almT,almG,almC,num_iter,weight);
00114 
00115     arr<double> temp(nlmax+1), pol(nlmax+1);
00116     if (nside_pixwin_in>0)
00117       {
00118       read_pixwin(datadir,nside_pixwin_in,temp,pol);
00119       for (int l=0; l<=nlmax; ++l)
00120         { temp[l] = 1/temp[l]; if (pol[l]!=0.) pol[l] = 1/pol[l]; }
00121       almT.ScaleL(temp); almG.ScaleL(pol); almC.ScaleL(pol);
00122       }
00123     if (nside_pixwin_out>0)
00124       {
00125       read_pixwin(datadir,nside_pixwin_out,temp,pol);
00126       almT.ScaleL(temp); almG.ScaleL(pol); almC.ScaleL(pol);
00127       }
00128 
00129     almT(0,0) += T(avg*sqrt(fourpi));
00130     mapT.SetNside(nside,RING);
00131     mapQ.SetNside(nside,RING);
00132     mapU.SetNside(nside,RING);
00133     alm2map_pol (almT,almG,almC,mapT,mapQ,mapU);
00134 
00135     write_Healpix_map_to_fits (outfile,mapT,mapQ,mapU,planckType<T>());
00136     }
00137   }
00138 
00139 } // unnamed namespace
00140 
00141 int udgrade_harmonic_cxx_module (int argc, const char **argv)
00142   {
00143   module_startup ("udgrade_harmonic_cxx", argc, argv);
00144   paramfile params (getParamsFromCmdline(argc,argv));
00145 
00146   bool dp = params.find<bool> ("double_precision",false);
00147   dp ? udgrade_harmonic_cxx<double>(params)
00148      : udgrade_harmonic_cxx<float> (params);
00149 
00150   return 0;
00151   }

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