Healpix C++  3.83
udgrade_harmonic_cxx_module.cc
1 /*
2  * This file is part of Healpix_cxx.
3  *
4  * Healpix_cxx is free software; you can redistribute it and/or modify
5  * it under the terms of the GNU General Public License as published by
6  * the Free Software Foundation; either version 2 of the License, or
7  * (at your option) any later version.
8  *
9  * Healpix_cxx is distributed in the hope that it will be useful,
10  * but WITHOUT ANY WARRANTY; without even the implied warranty of
11  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12  * GNU General Public License for more details.
13  *
14  * You should have received a copy of the GNU General Public License
15  * along with Healpix_cxx; if not, write to the Free Software
16  * Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
17  *
18  * For more information about HEALPix, see http://healpix.sourceforge.net
19  */
20 
21 /*
22  * Healpix_cxx is being developed at the Max-Planck-Institut fuer Astrophysik
23  * and financially supported by the Deutsches Zentrum fuer Luft- und Raumfahrt
24  * (DLR).
25  */
26 
27 /*
28  * Copyright (C) 2012-2017 Max-Planck-Society
29  * Author: Martin Reinecke
30  */
31 
32 #include "paramfile.h"
33 #include "alm.h"
34 #include "xcomplex.h"
35 #include "healpix_data_io.h"
36 #include "healpix_map.h"
37 #include "healpix_map_fitsio.h"
38 #include "alm_healpix_tools.h"
39 #include "weight_utils.h"
40 #include "levels_facilities.h"
41 #include "announce.h"
42 #include "lsconstants.h"
43 
44 using namespace std;
45 
46 namespace {
47 
48 template<typename T> void udgrade_harmonic_cxx (paramfile &params)
49  {
50  string infile = params.template find<string>("infile");
51  string outfile = params.template find<string>("outfile");
52  int nlmax = params.template find<int>("nlmax");
53  int nside = params.template find<int>("nside");
54  string pixwin_in = params.template find<string>("pixwin_in","");
55  string pixwin_out = params.template find<string>("pixwin_out","");
56  int num_iter = params.template find<int>("iter_order",0);
57  bool polarisation = params.template find<bool>("polarisation",false);
58  bool do_pwgt = params.param_present("pixelweights");
59  planck_assert(!(params.param_present("ringweights")&&do_pwgt),
60  "both pixel and ring weights specified");
61  planck_assert((num_iter==0)||(!do_pwgt),
62  "iterative analysis cannot be done in combination with pixel weights");
63 
64  if (!polarisation)
65  {
66  Healpix_Map<T> map;
67  read_Healpix_map_from_fits(infile,map,1,2);
68 
69  double avg=map.average();
70  map.Add(T(-avg));
71 
72  if (do_pwgt)
73  {
74  auto pwgt = read_fullweights_from_fits(
75  params.template find<string>("pixelweights"), map.Nside());
76  apply_fullweights(map,pwgt);
77  }
78 
79  arr<double> weight;
80  get_ring_weights (params,map.Nside(),weight);
81 
82  Alm<xcomplex<T> > alm(nlmax,nlmax);
83  if (map.Scheme()==NEST) map.swap_scheme();
84  map2alm_iter(map,alm,num_iter,weight);
85 
86  arr<double> temp(nlmax+1);
87  if (pixwin_in!="")
88  {
89  read_pixwin(pixwin_in,temp);
90  for (int l=0; l<=nlmax; ++l)
91  temp[l] = 1/temp[l];
92  alm.ScaleL (temp);
93  }
94  if (pixwin_out!="")
95  {
96  read_pixwin(pixwin_out,temp);
97  alm.ScaleL (temp);
98  }
99 
100  alm(0,0) += T(avg*sqrt(fourpi));
101  map.SetNside(nside,RING);
102  alm2map (alm,map);
103 
104  write_Healpix_map_to_fits (outfile,map,planckType<T>());
105  }
106  else
107  {
108  Healpix_Map<T> mapT, mapQ, mapU;
109  read_Healpix_map_from_fits(infile,mapT,mapQ,mapU,2);
110 
111  double avg=mapT.average();
112  mapT.Add(T(-avg));
113 
114  if (do_pwgt)
115  {
116  auto pwgt = read_fullweights_from_fits(
117  params.template find<string>("pixelweights"), mapT.Nside());
118  apply_fullweights(mapT,pwgt);
119  apply_fullweights(mapQ,pwgt);
120  apply_fullweights(mapU,pwgt);
121  }
122 
123  arr<double> weight;
124  get_ring_weights (params,mapT.Nside(),weight);
125 
126  Alm<xcomplex<T> > almT(nlmax,nlmax), almG(nlmax,nlmax), almC(nlmax,nlmax);
127  if (mapT.Scheme()==NEST) mapT.swap_scheme();
128  if (mapQ.Scheme()==NEST) mapQ.swap_scheme();
129  if (mapU.Scheme()==NEST) mapU.swap_scheme();
131  (mapT,mapQ,mapU,almT,almG,almC,num_iter,weight);
132 
133  arr<double> temp(nlmax+1), pol(nlmax+1);
134  if (pixwin_in!="")
135  {
136  read_pixwin(pixwin_in,temp,pol);
137  for (int l=0; l<=nlmax; ++l)
138  { temp[l] = 1/temp[l]; if (pol[l]!=0.) pol[l] = 1/pol[l]; }
139  almT.ScaleL(temp); almG.ScaleL(pol); almC.ScaleL(pol);
140  }
141  if (pixwin_out!="")
142  {
143  read_pixwin(pixwin_out,temp,pol);
144  almT.ScaleL(temp); almG.ScaleL(pol); almC.ScaleL(pol);
145  }
146 
147  almT(0,0) += T(avg*sqrt(fourpi));
148  mapT.SetNside(nside,RING);
149  mapQ.SetNside(nside,RING);
150  mapU.SetNside(nside,RING);
151  alm2map_pol (almT,almG,almC,mapT,mapQ,mapU);
152 
153  write_Healpix_map_to_fits (outfile,mapT,mapQ,mapU,planckType<T>());
154  }
155  }
156 
157 } // unnamed namespace
158 
159 int udgrade_harmonic_cxx_module (int argc, const char **argv)
160  {
161  module_startup ("udgrade_harmonic_cxx", argc, argv);
162  paramfile params (getParamsFromCmdline(argc,argv));
163 
164  bool dp = params.find<bool> ("double_precision",false);
165  dp ? udgrade_harmonic_cxx<double>(params)
166  : udgrade_harmonic_cxx<float> (params);
167 
168  return 0;
169  }
void SetNside(int nside, Healpix_Ordering_Scheme scheme)
Definition: healpix_map.h:83
I Nside() const
Definition: healpix_base.h:427
double average() const
Definition: healpix_map.h:244
Definition: alm.h:88
void map2alm_pol_iter(const Healpix_Map< T > &mapT, const Healpix_Map< T > &mapQ, const Healpix_Map< T > &mapU, Alm< xcomplex< T > > &almT, Alm< xcomplex< T > > &almG, Alm< xcomplex< T > > &almC, int num_iter, const arr< double > &weight)
void apply_fullweights(Healpix_Map< T > &map, const std::vector< double > &wgt)
void read_Healpix_map_from_fits(fitshandle &inp, Healpix_Map< T > &map, int colnum)
void alm2map_pol(const Alm< xcomplex< T > > &almT, const Alm< xcomplex< T > > &almG, const Alm< xcomplex< T > > &almC, Healpix_Map< T > &mapT, Healpix_Map< T > &mapQ, Healpix_Map< T > &mapU, bool add_map)
Healpix_Ordering_Scheme Scheme() const
Definition: healpix_base.h:431
void write_Healpix_map_to_fits(fitshandle &out, const Healpix_Map< T > &map, PDT datatype)
void map2alm_iter(const Healpix_Map< T > &map, Alm< xcomplex< T > > &alm, int num_iter, const arr< double > &weight)
void alm2map(const Alm< xcomplex< T > > &alm, Healpix_Map< T > &map, bool add_map)
void swap_scheme()
Definition: healpix_map.h:179
void Add(T val)
Definition: healpix_map.h:255

Generated on Wed Nov 13 2024 12:18:30 for Healpix C++