PURIFY
Next-generation radio interferometric imaging
Namespaces | Classes | Typedefs | Enumerations | Functions | Variables
purify Namespace Reference

Namespaces

 casa
 
 constant
 
 convol
 
 details
 
 distribute
 
 factory
 
 H5
 
 integration
 
 kernels
 
 logging
 
 measurementoperator
 
 operators
 
 pfitsio
 
 projection_kernels
 
 random_updater
 
 read_measurements
 
 utilities
 
 widefield
 
 wproj_utilities
 

Classes

class  IndexMapping
 
class  YamlParser
 

Typedefs

typedef std::complex< float > t_complexf
 
typedef Eigen::Triplet< t_complex > t_tripletList
 Root of the type hierarchy for triplet lists. More...
 
template<class T = t_real, class I = t_int>
using Sparse = Eigen::SparseMatrix< T, Eigen::RowMajor, I >
 A matrix of a given type. More...
 
template<class T = t_real, class I = t_int>
using SparseVector = Eigen::SparseVector< T, Eigen::RowMajor, I >
 

Enumerations

enum class  diff_func_type { L2Norm , L2Norm_with_CRR }
 
enum class  nondiff_func_type { L1Norm , Denoiser , RealIndicator }
 
enum class  stokes {
  I , Q , U , V ,
  XX , YY , XY , YX ,
  LL , RR , LR , RL ,
  P
}
 
enum class  dde_type { wkernel_radial , wkernel_2d }
 Types of DDEs in purify. More...
 

Functions

std::vector< t_int > all_to_all_send_sizes (const std::vector< t_int > &recv_sizes, const sopt::mpi::Communicator &comm)
 
std::string version ()
 Returns library version. More...
 
std::tuple< uint8_t, uint8_t, uint8_t > version_tuple ()
 Returns library version. More...
 
std::string gitref ()
 Returns library git reference, if known. More...
 
std::string default_logging_level ()
 Default logging level. More...
 
template<class T0 , class STORAGE_INDEX_TYPE = t_int>
std::set< STORAGE_INDEX_TYPE > non_empty_outers (Eigen::SparseMatrixBase< T0 > const &matrix)
 Indices of non empty outer indices. More...
 
template<class T0 >
Sparse< typename T0::Scalar > compress_outer (T0 const &matrix)
 Indices of non empty outer indices. More...
 
std::tuple< std::function< t_real(t_real)>, std::function< t_real(t_real)>, std::function< t_real(t_real)>, std::function< t_real(t_real)> > create_kernels (const kernels::kernel kernel_name_, const t_uint Ju_, const t_uint Jv_, const t_real imsizey_, const t_real imsizex_, const t_real oversample_ratio)
 
std::tuple< std::function< t_real(t_real)>, std::function< t_real(t_real)> > create_radial_ftkernel (const kernels::kernel kernel_name_, const t_uint Ju_, const t_real oversample_ratio)
 
std::vector< std::string > split (const std::string &s, char delim)
 splits string into vector More...
 
void mkdir_recursive (const std::string &path)
 recursively create directories when they do not exist More...
 
template<class T >
void mkdir_recursive (const T &path)
 boost wrapper More...
 
template<typename T >
void split (const std::string &s, char delim, T result)
 adds split string to container More...
 
template<typename T >
get_vector (const YAML::Node &node_map, const std::initializer_list< const char * > indicies)
 
template<class T >
std::enable_if< std::is_scalar< T >::value, std::vector< T > >::type read_data (const std::string &filename)
 read real values from data file More...
 
template<class T >
std::enable_if< std::is_same< t_complex, T >::value, std::vector< T > >::type read_data (const std::string &filename)
 read complex values from data file More...
 
std::vector< t_complex > read_data (const std::vector< t_real > &input)
 
std::string data_directory ()
 Holds data and such. More...
 
std::string models_directory ()
 Holds TF models. More...
 
std::string output_directory ()
 Where test outputs go. More...
 
std::string data_filename (std::string const &filename)
 Holds data and such. More...
 
std::string image_filename (std::string const &filename)
 Image filename. More...
 
std::string visibility_filename (std::string const &filename)
 Visibility filename. More...
 
std::string vla_filename (std::string const &filename)
 Specific vla data. More...
 
std::string atca_filename (std::string const &filename)
 Specific atca data. More...
 
std::string mwa_filename (std::string const &filename)
 Specific mwa data. More...
 
std::string gridding_filename (std::string const &filename)
 Some gridding regression data. More...
 
std::string degridding_filename (std::string const &filename)
 Some degridding regression data. More...
 
std::string output_filename (std::string const &filename)
 Test output file. More...
 
std::string ngc3256_ms ()
 

Variables

const std::map< std::string, diff_func_typediff_type_string
 
const std::map< std::string, nondiff_func_typenondiff_type_string
 
const std::map< stokes, t_int > stokes_int
 
const std::map< std::string, stokesstokes_string
 

Detailed Description

YamlParser definition for the purify project

Author
Roland Guichard
Version
1.0

YamlParser header for purify project

Author
Roland Guichard
Version
1.0 The YamlParser class definition.

It browses and searches through the config.yaml file in the data folder and sets the class variable members accordingly, to be used by the purify algorithm.

Parameters
filepathpath to the config file

Typedef Documentation

◆ Sparse

template<class T = t_real, class I = t_int>
using purify::Sparse = typedef Eigen::SparseMatrix<T, Eigen::RowMajor, I>

A matrix of a given type.

Operates as mathematical sparse matrix.

Definition at line 40 of file types.h.

◆ SparseVector

template<class T = t_real, class I = t_int>
using purify::SparseVector = typedef Eigen::SparseVector<T, Eigen::RowMajor, I>

Definition at line 42 of file types.h.

◆ t_complexf

typedef std::complex<float> purify::t_complexf

Definition at line 21 of file types.h.

◆ t_tripletList

typedef Eigen::Triplet<t_complex> purify::t_tripletList

Root of the type hierarchy for triplet lists.

Definition at line 23 of file types.h.

Enumeration Type Documentation

◆ dde_type

enum purify::dde_type
strong

Types of DDEs in purify.

Enumerator
wkernel_radial 
wkernel_2d 

Definition at line 59 of file types.h.

◆ diff_func_type

Enumerator
L2Norm 
L2Norm_with_CRR 

Definition at line 26 of file types.h.

◆ nondiff_func_type

Enumerator
L1Norm 
Denoiser 
RealIndicator 

Definition at line 31 of file types.h.

◆ stokes

enum purify::stokes
strong
Enumerator
XX 
YY 
XY 
YX 
LL 
RR 
LR 
RL 

Definition at line 44 of file types.h.

Function Documentation

◆ all_to_all_send_sizes()

std::vector<t_int> purify::all_to_all_send_sizes ( const std::vector< t_int > &  recv_sizes,
const sopt::mpi::Communicator &  comm 
)

Definition at line 5 of file AllToAllSparseVector.cc.

6  {
7  std::vector<t_int> send_sizes(recv_sizes.size(), 0);
8  for (t_int i = 0; i < comm.size(); i++) {
9  if (i == comm.rank())
10  send_sizes = comm.gather<t_int>(recv_sizes[i], i);
11  else
12  comm.gather<t_int>(recv_sizes[i], i);
13  }
14  return send_sizes;
15 }

◆ atca_filename()

std::string purify::atca_filename ( std::string const &  filename)
inline

Specific atca data.

Definition at line 32 of file directories.in.h.

32  {
33  return data_filename("atca/" + filename);
34 }
std::string data_filename(std::string const &filename)
Holds data and such.

References data_filename().

Referenced by TEST_CASE().

◆ compress_outer()

template<class T0 >
Sparse<typename T0::Scalar> purify::compress_outer ( T0 const &  matrix)

Indices of non empty outer indices.

Definition at line 71 of file IndexMapping.h.

71  {
72  static_assert(T0::IsRowMajor, "Not tested for col major");
73  auto const indices = non_empty_outers(matrix);
74 
75  SparseVector<typename T0::Index, typename T0::Index> mapping(matrix.cols());
76  mapping.reserve(indices.size());
77  t_int i(0);
78  for (auto const &index : indices) mapping.coeffRef(index) = i++;
79 
80  std::vector<typename Sparse<typename T0::Scalar>::Index> rows(matrix.rows() + 1, 0);
81  std::vector<typename Sparse<typename T0::Scalar>::Index> cols(matrix.nonZeros(), 0);
82  rows[matrix.rows()] = matrix.nonZeros();
83  t_int index = 0;
84  for (typename T0::Index k = 0; k < matrix.outerSize(); ++k) {
85  rows[k] = index;
86  for (typename T0::InnerIterator it(matrix, k); it; ++it) {
87  cols[index] = mapping.coeff(it.col());
88  index++;
89  }
90  }
91  return Eigen::MappedSparseMatrix<typename T0::Scalar, Eigen::RowMajor,
92  typename Sparse<typename T0::Scalar>::Index>(
93  matrix.rows(), indices.size(), matrix.nonZeros(), rows.data(), cols.data(),
94  const_cast<typename T0::Scalar *>(matrix.derived().valuePtr()));
95 }
std::set< STORAGE_INDEX_TYPE > non_empty_outers(Eigen::SparseMatrixBase< T0 > const &matrix)
Indices of non empty outer indices.
Definition: IndexMapping.h:61

References non_empty_outers().

Referenced by TEST_CASE().

◆ create_kernels()

std::tuple< std::function< t_real(t_real)>, std::function< t_real(t_real)>, std::function< t_real(t_real)>, std::function< t_real(t_real)> > purify::create_kernels ( const kernels::kernel  kernel_name_,
const t_uint  Ju_,
const t_uint  Jv_,
const t_real  imsizey_,
const t_real  imsizex_,
const t_real  oversample_ratio 
)

Definition at line 249 of file kernels.cc.

250  {
251  // PURIFY_MEDIUM_LOG("Kernel Name: {}", kernel_name_.c_str());
252  PURIFY_MEDIUM_LOG("Kernel Support: {} x {}", Ju_, Jv_);
253  const t_uint ftsizev_ = std::floor(imsizey_ * oversample_ratio);
254  const t_uint ftsizeu_ = std::floor(imsizex_ * oversample_ratio);
255  if ((kernel_name_ == kernels::kernel::pswf) and (Ju_ != 6 or Jv_ != 6)) {
256  PURIFY_ERROR("Error: Only a support of 6 is implemented for PSWFs.");
257  throw std::runtime_error("Incorrect input: PSWF requires a support of 6");
258  }
259  switch (kernel_name_) {
260  case kernels::kernel::kb: {
261  auto kbu = [=](const t_real x) { return kernels::kaiser_bessel(x, Ju_); };
262  auto kbv = [=](const t_real x) { return kernels::kaiser_bessel(x, Jv_); };
263  auto ftkbu = [=](const t_real x) { return kernels::ft_kaiser_bessel(x / ftsizeu_ - 0.5, Ju_); };
264  auto ftkbv = [=](const t_real x) { return kernels::ft_kaiser_bessel(x / ftsizev_ - 0.5, Jv_); };
265  return std::make_tuple(kbu, kbv, ftkbu, ftkbv);
266  break;
267  }
268  case kernels::kernel::kb_presample: {
269  if (Ju_ != Jv_)
270  throw std::runtime_error("Ju and Jv must be the same support size for presampling kernels.");
271  const auto samples = kernels::kernel_samples(
272  1e5, [&](const t_real x) { return kernels::kaiser_bessel(x, Ju_); });
273  auto kbu = [=](const t_real x) { return kernels::kernel_zero_interp(samples, x, Ju_); };
274  auto kbv = [=](const t_real x) { return kernels::kernel_zero_interp(samples, x, Jv_); };
275  auto ftkbu = [=](const t_real x) { return kernels::ft_kaiser_bessel(x / ftsizeu_ - 0.5, Ju_); };
276  auto ftkbv = [=](const t_real x) { return kernels::ft_kaiser_bessel(x / ftsizev_ - 0.5, Jv_); };
277  return std::make_tuple(kbu, kbv, ftkbu, ftkbv);
278  break;
279  }
280  case kernels::kernel::kbmin: {
281  const t_real kb_interp_alpha_Ju =
282  constant::pi * std::sqrt(Ju_ * Ju_ / (oversample_ratio * oversample_ratio) *
283  (oversample_ratio - 0.5) * (oversample_ratio - 0.5) -
284  0.8);
285  const t_real kb_interp_alpha_Jv =
286  constant::pi * std::sqrt(Jv_ * Jv_ / (oversample_ratio * oversample_ratio) *
287  (oversample_ratio - 0.5) * (oversample_ratio - 0.5) -
288  0.8);
289  auto kbu = [=](const t_real x) {
290  return kernels::kaiser_bessel_general(x, Ju_, kb_interp_alpha_Ju);
291  };
292  auto kbv = [=](const t_real x) {
293  return kernels::kaiser_bessel_general(x, Jv_, kb_interp_alpha_Jv);
294  };
295  auto ftkbu = [=](const t_real x) {
296  return kernels::ft_kaiser_bessel_general(x / ftsizeu_ - 0.5, Ju_, kb_interp_alpha_Ju);
297  };
298  auto ftkbv = [=](const t_real x) {
299  return kernels::ft_kaiser_bessel_general(x / ftsizev_ - 0.5, Jv_, kb_interp_alpha_Jv);
300  };
301  return std::make_tuple(kbu, kbv, ftkbu, ftkbv);
302  break;
303  }
304  case kernels::kernel::pswf: {
305  auto pswfu = [=](const t_real x) { return kernels::pswf(x, Ju_); };
306  auto pswfv = [=](const t_real x) { return kernels::pswf(x, Jv_); };
307  auto ftpswfu = [=](const t_real x) { return kernels::ft_pswf(x / ftsizeu_ - 0.5, Ju_); };
308  auto ftpswfv = [=](const t_real x) { return kernels::ft_pswf(x / ftsizev_ - 0.5, Jv_); };
309  return std::make_tuple(pswfu, pswfv, ftpswfu, ftpswfv);
310  break;
311  }
312  case kernels::kernel::gauss: {
313  auto gaussu = [=](const t_real x) { return kernels::gaussian(x, Ju_); };
314  auto gaussv = [=](const t_real x) { return kernels::gaussian(x, Jv_); };
315  auto ftgaussu = [=](const t_real x) { return kernels::ft_gaussian(x / ftsizeu_ - 0.5, Ju_); };
316  auto ftgaussv = [=](const t_real x) { return kernels::ft_gaussian(x / ftsizev_ - 0.5, Jv_); };
317  return std::make_tuple(gaussu, gaussv, ftgaussu, ftgaussv);
318  break;
319  }
320  case kernels::kernel::box: {
321  auto boxu = [=](const t_real x) { return kernels::pill_box(x, Ju_); };
322  auto boxv = [=](const t_real x) { return kernels::pill_box(x, Jv_); };
323  auto ftboxu = [=](const t_real x) { return kernels::ft_pill_box(x / ftsizeu_ - 0.5, Ju_); };
324  auto ftboxv = [=](const t_real x) { return kernels::ft_pill_box(x / ftsizev_ - 0.5, Jv_); };
325  return std::make_tuple(boxu, boxv, ftboxu, ftboxv);
326  break;
327  }
328  case kernels::kernel::gauss_alt: {
329  const t_real sigma = 1; // In units of radians, Rafael uses sigma = 2 * pi / ftsizeu_. However,
330  // this should be 1 in units of pixels.
331  auto gaussu = [=](const t_real x) { return kernels::gaussian_general(x, Ju_, sigma); };
332  auto gaussv = [=](const t_real x) { return kernels::gaussian_general(x, Jv_, sigma); };
333  auto ftgaussu = [=](const t_real x) {
334  return kernels::ft_gaussian_general(x / ftsizeu_ - 0.5, Ju_, sigma);
335  };
336  auto ftgaussv = [=](const t_real x) {
337  return kernels::ft_gaussian_general(x / ftsizev_ - 0.5, Jv_, sigma);
338  };
339  return std::make_tuple(gaussu, gaussv, ftgaussu, ftgaussv);
340  break;
341  }
342  default:
343  throw std::runtime_error("Did not choose valid kernel.");
344  }
345 }
#define PURIFY_ERROR(...)
\macro Something is definitely wrong, algorithm exits
Definition: logging.h:190
#define PURIFY_MEDIUM_LOG(...)
Medium priority message.
Definition: logging.h:205
const t_real pi
mathematical constant
Definition: types.h:70
t_real ft_kaiser_bessel(const t_real x, const t_real J)
Fourier transform of kaiser bessel kernel.
Definition: kernels.cc:40
t_real ft_gaussian_general(const t_real x, const t_real J, const t_real sigma)
Fourier transform of general Gaussian kernel.
Definition: kernels.cc:233
std::vector< t_real > kernel_samples(const t_int total_samples, const std::function< t_real(t_real)> kernelu)
Calculates samples of a kernel.
Definition: kernels.cc:144
t_real kaiser_bessel(const t_real x, const t_real J)
Kaiser-Bessel kernel.
Definition: kernels.cc:7
t_real kernel_zero_interp(const std::vector< t_real > &samples, const t_real x, const t_real J)
zeroth order interpolates from samples of kernel
Definition: kernels.cc:156
t_real pill_box(const t_real x, const t_real J)
Box car function for kernel.
Definition: kernels.cc:202
t_real ft_pswf(const t_real x, const t_real J)
Fourier transform of PSWF kernel.
Definition: kernels.cc:124
t_real gaussian_general(const t_real x, const t_real J, const t_real sigma)
Fourier transform of general Gaussian kernel.
Definition: kernels.cc:220
t_real kaiser_bessel_general(const t_real x, const t_real J, const t_real alpha)
More general Kaiser-Bessel kernel.
Definition: kernels.cc:15
t_real ft_pill_box(const t_real x, const t_real J)
Fourier transform of box car function, a Sinc function.
Definition: kernels.cc:211
t_real gaussian(const t_real x, const t_real J)
Gaussian kernel.
Definition: kernels.cc:49
t_real pswf(const t_real x, const t_real J)
PSWF kernel.
Definition: kernels.cc:107
t_real ft_gaussian(const t_real x, const t_real J)
Fourier transform of Gaussian kernel.
Definition: kernels.cc:60
t_real ft_kaiser_bessel_general(const t_real x, const t_real J, const t_real alpha)
Fourier transform of more general Kaiser-Bessel kernel.
Definition: kernels.cc:25

References purify::kernels::box, purify::kernels::ft_gaussian(), purify::kernels::ft_gaussian_general(), purify::kernels::ft_kaiser_bessel(), purify::kernels::ft_kaiser_bessel_general(), purify::kernels::ft_pill_box(), purify::kernels::ft_pswf(), purify::kernels::gauss, purify::kernels::gauss_alt, purify::kernels::gaussian(), purify::kernels::gaussian_general(), purify::kernels::kaiser_bessel(), purify::kernels::kaiser_bessel_general(), purify::kernels::kb, purify::kernels::kb_presample, purify::kernels::kbmin, purify::kernels::kernel_samples(), purify::kernels::kernel_zero_interp(), purify::constant::pi, purify::kernels::pill_box(), purify::kernels::pswf, PURIFY_ERROR, and PURIFY_MEDIUM_LOG.

Referenced by purify::operators::base_degrid_operator_2d(), GridOperatorFixture::SetUp(), and TEST_CASE().

◆ create_radial_ftkernel()

std::tuple< std::function< t_real(t_real)>, std::function< t_real(t_real)> > purify::create_radial_ftkernel ( const kernels::kernel  kernel_name_,
const t_uint  Ju_,
const t_real  oversample_ratio 
)

Definition at line 347 of file kernels.cc.

348  {
349  // PURIFY_MEDIUM_LOG("Kernel Name: {}", kernel_name_.c_str());
350  PURIFY_MEDIUM_LOG("Radial Kernel Support: {}", Ju_);
351  if ((kernel_name_ == kernels::kernel::pswf) and (Ju_ != 6)) {
352  PURIFY_ERROR("Error: Only a support of 6 is implemented for PSWFs.");
353  throw std::runtime_error("Incorrect input: PSWF requires a support of 6");
354  }
355  switch (kernel_name_) {
356  case kernels::kernel::kb: {
357  return std::make_tuple(
358  [=](const t_real x) { return kernels::ft_kaiser_bessel(x, static_cast<t_real>(Ju_)); },
359  [=](const t_real x) { return kernels::kaiser_bessel(x, static_cast<t_real>(Ju_)); });
360  break;
361  }
362  case kernels::kernel::kbmin: {
363  const t_real kb_interp_alpha_Ju =
364  constant::pi * std::sqrt(Ju_ * Ju_ / (oversample_ratio * oversample_ratio) *
365  (oversample_ratio - 0.5) * (oversample_ratio - 0.5) -
366  0.8);
367  auto kbuv = [=](const t_real x) {
368  return kernels::kaiser_bessel_general(x, Ju_, kb_interp_alpha_Ju);
369  };
370  auto ftkbuv = [=](const t_real x) {
371  return kernels::ft_kaiser_bessel_general(x, Ju_, kb_interp_alpha_Ju);
372  };
373  return std::make_tuple(ftkbuv, kbuv);
374  break;
375  }
376  case kernels::kernel::pswf: {
377  return std::make_tuple(
378  [=](const t_real x) { return kernels::ft_pswf(x, static_cast<t_real>(Ju_)); },
379  [=](const t_real x) { return kernels::pswf(x, Ju_); });
380  break;
381  }
382  case kernels::kernel::gauss: {
383  return std::make_tuple(
384  [=](const t_real x) { return kernels::ft_gaussian(x, static_cast<t_real>(Ju_)); },
385  [=](const t_real x) { return kernels::gaussian(x, Ju_); });
386  break;
387  }
388  case kernels::kernel::box: {
389  return std::make_tuple(
390  [=](const t_real x) { return kernels::ft_pill_box(x, static_cast<t_real>(Ju_)); },
391  [=](const t_real x) { return kernels::pill_box(x, Ju_); });
392  break;
393  }
394  default:
395  throw std::runtime_error("Did not choose valid radial kernel.");
396  }
397 }

References purify::kernels::box, purify::kernels::ft_gaussian(), purify::kernels::ft_kaiser_bessel(), purify::kernels::ft_kaiser_bessel_general(), purify::kernels::ft_pill_box(), purify::kernels::ft_pswf(), purify::kernels::gauss, purify::kernels::gaussian(), purify::kernels::kaiser_bessel(), purify::kernels::kaiser_bessel_general(), purify::kernels::kb, purify::kernels::kbmin, purify::constant::pi, purify::kernels::pill_box(), purify::kernels::pswf, PURIFY_ERROR, and PURIFY_MEDIUM_LOG.

Referenced by purify::operators::base_degrid_operator_2d().

◆ data_directory()

std::string purify::data_directory ( )
inline

Holds data and such.

Definition at line 9 of file directories.in.h.

9 { return "@CMAKE_INSTALL_PREFIX@/data"; }

Referenced by data_filename().

◆ data_filename()

std::string purify::data_filename ( std::string const &  filename)
inline

Holds data and such.

Definition at line 16 of file directories.in.h.

16  {
17  return data_directory() + "/" + filename;
18 }
std::string data_directory()
Holds data and such.
Definition: directories.in.h:9

References data_directory().

Referenced by atca_filename(), degridding_filename(), gridding_filename(), image_filename(), mwa_filename(), StochasticAlgoFixture::SetUp(), TEST_CASE(), visibility_filename(), and vla_filename().

◆ default_logging_level()

std::string purify::default_logging_level ( )
inline

Default logging level.

Definition at line 51 of file config.in.h.

51 { return "@PURIFY_TEST_LOG_LEVEL@"; }

Referenced by main().

◆ degridding_filename()

std::string purify::degridding_filename ( std::string const &  filename)
inline

Some degridding regression data.

Definition at line 44 of file directories.in.h.

44  {
45  return data_filename("expected/degridding/" + filename);
46 }

References data_filename().

◆ get_vector()

template<typename T >
T purify::get_vector ( const YAML::Node &  node_map,
const std::initializer_list< const char * >  indicies 
)

Definition at line 74 of file yaml-parser.cc.

74  {
75  YAML::Node node = YAML::Clone(node_map);
76  std::string faulty_variable;
77  for (const char* index : indicies) {
78  faulty_variable = std::string(index);
79  node = node[index];
80  if (!node.IsDefined()) {
81  throw std::runtime_error("The initialisation of " + faulty_variable + " is wrong.");
82  }
83  }
84  try {
85  T output;
86  for (int i = 0; i < node.size(); i++) output.push_back(node[i].as<typename T::value_type>());
87  return output;
88  } catch (std::exception e) {
89  throw std::runtime_error("There is a mismatch in the type conversion of " + faulty_variable);
90  }
91 }

◆ gitref()

std::string purify::gitref ( )
inline

Returns library git reference, if known.

Definition at line 49 of file config.in.h.

49 { return "@Purify_GITREF@"; }

◆ gridding_filename()

std::string purify::gridding_filename ( std::string const &  filename)
inline

Some gridding regression data.

Definition at line 40 of file directories.in.h.

40  {
41  return data_filename("expected/gridding/" + filename);
42 }

References data_filename().

◆ image_filename()

std::string purify::image_filename ( std::string const &  filename)
inline

Image filename.

Definition at line 20 of file directories.in.h.

20  {
21  return data_filename("images/" + filename);
22 }

References data_filename().

Referenced by main().

◆ mkdir_recursive() [1/2]

void purify::mkdir_recursive ( const std::string &  path)

recursively create directories when they do not exist

Definition at line 151 of file read_measurements.cc.

151  {
152  if (read_measurements::dir_exists(path)) return;
153  const auto folders = split(path, '/');
154  std::string current_path = "";
155  for (const auto f : folders) {
156  if (f == "") {
157  current_path += "/";
158  continue;
159  }
160  current_path += f;
161  if (not read_measurements::dir_exists(current_path)) {
162  const t_int status = mkdir(current_path.c_str(), ACCESSPERMS);
163  if (status != 0)
164  throw std::runtime_error("Error making recursive directory: " + current_path);
165  }
166  current_path += "/";
167  }
168  return;
169 }
bool dir_exists(const std::string &path)
check that directory path exists
std::vector< std::string > split(const std::string &s, char delim)
splits string into vector

References purify::read_measurements::dir_exists(), and split().

Referenced by main(), mkdir_recursive(), TEST_CASE(), and purify::YamlParser::writeOutput().

◆ mkdir_recursive() [2/2]

template<class T >
void purify::mkdir_recursive ( const T &  path)

boost wrapper

Definition at line 45 of file read_measurements.h.

45  {
46  mkdir_recursive(path.native());
47 }
void mkdir_recursive(const T &path)
boost wrapper

References mkdir_recursive().

◆ models_directory()

std::string purify::models_directory ( )
inline

Holds TF models.

Definition at line 11 of file directories.in.h.

11 { return "@CMAKE_INSTALL_PREFIX@/models"; }

◆ mwa_filename()

std::string purify::mwa_filename ( std::string const &  filename)
inline

Specific mwa data.

Definition at line 36 of file directories.in.h.

36  {
37  return data_filename("mwa/" + filename);
38 }

References data_filename().

Referenced by main(), and TEST_CASE().

◆ ngc3256_ms()

std::string purify::ngc3256_ms ( )
inline

Definition at line 53 of file directories.in.h.

53 { return "@NGC3256_MS@"; }

Referenced by main().

◆ non_empty_outers()

template<class T0 , class STORAGE_INDEX_TYPE = t_int>
std::set<STORAGE_INDEX_TYPE> purify::non_empty_outers ( Eigen::SparseMatrixBase< T0 > const &  matrix)

Indices of non empty outer indices.

Definition at line 61 of file IndexMapping.h.

61  {
62  std::set<STORAGE_INDEX_TYPE> result;
63  for (typename T0::Index k = 0; k < matrix.derived().outerSize(); ++k)
64  for (typename T0::InnerIterator it(matrix.derived(), k); it; ++it)
65  result.insert(static_cast<STORAGE_INDEX_TYPE>(it.col()));
66  return result;
67 }

Referenced by compress_outer(), and TEST_CASE().

◆ output_directory()

std::string purify::output_directory ( )
inline

Where test outputs go.

Definition at line 13 of file directories.in.h.

13 { return "@CMAKE_INSTALL_PREFIX@/outputs"; }

Referenced by output_filename().

◆ output_filename()

std::string purify::output_filename ( std::string const &  filename)
inline

Test output file.

Definition at line 49 of file directories.in.h.

49  {
50  return output_directory() + "/" + filename;
51 }
std::string output_directory()
Where test outputs go.

References output_directory().

Referenced by main(), padmm(), and TEST_CASE().

◆ read_data() [1/3]

template<class T >
std::enable_if<std::is_scalar<T>::value, std::vector<T> >::type purify::read_data ( const std::string &  filename)

read real values from data file

Definition at line 14 of file data.in.h.

15  {
16  std::ifstream data_file(filename);
17  if (!data_file) throw std::runtime_error("Test data is missing: " + filename);
18  std::string s;
19  // read data
20  std::vector<T> data;
21  for (std::string s; std::getline(data_file, s, ',');) {
22  data.push_back(std::stod(s));
23  }
24  return data;
25 }

◆ read_data() [2/3]

template<class T >
std::enable_if<std::is_same<t_complex, T>::value, std::vector<T> >::type purify::read_data ( const std::string &  filename)

read complex values from data file

Definition at line 28 of file data.in.h.

29  {
30  std::ifstream data_file(filename);
31  if (!data_file) throw std::runtime_error("Test data is missing: " + filename);
32  std::string real = "";
33  std::string imag = "";
34  // read data
35  std::vector<T> data;
36  for (std::string real; std::getline(data_file, real, ',');) {
37  real.erase(0, 1);
38  if (!std::getline(data_file, imag, ')')) break;
39  data.push_back(T(std::stod(real), std::stod(imag)));
40  real = "";
41  imag = "";
42  }
43  return data;
44 }

◆ read_data() [3/3]

std::vector<t_complex> purify::read_data ( const std::vector< t_real > &  input)

Definition at line 46 of file data.in.h.

46  {
47  std::vector<t_complex> output;
48  for (auto a = input.begin(); a != input.end(); a++) output.push_back(t_complex(*a, 0.));
49  return output;
50 }

◆ split() [1/2]

std::vector< std::string > purify::split ( const std::string &  s,
char  delim 
)

splits string into vector

Definition at line 146 of file read_measurements.cc.

146  {
147  std::vector<std::string> elems;
148  split(s, delim, std::back_inserter(elems));
149  return elems;
150 }

Referenced by mkdir_recursive().

◆ split() [2/2]

template<typename T >
void purify::split ( const std::string &  s,
char  delim,
result 
)

adds split string to container

Definition at line 50 of file read_measurements.h.

50  {
51  std::stringstream ss(s);
52  std::string item;
53  while (std::getline(ss, item, delim)) {
54  *(result++) = item;
55  }
56 };

◆ version()

std::string purify::version ( )
inline

Returns library version.

Definition at line 40 of file config.in.h.

40 { return "@Purify_VERSION@"; }

Referenced by main(), and TEST_CASE().

◆ version_tuple()

std::tuple<uint8_t, uint8_t, uint8_t> purify::version_tuple ( )
inline

Returns library version.

Definition at line 42 of file config.in.h.

42  {
43  // clang-format off
44  return std::tuple<uint8_t, uint8_t, uint8_t>(
45  @Purify_VERSION_MAJOR@, @Purify_VERSION_MINOR@, @Purify_VERSION_PATCH@);
46  // clang-format on
47 }

◆ visibility_filename()

std::string purify::visibility_filename ( std::string const &  filename)
inline

Visibility filename.

Definition at line 24 of file directories.in.h.

24  {
25  return data_filename("vis_" + filename);
26 }

References data_filename().

Referenced by b_utilities::random_measurements().

◆ vla_filename()

std::string purify::vla_filename ( std::string const &  filename)
inline

Specific vla data.

Definition at line 28 of file directories.in.h.

28  {
29  return data_filename("vla/" + filename);
30 }

References data_filename().

Referenced by main(), and TEST_CASE().

Variable Documentation

◆ diff_type_string

const std::map<std::string, diff_func_type> purify::diff_type_string
Initial value:
= {
{"l2", diff_func_type::L2Norm}, {"CRR", diff_func_type::L2Norm_with_CRR}}

Definition at line 27 of file types.h.

Referenced by purify::YamlParser::parseAndSetAlgorithmOptions().

◆ nondiff_type_string

const std::map<std::string, nondiff_func_type> purify::nondiff_type_string
Initial value:
= {
{"l1", nondiff_func_type::L1Norm},
{"denoiser", nondiff_func_type::Denoiser},
{"realIndicator", nondiff_func_type::RealIndicator}}

Definition at line 32 of file types.h.

Referenced by purify::YamlParser::parseAndSetAlgorithmOptions().

◆ stokes_int

const std::map<stokes, t_int> purify::stokes_int
Initial value:
= {{stokes::I, 1}, {stokes::Q, 2}, {stokes::U, 3},
{stokes::V, 4}, {stokes::RR, -1}, {stokes::LL, -2},
{stokes::RL, -3}, {stokes::LR, -4}, {stokes::XX, -5},
{stokes::YY, -6}, {stokes::XY, -7}, {stokes::YX, -8}}

Definition at line 45 of file types.h.

Referenced by purify::pfitsio::header_params::header_params(), and TEST_CASE().

◆ stokes_string

const std::map<std::string, stokes> purify::stokes_string
Initial value:
= {
{"I", stokes::I}, {"i", stokes::I}, {"Q", stokes::Q}, {"q", stokes::Q},
{"U", stokes::U}, {"u", stokes::U}, {"V", stokes::V}, {"v", stokes::V},
{"XX", stokes::XX}, {"xx", stokes::XX}, {"YY", stokes::YY}, {"yy", stokes::YY},
{"XY", stokes::XY}, {"xy", stokes::XY}, {"YX", stokes::YX}, {"yx", stokes::YX},
{"LL", stokes::LL}, {"ll", stokes::LL}, {"RR", stokes::RR}, {"rr", stokes::RR},
{"LR", stokes::LR}, {"lr", stokes::LR}, {"RL", stokes::RL}, {"rl", stokes::RL},
{"P", stokes::P}, {"p", stokes::P}}

Definition at line 49 of file types.h.

Referenced by purify::YamlParser::parseAndSetGeneralConfiguration(), and TEST_CASE().