DP3
Namespaces | Classes | Typedefs | Enumerations | Functions
dp3::ddecal Namespace Reference

Namespaces

 test
 

Classes

class  AmplitudeOnlyConstraint
 This class constrains the phases of the solution to be zero, but keeps the amplitude information. More...
 
class  AntennaConstraint
 This constraint averages the solutions of several groups of antennas, so that antennas within the same group have equal solutions. More...
 
class  AntennaIntervalConstraint
 
class  ApproximateTECConstraint
 
class  BdaSolverBuffer
 
class  Constraint
 This class is the base class for classes that implement a constraint on calibration solutions. Constraints are used to increase the converge of calibration by applying these inside the solving step. More...
 
struct  ConstraintResult
 
class  DiagonalLowRankSolver
 
class  DiagonalSolver
 
class  FaradayConstraint
 
class  FullJonesSolver
 
class  HybridSolver
 
class  IterativeDiagonalSolver
 
class  IterativeDiagonalSolverCuda
 
class  IterativeFullJonesSolver
 
class  IterativeScalarSolver
 
class  KernelSmoother
 Smooths a series of possibly irregularly gridded values by a given kernel. More...
 
class  KLFitter
 Creates KH base and fits screens from collection of PiercePoints. More...
 
class  LLSSolver
 
class  NormalEquationsSolver
 
class  PieceWisePhaseFitter
 
class  PiercePoint
 
class  PolarizationLeakageConstraint
 
class  QRSolver
 
class  RotationAndDiagonalConstraint
 
class  RotationConstraint
 
class  ScalarSolver
 
class  ScreenConstraint
 
class  ScreenFitter
 Class to perform screen fitting. More...
 
struct  Settings
 This struct parses the DDECal parset settings and stores them. More...
 
class  SmoothnessConstraint
 
class  SolutionResampler
 Class for resampling the solution into a square multi-dimensional array. More...
 
class  SolutionWriter
 
class  SolveData
 
class  SolverBase
 
class  SVDSolver
 
class  TecConstraint
 
class  TecOffsetDelayConstraint
 

Typedefs

using complex = std::complex< float >
 
using DuoSolveData = SolveData< aocommon::MC2x2FDiag >
 
using FullSolveData = SolveData< aocommon::MC2x2F >
 Stores all 4 polarizations of the data. More...
 
using UniSolveData = SolveData< std::complex< float > >
 
using SolutionTensor = xt::xtensor< std::complex< double >, 4 >
 
using SolutionSpan = aocommon::xt::Span< std::complex< double >, 4 >
 

Enumerations

enum class  LLSSolverType { QR , SVD , NORMAL_EQUATIONS }
 
enum class  SolverAlgorithm {
  kLowRank , kDirectionSolve , kDirectionIterative , kHybrid ,
  kLBFGS
}
 
enum class  SolverDataUse { kSingle , kDual , kFull }
 

Functions

void AssignAndWeight (std::vector< std::unique_ptr< base::DPBuffer >> &unweighted_buffers, const std::vector< std::string > &direction_names, std::vector< base::DPBuffer > &weighted_buffers, bool keep_unweighted_model_data, bool linear_weighting_mode)
 
std::vector< double > CalculateAntennaSmoothnessFactors (double ref_distance, const std::string &ref_antenna, const std::vector< std::array< double, 3 >> &antenna_positions, const std::vector< std::string > &antenna_factors, const std::vector< std::string > &antenna_names)
 
void cgels_ (char *trans, int *m, int *n, int *nrhs, complex *a, int *lda, complex *b, int *ldb, complex *work, int *lwork, int *info)
 
void cgelss_ (int *m, int *n, int *nrhs, complex *a, int *lda, complex *b, int *ldb, float *s, float *rcond, int *rank, complex *work, int *lwork, float *rwork, int *info)
 
void ConstrainDiagonal (std::array< std::complex< double >, 2 > &diagonal, base::CalType mode)
 
void cposv_ (char *uplo, int *n, int *nrhs, std::complex< float > *a, int *lda, std::complex< float > *b, int *ldb, int *info)
 
std::unique_ptr< SolverBaseCreateSolver (const Settings &settings, const std::vector< std::string > &station_names)
 
template<bool Add, typename MatrixType >
void DiagonalAddOrSubtractDirection (const typename SolveData< MatrixType >::ChannelBlockData &cb_data, std::vector< MatrixType > &v_residual, size_t direction, size_t n_solutions, const std::vector< std::complex< double >> &solutions, size_t max_threads)
 
float DominantEigenPair (const xt::xtensor< std::complex< float >, 2 > &matrix, xt::xtensor< std::complex< float >, 1 > &eigen_vector, size_t n_iterations)
 
float DominantEigenPairNear (const xt::xtensor< std::complex< float >, 2 > &matrix, xt::xtensor< std::complex< float >, 1 > &eigen_vector, size_t n_iterations)
 
std::vector< size_t > GetSolutionToDirectionVector (const std::vector< uint32_t > &sub_solutions_per_direction)
 
void InitializeSolverConstraints (SolverBase &solver, const Settings &settings, const std::vector< std::array< double, 3 >> &antenna_positions, const std::vector< std::string > &antenna_names, const std::vector< base::Direction > &source_positions, const std::vector< double > &frequencies)
 
std::vector< ConstraintResultMakeDiagonalResults (size_t n_antennas, size_t n_sub_solutions, size_t n_channels, base::CalType diagonal_solution_type)
 
void ShowConstraintSettings (std::ostream &output, const Settings &settings)
 
void StoreDiagonal (ConstraintResult *results, const std::array< std::complex< double >, 2 > &diagonal, size_t channel, size_t antenna, size_t subsolution, size_t n_channels, size_t n_sub_solutions, base::CalType diagonal_solution_type)
 
std::string ToString (SolverAlgorithm algorithm)
 
void Weigh (const base::DPBuffer::DataType &in, base::DPBuffer::DataType &out, const base::DPBuffer::WeightsType &weights)
 

Typedef Documentation

◆ complex

typedef std::complex< float > dp3::ddecal::complex

◆ DuoSolveData

using dp3::ddecal::DuoSolveData = typedef SolveData<aocommon::MC2x2FDiag>

Stores only the 2 diagonal values of the data (e.g. XX/YY). Because the word "diagonal" is extensively used for the solution type, the term "Duo" is used for this.

◆ FullSolveData

using dp3::ddecal::FullSolveData = typedef SolveData<aocommon::MC2x2F>

Stores all 4 polarizations of the data.

◆ SolutionSpan

using dp3::ddecal::SolutionSpan = typedef aocommon::xt::Span<std::complex<double>, 4>

◆ SolutionTensor

using dp3::ddecal::SolutionTensor = typedef xt::xtensor<std::complex<double>, 4>

Data structures for solutions. The dimensions are (nr. channel blocks, nr. antennas, nr. subsolutions, nr. polarizations). The number of subsolutions (third dimension) depends on the number of directions and the number of solutions per direction. Different directions may have different solution counts.

See also
dp3::ddecal::SolveData::ChannelBlockData::NSolutionsForDirection

◆ UniSolveData

using dp3::ddecal::UniSolveData = typedef SolveData<std::complex<float> >

Enumeration Type Documentation

◆ LLSSolverType

Linear least-squares solver to use.

Enumerator
QR 
SVD 
NORMAL_EQUATIONS 

◆ SolverAlgorithm

Enumerator
kLowRank 
kDirectionSolve 
kDirectionIterative 
kHybrid 
kLBFGS 

◆ SolverDataUse

Enumerator
kSingle 
kDual 
kFull 

Function Documentation

◆ AssignAndWeight()

void dp3::ddecal::AssignAndWeight ( std::vector< std::unique_ptr< base::DPBuffer >> &  unweighted_buffers,
const std::vector< std::string > &  direction_names,
std::vector< base::DPBuffer > &  weighted_buffers,
bool  keep_unweighted_model_data,
bool  linear_weighting_mode 
)

This function takes (unweighted) data and model data, as well as a weights array, weights these and writes the result.

Parameters
unweighted_buffersA vector with one buffer for each timestep. Each buffer should contain unweighted data and the corresponding weight values.
direction_namesA list with the names of the model data buffers in the DPBuffers.
weighted_buffersA vector with one buffer for each timestep. This function writes the weighted data to these buffers, using the same direction name. Existing data with that name is overwritten.
keep_unweighted_model_dataIf true, keep the model data in unweighted_buffers and create new model data buffers in weighted_buffers. If false, avoid creating new model data buffers by moving the model data from unweighted_buffers to weighted_buffers and weighing the data in-place.
linear_weighting_modeIf true, it has two effects: the data will be linearly weighted instead of weighting them by the sqrt of the weights, and the weights are stored in the buffer. Gradient descent and conjugate gradient-like methods require sqrt weighted data but do not need the weights, whereas a rank-based approach requires linear weighted data and needs to know the applied weights.

◆ CalculateAntennaSmoothnessFactors()

std::vector<double> dp3::ddecal::CalculateAntennaSmoothnessFactors ( double  ref_distance,
const std::string &  ref_antenna,
const std::vector< std::array< double, 3 >> &  antenna_positions,
const std::vector< std::string > &  antenna_factors,
const std::vector< std::string > &  antenna_names 
)

◆ cgels_()

void dp3::ddecal::cgels_ ( char *  trans,
int *  m,
int *  n,
int *  nrhs,
complex a,
int *  lda,
complex b,
int *  ldb,
complex work,
int *  lwork,
int *  info 
)

◆ cgelss_()

void dp3::ddecal::cgelss_ ( int *  m,
int *  n,
int *  nrhs,
complex a,
int *  lda,
complex b,
int *  ldb,
float *  s,
float *  rcond,
int *  rank,
complex work,
int *  lwork,
float *  rwork,
int *  info 
)

◆ ConstrainDiagonal()

void dp3::ddecal::ConstrainDiagonal ( std::array< std::complex< double >, 2 > &  diagonal,
base::CalType  mode 
)

Constrain the diagonal in accordance with the supplied mode. The mode should be a diagonal, scalar or rotational mode. In case of rotational, the diagonal is set to unity.

◆ cposv_()

void dp3::ddecal::cposv_ ( char *  uplo,
int *  n,
int *  nrhs,
std::complex< float > *  a,
int *  lda,
std::complex< float > *  b,
int *  ldb,
int *  info 
)

◆ CreateSolver()

std::unique_ptr<SolverBase> dp3::ddecal::CreateSolver ( const Settings settings,
const std::vector< std::string > &  station_names 
)

◆ DiagonalAddOrSubtractDirection()

template<bool Add, typename MatrixType >
void dp3::ddecal::DiagonalAddOrSubtractDirection ( const typename SolveData< MatrixType >::ChannelBlockData &  cb_data,
std::vector< MatrixType > &  v_residual,
size_t  direction,
size_t  n_solutions,
const std::vector< std::complex< double >> &  solutions,
size_t  max_threads 
)

◆ DominantEigenPair()

float dp3::ddecal::DominantEigenPair ( const xt::xtensor< std::complex< float >, 2 > &  matrix,
xt::xtensor< std::complex< float >, 1 > &  eigen_vector,
size_t  n_iterations 
)

Calculates the largest eigenvalue and corresponding eigen vector using the power iteration method. The eigen vector is calculated from the eigen vector using the Rayleigh quotient.

Parameters
[out]eigen_vectorUsed to store the calculated eigen vector.
n_iterationsNumber of power iterations. Typical values are 10-20.
Returns
Largest eigen value of the given matrix.

◆ DominantEigenPairNear()

float dp3::ddecal::DominantEigenPairNear ( const xt::xtensor< std::complex< float >, 2 > &  matrix,
xt::xtensor< std::complex< float >, 1 > &  eigen_vector,
size_t  n_iterations 
)

Like DominantEigenPair, but starts from an estimate of the eigen vector to speed up the calculation.

Parameters
[in,out]eigen_vectoron input, an estimate of the eigen vector corresponding to the dominant eigen value. On output, the calculated eigen vector.
n_iterationsNumber of power iterations. Typical values are 10-20.
Returns
Largest eigen value of the given matrix.

◆ GetSolutionToDirectionVector()

std::vector<size_t> dp3::ddecal::GetSolutionToDirectionVector ( const std::vector< uint32_t > &  sub_solutions_per_direction)
Returns
A vector that maps the solution index to the direction index it belongs to.

◆ InitializeSolverConstraints()

void dp3::ddecal::InitializeSolverConstraints ( SolverBase solver,
const Settings settings,
const std::vector< std::array< double, 3 >> &  antenna_positions,
const std::vector< std::string > &  antenna_names,
const std::vector< base::Direction > &  source_positions,
const std::vector< double > &  frequencies 
)

Initializes all constraints for a given solver.

While DPInfo can distinguish between used and unused antennas, the solvers work with the used antennas only. antenna_positions and antenna_names thus only contain used antennas. Their ordering should match the ordering while solving, when using Constraint::Apply.

Parameters
solverA solver that may have constraints.
settingsIncludes settings for the constraints.
antenna_positionsFor each antenna, the position.
antenna_namesFor each antenna, the name.
source_positionsFor each direction, the source position in ra/dec.

◆ MakeDiagonalResults()

std::vector<ConstraintResult> dp3::ddecal::MakeDiagonalResults ( size_t  n_antennas,
size_t  n_sub_solutions,
size_t  n_channels,
base::CalType  diagonal_solution_type 
)

◆ ShowConstraintSettings()

void dp3::ddecal::ShowConstraintSettings ( std::ostream &  output,
const Settings settings 
)

Writes the relevant constraints of the settings to the output.

◆ StoreDiagonal()

void dp3::ddecal::StoreDiagonal ( ConstraintResult results,
const std::array< std::complex< double >, 2 > &  diagonal,
size_t  channel,
size_t  antenna,
size_t  subsolution,
size_t  n_channels,
size_t  n_sub_solutions,
base::CalType  diagonal_solution_type 
)
inline

◆ ToString()

std::string dp3::ddecal::ToString ( SolverAlgorithm  algorithm)

◆ Weigh()

void dp3::ddecal::Weigh ( const base::DPBuffer::DataType in,
base::DPBuffer::DataType out,
const base::DPBuffer::WeightsType weights 
)

Weighs input data by multiplying it with weight values. All arguments must have the same shape.

This function is a manually simdized version of xt::noalias(out) = in * weights; (xt::noalias avoids an intermediate buffer for the result.)

On an Intel Xeon E5-2660, this function is ~3 x faster than xt::noalias(out) = in * weights; when keep_unweighted_model_data == false in AssignAndWeight.

Parameters
inInput data buffer.
outOutput data buffer. May be equal to 'in'.
weightsWeights buffer.