pyinterpol_ng.cc 7.75 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, int nthreads=0)
49
50
      : Interpolator<T>(makevec(slm, lmax, lmax),
                        makevec(blm, lmax, kmax),
51
                        separate, epsilon, nthreads) {}
52
    PyInterpolator(int64_t lmax, int64_t kmax, int64_t ncomp_, T epsilon, int nthreads=0)
53
      : Interpolator<T>(lmax, kmax, ncomp_, epsilon, nthreads) {}
54
    py::array pyinterpol(const py::array &ptg) const
Martin Reinecke's avatar
Martin Reinecke committed
55
56
      {
      auto ptg2 = to_mav<T,2>(ptg);
57
58
      auto res = make_Pyarr<T>({ptg2.shape(0),ncomp});
      auto res2 = to_mav<T,2>(res,true);
59
      interpol(ptg2, res2);
60
      return move(res);
Martin Reinecke's avatar
Martin Reinecke committed
61
62
      }

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

#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);
90
  return move(alm);
Martin Reinecke's avatar
Martin Reinecke committed
91
92
93
  }
#endif

94
95
96
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
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
----------
122
123
124
125
126
127
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.
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
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
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
146
147
148
ncomp : int
    the number of components which are going to input to `deinterpol`.
    Can be 1 or 3.
149
150
151
152
153
154
155
156
157
158
epsilon : float
    desired accuracy for the interpolation; a typical value is 1e-5
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
----------
159
ptg : numpy.ndarray((N, 3), dtype=numpy.float64)
160
161
162
163
164
165
166
    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
-------
167
numpy.array((N, n2), dtype=numpy.float64)
168
    the interpolated values
169
170
    n2 is either 1 (if separate=True was used in the constructor) or the
    second dimension of the input slm and blm arrays (otherwise)
171
172
173
174
175
176
177
178
179
180
181
182
183
184

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
----------
185
ptg : numpy.ndarray((N,3), dtype=numpy.float64)
186
187
188
189
    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]
190
data : numpy.ndarray((N, n2), dtype=numpy.float64)
191
    the interpolated values
192
    n2 must match the `ncomp` value specified in the constructor.
193
194
195
196
197
198
199
200
201
202

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
203
interpolation
204
205
206

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

Returns
-------
214
215
numpy.array(nalm_sky, nbeam), dtype=numpy.complex):
    spherical harmonic coefficients of the sky with lmax defined
216
217
218
219
220
221
222
223
    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
224
225
226
227
228
229
} // unnamed namespace

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

230
231
  m.doc() = pyinterpol_ng_DS;

232
233
  py::class_<PyInterpolator<fptype>> (m, "PyInterpolator", pyinterpolator_DS)
    .def(py::init<const py::array &, const py::array &, bool, int64_t, int64_t, fptype, int>(),
234
      initnormal_DS, "sky"_a, "beam"_a, "separate"_a, "lmax"_a, "kmax"_a, "epsilon"_a,
235
      "nthreads"_a)
236
    .def(py::init<int64_t, int64_t, int64_t, fptype, int>(), initadjoint_DS,
237
      "lmax"_a, "kmax"_a, "ncomp"_a, "epsilon"_a, "nthreads"_a)
238
239
240
    .def ("interpol", &PyInterpolator<fptype>::pyinterpol, interpol_DS, "ptg"_a)
    .def ("deinterpol", &PyInterpolator<fptype>::pydeinterpol, deinterpol_DS, "ptg"_a, "data"_a)
    .def ("getSlm", &PyInterpolator<fptype>::pygetSlm, getSlm_DS, "beam"_a);
Martin Reinecke's avatar
Martin Reinecke committed
241
#if 1
242
  m.def("rotate_alm", &pyrotate_alm<fptype>, "alm"_a, "lmax"_a, "psi"_a, "theta"_a,
Martin Reinecke's avatar
Martin Reinecke committed
243
244
245
    "phi"_a);
#endif
  }