Common Pipeline Library Reference 7.3.2
|
Functions | |
cpl_error_code | cpl_bivector_copy (cpl_bivector *self, const cpl_bivector *other) |
Copy contents of a bivector into another bivector. | |
void | cpl_bivector_delete (cpl_bivector *f) |
Delete a cpl_bivector. | |
void | cpl_bivector_dump (const cpl_bivector *f, FILE *stream) |
Dump a cpl_bivector as ASCII to a stream. | |
cpl_bivector * | cpl_bivector_duplicate (const cpl_bivector *in) |
Duplicate a cpl_bivector. | |
cpl_size | cpl_bivector_get_size (const cpl_bivector *in) |
Get the size of the cpl_bivector. | |
cpl_vector * | cpl_bivector_get_x (cpl_bivector *in) |
Get a pointer to the x vector of the cpl_bivector. | |
const cpl_vector * | cpl_bivector_get_x_const (const cpl_bivector *in) |
Get a pointer to the x vector of the cpl_bivector. | |
double * | cpl_bivector_get_x_data (cpl_bivector *in) |
Get a pointer to the x data part of the cpl_bivector. | |
const double * | cpl_bivector_get_x_data_const (const cpl_bivector *in) |
Get a pointer to the x data part of the cpl_bivector. | |
cpl_vector * | cpl_bivector_get_y (cpl_bivector *in) |
Get a pointer to the y vector of the cpl_bivector. | |
const cpl_vector * | cpl_bivector_get_y_const (const cpl_bivector *in) |
Get a pointer to the y vector of the cpl_bivector. | |
double * | cpl_bivector_get_y_data (cpl_bivector *in) |
Get a pointer to the y data part of the cpl_bivector. | |
const double * | cpl_bivector_get_y_data_const (const cpl_bivector *in) |
Get a pointer to the y data part of the cpl_bivector. | |
cpl_error_code | cpl_bivector_interpolate_linear (cpl_bivector *fout, const cpl_bivector *fref) |
Linear interpolation of a 1d-function. | |
cpl_bivector * | cpl_bivector_new (cpl_size n) |
Create a new cpl_bivector. | |
cpl_bivector * | cpl_bivector_read (const char *filename) |
Read a list of values from an ASCII file and create a cpl_bivector. | |
cpl_error_code | cpl_bivector_sort (cpl_bivector *self, const cpl_bivector *other, cpl_sort_direction dir, cpl_sort_mode mode) |
Sort a cpl_bivector. | |
void | cpl_bivector_unwrap_vectors (cpl_bivector *f) |
Free memory associated to a cpl_bivector, excluding the two vectors. | |
cpl_bivector * | cpl_bivector_wrap_vectors (cpl_vector *x, cpl_vector *y) |
Create a new cpl_bivector from two cpl_vectors. | |
This module provides functions to handle cpl_bivector.
A cpl_bivector is composed of two vectors of the same size. It can be used to store 1d functions, with the x and y positions of the samples, offsets in x and y or simply positions in an image.
This module provides among other things functions for interpolation and for sorting one vector according to another.
cpl_error_code cpl_bivector_copy | ( | cpl_bivector * | self, |
const cpl_bivector * | other | ||
) |
Copy contents of a bivector into another bivector.
self | destination cpl_vector |
other | source cpl_vector |
Possible _cpl_error_code_ set in this function:
References cpl_bivector_get_x(), cpl_bivector_get_x_const(), cpl_bivector_get_y(), cpl_bivector_get_y_const(), CPL_ERROR_NONE, and cpl_vector_copy().
void cpl_bivector_delete | ( | cpl_bivector * | f | ) |
Delete a cpl_bivector.
f | cpl_bivector to delete |
This function deletes a bivector. If the input bivector f is NULL
, nothing is done, and no error is set.
References cpl_free(), and cpl_vector_delete().
Referenced by cpl_apertures_get_fwhm(), cpl_flux_get_noise_ring(), cpl_geom_img_offset_combine(), cpl_image_fit_gaussian(), cpl_wlcalib_slitmodel_delete(), and cpl_wlcalib_slitmodel_set_catalog().
void cpl_bivector_dump | ( | const cpl_bivector * | f, |
FILE * | stream | ||
) |
Dump a cpl_bivector as ASCII to a stream.
f | Input cpl_bivector to dump or NULL |
stream | Output stream, accepts stdout or stderr or NULL |
Comment lines start with the hash character.
stream may be NULL in which case stdout
is used.
References cpl_bivector_get_size(), cpl_vector_get_data_const(), and cpl_vector_get_size().
cpl_bivector * cpl_bivector_duplicate | ( | const cpl_bivector * | in | ) |
Duplicate a cpl_bivector.
in | cpl_bivector to duplicate |
The returned object must be deallocated using cpl_bivector_delete()
Possible _cpl_error_code_ set in this function:
References cpl_ensure, CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_NULL_INPUT, cpl_malloc(), cpl_vector_duplicate(), and cpl_vector_get_size().
cpl_size cpl_bivector_get_size | ( | const cpl_bivector * | in | ) |
Get the size of the cpl_bivector.
in | the input bivector |
Possible _cpl_error_code_ set in this function:
References cpl_ensure, CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_NULL_INPUT, and cpl_vector_get_size().
Referenced by cpl_bivector_dump(), cpl_bivector_interpolate_linear(), cpl_bivector_sort(), cpl_geom_img_offset_combine(), cpl_geom_img_offset_fine(), cpl_geom_img_offset_saa(), cpl_plot_bivector(), and cpl_plot_bivectors().
cpl_vector * cpl_bivector_get_x | ( | cpl_bivector * | in | ) |
Get a pointer to the x vector of the cpl_bivector.
in | a cpl_bivector |
The returned pointer refers to an already created cpl_vector.
Possible _cpl_error_code_ set in this function:
References cpl_ensure, and CPL_ERROR_NULL_INPUT.
Referenced by cpl_bivector_copy(), and cpl_geom_img_offset_combine().
const cpl_vector * cpl_bivector_get_x_const | ( | const cpl_bivector * | in | ) |
Get a pointer to the x vector of the cpl_bivector.
in | a cpl_bivector |
References cpl_ensure, and CPL_ERROR_NULL_INPUT.
Referenced by cpl_bivector_copy(), cpl_bivector_interpolate_linear(), cpl_geom_img_offset_saa(), and cpl_polynomial_fit_2d_create().
double * cpl_bivector_get_x_data | ( | cpl_bivector * | in | ) |
Get a pointer to the x data part of the cpl_bivector.
in | a cpl_bivector |
Possible _cpl_error_code_ set in this function:
References cpl_ensure, CPL_ERROR_NULL_INPUT, and cpl_vector_get_data().
Referenced by cpl_apertures_get_fwhm(), cpl_bivector_interpolate_linear(), cpl_bivector_sort(), cpl_geom_img_offset_combine(), cpl_geom_img_offset_fine(), cpl_image_fit_gaussian(), and cpl_image_iqe().
const double * cpl_bivector_get_x_data_const | ( | const cpl_bivector * | in | ) |
Get a pointer to the x data part of the cpl_bivector.
in | a cpl_bivector |
References cpl_ensure, CPL_ERROR_NULL_INPUT, and cpl_vector_get_data_const().
Referenced by cpl_bivector_interpolate_linear(), cpl_bivector_sort(), cpl_flux_get_noise_ring(), cpl_geom_img_offset_fine(), cpl_plot_bivector(), and cpl_plot_bivectors().
cpl_vector * cpl_bivector_get_y | ( | cpl_bivector * | in | ) |
Get a pointer to the y vector of the cpl_bivector.
in | a cpl_bivector |
The returned pointer refers to an already created cpl_vector.
Possible _cpl_error_code_ set in this function:
References cpl_ensure, and CPL_ERROR_NULL_INPUT.
Referenced by cpl_bivector_copy(), and cpl_geom_img_offset_combine().
const cpl_vector * cpl_bivector_get_y_const | ( | const cpl_bivector * | in | ) |
Get a pointer to the y vector of the cpl_bivector.
in | a cpl_bivector |
References cpl_ensure, and CPL_ERROR_NULL_INPUT.
Referenced by cpl_bivector_copy(), cpl_geom_img_offset_saa(), and cpl_polynomial_fit_2d_create().
double * cpl_bivector_get_y_data | ( | cpl_bivector * | in | ) |
Get a pointer to the y data part of the cpl_bivector.
in | a cpl_bivector |
Possible _cpl_error_code_ set in this function:
References cpl_ensure, CPL_ERROR_NULL_INPUT, and cpl_vector_get_data().
Referenced by cpl_apertures_get_fwhm(), cpl_bivector_interpolate_linear(), cpl_bivector_sort(), cpl_geom_img_offset_combine(), cpl_geom_img_offset_fine(), and cpl_image_iqe().
const double * cpl_bivector_get_y_data_const | ( | const cpl_bivector * | in | ) |
Get a pointer to the y data part of the cpl_bivector.
in | a cpl_bivector |
References cpl_ensure, CPL_ERROR_NULL_INPUT, and cpl_vector_get_data_const().
Referenced by cpl_bivector_interpolate_linear(), cpl_bivector_sort(), cpl_flux_get_noise_ring(), cpl_geom_img_offset_fine(), cpl_plot_bivector(), and cpl_plot_bivectors().
cpl_error_code cpl_bivector_interpolate_linear | ( | cpl_bivector * | fout, |
const cpl_bivector * | fref | ||
) |
Linear interpolation of a 1d-function.
fout | Preallocated with X-vector set, to hold interpolation in Y |
fref | Reference 1d-function |
fref must have both its abscissa and ordinate defined. fout must have its abscissa defined and its ordinate allocated.
The linear interpolation will be done from the values in fref to the abscissa points in fout.
For each abscissa point in fout, fref must either have two neigboring abscissa points such that xref_i < xout_j < xref{i+1}, or a single identical abscissa point, such that xref_i == xout_j.
This is ensured by monotonely growing abscissa points in both fout and fref (and by min(xref) <= min(xout) and max(xout) < max(xref)).
However, for efficiency reasons (since fref can be very long) the monotonicity is only verified to the extent necessary to actually perform the interpolation.
This input requirement implies that extrapolation is not allowed.
Possible _cpl_error_code_ set in this function:
References cpl_bivector_get_size(), cpl_bivector_get_x_const(), cpl_bivector_get_x_data(), cpl_bivector_get_x_data_const(), cpl_bivector_get_y_data(), cpl_bivector_get_y_data_const(), cpl_ensure_code, CPL_ERROR_DATA_NOT_FOUND, cpl_error_get_code(), CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_NONE, and cpl_vector_find().
cpl_bivector * cpl_bivector_new | ( | cpl_size | n | ) |
Create a new cpl_bivector.
n | Positive number of points |
The returned object must be deallocated using cpl_bivector_delete() or cpl_bivector_unwrap_vectors(), provided the two cpl_vectors are deallocated separately.
Possible _cpl_error_code_ set in this function:
References cpl_ensure, CPL_ERROR_ILLEGAL_INPUT, cpl_malloc(), and cpl_vector_new().
Referenced by cpl_apertures_get_fwhm(), cpl_geom_img_offset_combine(), cpl_geom_img_offset_fine(), and cpl_image_iqe().
cpl_bivector * cpl_bivector_read | ( | const char * | filename | ) |
Read a list of values from an ASCII file and create a cpl_bivector.
filename | Name of the input ASCII file |
The input ASCII file must contain two values per line.
The returned object must be deallocated using cpl_bivector_delete() Two columns of numbers are expected in the input file.
In addition to normal files, FIFO (see man mknod) are also supported.
Possible _cpl_error_code_ set in this function:
References cpl_bivector_wrap_vectors(), cpl_ensure, CPL_ERROR_BAD_FILE_FORMAT, CPL_ERROR_FILE_IO, CPL_ERROR_NULL_INPUT, cpl_vector_delete(), cpl_vector_new(), cpl_vector_set(), and cpl_vector_set_size().
cpl_error_code cpl_bivector_sort | ( | cpl_bivector * | self, |
const cpl_bivector * | other, | ||
cpl_sort_direction | dir, | ||
cpl_sort_mode | mode | ||
) |
Sort a cpl_bivector.
self | cpl_bivector to hold sorted result |
other | Input cpl_bivector to sort, may equal self |
dir | CPL_SORT_ASCENDING or CPL_SORT_DESCENDING |
mode | CPL_SORT_BY_X or CPL_SORT_BY_Y |
The values in the input are sorted according to direction and mode, and the result is placed self which must be of the same size as other.
As for qsort(): If two members compare as equal, their order in the sorted array is undefined.
In place sorting is supported.
Possible _cpl_error_code_ set in this function:
References cpl_bivector_get_size(), cpl_bivector_get_x_data(), cpl_bivector_get_x_data_const(), cpl_bivector_get_y_data(), cpl_bivector_get_y_data_const(), cpl_ensure_code, CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_INCOMPATIBLE_INPUT, CPL_ERROR_NONE, CPL_ERROR_NULL_INPUT, and CPL_ERROR_UNSUPPORTED_MODE.
void cpl_bivector_unwrap_vectors | ( | cpl_bivector * | f | ) |
Free memory associated to a cpl_bivector, excluding the two vectors.
f | cpl_bivector to delete |
References cpl_free().
Referenced by cpl_plot_column(), and cpl_polynomial_fit().
cpl_bivector * cpl_bivector_wrap_vectors | ( | cpl_vector * | x, |
cpl_vector * | y | ||
) |
Create a new cpl_bivector from two cpl_vectors.
x | the x cpl_vector |
y | the y cpl_vector |
The returned object must be deallocated using cpl_bivector_delete() or with cpl_bivector_unwrap_vectors(), provided the two cpl_vectors are deallocated separately.
Possible _cpl_error_code_ set in this function:
References cpl_ensure, CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_NULL_INPUT, cpl_malloc(), and cpl_vector_get_size().
Referenced by cpl_bivector_read(), cpl_plot_column(), cpl_polynomial_fit(), and cpl_ppm_match_positions().