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-2013 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     double avg=map.average();
00070     map.Add(T(-avg));
00071 
00072     arr<double> weight;
00073     get_ring_weights (params,map.Nside(),weight);
00074 
00075     Alm<xcomplex<T> > alm(nlmax,nlmax);
00076     if (map.Scheme()==NEST) map.swap_scheme();
00077     map2alm_iter(map,alm,num_iter,weight);
00078 
00079     arr<double> temp(nlmax+1);
00080     if (nside_pixwin_in>0)
00081       {
00082       read_pixwin(datadir,nside_pixwin_in,temp);
00083       for (int l=0; l<=nlmax; ++l)
00084         temp[l] = 1/temp[l];
00085       alm.ScaleL (temp);
00086       }
00087     if (nside_pixwin_out>0)
00088       {
00089       read_pixwin(datadir,nside_pixwin_out,temp);
00090       alm.ScaleL (temp);
00091       }
00092 
00093     alm(0,0) += T(avg*sqrt(fourpi));
00094     map.SetNside(nside,RING);
00095     alm2map (alm,map);
00096 
00097     write_Healpix_map_to_fits (outfile,map,planckType<T>());
00098     }
00099   else
00100     {
00101     Healpix_Map<T> mapT, mapQ, mapU;
00102     read_Healpix_map_from_fits(infile,mapT,mapQ,mapU,2);
00103 
00104     double avg=mapT.average();
00105     mapT.Add(T(-avg));
00106 
00107     arr<double> weight;
00108     get_ring_weights (params,mapT.Nside(),weight);
00109 
00110     Alm<xcomplex<T> > almT(nlmax,nlmax), almG(nlmax,nlmax), almC(nlmax,nlmax);
00111     if (mapT.Scheme()==NEST) mapT.swap_scheme();
00112     if (mapQ.Scheme()==NEST) mapQ.swap_scheme();
00113     if (mapU.Scheme()==NEST) mapU.swap_scheme();
00114     map2alm_pol_iter
00115       (mapT,mapQ,mapU,almT,almG,almC,num_iter,weight);
00116 
00117     arr<double> temp(nlmax+1), pol(nlmax+1);
00118     if (nside_pixwin_in>0)
00119       {
00120       read_pixwin(datadir,nside_pixwin_in,temp,pol);
00121       for (int l=0; l<=nlmax; ++l)
00122         { temp[l] = 1/temp[l]; if (pol[l]!=0.) pol[l] = 1/pol[l]; }
00123       almT.ScaleL(temp); almG.ScaleL(pol); almC.ScaleL(pol);
00124       }
00125     if (nside_pixwin_out>0)
00126       {
00127       read_pixwin(datadir,nside_pixwin_out,temp,pol);
00128       almT.ScaleL(temp); almG.ScaleL(pol); almC.ScaleL(pol);
00129       }
00130 
00131     almT(0,0) += T(avg*sqrt(fourpi));
00132     mapT.SetNside(nside,RING);
00133     mapQ.SetNside(nside,RING);
00134     mapU.SetNside(nside,RING);
00135     alm2map_pol (almT,almG,almC,mapT,mapQ,mapU);
00136 
00137     write_Healpix_map_to_fits (outfile,mapT,mapQ,mapU,planckType<T>());
00138     }
00139   }
00140 
00141 } // unnamed namespace
00142 
00143 int udgrade_harmonic_cxx_module (int argc, const char **argv)
00144   {
00145   module_startup ("udgrade_harmonic_cxx", argc, argv);
00146   paramfile params (getParamsFromCmdline(argc,argv));
00147 
00148   bool dp = params.find<bool> ("double_precision",false);
00149   dp ? udgrade_harmonic_cxx<double>(params)
00150      : udgrade_harmonic_cxx<float> (params);
00151 
00152   return 0;
00153   }

Generated on Thu Oct 8 14:48:52 2015 for Healpix C++