totalconvolve.cc 8.87 KB
Newer Older
Martin Reinecke's avatar
Martin Reinecke committed
1
2
3
4
5
6
7
/*
 *  Copyright (C) 2020 Max-Planck-Society
 *  Author: Martin Reinecke
 */

#include <pybind11/pybind11.h>
#include <pybind11/numpy.h>
8
#include "python/totalconvolve.h"
Martin Reinecke's avatar
Martin Reinecke committed
9

10
11
namespace mr {

Martin Reinecke's avatar
Martin Reinecke committed
12
namespace detail_pymodule_totalconvolve {
13

Martin Reinecke's avatar
Martin Reinecke committed
14
15
16
17
18
19
20
21
22
using namespace std;

namespace py = pybind11;

template<typename T> class PyInterpolator: public Interpolator<T>
  {
  protected:
    using Interpolator<T>::lmax;
    using Interpolator<T>::kmax;
23
    using Interpolator<T>::ncomp;
24
25
26
    using Interpolator<T>::interpol;
    using Interpolator<T>::deinterpol;
    using Interpolator<T>::getSlm;
Martin Reinecke's avatar
Martin Reinecke committed
27

28
vector<Alm<complex<T>>> makevec(const py::array &inp, int64_t lmax, int64_t kmax)
29
  {
30
  auto inp2 = to_mav<complex<T>,2>(inp);
31
32
33
34
35
  vector<Alm<complex<T>>> res;
  for (size_t i=0; i<inp2.shape(1); ++i)
    res.push_back(Alm<complex<T>>(inp2.template subarray<1>({0,i},{inp2.shape(0),0}),lmax, kmax));
  return res;
  }
36
37
38
39
40
41
42
43
44
void makevec_v(py::array &inp, int64_t lmax, int64_t kmax, vector<Alm<complex<T>>> &res)
  {
  auto inp2 = to_mav<complex<T>,2>(inp, true);
  for (size_t i=0; i<inp2.shape(1); ++i)
    {
    auto xtmp = inp2.template subarray<1>({0,i},{inp2.shape(0),0});
    res.emplace_back(xtmp, lmax, kmax);
    }
  }
Martin Reinecke's avatar
Martin Reinecke committed
45
  public:
46
    PyInterpolator(const py::array &slm, const py::array &blm,
47
      bool separate, int64_t lmax, int64_t kmax, T epsilon, T ofactor, int nthreads)
48
49
      : Interpolator<T>(makevec(slm, lmax, lmax),
                        makevec(blm, lmax, kmax),
50
51
52
53
54
55
                        separate, epsilon, ofactor, nthreads) {}
    PyInterpolator(int64_t lmax, int64_t kmax, int64_t ncomp_, T epsilon, T ofactor, int nthreads)
      : Interpolator<T>(lmax, kmax, ncomp_, epsilon, ofactor, nthreads) {}

    using Interpolator<T>::support;

56
    py::array pyinterpol(const py::array &ptg) const
Martin Reinecke's avatar
Martin Reinecke committed
57
58
      {
      auto ptg2 = to_mav<T,2>(ptg);
59
60
      auto res = make_Pyarr<T>({ptg2.shape(0),ncomp});
      auto res2 = to_mav<T,2>(res,true);
61
      interpol(ptg2, res2);
62
      return move(res);
Martin Reinecke's avatar
Martin Reinecke committed
63
64
      }

65
    void pydeinterpol(const py::array &ptg, const py::array &data)
Martin Reinecke's avatar
Martin Reinecke committed
66
67
      {
      auto ptg2 = to_mav<T,2>(ptg);
68
      auto data2 = to_mav<T,2>(data);
69
      deinterpol(ptg2, data2);
Martin Reinecke's avatar
Martin Reinecke committed
70
      }
71
    py::array pygetSlm(const py::array &blm_)
Martin Reinecke's avatar
Martin Reinecke committed
72
      {
73
74
      auto blm = makevec(blm_, lmax, kmax);
      auto res = make_Pyarr<complex<T>>({Alm_Base::Num_Alms(lmax, lmax),blm.size()});
75
      vector<Alm<complex<T>>> slm;
76
      makevec_v(res, lmax, lmax, slm);
77
      getSlm(blm, slm);
78
      return move(res);
Martin Reinecke's avatar
Martin Reinecke committed
79
80
81
      }
  };

Martin Reinecke's avatar
Martin Reinecke committed
82
constexpr const char *totalconvolve_DS = R"""(
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
Python interface for total convolution/interpolation library

All arrays containing spherical harmonic coefficients are assumed to have the
following format:
- values for m=0, l going from 0 to lmax
  (these values must have an imaginary part of zero)
- values for m=1, l going from 1 to lmax
  (these values can be fully complex)
- values for m=2, l going from 2 to lmax
- ...
- values for m=mmax, l going from mmax to lmax 

Error conditions are reported by raising exceptions.
)""";

constexpr const char *pyinterpolator_DS = R"""(
Class encapsulating the convolution/interpolation functionality

The class can be configured for interpolation or for adjoint interpolation, by
means of two different constructors.
)""";

constexpr const char *initnormal_DS = R"""(
Constructor for interpolation mode

Parameters
----------
110
111
112
113
114
115
sky : numpy.ndarray((nalm_sky, ncomp), dtype=numpy.complex)
    spherical harmonic coefficients of the sky. ncomp can be 1 or 3.
beam : numpy.ndarray((nalm_beam, ncomp), dtype=numpy.complex)
    spherical harmonic coefficients of the beam. ncomp can be 1 or 3
separate : bool
    whether contributions of individual components should be added together.
116
117
118
119
120
121
lmax : int
    maximum l in the coefficient arays
kmax : int
    maximum azimuthal moment in the beam coefficients
epsilon : float
    desired accuracy for the interpolation; a typical value is 1e-5
122
123
124
125
126
ofactor : float
    oversampling factor to be used for the interpolation grids.
    Should be in the range [1.2; 2], a typical value is 1.5
    Increasing this factor makes (adjoint) convolution slower and
    increases memory consumption, but speeds up interpolation/deinterpolation.
127
128
129
130
131
132
133
134
135
136
137
138
nthreads : the number of threads to use for computation
)""";

constexpr const char *initadjoint_DS = R"""(
Constructor for adjoint interpolation mode

Parameters
----------
lmax : int
    maximum l in the coefficient arays
kmax : int
    maximum azimuthal moment in the beam coefficients
139
140
141
ncomp : int
    the number of components which are going to input to `deinterpol`.
    Can be 1 or 3.
142
143
epsilon : float
    desired accuracy for the interpolation; a typical value is 1e-5
144
145
146
147
148
ofactor : float
    oversampling factor to be used for the interpolation grids.
    Should be in the range [1.2; 2], a typical value is 1.5
    Increasing this factor makes (adjoint) convolution slower and
    increases memory consumption, but speeds up interpolation/deinterpolation.
149
150
151
152
153
154
155
156
nthreads : the number of threads to use for computation
)""";

constexpr const char *interpol_DS = R"""(
Computes the interpolated values for a given set of angle triplets

Parameters
----------
157
ptg : numpy.ndarray((N, 3), dtype=numpy.float64)
158
159
160
161
162
163
164
    theta, phi and psi angles (in radian) for N pointings
    theta must be in the range [0; pi]
    phi must be in the range [0; 2pi]
    psi should be in the range [-2pi; 2pi]

Returns
-------
165
numpy.array((N, n2), dtype=numpy.float64)
166
    the interpolated values
167
168
    n2 is either 1 (if separate=True was used in the constructor) or the
    second dimension of the input slm and blm arrays (otherwise)
169
170
171
172
173
174
175
176
177
178
179
180
181
182

Notes
-----
    - Can only be called in "normal" (i.e. not adjoint) mode
    - repeated calls to this method are fine, but for good performance the
      number of pointings passed per call should be as large as possible.
)""";

constexpr const char *deinterpol_DS = R"""(
Takes a set of angle triplets and interpolated values and spreads them onto the
data cube.

Parameters
----------
183
ptg : numpy.ndarray((N,3), dtype=numpy.float64)
184
185
186
187
    theta, phi and psi angles (in radian) for N pointings
    theta must be in the range [0; pi]
    phi must be in the range [0; 2pi]
    psi should be in the range [-2pi; 2pi]
188
data : numpy.ndarray((N, n2), dtype=numpy.float64)
189
    the interpolated values
190
    n2 must match the `ncomp` value specified in the constructor.
191
192
193
194
195
196
197
198
199
200

Notes
-----
    - Can only be called in adjoint mode
    - repeated calls to this method are fine, but for good performance the
      number of pointings passed per call should be as large as possible.
)""";

constexpr const char *getSlm_DS = R"""(
Returns a set of sky spherical hamonic coefficients resulting from adjoint
201
interpolation
202
203
204

Parameters
----------
205
beam : numpy.array(nalm_beam, nbeam), dtype=numpy.complex)
206
207
    spherical harmonic coefficients of the beam with lmax and kmax defined
    in the constructor call
208
    nbeam must match the ncomp specified in the constructor, unless ncomp was 1.
209
210
211

Returns
-------
212
213
numpy.array(nalm_sky, nbeam), dtype=numpy.complex):
    spherical harmonic coefficients of the sky with lmax defined
214
215
216
217
218
219
220
221
    in the constructor call

Notes
-----
    - Can only be called in adjoint mode
    - must be the last call to the object
)""";

Martin Reinecke's avatar
Martin Reinecke committed
222
void add_totalconvolve(py::module &msup)
Martin Reinecke's avatar
Martin Reinecke committed
223
224
  {
  using namespace pybind11::literals;
Martin Reinecke's avatar
Martin Reinecke committed
225
  auto m = msup.def_submodule("totalconvolve");
Martin Reinecke's avatar
Martin Reinecke committed
226

Martin Reinecke's avatar
Martin Reinecke committed
227
  m.doc() = totalconvolve_DS;
228

Martin Reinecke's avatar
Martin Reinecke committed
229
  using inter_d = PyInterpolator<double>;
Martin Reinecke's avatar
Martin Reinecke committed
230
  py::class_<inter_d> (m, "PyInterpolator", py::module_local(), pyinterpolator_DS)
231
232
    .def(py::init<const py::array &, const py::array &, bool, int64_t, int64_t, double, double, int>(),
      initnormal_DS, "sky"_a, "beam"_a, "separate"_a, "lmax"_a, "kmax"_a, "epsilon"_a, "ofactor"_a=1.5,
233
      "nthreads"_a=0)
234
235
    .def(py::init<int64_t, int64_t, int64_t, double, double, int>(), initadjoint_DS,
      "lmax"_a, "kmax"_a, "ncomp"_a, "epsilon"_a, "ofactor"_a=1.5, "nthreads"_a=0)
Martin Reinecke's avatar
Martin Reinecke committed
236
237
238
239
    .def ("interpol", &inter_d::pyinterpol, interpol_DS, "ptg"_a)
    .def ("deinterpol", &inter_d::pydeinterpol, deinterpol_DS, "ptg"_a, "data"_a)
    .def ("getSlm", &inter_d::pygetSlm, getSlm_DS, "beam"_a)
    .def ("support", &inter_d::support);
240
  using inter_f = PyInterpolator<float>;
Martin Reinecke's avatar
Martin Reinecke committed
241
  py::class_<inter_f> (m, "PyInterpolator_f", py::module_local(), pyinterpolator_DS)
242
243
244
245
246
247
248
249
250
    .def(py::init<const py::array &, const py::array &, bool, int64_t, int64_t, float, float, int>(),
      initnormal_DS, "sky"_a, "beam"_a, "separate"_a, "lmax"_a, "kmax"_a, "epsilon"_a, "ofactor"_a=1.5f,
      "nthreads"_a=0)
    .def(py::init<int64_t, int64_t, int64_t, float, float, int>(), initadjoint_DS,
      "lmax"_a, "kmax"_a, "ncomp"_a, "epsilon"_a, "ofactor"_a=1.5f, "nthreads"_a=0)
    .def ("interpol", &inter_f::pyinterpol, interpol_DS, "ptg"_a)
    .def ("deinterpol", &inter_f::pydeinterpol, deinterpol_DS, "ptg"_a, "data"_a)
    .def ("getSlm", &inter_f::pygetSlm, getSlm_DS, "beam"_a)
    .def ("support", &inter_f::support);
251
  m.def("epsilon_guess", &epsilon_guess, "support"_a, "ofactor"_a);
Martin Reinecke's avatar
Martin Reinecke committed
252
  }
253
254
255

}

Martin Reinecke's avatar
Martin Reinecke committed
256
using detail_pymodule_totalconvolve::add_totalconvolve;
257
258

}