PURIFY
Next-generation radio interferometric imaging
utilities.cc
Go to the documentation of this file.
1 #include <fstream>
2 #include <sstream>
3 #include <benchmarks/utilities.h>
4 #include "purify/directories.h"
5 #include "purify/distribute.h"
6 #include "purify/logging.h"
7 #include "purify/mpi_utilities.h"
8 #include "purify/operators.h"
9 #include "purify/pfitsio.h"
10 #include <sopt/linear_transform.h>
11 
12 using namespace purify;
13 
14 namespace b_utilities {
15 
16 void Arguments(benchmark::internal::Benchmark *b) {
17  int im_size_max = 4096; // 4096
18  int uv_size_max = 10000000; // 1M, 10M, 100M
19  int kernel_max = 16; // 16
20  for (int i = 128; i <= im_size_max; i *= 2)
21  for (int j = 1000000; j <= uv_size_max; j *= 10)
22  for (int k = 2; k <= kernel_max; k *= 2)
23  if (k * k < i) b->Args({i, j, k});
24 }
25 
26 double duration(std::chrono::high_resolution_clock::time_point start,
27  std::chrono::high_resolution_clock::time_point end) {
28  auto elapsed_seconds = std::chrono::duration_cast<std::chrono::duration<double>>(end - start);
29  return elapsed_seconds.count();
30 }
31 
32 bool updateImage(t_uint newSize, Image<t_complex> &image, t_uint &sizex, t_uint &sizey) {
33  if (sizex == newSize) {
34  return false;
35  }
36  image = Image<t_complex>::Random(newSize, newSize);
37  sizex = image.cols();
38  sizey = image.rows();
39  t_real const max = image.array().abs().maxCoeff();
40  image = image * 1. / max;
41  return true;
42 }
43 
44 bool updateEmptyImage(t_uint newSize, Vector<t_complex> &image, t_uint &sizex, t_uint &sizey) {
45  if (sizex == newSize) {
46  return false;
47  }
48  image.resize(newSize * newSize);
49  sizex = newSize;
50  sizey = newSize;
51  return true;
52 }
53 
54 bool updateMeasurements(t_uint newSize, utilities::vis_params &data) {
55  if (data.vis.size() == newSize) {
56  return false;
57  }
58  data = b_utilities::random_measurements(newSize);
59  return true;
60 }
61 
62 bool updateMeasurements(t_uint newSize, utilities::vis_params &data, t_real &epsilon, bool newImage,
63  Image<t_complex> &image) {
64  if (data.vis.size() == newSize && !newImage) {
65  return false;
66  }
67  const t_real FoV = 1; // deg
68  const t_real cellsize = FoV / image.size() * 60. * 60.;
69  std::tuple<utilities::vis_params, t_real> temp =
70  b_utilities::dirty_measurements(image, newSize, 30., cellsize);
71  data = std::get<0>(temp);
72  epsilon = utilities::calculate_l2_radius(data.vis.size(), std::get<1>(temp));
73 
74  return true;
75 }
76 
77 std::tuple<utilities::vis_params, t_real> dirty_measurements(
78  Image<t_complex> const &ground_truth_image, t_uint number_of_vis, t_real snr,
79  const t_real &cellsize) {
80  auto uv_data = random_measurements(number_of_vis);
81  // creating operator to generate measurements
82  auto measurement_op = measurementoperator::init_degrid_operator_2d<Vector<t_complex>>(
83  uv_data, ground_truth_image.rows(), ground_truth_image.cols(), cellsize, cellsize, 2,
84  kernels::kernel::kb, 8, 8, false);
85  // Generates measurements from image
86  uv_data.vis = (*measurement_op) *
87  Image<t_complex>::Map(ground_truth_image.data(), ground_truth_image.size(), 1);
88 
89  // working out value of signal given SNR
90  auto const sigma = utilities::SNR_to_standard_deviation(uv_data.vis, snr);
91  // adding noise to visibilities
92  uv_data.vis = utilities::add_noise(uv_data.vis, 0., sigma);
93  return std::make_tuple(uv_data, sigma);
94 }
95 
96 utilities::vis_params random_measurements(t_int size, const t_real max_w, const t_int id,
97  const bool cache_visibilities) {
98  utilities::vis_params uv_data;
99 
100  std::stringstream filename;
101  filename << "random_" << size << "_";
102  filename << std::to_string(id) << ".vis";
103  std::string const vis_file = visibility_filename(filename.str());
104  std::ifstream vis_file_str(vis_file);
105 
106  if (cache_visibilities and vis_file_str.good()) {
107  PURIFY_INFO("Reading random visibilities from file {}", vis_file);
108  uv_data = utilities::read_visibility(vis_file, true);
110  } else {
111  PURIFY_INFO("Generating random visibilities");
112  t_real const sigma_m = constant::pi / 3;
113  uv_data = utilities::random_sample_density(size, 0, sigma_m, max_w);
115  if (cache_visibilities) {
116  utilities::write_visibility(uv_data, vis_file, true);
117  }
118  }
119  return uv_data;
120 }
121 
122 #ifdef PURIFY_MPI
123 utilities::vis_params random_measurements(t_int size, sopt::mpi::Communicator const &comm) {
124  if (comm.is_root()) {
125  // Generate random measurements
126  auto uv_data = random_measurements(size);
127  if (comm.size() == 1) return uv_data;
128 
129  // Distribute them
130  auto const order = distribute::distribute_measurements(uv_data, comm, distribute::plan::radial);
131  uv_data = utilities::regroup_and_scatter(uv_data, order, comm);
132  return uv_data;
133  }
134  return utilities::scatter_visibilities(comm);
135 }
136 
137 bool updateMeasurements(t_uint newSize, utilities::vis_params &data,
138  sopt::mpi::Communicator &comm) {
139  if (data.vis.size() == newSize) {
140  return false;
141  }
142  comm = sopt::mpi::Communicator::World();
143  data = b_utilities::random_measurements(newSize, comm);
144  return true;
145 }
146 
147 bool updateMeasurements(t_uint newSize, utilities::vis_params &data, t_real &epsilon, bool newImage,
148  Image<t_complex> &image, sopt::mpi::Communicator &comm) {
149  if (data.vis.size() == newSize && !newImage) {
150  return false;
151  }
152  comm = sopt::mpi::Communicator::World();
153  const t_real FoV = 1; // deg
154  const t_real cellsize = FoV / image.size() * 60. * 60.;
155  std::tuple<utilities::vis_params, t_real> temp =
156  b_utilities::dirty_measurements(image, newSize, 30., cellsize, comm);
157  data = std::get<0>(temp);
158  epsilon = utilities::calculate_l2_radius(data.vis.size(), std::get<1>(temp));
159 
160  return true;
161 }
162 double duration(std::chrono::high_resolution_clock::time_point start,
163  std::chrono::high_resolution_clock::time_point end,
164  sopt::mpi::Communicator const &comm) {
165  auto elapsed_seconds = duration(start, end);
166  // Now get the max time across all procs: the slowest processor is the one that is
167  // holding back the others in the benchmark.
168  return comm.all_reduce(elapsed_seconds, MPI_MAX);
169 }
170 std::tuple<utilities::vis_params, t_real> dirty_measurements(
171  Image<t_complex> const &ground_truth_image, t_uint number_of_vis, t_real snr,
172  const t_real &cellsize, sopt::mpi::Communicator const &comm) {
173  if (comm.size() == 1) return dirty_measurements(ground_truth_image, number_of_vis, snr, cellsize);
174  if (comm.is_root()) {
175  auto result = dirty_measurements(ground_truth_image, number_of_vis, snr, cellsize);
176  comm.broadcast(std::get<1>(result));
177  auto const order =
179  std::get<0>(result) = utilities::regroup_and_scatter(std::get<0>(result), order, comm);
180  return result;
181  }
182  auto const sigma = comm.broadcast<t_real>();
183  return std::make_tuple(utilities::scatter_visibilities(comm), sigma);
184 }
185 void update_comm(sopt::mpi::Communicator &comm) { comm = sopt::mpi::Communicator::World(); }
186 #endif
187 } // namespace b_utilities
#define PURIFY_INFO(...)
Definition: logging.h:195
bool updateImage(t_uint newSize, Image< t_complex > &image, t_uint &sizex, t_uint &sizey)
Definition: utilities.cc:32
std::tuple< utilities::vis_params, t_real > dirty_measurements(Image< t_complex > const &ground_truth_image, t_uint number_of_vis, t_real snr, const t_real &cellsize)
Definition: utilities.cc:77
bool updateEmptyImage(t_uint newSize, Vector< t_complex > &image, t_uint &sizex, t_uint &sizey)
Definition: utilities.cc:44
void Arguments(benchmark::internal::Benchmark *b)
Definition: utilities.cc:16
bool updateMeasurements(t_uint newSize, utilities::vis_params &data, t_real &epsilon, bool newImage, Image< t_complex > &image)
Definition: utilities.cc:62
double duration(std::chrono::high_resolution_clock::time_point start, std::chrono::high_resolution_clock::time_point end)
Definition: utilities.cc:26
utilities::vis_params random_measurements(t_int size, const t_real max_w, const t_int id, const bool cache_visibilities)
Definition: utilities.cc:96
const t_real pi
mathematical constant
Definition: types.h:70
std::vector< t_int > distribute_measurements(Vector< t_real > const &u, Vector< t_real > const &v, Vector< t_real > const &w, t_int const number_of_nodes, distribute::plan const distribution_plan, t_int const &grid_size)
Distribute visiblities into groups.
Definition: distribute.cc:6
void write_visibility(const utilities::vis_params &uv_vis, const std::string &file_name, const bool w_term)
Writes visibilities to txt.
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
vis_params scatter_visibilities(vis_params const &params, std::vector< t_int > const &sizes, sopt::mpi::Communicator const &comm)
vis_params regroup_and_scatter(vis_params const &params, std::vector< t_int > const &groups, sopt::mpi::Communicator const &comm)
utilities::vis_params read_visibility(const std::vector< std::string > &names, const bool w_term)
Read visibility files from name of vector.
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.
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
std::string visibility_filename(std::string const &filename)
Visibility filename.
Vector< t_complex > vis
Definition: uvw_utilities.h:22