PURIFY
Next-generation radio interferometric imaging
Functions
purify::widefield Namespace Reference

Functions

t_int w_support (const t_real w, const t_real du, const t_int min, const t_int max)
 estimate support size of w given u resolution du More...
 
t_real pixel_to_lambda (const t_real cell, const t_uint imsize, const t_real oversample_ratio)
 return factors to convert between arcsecond pixel size image space and lambda for uv space More...
 
t_real estimate_cell_size (const t_real max_u, const t_uint imsize, const t_real oversample_ratio)
 return cell size from the bandwidth More...
 
t_real fov_cosine (t_real const cell, t_uint const imsize)
 Work out max L and M directional cosines from image parameters. More...
 
t_real equivalent_miriad_cell_size (const t_real cell, const t_uint imsize, const t_real oversample_ratio)
 for a given purify cell size in arcsec provide the equivalent miriad cell size in arcsec More...
 
Matrix< t_complex > generate_chirp (const t_real w_rate, const t_real cell_x, const t_real cell_y, const t_uint x_size, const t_uint y_size)
 Generates image of chirp. More...
 
Matrix< t_complex > estimate_sample_density (const Vector< t_real > &u, const Vector< t_real > &v, const t_real cellx, const t_real celly, const t_uint imsizex, const t_uint imsizey, const t_real oversample_ratio, const t_real scale)
 estimate sample desity grid for a given field of view More...
 
Vector< t_complex > sample_density_weights (const Vector< t_real > &u, const Vector< t_real > &v, const t_real cellx, const t_real celly, const t_uint imsizex, const t_uint imsizey, const t_real oversample_ratio, const t_real scale)
 create sample density weights for a given field of view, uniform weighting More...
 
template<class DDE >
Matrix< t_complex > generate_dde (const DDE &dde, const t_real cell_x, const t_real cell_y, const t_uint x_size, const t_uint y_size, const t_real stop_gap)
 Generate image of DDE for aw-stacking. More...
 
template<class DDE >
Matrix< t_complex > generate_chirp (const DDE &dde, const t_real w_rate, const t_real cell_x, const t_real cell_y, const t_uint x_size, const t_uint y_size, const t_real stop_gap=0.1)
 generates image of chirp and DDE More...
 

Function Documentation

◆ equivalent_miriad_cell_size()

t_real purify::widefield::equivalent_miriad_cell_size ( const t_real  cell,
const t_uint  imsize,
const t_real  oversample_ratio 
)

for a given purify cell size in arcsec provide the equivalent miriad cell size in arcsec

Definition at line 28 of file wide_field_utilities.cc.

29  {
30  const t_real du = pixel_to_lambda(cell, imsize, oversample_ratio);
31  return 60. * 60. * 180. / (static_cast<t_real>(imsize) * oversample_ratio * du * constant::pi);
32 }
const t_real pi
mathematical constant
Definition: types.h:70
t_real pixel_to_lambda(const t_real cell, const t_uint imsize, const t_real oversample_ratio)
return factors to convert between arcsecond pixel size image space and lambda for uv space

References purify::constant::pi, and pixel_to_lambda().

Referenced by getInputData(), and TEST_CASE().

◆ estimate_cell_size()

t_real purify::widefield::estimate_cell_size ( const t_real  max_u,
const t_uint  imsize,
const t_real  oversample_ratio 
)

return cell size from the bandwidth

Definition at line 15 of file wide_field_utilities.cc.

15  {
16  return (2. / static_cast<t_real>(imsize)) *
17  std::asin(static_cast<t_real>(imsize) / (4 * oversample_ratio * max_u)) * 60. * 60. *
18  180. / constant::pi;
19 }

References purify::constant::pi.

Referenced by getInputData().

◆ estimate_sample_density()

Matrix< t_complex > purify::widefield::estimate_sample_density ( const Vector< t_real > &  u,
const Vector< t_real > &  v,
const t_real  cellx,
const t_real  celly,
const t_uint  imsizex,
const t_uint  imsizey,
const t_real  oversample_ratio,
const t_real  scale 
)

estimate sample desity grid for a given field of view

Definition at line 39 of file wide_field_utilities.cc.

42  {
43  const t_int ftsizeu = std::floor(oversample_ratio * imsizex);
44  const t_int ftsizev = std::floor(oversample_ratio * imsizex);
45  Matrix<t_complex> sample_density = Matrix<t_complex>::Zero(ftsizev, ftsizeu);
46  const t_real du = pixel_to_lambda(cellx, imsizex, oversample_ratio) * scale;
47  const t_real dv = pixel_to_lambda(celly, imsizey, oversample_ratio) * scale;
48  for (t_int i = 0; i < u.size(); ++i) {
49  const t_int q = utilities::mod(floor(u(i) * du), ftsizeu);
50  const t_int p = utilities::mod(floor(v(i) * dv), ftsizev);
51  sample_density(p, q) += 1.;
52  }
53  return sample_density;
54 }
const std::vector< t_real > u
data for u coordinate
Definition: operators.cc:18
const std::vector< t_real > v
data for v coordinate
Definition: operators.cc:20
t_real mod(const t_real &x, const t_real &y)
Mod function modified to wrap circularly for negative numbers.
Definition: utilities.cc:41

References purify::utilities::mod(), pixel_to_lambda(), operators_test::u, and operators_test::v.

Referenced by sample_density_weights().

◆ fov_cosine()

t_real purify::widefield::fov_cosine ( t_real const  cell,
t_uint const  imsize 
)

Work out max L and M directional cosines from image parameters.

Definition at line 21 of file wide_field_utilities.cc.

21  {
22  const t_real theta_FoV_L = imsize * cell;
23  const t_real extra_theta = (theta_FoV_L > 180. * 60. * 60.) ? theta_FoV_L - 180. * 60. * 60. : 0.;
24  return 2 * std::sin(constant::pi / 180. * (theta_FoV_L - extra_theta) / (60. * 60.) * 0.5) +
25  2 * std::sin(constant::pi / 180. * extra_theta / (60. * 60.) * 0.5);
26 }

References purify::constant::pi.

Referenced by generate_dde(), and pixel_to_lambda().

◆ generate_chirp() [1/2]

template<class DDE >
Matrix<t_complex> purify::widefield::generate_chirp ( const DDE &  dde,
const t_real  w_rate,
const t_real  cell_x,
const t_real  cell_y,
const t_uint  x_size,
const t_uint  y_size,
const t_real  stop_gap = 0.1 
)

generates image of chirp and DDE

Definition at line 64 of file wide_field_utilities.h.

66  {
67  const t_real nz = y_size * x_size;
68  const t_complex I(0, 1);
69  const auto chirp = [=](const t_real y, const t_real x) {
70  return dde(y, x) *
71  (std::exp(-2 * constant::pi * I * w_rate * (std::sqrt(1 - x * x - y * y) - 1))) /
72  std::sqrt(1 - x * x - y * y) / nz;
73  };
74  return generate_dde(chirp, cell_x, cell_y, x_size, y_size, stop_gap);
75 }
Matrix< t_complex > generate_dde(const DDE &dde, const t_real cell_x, const t_real cell_y, const t_uint x_size, const t_uint y_size, const t_real stop_gap)
Generate image of DDE for aw-stacking.

References generate_dde(), purify::I, and purify::constant::pi.

◆ generate_chirp() [2/2]

Matrix< t_complex > purify::widefield::generate_chirp ( const t_real  w_rate,
const t_real  cell_x,
const t_real  cell_y,
const t_uint  x_size,
const t_uint  y_size 
)

Generates image of chirp.

Definition at line 34 of file wide_field_utilities.cc.

35  {
36  return generate_chirp([](t_real, t_real) { return 1.; }, w_rate, cell_x, cell_y, x_size, y_size);
37 }
Matrix< t_complex > generate_chirp(const t_real w_rate, const t_real cell_x, const t_real cell_y, const t_uint x_size, const t_uint y_size)
Generates image of chirp.

Referenced by purify::details::init_correction2d(), and purify::details::init_correction_radial_2d().

◆ generate_dde()

template<class DDE >
Matrix<t_complex> purify::widefield::generate_dde ( const DDE &  dde,
const t_real  cell_x,
const t_real  cell_y,
const t_uint  x_size,
const t_uint  y_size,
const t_real  stop_gap 
)

Generate image of DDE for aw-stacking.

Definition at line 43 of file wide_field_utilities.h.

44  {
45  assert(stop_gap <= 1);
46  const t_real L = fov_cosine(cell_x, x_size);
47  const t_real M = fov_cosine(cell_y, y_size);
48 
49  const t_real delt_x = L / x_size;
50  const t_real delt_y = M / y_size;
51  Image<t_complex> output = Image<t_complex>::Zero(y_size, x_size);
52 
53  for (t_int l = 0; l < x_size; ++l)
54  for (t_int m = 0; m < y_size; ++m) {
55  const t_real x = (l - x_size * 0.5) * delt_x;
56  const t_real y = (m - y_size * 0.5) * delt_y;
57  output(m, l) = ((x * x + y * y) < 1 - stop_gap) ? dde(y, x) : 0.;
58  }
59 
60  return output;
61 };
t_real fov_cosine(t_real const cell, t_uint const imsize)
Work out max L and M directional cosines from image parameters.

References fov_cosine().

Referenced by generate_chirp(), and TEST_CASE().

◆ pixel_to_lambda()

t_real purify::widefield::pixel_to_lambda ( const t_real  cell,
const t_uint  imsize,
const t_real  oversample_ratio 
)

return factors to convert between arcsecond pixel size image space and lambda for uv space

Definition at line 11 of file wide_field_utilities.cc.

11  {
12  return 1. / (oversample_ratio * fov_cosine(cell, imsize));
13 }

References fov_cosine().

Referenced by purify::utilities::convert_to_pixels(), equivalent_miriad_cell_size(), estimate_sample_density(), getInputData(), purify::details::init_gridding_matrix_2d(), main(), sample_density_weights(), DegridOperatorCtorFixturePar::SetUp(), and TEST_CASE().

◆ sample_density_weights()

Vector< t_complex > purify::widefield::sample_density_weights ( const Vector< t_real > &  u,
const Vector< t_real > &  v,
const t_real  cellx,
const t_real  celly,
const t_uint  imsizex,
const t_uint  imsizey,
const t_real  oversample_ratio,
const t_real  scale 
)

create sample density weights for a given field of view, uniform weighting

Definition at line 56 of file wide_field_utilities.cc.

59  {
60  Vector<t_complex> weights = Vector<t_complex>::Zero(u.size());
61  const t_int ftsizeu = std::floor(oversample_ratio * imsizex);
62  const t_int ftsizev = std::floor(oversample_ratio * imsizex);
63  const Matrix<t_complex> sample_density =
64  estimate_sample_density(u, v, cellx, celly, imsizex, imsizey, oversample_ratio, scale);
65  const t_real du = pixel_to_lambda(cellx, imsizex, oversample_ratio) * scale;
66  const t_real dv = pixel_to_lambda(celly, imsizey, oversample_ratio) * scale;
67  for (t_int i = 0; i < u.size(); ++i) {
68  const t_int q = utilities::mod(floor(u(i) * du), ftsizeu);
69  const t_int p = utilities::mod(floor(v(i) * dv), ftsizev);
70  weights(i) = 1. / sample_density(p, q);
71  }
72  const t_real max_weight = weights.array().cwiseAbs().maxCoeff();
73  weights /= max_weight;
74  return weights;
75 }
Matrix< t_complex > estimate_sample_density(const Vector< t_real > &u, const Vector< t_real > &v, const t_real cellx, const t_real celly, const t_uint imsizex, const t_uint imsizey, const t_real oversample_ratio, const t_real scale)
estimate sample desity grid for a given field of view

References estimate_sample_density(), purify::utilities::mod(), pixel_to_lambda(), operators_test::u, and operators_test::v.

Referenced by main(), and TEST_CASE().

◆ w_support()

t_int purify::widefield::w_support ( const t_real  w,
const t_real  du,
const t_int  min,
const t_int  max 
)

estimate support size of w given u resolution du

Definition at line 7 of file wide_field_utilities.cc.

7  {
8  return std::min<int>(std::max<int>(static_cast<int>(2 * std::floor(std::abs(w / du))), min), max);
9 }

Referenced by purify::details::init_gridding_matrix_2d(), main(), TEST_CASE(), and purify::utilities::w_stacking_with_all_to_all().