Common Pipeline Library Reference 7.3.2
|
Functions | |
cpl_error_code | cpl_vector_add (cpl_vector *v1, const cpl_vector *v2) |
Add a cpl_vector to another. | |
cpl_error_code | cpl_vector_add_scalar (cpl_vector *v, double addend) |
Elementwise addition of a scalar to a vector. | |
cpl_error_code | cpl_vector_convolve_symmetric (cpl_vector *smoothed, const cpl_vector *conv_kernel) |
Convolve a 1d-signal with a symmetric 1D-signal. | |
cpl_error_code | cpl_vector_copy (cpl_vector *destination, const cpl_vector *source) |
This function copies contents of a vector into another vector. | |
cpl_size | cpl_vector_correlate (cpl_vector *vxc, const cpl_vector *v1, const cpl_vector *v2) |
Cross-correlation of two vectors. | |
cpl_error_code | cpl_vector_cycle (cpl_vector *self, const cpl_vector *other, double shift) |
Perform a cyclic shift to the right of the elements of a vector. | |
void | cpl_vector_delete (cpl_vector *v) |
Delete a cpl_vector. | |
cpl_error_code | cpl_vector_divide (cpl_vector *v1, const cpl_vector *v2) |
Divide two vectors element-wise. | |
cpl_error_code | cpl_vector_divide_scalar (cpl_vector *v, double divisor) |
Elementwise division of a vector with a scalar. | |
void | cpl_vector_dump (const cpl_vector *v, FILE *stream) |
Dump a cpl_vector as ASCII to a stream. | |
cpl_vector * | cpl_vector_duplicate (const cpl_vector *v) |
This function duplicates an existing vector and allocates memory. | |
cpl_error_code | cpl_vector_exponential (cpl_vector *v, double base) |
Compute the exponential of all vector elements. | |
cpl_vector * | cpl_vector_extract (const cpl_vector *v, cpl_size istart, cpl_size istop, cpl_size istep) |
Extract a sub_vector from a vector. | |
cpl_error_code | cpl_vector_fill (cpl_vector *v, double val) |
Fill a cpl_vector. | |
cpl_error_code | cpl_vector_fill_kernel_profile (cpl_vector *profile, cpl_kernel type, double radius) |
Fill a vector with a kernel profile. | |
cpl_vector * | cpl_vector_filter_lowpass_create (const cpl_vector *v, cpl_lowpass filter_type, cpl_size hw) |
Apply a low-pass filter to a cpl_vector. | |
cpl_vector * | cpl_vector_filter_median_create (const cpl_vector *self, cpl_size hw) |
Apply a 1D median filter of given half-width to a cpl_vector. | |
cpl_size | cpl_vector_find (const cpl_vector *sorted, double key) |
In a sorted vector find the element closest to the given value. | |
cpl_error_code | cpl_vector_fit_gaussian (const cpl_vector *x, const cpl_vector *sigma_x, const cpl_vector *y, const cpl_vector *sigma_y, cpl_fit_mode fit_pars, double *x0, double *sigma, double *area, double *offset, double *mse, double *red_chisq, cpl_matrix **covariance) |
Apply a 1d gaussian fit. | |
double | cpl_vector_get (const cpl_vector *in, cpl_size idx) |
Get an element of the vector. | |
double * | cpl_vector_get_data (cpl_vector *in) |
Get a pointer to the data part of the vector. | |
const double * | cpl_vector_get_data_const (const cpl_vector *in) |
Get a pointer to the data part of the vector. | |
double | cpl_vector_get_max (const cpl_vector *v) |
Get the maximum of the cpl_vector. | |
cpl_size | cpl_vector_get_maxpos (const cpl_vector *self) |
Get the index of the maximum element of the cpl_vector. | |
double | cpl_vector_get_mean (const cpl_vector *v) |
Compute the mean value of vector elements. | |
double | cpl_vector_get_median (cpl_vector *v) |
Compute the median of the elements of a vector. | |
double | cpl_vector_get_median_const (const cpl_vector *v) |
Compute the median of the elements of a vector. | |
double | cpl_vector_get_min (const cpl_vector *v) |
Get the minimum of the cpl_vector. | |
cpl_size | cpl_vector_get_minpos (const cpl_vector *self) |
Get the index of the minimum element of the cpl_vector. | |
cpl_size | cpl_vector_get_size (const cpl_vector *in) |
Get the size of the vector. | |
double | cpl_vector_get_stdev (const cpl_vector *v) |
Compute the bias-corrected standard deviation of a vectors elements. | |
double | cpl_vector_get_sum (const cpl_vector *v) |
Get the sum of the elements of the cpl_vector. | |
cpl_vector * | cpl_vector_load (const char *filename, cpl_size xtnum) |
Load a list of values from a FITS file. | |
cpl_error_code | cpl_vector_logarithm (cpl_vector *v, double base) |
Compute the element-wise logarithm. | |
cpl_error_code | cpl_vector_multiply (cpl_vector *v1, const cpl_vector *v2) |
Multiply two vectors component-wise. | |
cpl_error_code | cpl_vector_multiply_scalar (cpl_vector *v, double factor) |
Elementwise multiplication of a vector with a scalar. | |
cpl_vector * | cpl_vector_new (cpl_size n) |
Create a new cpl_vector. | |
cpl_vector * | cpl_vector_new_lss_kernel (double slitw, double fwhm) |
Create Right Half of a symmetric smoothing kernel for LSS. | |
cpl_error_code | cpl_vector_power (cpl_vector *v, double exponent) |
Compute the power of all vector elements. | |
double | cpl_vector_product (const cpl_vector *v1, const cpl_vector *v2) |
Compute the vector dot product. | |
cpl_vector * | cpl_vector_read (const char *filename) |
Read a list of values from an ASCII file and create a cpl_vector. | |
cpl_error_code | cpl_vector_save (const cpl_vector *self, const char *filename, cpl_type type, const cpl_propertylist *plist, unsigned mode) |
Save a vector to a FITS file. | |
cpl_error_code | cpl_vector_set (cpl_vector *in, cpl_size idx, double value) |
Set an element of the vector. | |
cpl_error_code | cpl_vector_set_size (cpl_vector *in, cpl_size newsize) |
Resize the vector. | |
cpl_error_code | cpl_vector_sort (cpl_vector *self, cpl_sort_direction dir) |
Sort a cpl_vector. | |
cpl_error_code | cpl_vector_sqrt (cpl_vector *v) |
Compute the sqrt of a cpl_vector. | |
cpl_error_code | cpl_vector_subtract (cpl_vector *v1, const cpl_vector *v2) |
Subtract a cpl_vector from another. | |
cpl_error_code | cpl_vector_subtract_scalar (cpl_vector *v, double subtrahend) |
Elementwise subtraction of a scalar from a vector. | |
void * | cpl_vector_unwrap (cpl_vector *v) |
Delete a cpl_vector except the data array. | |
cpl_vector * | cpl_vector_wrap (cpl_size n, double *data) |
Create a cpl_vector from existing data. | |
This module provides functions to handle cpl_vector.
A cpl_vector is an object containing a list of values (as doubles) and the (always positive) number of values. The functionalities provided here are simple ones like sorting, statistics, or simple operations. The cpl_bivector object is composed of two of these vectors. No special provisions are made to handle special values like NaN or Inf, for data with such elements, the cpl_array object may be preferable.
cpl_error_code cpl_vector_add | ( | cpl_vector * | v1, |
const cpl_vector * | v2 | ||
) |
Add a cpl_vector to another.
v1 | First cpl_vector (modified) |
v2 | Second cpl_vector |
The second vector is added to the first one. The input first vector is modified.
The input vectors must have the same size.
Possible _cpl_error_code_ set in this function:
References cpl_ensure_code, CPL_ERROR_INCOMPATIBLE_INPUT, CPL_ERROR_NONE, and CPL_ERROR_NULL_INPUT.
cpl_error_code cpl_vector_add_scalar | ( | cpl_vector * | v, |
double | addend | ||
) |
Elementwise addition of a scalar to a vector.
v | cpl_vector to modify |
addend | Number to add |
Add a number to each element of the cpl_vector.
Possible _cpl_error_code_ set in this function:
References cpl_ensure_code, CPL_ERROR_NONE, and CPL_ERROR_NULL_INPUT.
cpl_error_code cpl_vector_convolve_symmetric | ( | cpl_vector * | smoothed, |
const cpl_vector * | conv_kernel | ||
) |
Convolve a 1d-signal with a symmetric 1D-signal.
smoothed | Preallocated vector to be smoothed in place |
conv_kernel | Vector with symmetric convolution function |
References cpl_ensure_code, CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_NONE, CPL_ERROR_NULL_INPUT, cpl_vector_delete(), cpl_vector_duplicate(), cpl_vector_get_data(), cpl_vector_get_data_const(), and cpl_vector_get_size().
cpl_error_code cpl_vector_copy | ( | cpl_vector * | destination, |
const cpl_vector * | source | ||
) |
This function copies contents of a vector into another vector.
destination | destination cpl_vector |
source | source cpl_vector |
Possible _cpl_error_code_ set in this function:
References cpl_ensure_code, cpl_error_get_code(), CPL_ERROR_NONE, CPL_ERROR_NULL_INPUT, and cpl_vector_set_size().
Referenced by cpl_bivector_copy(), and cpl_vector_duplicate().
cpl_size cpl_vector_correlate | ( | cpl_vector * | vxc, |
const cpl_vector * | v1, | ||
const cpl_vector * | v2 | ||
) |
Cross-correlation of two vectors.
vxc | Odd-sized vector with the computed cross-correlations |
v1 | 1st vector to correlate |
v2 | 2nd vector to correlate |
vxc must have an odd number of elements, 2*half_search+1, where half_search is the half-size of the search domain.
The length of v2 may not exceed that of v1. If the difference in length between v1 and v2 is less than half_search then this difference must be even (if the difference is odd resampling of v2 may be useful).
The cross-correlation is computed with shifts ranging from -half_search to half_search.
On succesful return element i (starting with 0) of vxc contains the cross- correlation at offset i-half_search. On error vxc is unmodified.
The cross-correlation is in fact the dot-product of two unit-vectors and ranges therefore from -1 to 1.
The cross-correlation is, in absence of rounding errors, commutative only for equal-sized vectors, i.e. changing the order of v1 and v2 will move element j in vxc to 2*half_search - j and thus change the return value from i to 2*half_search - i.
If, in absence of rounding errors, more than one shift would give the maximum cross-correlation, rounding errors may cause any one of those shifts to be returned. If rounding errors have no effect the index corresponding to the shift with the smallest absolute value is returned (with preference given to the smaller of two indices that correspond to the same absolute shift).
If v1 is longer than v2, the first element in v1 used for the resulting cross-correlation is max(0,shift + (cpl_vector_get_size(v1)-cpl_vector_get_size(v2))/2).
Cross-correlation with half_search == 0 requires about 8n FLOPs, where n = cpl_vector_get_size(v2). Each increase of half_search by 1 requires about 4n FLOPs more, when all of v2's elements can be cross-correlated, otherwise the extra cost is about 4m, where m is the number of elements in v2 that can be cross-correlated, n - half_search <= m < n.
In case of error, the _cpl_error_code_ code is set, and the returned delta and cross-correlation is undefined.
Example of 1D-wavelength calibration (error handling omitted for brevity):
Possible _cpl_error_code_ set in this function:
References cpl_ensure, CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_NULL_INPUT, cpl_vector_get_size(), and cpl_vector_set().
cpl_error_code cpl_vector_cycle | ( | cpl_vector * | self, |
const cpl_vector * | other, | ||
double | shift | ||
) |
Perform a cyclic shift to the right of the elements of a vector.
self | Vector to hold the shifted result |
other | Vector to read from, or NULL for in-place shift of self |
shift | The number of positions to cyclic right-shift |
A shift of +1 will move the last element to the first, a shift of -1 will move the first element to the last, a zero-shift will perform a copy (or do nothing in case of an in-place operation).
A non-integer shift will perform the shift in the Fourier domain. Large discontinuities in the vector to shift will thus lead to FFT artifacts around each discontinuity.
Possible _cpl_error_code_ set in this function:
References cpl_ensure_code, CPL_ERROR_INCOMPATIBLE_INPUT, CPL_ERROR_NONE, CPL_ERROR_NULL_INPUT, cpl_free(), cpl_image_unwrap(), cpl_image_wrap_double(), cpl_malloc(), cpl_vector_get_data(), and cpl_vector_get_size().
void cpl_vector_delete | ( | cpl_vector * | v | ) |
Delete a cpl_vector.
v | cpl_vector to delete |
If the vector v is NULL
, nothing is done and no error is set.
References cpl_free().
Referenced by cpl_bivector_delete(), cpl_bivector_read(), cpl_fit_image_gaussian(), cpl_flux_get_noise_ring(), cpl_geom_img_offset_combine(), cpl_geom_img_offset_saa(), cpl_image_fill_jacobian_polynomial(), cpl_image_get_fwhm(), cpl_image_warp_polynomial(), cpl_matrix_solve_svd(), cpl_matrix_solve_svd_threshold(), cpl_plot_column(), cpl_vector_convolve_symmetric(), cpl_vector_fill_polynomial_fit_residual(), cpl_vector_filter_lowpass_create(), cpl_vector_fit_gaussian(), cpl_vector_new_from_image_column(), cpl_vector_new_from_image_row(), cpl_vector_new_lss_kernel(), cpl_vector_read(), cpl_wlcalib_find_best_1d(), and cpl_wlcalib_slitmodel_delete().
cpl_error_code cpl_vector_divide | ( | cpl_vector * | v1, |
const cpl_vector * | v2 | ||
) |
Divide two vectors element-wise.
v1 | First cpl_vector |
v2 | Second cpl_vector |
If an element in v2 is zero v1 is not modified and CPL_ERROR_DIVISION_BY_ZERO is returned.
Possible _cpl_error_code_ set in this function:
References cpl_ensure_code, CPL_ERROR_DIVISION_BY_ZERO, CPL_ERROR_INCOMPATIBLE_INPUT, CPL_ERROR_NONE, and CPL_ERROR_NULL_INPUT.
cpl_error_code cpl_vector_divide_scalar | ( | cpl_vector * | v, |
double | divisor | ||
) |
Elementwise division of a vector with a scalar.
v | cpl_vector to modify |
divisor | Non-zero number to divide with |
Divide each element of the cpl_vector with a number.
Possible _cpl_error_code_ set in this function:
References cpl_ensure_code, CPL_ERROR_DIVISION_BY_ZERO, CPL_ERROR_NONE, and CPL_ERROR_NULL_INPUT.
void cpl_vector_dump | ( | const cpl_vector * | v, |
FILE * | stream | ||
) |
Dump a cpl_vector as ASCII to a stream.
v | Input cpl_vector to dump |
stream | Output stream, accepts stdout or stderr |
Each element is preceded by its index number (starting with 1!) and written on a single line.
Comment lines start with the hash character.
stream may be NULL in which case stdout
is used.
References CPL_SIZE_FORMAT.
cpl_vector * cpl_vector_duplicate | ( | const cpl_vector * | v | ) |
This function duplicates an existing vector and allocates memory.
v | the input cpl_vector |
The returned object must be deallocated using cpl_vector_delete()
Possible _cpl_error_code_ set in this function:
References cpl_ensure, CPL_ERROR_NULL_INPUT, cpl_vector_copy(), and cpl_vector_new().
Referenced by cpl_bivector_duplicate(), cpl_geom_img_offset_saa(), cpl_plot_column(), cpl_vector_convolve_symmetric(), cpl_vector_filter_median_create(), and cpl_vector_fit_gaussian().
cpl_error_code cpl_vector_exponential | ( | cpl_vector * | v, |
double | base | ||
) |
Compute the exponential of all vector elements.
v | Target cpl_vector. |
base | Exponential base. |
If the base is zero all vector elements must be positive and if the base is negative all vector elements must be integer, otherwise a cpl_error_code is returned and the vector is unmodified.
Possible _cpl_error_code_ set in this function:
References cpl_ensure_code, CPL_ERROR_DIVISION_BY_ZERO, CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_NONE, and CPL_ERROR_NULL_INPUT.
cpl_vector * cpl_vector_extract | ( | const cpl_vector * | v, |
cpl_size | istart, | ||
cpl_size | istop, | ||
cpl_size | istep | ||
) |
Extract a sub_vector from a vector.
v | Input cpl_vector |
istart | Start index (from 0 to number of elements - 1) |
istop | Stop index (from 0 to number of elements - 1) |
istep | Extract every step element |
The returned object must be deallocated using cpl_vector_delete()
FIXME: Currently istop must be greater than istart. FIXME: Currently istep must equal 1.
Possible _cpl_error_code_ set in this function:
References cpl_ensure, CPL_ERROR_ACCESS_OUT_OF_RANGE, CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_NULL_INPUT, and cpl_vector_new().
cpl_error_code cpl_vector_fill | ( | cpl_vector * | v, |
double | val | ||
) |
Fill a cpl_vector.
v | cpl_vector to be filled with the value val |
val | Value used to fill the cpl_vector |
Input vector is modified
Possible _cpl_error_code_ set in this function:
References cpl_ensure_code, CPL_ERROR_NONE, and CPL_ERROR_NULL_INPUT.
Referenced by cpl_wlcalib_find_best_1d().
cpl_error_code cpl_vector_fill_kernel_profile | ( | cpl_vector * | profile, |
cpl_kernel | type, | ||
double | radius | ||
) |
Fill a vector with a kernel profile.
profile | cpl_vector to be filled |
type | Type of kernel profile |
radius | Radius of the profile in pixels |
A number of predefined kernel profiles are available:
Possible _cpl_error_code_ set in this function:
References cpl_ensure_code, cpl_error_get_code(), CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_INVALID_TYPE, CPL_ERROR_NONE, and cpl_vector_get_size().
Referenced by cpl_geom_img_offset_saa().
cpl_vector * cpl_vector_filter_lowpass_create | ( | const cpl_vector * | v, |
cpl_lowpass | filter_type, | ||
cpl_size | hw | ||
) |
Apply a low-pass filter to a cpl_vector.
v | cpl_vector |
filter_type | Type of filter to use |
hw | Filter half-width |
This type of low-pass filtering consists in a convolution with a given kernel. The chosen filter type determines the kind of kernel to apply for convolution. Supported kernels are CPL_LOWPASS_LINEAR and CPL_LOWPASS_GAUSSIAN.
In the case of CPL_LOWPASS_GAUSSIAN, the gaussian sigma used is 1/sqrt(2). As this function is not meant to be general and cover all possible cases, this sigma is hardcoded and cannot be changed.
The returned smooth cpl_vector must be deallocated using cpl_vector_delete(). The returned signal has exactly as many samples as the input signal.
Possible _cpl_error_code_ set in this function:
References cpl_ensure, CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_NULL_INPUT, cpl_vector_delete(), and cpl_vector_new().
cpl_vector * cpl_vector_filter_median_create | ( | const cpl_vector * | self, |
cpl_size | hw | ||
) |
Apply a 1D median filter of given half-width to a cpl_vector.
self | Input vector to be filtered |
hw | Filter half-width |
This function applies a median smoothing to a cpl_vector and returns a newly allocated cpl_vector containing a median-smoothed version of the input. The returned cpl_vector must be deallocated using cpl_vector_delete().
The returned cpl_vector has exactly as many samples as the input one. The outermost hw values are copies of the input, each of the others is set to the median of its surrounding 1 + 2 * hw values.
For historical reasons twice the half-width is allowed to equal the vector length, although in this case the returned vector is simply a duplicate of the input one.
If different processing of the outer values is needed or if a more general kernel is needed, then cpl_image_filter_mask() can be called instead with CPL_FILTER_MEDIAN and the 1D-image input wrapped around self.
Possible _cpl_error_code_ set in this function:
References CPL_BORDER_COPY, cpl_ensure, CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_NONE, CPL_ERROR_NULL_INPUT, CPL_FILTER_MEDIAN, cpl_image_filter_mask(), cpl_image_unwrap(), cpl_image_wrap_double(), cpl_mask_delete(), cpl_mask_new(), cpl_mask_not(), cpl_vector_duplicate(), and cpl_vector_new().
cpl_size cpl_vector_find | ( | const cpl_vector * | sorted, |
double | key | ||
) |
In a sorted vector find the element closest to the given value.
sorted | CPL vector sorted using CPL_SORT_ASCENDING |
key | Value to find |
Bisection is used to find the element.
If two (neighboring) elements with different values both minimize fabs(sorted[index] - key) the index of the larger element is returned.
If the vector contains identical elements that minimize fabs(sorted[index] - key) then it is undefined which element has its index returned.
Possible _cpl_error_code_ set in this function:
References cpl_ensure, CPL_ERROR_ILLEGAL_INPUT, and CPL_ERROR_NULL_INPUT.
Referenced by cpl_bivector_interpolate_linear().
cpl_error_code cpl_vector_fit_gaussian | ( | const cpl_vector * | x, |
const cpl_vector * | sigma_x, | ||
const cpl_vector * | y, | ||
const cpl_vector * | sigma_y, | ||
cpl_fit_mode | fit_pars, | ||
double * | x0, | ||
double * | sigma, | ||
double * | area, | ||
double * | offset, | ||
double * | mse, | ||
double * | red_chisq, | ||
cpl_matrix ** | covariance | ||
) |
Apply a 1d gaussian fit.
x | Positions to fit |
sigma_x | Uncertainty (one sigma, gaussian errors assumed) assosiated with x. Taking into account the uncertainty of the independent variable is currently unsupported, and this parameter must therefore be set to NULL. |
y | Values to fit |
sigma_y | Uncertainty (one sigma, gaussian errors assumed) associated with y. If NULL, constant uncertainties are assumed. |
fit_pars | Specifies which parameters participate in the fit (any other parameters will be held constant). Possible values are CPL_FIT_CENTROID, CPL_FIT_STDEV, CPL_FIT_AREA, CPL_FIT_OFFSET and any bitwise combination of these. As a shorthand for including all four parameters in the fit, use CPL_FIT_ALL. |
x0 | (output) Center of best fit gaussian. If CPL_FIT_CENTROID is not set, this is also an input parameter. |
sigma | (output) Width of best fit gaussian. A positive number on success. If CPL_FIT_STDEV is not set, this is also an input parameter. |
area | (output) Area of gaussian. A positive number on succes. If CPL_FIT_AREA is not set, this is also an input parameter. |
offset | (output) Fitted background level. If CPL_FIT_OFFSET is not set, this is also an input parameter. |
mse | (output) If non-NULL, the mean squared error of the best fit is returned. |
red_chisq | (output) If non-NULL, the reduced chi square of the best fit is returned. This requires the noise vector to be specified. |
covariance | (output) If non-NULL, the formal covariance matrix of the best fit is returned. This requires sigma_y to be specified. The order of fit parameters in the covariance matrix is defined as (x0, sigma, area, offset), for example the (3,3) element of the matrix (counting from zero) is the variance of the fitted offset. The matrix must be deallocated by calling cpl_matrix_delete() . On error, NULL is returned. |
This function fits to the input vectors a 1d gaussian function of the form
f(x) = area / sqrt(2 pi sigma^2) * exp( -(x - x0)^2/(2 sigma^2)) + offset
(area > 0) by minimizing chi^2 using a Levenberg-Marquardt algorithm.
The values to fit are read from the input vector x. Optionally, a vector sigma_x (of same size as x) may be specified.
Optionally, the mean squared error, the reduced chi square and the covariance matrix of the best fit are computed . Set corresponding parameter to NULL to ignore.
If the covariance matrix is requested and successfully computed, the diagonal elements (the variances) are guaranteed to be positive.
Occasionally, the Levenberg-Marquardt algorithm fails to converge to a set of sensible parameters. In this case (and only in this case), a CPL_ERROR_CONTINUE is set. To allow the caller to recover from this particular error, the parameters x0, sigma, area and offset will on output contain estimates of the best fit parameters, specifically estimated as the median position, the median of the absolute residuals multiplied by 1.4828, the minimum flux value and the maximum flux difference multiplied by sqrt(2 pi sigma^2), respectively. In this case, mse, red_chisq and covariance are not computed. Note that the variance of x0 (the (0,0) element of the covariance matrix) in this case can be estimated by sigma^2 / area .
A CPL_ERROR_SINGULAR_MATRIX occurs if the covariance matrix cannot be computed. In that case all other output parameters are valid.
Current limitations
Possible _cpl_error_code_ set in this function:
References cpl_ensure_code, CPL_ERROR_CONTINUE, CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_ILLEGAL_OUTPUT, CPL_ERROR_INCOMPATIBLE_INPUT, CPL_ERROR_INVALID_TYPE, CPL_ERROR_NONE, CPL_ERROR_NULL_INPUT, CPL_ERROR_SINGULAR_MATRIX, CPL_ERROR_UNSUPPORTED_MODE, cpl_free(), cpl_malloc(), CPL_MATH_SQRT2PI, cpl_matrix_delete(), cpl_matrix_get_mean(), cpl_matrix_unwrap(), cpl_matrix_wrap(), cpl_vector_delete(), cpl_vector_duplicate(), cpl_vector_get(), cpl_vector_get_data(), cpl_vector_get_data_const(), cpl_vector_get_max(), cpl_vector_get_min(), cpl_vector_get_size(), cpl_vector_new(), and cpl_vector_set().
double cpl_vector_get | ( | const cpl_vector * | in, |
cpl_size | idx | ||
) |
Get an element of the vector.
in | the input vector |
idx | the index of the element (0 to nelem-1) |
In case of error, the _cpl_error_code_ code is set, and the returned double is undefined.
Possible _cpl_error_code_ set in this function:
References cpl_ensure, CPL_ERROR_ACCESS_OUT_OF_RANGE, CPL_ERROR_ILLEGAL_INPUT, and CPL_ERROR_NULL_INPUT.
Referenced by cpl_apertures_extract(), cpl_imagelist_erase(), cpl_vector_fit_gaussian(), and cpl_wlcalib_find_best_1d().
double * cpl_vector_get_data | ( | cpl_vector * | in | ) |
Get a pointer to the data part of the vector.
in | the input vector |
The returned pointer refers to already allocated data.
Possible _cpl_error_code_ set in this function:
References cpl_ensure, and CPL_ERROR_NULL_INPUT.
Referenced by cpl_bivector_get_x_data(), cpl_bivector_get_y_data(), cpl_geom_img_offset_fine(), cpl_image_fill_jacobian_polynomial(), cpl_image_warp_polynomial(), cpl_photom_fill_blackbody(), cpl_vector_convolve_symmetric(), cpl_vector_cycle(), cpl_vector_fit_gaussian(), cpl_vector_new_from_image_column(), cpl_vector_new_from_image_row(), and cpl_wlcalib_find_best_1d().
const double * cpl_vector_get_data_const | ( | const cpl_vector * | in | ) |
Get a pointer to the data part of the vector.
in | the input vector |
References cpl_ensure, and CPL_ERROR_NULL_INPUT.
Referenced by cpl_bivector_dump(), cpl_bivector_get_x_data_const(), cpl_bivector_get_y_data_const(), cpl_geom_img_offset_combine(), cpl_geom_img_offset_saa(), cpl_image_get_interpolated(), cpl_image_warp(), cpl_image_warp_polynomial(), cpl_matrix_solve_svd(), cpl_matrix_solve_svd_threshold(), cpl_photom_fill_blackbody(), cpl_plot_vector(), cpl_plot_vectors(), cpl_polynomial_fit_2d_create(), cpl_ppm_match_positions(), cpl_vector_convolve_symmetric(), cpl_vector_fit_gaussian(), cpl_vector_product(), and cpl_wlcalib_find_best_1d().
double cpl_vector_get_max | ( | const cpl_vector * | v | ) |
Get the maximum of the cpl_vector.
v | const cpl_vector |
References cpl_vector_get_maxpos().
Referenced by cpl_geom_img_offset_saa(), cpl_matrix_solve_svd_threshold(), and cpl_vector_fit_gaussian().
cpl_size cpl_vector_get_maxpos | ( | const cpl_vector * | self | ) |
Get the index of the maximum element of the cpl_vector.
self | The vector to process |
Possible _cpl_error_code_ set in this function:
References cpl_ensure, and CPL_ERROR_NULL_INPUT.
Referenced by cpl_vector_get_max().
double cpl_vector_get_mean | ( | const cpl_vector * | v | ) |
Compute the mean value of vector elements.
v | Input const cpl_vector |
References cpl_ensure, and CPL_ERROR_NULL_INPUT.
double cpl_vector_get_median | ( | cpl_vector * | v | ) |
Compute the median of the elements of a vector.
v | Input cpl_vector |
References cpl_ensure, and CPL_ERROR_NULL_INPUT.
Referenced by cpl_flux_get_noise_ring().
double cpl_vector_get_median_const | ( | const cpl_vector * | v | ) |
Compute the median of the elements of a vector.
v | Input const cpl_vector |
For a finite population or sample, the median is the middle value of an odd number of values (arranged in ascending order) or any value between the two middle values of an even number of values. The criteria used for an even number of values in the input array is to choose the mean between the two middle values. Note that in this case, the median might not be a value of the input array. Also, note that in the case of integer data types, the result will be converted to an integer. Consider to transform your int array to float if that is not the desired behavior.
References cpl_ensure, CPL_ERROR_NULL_INPUT, cpl_free(), and cpl_malloc().
double cpl_vector_get_min | ( | const cpl_vector * | v | ) |
Get the minimum of the cpl_vector.
v | const cpl_vector |
In case of error, the _cpl_error_code_ code is set, and the returned double is undefined.
Possible _cpl_error_code_ set in this function:
References cpl_vector_get_minpos().
Referenced by cpl_geom_img_offset_saa(), and cpl_vector_fit_gaussian().
cpl_size cpl_vector_get_minpos | ( | const cpl_vector * | self | ) |
Get the index of the minimum element of the cpl_vector.
self | The vector to process |
Possible _cpl_error_code_ set in this function:
References cpl_ensure, and CPL_ERROR_NULL_INPUT.
Referenced by cpl_vector_get_min().
cpl_size cpl_vector_get_size | ( | const cpl_vector * | in | ) |
Get the size of the vector.
in | the input vector |
Possible _cpl_error_code_ set in this function:
References cpl_ensure, and CPL_ERROR_NULL_INPUT.
Referenced by cpl_apertures_extract(), cpl_bivector_dump(), cpl_bivector_duplicate(), cpl_bivector_get_size(), cpl_bivector_wrap_vectors(), cpl_fit_imagelist_polynomial_window(), cpl_geom_img_offset_fine(), cpl_image_get_interpolated(), cpl_image_warp(), cpl_image_warp_polynomial(), cpl_imagelist_erase(), cpl_photom_fill_blackbody(), cpl_plot_column(), cpl_plot_vector(), cpl_plot_vectors(), cpl_polynomial_fit(), cpl_polynomial_fit_2d_create(), cpl_ppm_match_positions(), cpl_test_get_bytes_vector(), cpl_vector_convolve_symmetric(), cpl_vector_correlate(), cpl_vector_cycle(), cpl_vector_fill_kernel_profile(), cpl_vector_fill_polynomial_fit_residual(), cpl_vector_fit_gaussian(), cpl_vector_product(), and cpl_wlcalib_find_best_1d().
double cpl_vector_get_stdev | ( | const cpl_vector * | v | ) |
Compute the bias-corrected standard deviation of a vectors elements.
v | Input const cpl_vector |
S(n-1) = sqrt((1/n-1) sum((xi-mean)^2) (i=1 -> n)
The length of v must be at least 2.
References cpl_ensure, CPL_ERROR_ILLEGAL_INPUT, and CPL_ERROR_NULL_INPUT.
Referenced by cpl_flux_get_noise_ring().
double cpl_vector_get_sum | ( | const cpl_vector * | v | ) |
Get the sum of the elements of the cpl_vector.
v | const cpl_vector |
References cpl_ensure, and CPL_ERROR_NULL_INPUT.
cpl_vector * cpl_vector_load | ( | const char * | filename, |
cpl_size | xtnum | ||
) |
Load a list of values from a FITS file.
filename | Name of the input file |
xtnum | Extension number in the file (0 for primary HDU) |
This function loads a vector from a FITS file (NAXIS=1), using cfitsio. The returned image has to be deallocated with cpl_vector_delete().
'xtnum' specifies from which extension the vector should be loaded. This could be 0 for the main data section or any number between 1 and N, where N is the number of extensions present in the file.
Possible _cpl_error_code_ set in this function:
References cpl_ensure, CPL_ERROR_BAD_FILE_FORMAT, CPL_ERROR_FILE_IO, CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_NULL_INPUT, CPL_ERROR_UNSUPPORTED_MODE, cpl_free(), cpl_malloc(), CPL_SIZE_FORMAT, and cpl_vector_wrap().
cpl_error_code cpl_vector_logarithm | ( | cpl_vector * | v, |
double | base | ||
) |
Compute the element-wise logarithm.
v | cpl_vector to modify. |
base | Logarithm base. |
The base and all the vector elements must be positive and the base must be different from 1, or a cpl_error_code will be returned and the vector will be left unmodified.
Possible _cpl_error_code_ set in this function:
References cpl_ensure_code, CPL_ERROR_DIVISION_BY_ZERO, CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_NONE, and CPL_ERROR_NULL_INPUT.
cpl_error_code cpl_vector_multiply | ( | cpl_vector * | v1, |
const cpl_vector * | v2 | ||
) |
Multiply two vectors component-wise.
v1 | First cpl_vector |
v2 | Second cpl_vector |
References cpl_ensure_code, CPL_ERROR_INCOMPATIBLE_INPUT, CPL_ERROR_NONE, and CPL_ERROR_NULL_INPUT.
cpl_error_code cpl_vector_multiply_scalar | ( | cpl_vector * | v, |
double | factor | ||
) |
Elementwise multiplication of a vector with a scalar.
v | cpl_vector to modify |
factor | Number to multiply with |
Multiply each element of the cpl_vector with a number.
Possible _cpl_error_code_ set in this function:
References cpl_ensure_code, CPL_ERROR_NONE, and CPL_ERROR_NULL_INPUT.
cpl_vector * cpl_vector_new | ( | cpl_size | n | ) |
Create a new cpl_vector.
n | Number of element of the cpl_vector |
The returned object must be deallocated using cpl_vector_delete(). There is no default values assigned to the created object, they are undefined until they are set.
Possible _cpl_error_code_ set in this function:
References CPL_ERROR_ILLEGAL_INPUT, cpl_malloc(), and CPL_SIZE_FORMAT.
Referenced by cpl_bivector_new(), cpl_bivector_read(), cpl_flux_get_noise_ring(), cpl_geom_img_offset_combine(), cpl_geom_img_offset_saa(), cpl_image_fill_jacobian_polynomial(), cpl_image_warp_polynomial(), cpl_matrix_solve_svd(), cpl_matrix_solve_svd_threshold(), cpl_vector_duplicate(), cpl_vector_extract(), cpl_vector_filter_lowpass_create(), cpl_vector_filter_median_create(), cpl_vector_fit_gaussian(), cpl_vector_new_from_image_column(), cpl_vector_new_from_image_row(), cpl_vector_new_lss_kernel(), cpl_vector_read(), and cpl_wlcalib_find_best_1d().
cpl_vector * cpl_vector_new_lss_kernel | ( | double | slitw, |
double | fwhm | ||
) |
Create Right Half of a symmetric smoothing kernel for LSS.
slitw | The slit width [pixel] |
fwhm | The spectral FWHM [pixel] |
References CPL_MATH_SIG_FWHM, cpl_vector_delete(), and cpl_vector_new().
cpl_error_code cpl_vector_power | ( | cpl_vector * | v, |
double | exponent | ||
) |
Compute the power of all vector elements.
v | Target cpl_vector. |
exponent | Constant exponent. |
If the exponent is negative all vector elements must be non-zero and if the exponent is non-integer all vector elements must be non-negative, otherwise a cpl_error_code is returned and the vector is unmodified.
Following the behaviour of C99 pow() function, this function sets 0^0 = 1.
Possible _cpl_error_code_ set in this function:
References cpl_ensure_code, CPL_ERROR_DIVISION_BY_ZERO, CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_NONE, and CPL_ERROR_NULL_INPUT.
double cpl_vector_product | ( | const cpl_vector * | v1, |
const cpl_vector * | v2 | ||
) |
Compute the vector dot product.
v1 | One vector |
v2 | Another vector of the same size |
The same vector may be passed twice, in which case the square of its 2-norm is computed.
Possible _cpl_error_code_ set in this function:
References cpl_ensure, CPL_ERROR_INCOMPATIBLE_INPUT, CPL_ERROR_NULL_INPUT, cpl_vector_get_data_const(), and cpl_vector_get_size().
Referenced by cpl_vector_fill_polynomial_fit_residual().
cpl_vector * cpl_vector_read | ( | const char * | filename | ) |
Read a list of values from an ASCII file and create a cpl_vector.
filename | Name of the input ASCII file |
Parse an input ASCII file values and create a cpl_vector from it Lines beginning with a hash are ignored, blank lines also. In valid lines the value is preceeded by an integer, which is ignored.
The returned object must be deallocated using cpl_vector_delete()
In addition to normal files, FIFO (see man mknod) are also supported.
Possible _cpl_error_code_ set in this function:
References 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_vector_save | ( | const cpl_vector * | self, |
const char * | filename, | ||
cpl_type | type, | ||
const cpl_propertylist * | plist, | ||
unsigned | mode | ||
) |
Save a vector to a FITS file.
self | Vector to write to disk or NULL |
filename | Name of the file to write |
type | The type used to represent the data in the file |
plist | Property list for the output header or NULL |
mode | The desired output options (combined with bitwise or) |
This function saves a vector to a FITS file (NAXIS=1), using cfitsio. If a property list is provided, it is written to the named file before the pixels are written. If the image is not provided, the created file will only contain the primary header. This can be useful to create multiple extension files.
The type used in the file can be one of: CPL_TYPE_UCHAR (8 bit unsigned), CPL_TYPE_SHORT (16 bit signed), CPL_TYPE_USHORT (16 bit unsigned), CPL_TYPE_INT (32 bit signed), CPL_TYPE_FLOAT (32 bit floating point), or CPL_TYPE_DOUBLE (64 bit floating point). Use CPL_TYPE_DOUBLE when no loss of information is required.
Supported output modes are CPL_IO_CREATE (create a new file) and CPL_IO_EXTEND (append to an existing file)
If you are in append mode, make sure that the file has writing permissions. You may have problems if you create a file in your application and append something to it with the umask set to 222. In this case, the file created by your application would not be writable, and the append would fail.
Possible _cpl_error_code_ set in this function:
References cpl_ensure_code, CPL_ERROR_FILE_IO, CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_NONE, CPL_ERROR_NULL_INPUT, CPL_ERROR_UNSUPPORTED_MODE, cpl_free(), CPL_IO_CREATE, CPL_IO_EXTEND, cpl_propertylist_save(), and cpl_sprintf().
cpl_error_code cpl_vector_set | ( | cpl_vector * | in, |
cpl_size | idx, | ||
double | value | ||
) |
Set an element of the vector.
in | the input vector |
idx | the index of the element (0 to nelem-1) |
value | the value to set in the vector |
Possible _cpl_error_code_ set in this function:
References cpl_ensure_code, CPL_ERROR_ACCESS_OUT_OF_RANGE, CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_NONE, and CPL_ERROR_NULL_INPUT.
Referenced by cpl_bivector_read(), cpl_flux_get_noise_ring(), cpl_plot_column(), cpl_vector_correlate(), cpl_vector_fit_gaussian(), cpl_vector_read(), and cpl_wlcalib_find_best_1d().
cpl_error_code cpl_vector_set_size | ( | cpl_vector * | in, |
cpl_size | newsize | ||
) |
Resize the vector.
in | The vector to be resized |
newsize | The new (positive) number of elements in the vector |
cpl_vector_get_data()
should be discarded. If the vector was created with cpl_vector_wrap() the argument pointer to that call should be discarded as well.Possible _cpl_error_code_ set in this function:
References cpl_ensure_code, CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_NONE, CPL_ERROR_NULL_INPUT, and cpl_realloc().
Referenced by cpl_bivector_read(), cpl_geom_img_offset_combine(), cpl_vector_copy(), cpl_vector_fill_polynomial_fit_residual(), and cpl_vector_read().
cpl_error_code cpl_vector_sort | ( | cpl_vector * | self, |
cpl_sort_direction | dir | ||
) |
Sort a cpl_vector.
self | cpl_vector to sort in place |
dir | CPL_SORT_ASCENDING or CPL_SORT_DESCENDING |
The input cpl_vector is modified to sort its values in either ascending (CPL_SORT_ASCENDING) or descending (CPL_SORT_DESCENDING) order.
Possible _cpl_error_code_ set in this function:
References cpl_ensure_code, CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_NONE, and CPL_ERROR_NULL_INPUT.
Referenced by cpl_geom_img_offset_saa().
cpl_error_code cpl_vector_sqrt | ( | cpl_vector * | v | ) |
Compute the sqrt of a cpl_vector.
v | cpl_vector |
The sqrt of the data is computed. The input cpl_vector is modified
If an element in v is negative v is not modified and CPL_ERROR_ILLEGAL_INPUT is returned.
Possible _cpl_error_code_ set in this function:
References cpl_ensure_code, CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_NONE, and CPL_ERROR_NULL_INPUT.
cpl_error_code cpl_vector_subtract | ( | cpl_vector * | v1, |
const cpl_vector * | v2 | ||
) |
Subtract a cpl_vector from another.
v1 | First cpl_vector (modified) |
v2 | Second cpl_vector |
References cpl_ensure_code, CPL_ERROR_INCOMPATIBLE_INPUT, CPL_ERROR_NONE, and CPL_ERROR_NULL_INPUT.
cpl_error_code cpl_vector_subtract_scalar | ( | cpl_vector * | v, |
double | subtrahend | ||
) |
Elementwise subtraction of a scalar from a vector.
v | cpl_vector to modify |
subtrahend | Number to subtract |
Subtract a number from each element of the cpl_vector.
Possible _cpl_error_code_ set in this function:
References cpl_ensure_code, CPL_ERROR_NONE, and CPL_ERROR_NULL_INPUT.
void * cpl_vector_unwrap | ( | cpl_vector * | v | ) |
Delete a cpl_vector except the data array.
v | cpl_vector to delete |
References cpl_free().
Referenced by cpl_fit_image_gaussian(), cpl_flux_get_noise_ring(), cpl_plot_column(), cpl_polynomial_fit(), cpl_polynomial_fit_2d_create(), and cpl_wlcalib_find_best_1d().
cpl_vector * cpl_vector_wrap | ( | cpl_size | n, |
double * | data | ||
) |
Create a cpl_vector from existing data.
n | Number of elements in the vector |
data | Pointer to array of n doubles |
Possible _cpl_error_code_ set in this function:
References CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_NULL_INPUT, cpl_malloc(), and CPL_SIZE_FORMAT.
Referenced by cpl_fit_image_gaussian(), cpl_flux_get_noise_ring(), cpl_plot_column(), cpl_polynomial_fit(), cpl_polynomial_fit_2d_create(), cpl_ppm_match_positions(), cpl_vector_fill_polynomial_fit_residual(), cpl_vector_load(), and cpl_wlcalib_find_best_1d().