/* * This file is part of Healpix_cxx. * * Healpix_cxx 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. * * Healpix_cxx 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 Healpix_cxx; if not, write to the Free Software * Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA * * For more information about HEALPix, see http://healpix.jpl.nasa.gov */ /* * Healpix_cxx is being developed at the Max-Planck-Institut fuer Astrophysik * and financially supported by the Deutsches Zentrum fuer Luft- und Raumfahrt * (DLR). */ /* * Copyright (C) 2003-2010 Max-Planck-Society * Author: Martin Reinecke */ #include "xcomplex.h" #include "cxxutils.h" #include "paramfile.h" #include "healpix_data_io.h" #include "alm.h" #include "alm_fitsio.h" #include "alm_powspec_tools.h" #include "fitshandle.h" #include "levels_facilities.h" #include "lsconstants.h" using namespace std; namespace { template void mult_alm (paramfile ¶ms) { string infile = params.template find("infile"); string outfile = params.template find("outfile"); int nside_pixwin_in = params.template find("nside_pixwin_in",0); planck_assert (nside_pixwin_in>=0,"nside_pixwin_in must be >= 0"); int nside_pixwin_out = params.template find("nside_pixwin_out",0); planck_assert (nside_pixwin_out>=0,"nside_pixwin_out must be >= 0"); double fwhm_in = arcmin2rad*params.template find("fwhm_arcmin_in",0); planck_assert (fwhm_in>=0,"fwhm_arcmin_in must be >= 0"); double fwhm_out = arcmin2rad*params.template find("fwhm_arcmin_out",0); planck_assert (fwhm_out>=0,"fwhm_arcmin_out must be >= 0"); string datadir; if ((nside_pixwin_in>0) || (nside_pixwin_out>0)) datadir = params.template find("healpix_data"); bool polarisation = params.template find("polarisation"); if (!polarisation) { int nlmax, nmmax; get_almsize(infile, nlmax, nmmax, 2); Alm > alm; read_Alm_from_fits(infile,alm,nlmax,nmmax,2); if (fwhm_in>0) smoothWithGauss (alm, -fwhm_in); arr temp(nlmax+1); if (nside_pixwin_in>0) { read_pixwin(datadir,nside_pixwin_in,temp); for (int l=0; l<=nlmax; ++l) temp[l] = 1/temp[l]; alm.ScaleL (temp); } if (nside_pixwin_out>0) { read_pixwin(datadir,nside_pixwin_out,temp); alm.ScaleL (temp); } if (fwhm_out>0) smoothWithGauss (alm, fwhm_out); write_Alm_to_fits (outfile,alm,nlmax,nmmax,planckType()); } else { int nlmax, nmmax; get_almsize_pol(infile, nlmax, nmmax); Alm > almT, almG, almC; read_Alm_from_fits(infile,almT,almG,almC,nlmax,nmmax,2); if (fwhm_in>0) smoothWithGauss (almT, almG, almC, -fwhm_in); arr temp(nlmax+1), pol(nlmax+1); if (nside_pixwin_in>0) { read_pixwin(datadir,nside_pixwin_in,temp,pol); for (int l=0; l<=nlmax; ++l) { temp[l] = 1/temp[l]; pol[l] = 1/pol[l]; } almT.ScaleL(temp); almG.ScaleL(pol); almC.ScaleL(pol); } if (nside_pixwin_out>0) { read_pixwin(datadir,nside_pixwin_out,temp,pol); almT.ScaleL(temp); almG.ScaleL(pol); almC.ScaleL(pol); } if (fwhm_out>0) smoothWithGauss (almT, almG, almC, fwhm_out); write_Alm_to_fits (outfile,almT,almG,almC,nlmax,nmmax,planckType()); } } } // unnamed namespace int mult_alm_module (int argc, const char **argv) { module_startup ("mult_alm", argc, argv, 2, ""); paramfile params (argv[1]); bool dp = params.find ("double_precision",false); dp ? mult_alm(params) : mult_alm(params); return 0; }