PURIFY
Next-generation radio interferometric imaging
Functions
padmm_random_coverage.cc File Reference
#include <array>
#include <memory>
#include <random>
#include <boost/math/special_functions/erf.hpp>
#include "purify/directories.h"
#include "purify/logging.h"
#include "purify/operators.h"
#include "purify/wproj_operators.h"
#include <sopt/credible_region.h>
#include <sopt/imaging_padmm.h>
#include <sopt/power_method.h>
#include <sopt/relative_variation.h>
#include <sopt/utilities.h>
#include <sopt/wavelets.h>
#include <sopt/wavelets/sara.h>
#include "purify/types.h"
#include "purify/cimg.h"
#include "purify/pfitsio.h"
#include "purify/utilities.h"
+ Include dependency graph for padmm_random_coverage.cc:

Go to the source code of this file.

Functions

void padmm (const std::string &name, const Image< t_complex > &M31, const std::string &kernel, const t_int J, const utilities::vis_params &uv_data, const t_real sigma, const std::tuple< bool, t_real > &w_term)
 
int main (int, char **)
 

Function Documentation

◆ main()

int main ( int  ,
char **   
)

Definition at line 134 of file padmm_random_coverage.cc.

134  {
135  // sopt::logging::set_level("debug");
136  // purify::logging::set_level("debug");
137  const std::string &name = "M31";
138  const t_real FoV = 15; // deg
139  const t_real max_w = 15.; // lambda
140  const t_real snr = 30;
141  const bool w_term = false;
142  const std::string kernel = "kb";
143  std::string const fitsfile = image_filename(name + ".fits");
144  auto M31 = pfitsio::read2d(fitsfile);
145  const t_real cellsize = FoV / M31.cols() * 60. * 60.;
146  std::string const inputfile = output_filename(name + "_" + "input.fits");
147 
148  t_real const max = M31.array().abs().maxCoeff();
149  M31 = M31 * 1. / max;
150  pfitsio::write2d(M31.real(), inputfile);
151 
152  t_int const number_of_pxiels = M31.size();
153  t_int const number_of_vis = std::floor(number_of_pxiels * 2);
154  // Generating random uv(w) coverage
155  t_real const sigma_m = constant::pi / 3;
156  auto uv_data = utilities::random_sample_density(number_of_vis, 0, sigma_m, max_w);
157  uv_data.units = utilities::vis_units::radians;
158  PURIFY_MEDIUM_LOG("Number of measurements / number of pixels: {}",
159  uv_data.u.size() * 1. / number_of_pxiels);
160 
161 #ifndef PURIFY_GPU
162  auto const sky_measurements = std::get<2>(sopt::algorithm::normalise_operator<Vector<t_complex>>(
164  uv_data, M31.rows(), M31.cols(), cellsize, cellsize, 2,
165  kernels::kernel_from_string.at("kb"), 8, 30, true, 1e-6, 1e-6, dde_type::wkernel_radial),
166  1000, 0.0001, Vector<t_complex>::Random(M31.size())));
167 #else
168  auto const sky_measurements = std::get<2>(sopt::algorithm::normalise_operator<Vector<t_complex>>(
170  uv_data, M31.rows(), M31.cols(), cellsize, cellsize, 2,
171  kernels::kernel_from_string.at("kb"), 8, 8, w_term),
172  1000, 0.0001, Vector<t_complex>::Random(M31.size())));
173 
174 #endif
175  uv_data.vis = (*sky_measurements) * Image<t_complex>::Map(M31.data(), M31.size(), 1);
176  Vector<t_complex> const y0 = uv_data.vis;
177  // working out value of signal given SNR of 30
178  t_real const sigma = utilities::SNR_to_standard_deviation(y0, snr);
179 
180  std::cout << std::setprecision(13) << sigma << std::endl;
181  // adding noise to visibilities
182  uv_data.vis = utilities::add_noise(y0, 0., sigma);
183  padmm(name + "30", M31, kernel, 4, uv_data, sigma, std::make_tuple(w_term, cellsize));
184  return 0;
185 }
#define PURIFY_MEDIUM_LOG(...)
Medium priority message.
Definition: logging.h:205
const t_real pi
mathematical constant
Definition: types.h:70
const std::map< std::string, kernel > kernel_from_string
Definition: kernels.h:16
std::shared_ptr< sopt::LinearTransform< T > > init_degrid_operator_2d(const Vector< t_real > &u, const Vector< t_real > &v, const Vector< t_real > &w, const Vector< t_complex > &weights, const t_uint &imsizey, const t_uint &imsizex, const t_real &oversample_ratio=2, const kernels::kernel kernel=kernels::kernel::kb, const t_uint Ju=4, const t_uint Jv=4, const bool w_stacking=false, const t_real &cellx=1, const t_real &celly=1)
Returns linear transform that is the standard degridding operator.
Definition: operators.h:608
void write2d(const Image< t_real > &eigen_image, const pfitsio::header_params &header, const bool &overwrite)
Write image to fits file.
Definition: pfitsio.cc:30
Image< t_complex > read2d(const std::string &fits_name)
Read image from fits file.
Definition: pfitsio.cc:109
t_real SNR_to_standard_deviation(const Vector< t_complex > &y0, const t_real &SNR)
Converts SNR to RMS noise.
Definition: utilities.cc:101
Vector< t_complex > add_noise(const Vector< t_complex > &y0, const t_complex &mean, const t_real &standard_deviation)
Add guassian noise to vector.
Definition: utilities.cc:113
utilities::vis_params random_sample_density(const t_int vis_num, const t_real mean, const t_real standard_deviation, const t_real rms_w)
Generates a random visibility coverage.
std::string output_filename(std::string const &filename)
Test output file.
std::string image_filename(std::string const &filename)
Image filename.
void padmm(const std::string &name, const Image< t_complex > &M31, const std::string &kernel, const t_int J, const utilities::vis_params &uv_data, const t_real sigma, const std::tuple< bool, t_real > &w_term)

References purify::utilities::add_noise(), purify::image_filename(), purify::measurementoperator::init_degrid_operator_2d(), purify::kernels::kernel_from_string, purify::output_filename(), padmm(), purify::constant::pi, PURIFY_MEDIUM_LOG, purify::utilities::radians, purify::utilities::random_sample_density(), purify::pfitsio::read2d(), purify::utilities::SNR_to_standard_deviation(), purify::wkernel_radial, and purify::pfitsio::write2d().

◆ padmm()

void padmm ( const std::string &  name,
const Image< t_complex > &  M31,
const std::string &  kernel,
const t_int  J,
const utilities::vis_params uv_data,
const t_real  sigma,
const std::tuple< bool, t_real > &  w_term 
)

Definition at line 25 of file padmm_random_coverage.cc.

27  {
28  std::string const outfile = output_filename(name + "_" + kernel + ".tiff");
29  std::string const outfile_fits = output_filename(name + "_" + kernel + "_solution.fits");
30  std::string const residual_fits = output_filename(name + "_" + kernel + "_residual.fits");
31  std::string const dirty_image = output_filename(name + "_" + kernel + "_dirty.tiff");
32  std::string const dirty_image_fits = output_filename(name + "_" + kernel + "_dirty.fits");
33 
34  t_real const over_sample = 2;
35  t_uint const imsizey = M31.rows();
36  t_uint const imsizex = M31.cols();
37  utilities::write_visibility(uv_data, output_filename("input_data.vis"));
38 #ifndef PURIFY_GPU
39  auto const measurements_transform =
40  std::get<2>(sopt::algorithm::normalise_operator<Vector<t_complex>>(
42  uv_data, imsizey, imsizex, std::get<1>(w_term), std::get<1>(w_term), over_sample,
43  kernels::kernel_from_string.at(kernel), J, 30, true, 1e-6, 1e-6,
44  dde_type::wkernel_radial),
45  1000, 1e-4, Vector<t_complex>::Random(imsizex * imsizey)));
46 #else
47  af::setDevice(0);
48  auto const measurements_transform =
49  std::get<2>(sopt::algorithm::normalise_operator<Vector<t_complex>>(
51  uv_data, imsizey, imsizex, std::get<1>(w_term), std::get<1>(w_term), over_sample,
52  kernels::kernel_from_string.at(kernel), J, J, std::get<0>(w_term)),
53  1000, 1e-4, Vector<t_complex>::Random(imsizex * imsizey)));
54 #endif
55  sopt::wavelets::SARA const sara{
56  std::make_tuple("Dirac", 3u), std::make_tuple("DB1", 3u), std::make_tuple("DB2", 3u),
57  std::make_tuple("DB3", 3u), std::make_tuple("DB4", 3u), std::make_tuple("DB5", 3u),
58  std::make_tuple("DB6", 3u), std::make_tuple("DB7", 3u), std::make_tuple("DB8", 3u)};
59 
60  auto const Psi = sopt::linear_transform<t_complex>(sara, imsizey, imsizex);
61  Vector<> dimage = (measurements_transform->adjoint() * uv_data.vis).real();
62  assert(dimage.size() == M31.size());
63  Vector<t_complex> initial_estimate = Vector<t_complex>::Zero(dimage.size());
64  sopt::utilities::write_tiff(Image<t_real>::Map(dimage.data(), imsizey, imsizex), dirty_image);
65  pfitsio::write2d(Image<t_real>::Map(dimage.data(), imsizey, imsizex), dirty_image_fits);
66 
67  auto const epsilon = utilities::calculate_l2_radius(uv_data.vis.size(), sigma);
68  PURIFY_HIGH_LOG("Using epsilon of {}", epsilon);
69 #ifdef PURIFY_CImg
70  auto const canvas =
71  std::make_shared<CDisplay>(cimg::make_display(Vector<t_real>::Zero(1024 * 512), 1024, 512));
72  canvas->resize(true);
73  auto const show_image = [&, measurements_transform](const Vector<t_complex> &x) -> bool {
74  if (!canvas->is_closed()) {
75  const Vector<t_complex> res =
76  (measurements_transform->adjoint() * (uv_data.vis - (*measurements_transform * x)));
77  const auto img1 = cimg::make_image(x.real().eval(), imsizey, imsizex)
78  .get_normalize(0, 1)
79  .get_resize(512, 512);
80  const auto img2 = cimg::make_image(res.real().eval(), imsizey, imsizex)
81  .get_normalize(0, 1)
82  .get_resize(512, 512);
83  const auto results = CImageList<t_real>(img1, img2);
84  canvas->display(results);
85  canvas->resize(true);
86  }
87  return true;
88  };
89 #endif
90  auto const padmm =
91  sopt::algorithm::ImagingProximalADMM<t_complex>(uv_data.vis)
92  .itermax(500)
93  .regulariser_strength(
94  (Psi.adjoint() * (measurements_transform->adjoint() * uv_data.vis).eval())
95  .cwiseAbs()
96  .maxCoeff() *
97  1e-3)
98  .relative_variation(1e-3)
99  .l2ball_proximal_epsilon(epsilon)
100  .tight_frame(false)
101  .l1_proximal_tolerance(1e-2)
102  .l1_proximal_nu(1.)
103  .l1_proximal_itermax(50)
104  .l1_proximal_positivity_constraint(true)
105  .l1_proximal_real_constraint(true)
106  .residual_convergence(epsilon)
107  .lagrange_update_scale(0.9)
108 #ifdef PURIFY_CImg
109  .is_converged(show_image)
110 #endif
111  .Psi(Psi)
112  .Phi(*measurements_transform);
113 
114  auto const diagnostic = padmm();
115  assert(diagnostic.x.size() == M31.size());
116  Image<t_complex> image = Image<t_complex>::Map(diagnostic.x.data(), imsizey, imsizex);
117  pfitsio::write2d(image.real(), outfile_fits);
118  Vector<t_complex> residuals = measurements_transform->adjoint() *
119  (uv_data.vis - ((*measurements_transform) * diagnostic.x));
120  Image<t_complex> residual_image = Image<t_complex>::Map(residuals.data(), imsizey, imsizex);
121  pfitsio::write2d(residual_image.real(), residual_fits);
122 #ifdef PURIFY_CImg
123  const auto results = CImageList<t_real>(
124  cimg::make_image(diagnostic.x.real().eval(), imsizey, imsizex).get_resize(512, 512),
125  cimg::make_image(residuals.real().eval(), imsizey, imsizex).get_resize(512, 512));
126  canvas->display(results);
127  cimg::make_image(residuals.real().eval(), imsizey, imsizex)
128  .histogram(256)
129  .display_graph("Residual Histogram", 2);
130  while (!canvas->is_closed()) canvas->wait();
131 #endif
132 }
#define PURIFY_HIGH_LOG(...)
High priority message.
Definition: logging.h:203
void write_visibility(const utilities::vis_params &uv_vis, const std::string &h5name, const bool w_term, const size_t chunksize=0)
Write an HDF5 file with u,v visibilities from a vis_params object.
Definition: h5reader.h:244
t_real calculate_l2_radius(const t_uint y_size, const t_real &sigma, const t_real &n_sigma, const std::string distirbution)
A function that calculates the l2 ball radius for sopt.
Definition: utilities.cc:75
Vector< t_complex > vis
Definition: uvw_utilities.h:22

References purify::utilities::calculate_l2_radius(), purify::measurementoperator::init_degrid_operator_2d(), purify::kernels::kernel_from_string, purify::output_filename(), padmm(), PURIFY_HIGH_LOG, purify::utilities::vis_params::vis, purify::wkernel_radial, purify::pfitsio::write2d(), and purify::utilities::write_visibility().

Referenced by main(), padmm(), purify::factory::padmm_factory(), padmm_factory(), and TEST_CASE().