pyinterpol_ng.cc 9.32 KB
Newer Older
Martin Reinecke's avatar
Martin Reinecke committed
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
/*
 *  Copyright (C) 2020 Max-Planck-Society
 *  Author: Martin Reinecke
 */

#include <pybind11/pybind11.h>
#include <pybind11/numpy.h>
#include "interpol_ng.h"

using namespace std;
using namespace mr;

namespace py = pybind11;

namespace {

Martin Reinecke's avatar
Martin Reinecke committed
17
using fptype = double;
18

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

29
vector<Alm<complex<T>>> makevec(const py::array &inp, int64_t lmax, int64_t kmax)
30
  {
31
  auto inp2 = to_mav<complex<T>,2>(inp);
32
33
34
35
36
  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;
  }
37
38
39
40
41
42
43
44
45
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
46
  public:
47
    PyInterpolator(const py::array &slm, const py::array &blm,
48
      bool separate, int64_t lmax, int64_t kmax, T epsilon, T ofactor, int nthreads)
49
50
      : Interpolator<T>(makevec(slm, lmax, lmax),
                        makevec(blm, lmax, kmax),
51
52
53
54
55
56
                        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;

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

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

#if 1
template<typename T> py::array pyrotate_alm(const py::array &alm_, int64_t lmax,
  double psi, double theta, double phi)
  {
  auto a1 = to_mav<complex<T>,1>(alm_);
  auto alm = make_Pyarr<complex<T>>({a1.shape(0)});
  auto a2 = to_mav<complex<T>,1>(alm,true);
  for (size_t i=0; i<a1.shape(0); ++i) a2.v(i)=a1(i);
  auto tmp = Alm<complex<T>>(a2,lmax,lmax);
  rotate_alm(tmp, psi, theta, phi);
93
  return move(alm);
Martin Reinecke's avatar
Martin Reinecke committed
94
95
96
  }
#endif

97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
constexpr const char *pyinterpol_ng_DS = R"""(
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
----------
125
126
127
128
129
130
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.
131
132
133
134
135
136
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
137
138
139
140
141
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.
142
143
144
145
146
147
148
149
150
151
152
153
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
154
155
156
ncomp : int
    the number of components which are going to input to `deinterpol`.
    Can be 1 or 3.
157
158
epsilon : float
    desired accuracy for the interpolation; a typical value is 1e-5
159
160
161
162
163
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.
164
165
166
167
168
169
170
171
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
----------
172
ptg : numpy.ndarray((N, 3), dtype=numpy.float64)
173
174
175
176
177
178
179
    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
-------
180
numpy.array((N, n2), dtype=numpy.float64)
181
    the interpolated values
182
183
    n2 is either 1 (if separate=True was used in the constructor) or the
    second dimension of the input slm and blm arrays (otherwise)
184
185
186
187
188
189
190
191
192
193
194
195
196
197

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
----------
198
ptg : numpy.ndarray((N,3), dtype=numpy.float64)
199
200
201
202
    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]
203
data : numpy.ndarray((N, n2), dtype=numpy.float64)
204
    the interpolated values
205
    n2 must match the `ncomp` value specified in the constructor.
206
207
208
209
210
211
212
213
214
215

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
216
interpolation
217
218
219

Parameters
----------
220
beam : numpy.array(nalm_beam, nbeam), dtype=numpy.complex)
221
222
    spherical harmonic coefficients of the beam with lmax and kmax defined
    in the constructor call
223
    nbeam must match the ncomp specified in the constructor, unless ncomp was 1.
224
225
226

Returns
-------
227
228
numpy.array(nalm_sky, nbeam), dtype=numpy.complex):
    spherical harmonic coefficients of the sky with lmax defined
229
230
231
232
233
234
235
236
    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
237
238
239
240
241
242
} // unnamed namespace

PYBIND11_MODULE(pyinterpol_ng, m)
  {
  using namespace pybind11::literals;

243
244
  m.doc() = pyinterpol_ng_DS;

Martin Reinecke's avatar
Martin Reinecke committed
245
246
  using inter_d = PyInterpolator<double>;
  py::class_<inter_d> (m, "PyInterpolator", pyinterpolator_DS)
247
248
249
250
251
    .def(py::init<const py::array &, const py::array &, bool, int64_t, int64_t, fptype, fptype, int>(),
      initnormal_DS, "sky"_a, "beam"_a, "separate"_a, "lmax"_a, "kmax"_a, "epsilon"_a, "ofactor"_a=fptype(1.5),
      "nthreads"_a=0)
    .def(py::init<int64_t, int64_t, int64_t, fptype, fptype, int>(), initadjoint_DS,
      "lmax"_a, "kmax"_a, "ncomp"_a, "epsilon"_a, "ofactor"_a=fptype(1.5), "nthreads"_a=0)
Martin Reinecke's avatar
Martin Reinecke committed
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
    .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);
//   using inter_f = PyInterpolator<float>;
//   py::class_<inter_f> (m, "PyInterpolator_f", pyinterpolator_DS)
//     .def(py::init<const py::array &, const py::array &, bool, int64_t, int64_t, fptype, fptype, int>(),
//       initnormal_DS, "sky"_a, "beam"_a, "separate"_a, "lmax"_a, "kmax"_a, "epsilon"_a, "ofactor"_a=fptype(1.5),
//       "nthreads"_a=0)
//     .def(py::init<int64_t, int64_t, int64_t, fptype, fptype, int>(), initadjoint_DS,
//       "lmax"_a, "kmax"_a, "ncomp"_a, "epsilon"_a, "ofactor"_a=fptype(1.5), "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);
Martin Reinecke's avatar
Martin Reinecke committed
267
#if 1
268
  m.def("rotate_alm", &pyrotate_alm<fptype>, "alm"_a, "lmax"_a, "psi"_a, "theta"_a,
Martin Reinecke's avatar
Martin Reinecke committed
269
270
    "phi"_a);
#endif
271
  m.def("epsilon_guess", &epsilon_guess, "support"_a, "ofactor"_a);
Martin Reinecke's avatar
Martin Reinecke committed
272
  }