File: measurement_operator.cc

package info (click to toggle)
purify 5.0.1-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 186,836 kB
  • sloc: cpp: 17,731; python: 510; xml: 182; makefile: 7; sh: 6
file content (314 lines) | stat: -rw-r--r-- 14,501 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
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
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
#include <iomanip>
#include "catch2/catch_all.hpp"

#include "purify/directories.h"
#include "purify/kernels.h"
#include "purify/operators.h"
#include "purify/pfitsio.h"
#include "purify/test_data.h"
#include "purify/utilities.h"
#include "purify/wproj_operators.h"
#include <sopt/power_method.h>

using namespace purify;
using Catch::Approx;

TEST_CASE("regression_degrid") {
  const t_int imsizex = 256;
  const t_int imsizey = 256;
  Vector<t_real> u = Vector<t_real>::Random(10) * imsizex / 2;
  Vector<t_real> v = Vector<t_real>::Random(10) * imsizey / 2;
  Vector<t_complex> input = Vector<t_complex>::Zero(imsizex * imsizey);
  input(static_cast<t_int>(imsizex * 0.5 + imsizey * 0.5 * imsizex)) = 1.;
  const t_uint M = u.size();
  const t_real oversample_ratio = 2;
  const t_int ftsizev = std::floor(imsizey * oversample_ratio);
  const t_int ftsizeu = std::floor(imsizex * oversample_ratio);
  const t_uint Ju = 8;
  const t_uint Jv = 8;

  const Vector<t_complex> y = Vector<t_complex>::Ones(u.size());
  CAPTURE(u);
  CAPTURE(v);

  SECTION("kb") {
    const kernels::kernel kernel = kernels::kernel::kb;
    const std::shared_ptr<sopt::LinearTransform<Vector<t_complex>>> measure_op =
        measurementoperator::init_degrid_operator_2d<Vector<t_complex>>(
            u, v, Vector<t_real>::Zero(M), Vector<t_complex>::Ones(M), imsizey, imsizex,
            oversample_ratio, kernel, Ju, Jv);
    const Vector<t_complex> y_test = *measure_op * input;
    CAPTURE(y_test.array() / y.array());
    const t_real max_test = y_test.cwiseAbs().mean();
    CAPTURE(y_test / max_test);
    CAPTURE(y);
    CHECK((y_test / max_test).isApprox((y), 1e-6));
  }
  SECTION("pswf") {
    const kernels::kernel kernel = kernels::kernel::pswf;
    const auto measure_op = measurementoperator::init_degrid_operator_2d<Vector<t_complex>>(
        u, v, Vector<t_real>::Zero(M), Vector<t_complex>::Ones(M), imsizey, imsizex,
        oversample_ratio, kernel, 6, 6);
    const Vector<t_complex> y_test = *measure_op * input;
    CAPTURE(y_test.array() / y.array());
    const t_real max_test = y_test.cwiseAbs().mean();
    CAPTURE(y_test / max_test);
    CAPTURE(y);
    CHECK((y_test / max_test).isApprox((y), 1e-3));
  }
  SECTION("gauss") {
    const kernels::kernel kernel = kernels::kernel::gauss;
    const auto measure_op = measurementoperator::init_degrid_operator_2d<Vector<t_complex>>(
        u, v, Vector<t_real>::Zero(M), Vector<t_complex>::Ones(M), imsizey, imsizex,
        oversample_ratio, kernel, 10, 10);
    const Vector<t_complex> y_test = *measure_op * input;
    CAPTURE(y_test.array() / y.array());
    const t_real max_test = y_test.cwiseAbs().mean();
    CAPTURE(y_test / max_test);
    CAPTURE(y);
    CHECK((y_test / max_test).isApprox((y), 1e-4));
  }
}

TEST_CASE("flux units") {
  const t_real oversample_ratio = 2;
  Vector<t_real> u = Vector<t_real>::Random(10);
  Vector<t_real> v = Vector<t_real>::Random(10);
  const t_uint M = u.size();
  const Vector<t_complex> y = Vector<t_complex>::Ones(u.size());
  SECTION("kb") {
    const kernels::kernel kernel = kernels::kernel::kb;
    for (auto& J : {4, 5, 6, 7, 8}) {
      for (auto& imsize : {128, 256, 512}) {
        const std::shared_ptr<sopt::LinearTransform<Vector<t_complex>>> measure_op =
            measurementoperator::init_degrid_operator_2d<Vector<t_complex>>(
                u * imsize / 2, v * imsize / 2, Vector<t_real>::Zero(M), Vector<t_complex>::Ones(M),
                imsize, imsize, oversample_ratio, kernel, J, J);
        Vector<t_complex> input = Vector<t_complex>::Zero(imsize * imsize);
        input(static_cast<t_int>(imsize * 0.5 + imsize * 0.5 * imsize)) = 1.;
        const Vector<t_complex> y_test = (*measure_op * input).eval();
        CAPTURE(y_test.cwiseAbs().mean());
        CAPTURE(y);
        CAPTURE(y_test);
        CAPTURE(J);
        CAPTURE(imsize);
        CHECK(y_test.isApprox(y, 1e-3));
        const Vector<t_complex> psf = (measure_op->adjoint() * y) * 1. / M;
        CHECK(std::real(psf(static_cast<t_int>(imsize * 0.5 + imsize * 0.5 * imsize))) ==
              Approx(1.).margin(0.001));
      }
    }
  }
  SECTION("pswf") {
    const kernels::kernel kernel = kernels::kernel::pswf;
    for (auto& J : {4, 5, 6, 7, 8}) {
      for (auto& imsize : {128, 256, 512}) {
        if (J != 6)
          CHECK_THROWS(measurementoperator::init_degrid_operator_2d<Vector<t_complex>>(
              u * imsize / 2, v * imsize / 2, Vector<t_real>::Zero(M), Vector<t_complex>::Ones(M),
              imsize, imsize, oversample_ratio, kernel, J, J));
        else {
          const std::shared_ptr<sopt::LinearTransform<Vector<t_complex>>> measure_op =
              measurementoperator::init_degrid_operator_2d<Vector<t_complex>>(
                  u * imsize / 2, v * imsize / 2, Vector<t_real>::Zero(M),
                  Vector<t_complex>::Ones(M), imsize, imsize, oversample_ratio, kernel, J, J);
          Vector<t_complex> input = Vector<t_complex>::Zero(imsize * imsize);
          input(static_cast<t_int>(imsize * 0.5 + imsize * 0.5 * imsize)) = 1.;
          const Vector<t_complex> y_test = (*measure_op * input).eval();
          CAPTURE(y_test.cwiseAbs().mean());
          CAPTURE(y);
          CAPTURE(y_test);
          CAPTURE(J);
          CAPTURE(imsize);
          CHECK(y_test.isApprox(y, 1e-3));
          const Vector<t_complex> psf = (measure_op->adjoint() * y) * 1. / M;
          CHECK(std::real(psf(static_cast<t_int>(imsize * 0.5 + imsize * 0.5 * imsize))) ==
                Approx(1.).margin(0.001));
        }
      }
    }
  }
  SECTION("gauss") {
    const kernels::kernel kernel = kernels::kernel::gauss;
    for (auto& J : {4, 5, 6, 7, 8}) {
      for (auto& imsize : {128, 256, 512}) {
        const std::shared_ptr<sopt::LinearTransform<Vector<t_complex>>> measure_op =
            measurementoperator::init_degrid_operator_2d<Vector<t_complex>>(
                u * imsize / 2, v * imsize / 2, Vector<t_real>::Zero(M), Vector<t_complex>::Ones(M),
                imsize, imsize, oversample_ratio, kernel, J, J);
        Vector<t_complex> input = Vector<t_complex>::Zero(imsize * imsize);
        input(static_cast<t_int>(imsize * 0.5 + imsize * 0.5 * imsize)) = 1.;
        const Vector<t_complex> y_test = (*measure_op * input).eval();
        CAPTURE(y_test.cwiseAbs().mean());
        CAPTURE(y);
        CAPTURE(y_test);
        CAPTURE(J);
        CAPTURE(imsize);
        CHECK(y_test.isApprox(y, 1e-2));
        const Vector<t_complex> psf = (measure_op->adjoint() * y) * 1. / M;
        CHECK(std::real(psf(static_cast<t_int>(imsize * 0.5 + imsize * 0.5 * imsize))) ==
              Approx(1.).margin(0.01));
      }
    }
  }
  SECTION("box") {
    const kernels::kernel kernel = kernels::kernel::box;
    for (auto& J : {4, 5, 6, 7, 8}) {
      for (auto& imsize : {128, 256, 512}) {
        const std::shared_ptr<sopt::LinearTransform<Vector<t_complex>>> measure_op =
            measurementoperator::init_degrid_operator_2d<Vector<t_complex>>(
                u * imsize / 2, v * imsize / 2, Vector<t_real>::Zero(M), Vector<t_complex>::Ones(M),
                imsize, imsize, oversample_ratio, kernel, J, J);
        Vector<t_complex> input = Vector<t_complex>::Zero(imsize * imsize);
        input(static_cast<t_int>(imsize * 0.5 + imsize * 0.5 * imsize)) = 1.;
        const Vector<t_complex> y_test = (*measure_op * input).eval();
        CAPTURE(y_test.cwiseAbs().mean());
        CAPTURE(y);
        CAPTURE(y_test);
        CAPTURE(J);
        CAPTURE(imsize);
        CHECK(y_test.isApprox(y, 1e-3));
        const Vector<t_complex> psf = (measure_op->adjoint() * y) * 1. / M;
        CHECK(std::real(psf(static_cast<t_int>(imsize * 0.5 + imsize * 0.5 * imsize))) ==
              Approx(1.).margin(0.001));
      }
    }
  }
  SECTION("wproj kb") {
    const kernels::kernel kernel = kernels::kernel::kb;
    for (auto& J : {4, 5, 6, 7, 8}) {
      for (auto& imsize : {128, 256, 512}) {
        const std::shared_ptr<sopt::LinearTransform<Vector<t_complex>>> measure_op =
            measurementoperator::init_degrid_operator_2d<Vector<t_complex>>(
                u * imsize / 2, v * imsize / 2, Vector<t_real>::Zero(M), Vector<t_complex>::Ones(M),
                imsize, imsize, oversample_ratio, kernel, J, J, false, 1e-6, 1e-6);
        Vector<t_complex> input = Vector<t_complex>::Zero(imsize * imsize);
        input(static_cast<t_int>(imsize * 0.5 + imsize * 0.5 * imsize)) = 1.;
        const Vector<t_complex> y_test = (*measure_op * input).eval();
        CAPTURE(y_test.cwiseAbs().mean());
        CAPTURE(y);
        CAPTURE(y_test);
        CAPTURE(J);
        CAPTURE(imsize);
        CHECK(y_test.isApprox(y, 1e-3));
        const Vector<t_complex> psf = (measure_op->adjoint() * y) * 1. / M;
        CHECK(std::real(psf(static_cast<t_int>(imsize * 0.5 + imsize * 0.5 * imsize))) ==
              Approx(1.).margin(0.001));
      }
    }
  }
  SECTION("wproj pswf") {
    const kernels::kernel kernel = kernels::kernel::pswf;
    for (auto& J : {4, 5, 6, 7, 8}) {
      for (auto& imsize : {128, 256, 512}) {
        if (J != 6)
          CHECK_THROWS(measurementoperator::init_degrid_operator_2d<Vector<t_complex>>(
              u * imsize / 2, v * imsize / 2, Vector<t_real>::Zero(M), Vector<t_complex>::Ones(M),
              imsize, imsize, oversample_ratio, kernel, J, J, false, 1e-6, 1e-6));
        else {
          const std::shared_ptr<sopt::LinearTransform<Vector<t_complex>>> measure_op =
              measurementoperator::init_degrid_operator_2d<Vector<t_complex>>(
                  u * imsize / 2, v * imsize / 2, Vector<t_real>::Zero(M),
                  Vector<t_complex>::Ones(M), imsize, imsize, oversample_ratio, kernel, J, J, false,
                  1e-6, 1e-6);
          Vector<t_complex> input = Vector<t_complex>::Zero(imsize * imsize);
          input(static_cast<t_int>(imsize * 0.5 + imsize * 0.5 * imsize)) = 1.;
          const Vector<t_complex> y_test = (*measure_op * input).eval();
          CAPTURE(y_test.cwiseAbs().mean());
          CAPTURE(y);
          CAPTURE(y_test);
          CAPTURE(J);
          CAPTURE(imsize);
          CHECK(y_test.isApprox(y, 1e-3));
          const Vector<t_complex> psf = (measure_op->adjoint() * y) * 1. / M;
          CHECK(std::real(psf(static_cast<t_int>(imsize * 0.5 + imsize * 0.5 * imsize))) ==
                Approx(1.).margin(0.001));
        }
      }
    }
  }
  SECTION("wproj gauss") {
    const kernels::kernel kernel = kernels::kernel::gauss;
    for (auto& J : {4, 5, 6, 7, 8}) {
      for (auto& imsize : {128, 256, 512}) {
        const std::shared_ptr<sopt::LinearTransform<Vector<t_complex>>> measure_op =
            measurementoperator::init_degrid_operator_2d<Vector<t_complex>>(
                u * imsize / 2, v * imsize / 2, Vector<t_real>::Zero(M), Vector<t_complex>::Ones(M),
                imsize, imsize, oversample_ratio, kernel, J, J, false, 1e-6, 1e-6);
        Vector<t_complex> input = Vector<t_complex>::Zero(imsize * imsize);
        input(static_cast<t_int>(imsize * 0.5 + imsize * 0.5 * imsize)) = 1.;
        const Vector<t_complex> y_test = (*measure_op * input).eval();
        CAPTURE(y_test.cwiseAbs().mean());
        CAPTURE(y);
        CAPTURE(y_test);
        CAPTURE(J);
        CAPTURE(imsize);
        CHECK(y_test.isApprox(y, 1e-2));
        const Vector<t_complex> psf = (measure_op->adjoint() * y) * 1. / M;
        CHECK(std::real(psf(static_cast<t_int>(imsize * 0.5 + imsize * 0.5 * imsize))) ==
              Approx(1.).margin(0.01));
      }
    }
  }
  SECTION("wproj box") {
    const kernels::kernel kernel = kernels::kernel::box;
    for (auto& J : {4, 5, 6, 7, 8}) {
      for (auto& imsize : {128, 256, 512}) {
        const std::shared_ptr<sopt::LinearTransform<Vector<t_complex>>> measure_op =
            measurementoperator::init_degrid_operator_2d<Vector<t_complex>>(
                u * imsize / 2, v * imsize / 2, Vector<t_real>::Zero(M), Vector<t_complex>::Ones(M),
                imsize, imsize, oversample_ratio, kernel, J, J, false, 1e-6, 1e-6);
        Vector<t_complex> input = Vector<t_complex>::Zero(imsize * imsize);
        input(static_cast<t_int>(imsize * 0.5 + imsize * 0.5 * imsize)) = 1.;
        const Vector<t_complex> y_test = (*measure_op * input).eval();
        CAPTURE(y_test.cwiseAbs().mean());
        CAPTURE(y);
        CAPTURE(y_test);
        CAPTURE(J);
        CAPTURE(imsize);
        CHECK(y_test.isApprox(y, 1e-2));
        const Vector<t_complex> psf = (measure_op->adjoint() * y) * 1. / M;
        CHECK(std::real(psf(static_cast<t_int>(imsize * 0.5 + imsize * 0.5 * imsize))) ==
              Approx(1.).margin(0.001));
      }
    }
  }
}
TEST_CASE("normed operator") {
  const t_real oversample_ratio = 2;
  const t_int power_iters = 1000;
  const t_real power_tol = 1e-4;
  Vector<t_real> u = Vector<t_real>::Random(10);
  Vector<t_real> v = Vector<t_real>::Random(10);
  const t_uint M = u.size();
  const Vector<t_complex> y = Vector<t_complex>::Ones(u.size());
  SECTION("kb") {
    const kernels::kernel kernel = kernels::kernel::kb;
    for (auto& J : {4, 5, 6, 7, 8}) {
      for (auto& imsize : {128, 256, 512}) {
        const Vector<t_complex> init = Vector<t_complex>::Ones(imsize * imsize);
        auto power_method_result = sopt::algorithm::normalise_operator<Vector<t_complex>>(
            measurementoperator::init_degrid_operator_2d<Vector<t_complex>>(
                u * imsize / 2, v * imsize / 2, Vector<t_real>::Zero(M), Vector<t_complex>::Ones(M),
                imsize, imsize, oversample_ratio, kernel, J, J),
            power_iters, power_tol, init);
        const t_real norm = std::get<0>(power_method_result);
        const std::shared_ptr<sopt::LinearTransform<Vector<t_complex>>> measure_op =
            std::get<2>(power_method_result);

        Vector<t_complex> input = Vector<t_complex>::Zero(imsize * imsize);
        input(static_cast<t_int>(imsize * 0.5 + imsize * 0.5 * imsize)) = 1.;
        const Vector<t_complex> y_test = (*measure_op * input).eval();
        CAPTURE(y_test.cwiseAbs().mean());
        CAPTURE(y);
        CAPTURE(y_test);
        CAPTURE(J);
        CAPTURE(imsize);
        CHECK((y_test * norm).isApprox(y, 1e-3));
        const Vector<t_complex> psf = (measure_op->adjoint() * y) * 1. / M * norm;
        CHECK(std::real(psf(static_cast<t_int>(imsize * 0.5 + imsize * 0.5 * imsize))) ==
              Approx(1.).margin(0.001));
      }
    }
  }
}