LibSakura
5.2
|
Sakura main header file. More...
#include <stddef.h>
#include <stdbool.h>
#include <stdint.h>
#include <sys/types.h>
#include <libsakura/config.h>
Go to the source code of this file.
Classes | |
struct | sakura_StatisticsResultFloat |
A structure in which the result of sakura_ComputeStatisticsFloat and sakura_ComputeAccurateStatisticsFloat is stored. More... | |
Macros | |
#define | LIBSAKURA_WARN_UNUSED_RESULT /* Don't ignore result value */ |
#define | LIBSAKURA_NOEXCEPT /* noexcept */ |
Typedefs | |
typedef void *(* | sakura_UserAllocator) (size_t size) |
A type of the allocator function used by Sakura Library. More... | |
typedef void(* | sakura_UserDeallocator) (void *pointer) |
A type of the deallocator function used by Sakura Library. More... | |
Enumerations | |
enum | sakura_Status { sakura_Status_kOK = 0 , sakura_Status_kNG = 1 , sakura_Status_kInvalidArgument = 2 , sakura_Status_kNoMemory = 3 , sakura_Status_kUnknownError = 99 } |
A result of function call. More... | |
enum | sakura_InterpolationMethod { sakura_InterpolationMethod_kNearest , sakura_InterpolationMethod_kLinear , sakura_InterpolationMethod_kPolynomial , sakura_InterpolationMethod_kSpline , sakura_InterpolationMethod_kNumElements } |
Enumerations to define interpolation types. More... | |
enum | sakura_LSQFitStatus { sakura_LSQFitStatus_kOK = 0 , sakura_LSQFitStatus_kNG = 1 , sakura_LSQFitStatus_kNotEnoughData = 2 , sakura_LSQFitStatus_kNumElements } |
Enumerations to define least-square fitting specific error code. More... | |
enum | sakura_LSQFitType { sakura_LSQFitType_kPolynomial , sakura_LSQFitType_kChebyshev , sakura_LSQFitType_kNumElements } |
Enumerations to define type of least-square fitting models. More... | |
Functions | |
sakura_Status | sakura_Initialize (sakura_UserAllocator allocator, sakura_UserDeallocator deallocator) |
Initializes Sakura Library. More... | |
void | sakura_CleanUp () |
Cleans up Sakura Library. More... | |
bool | sakura_IsAligned (void const *ptr) |
Checks if ptr points the aligned address Sakura Library requires. More... | |
size_t | sakura_GetAlignment () |
Returns an alignment that Sakura Library expects for arrays on which vector operations are performed. More... | |
void * | sakura_AlignAny (size_t size_of_arena, void *arena, size_t size_required) |
Returns an aligned address close to arena by adding 0 or minimum offset. More... | |
float * | sakura_AlignFloat (size_t elements_in_arena, float *arena, size_t elements_required) |
Returns an aligned address close to arena by adding 0 or minimum offset. More... | |
double * | sakura_AlignDouble (size_t elements_in_arena, double *arena, size_t elements_required) |
Returns an aligned address close to arena by adding 0 or minimum offset. More... | |
sakura_Status | sakura_ComputeStatisticsFloat (size_t num_data, float const data[], bool const is_valid[], sakura_StatisticsResultFloat *result) |
Computes statistics. More... | |
sakura_Status | sakura_ComputeAccurateStatisticsFloat (size_t num_data, float const data[], bool const is_valid[], sakura_StatisticsResultFloat *result) |
Computes statistics. More... | |
sakura_Status | sakura_SortValidValuesDenselyFloat (size_t num_data, bool const is_valid[], float data[], size_t *new_num_data) |
Sorts only valid data in ascending order. More... | |
sakura_Status | sakura_ComputeMedianAbsoluteDeviationFloat (size_t num_data, float const data[], float new_data[]) |
Computes median absolute deviation. More... | |
sakura_Status | sakura_GridConvolvingFloat (size_t num_spectra, size_t start_spectrum, size_t end_spectrum, bool const spectrum_mask[], double const x[], double const y[], size_t support, size_t sampling, size_t num_polarizations, uint32_t const polarization_map[], size_t num_channels, uint32_t const channel_map[], bool const mask[], float const value[], float const weight[], bool weight_only, size_t num_convolution_table, float const convolution_table[], size_t num_polarizations_for_grid, size_t num_channels_for_grid, size_t width, size_t height, double weight_sum[], float weight_of_grid[], float grid[]) |
Grids data with convolution. More... | |
sakura_Status | sakura_SetTrueIfInRangesInclusiveFloat (size_t num_data, float const data[], size_t num_condition, float const lower_bounds[], float const upper_bounds[], bool result[]) |
Sets true if the values in an input array (data ) are in any of specified range (inclusive). More... | |
sakura_Status | sakura_SetTrueIfInRangesInclusiveInt (size_t num_data, int const data[], size_t num_condition, int const lower_bounds[], int const upper_bounds[], bool result[]) |
Sets true if the values in an input array (data ) are in any of specified range (inclusive). More... | |
sakura_Status | sakura_SetTrueIfInRangesExclusiveFloat (size_t num_data, float const data[], size_t num_condition, float const lower_bounds[], float const upper_bounds[], bool result[]) |
Sets true if the values in an input array (data ) are in any of specified range (exclusive). More... | |
sakura_Status | sakura_SetTrueIfInRangesExclusiveInt (size_t num_data, int const data[], size_t num_condition, int const lower_bounds[], int const upper_bounds[], bool result[]) |
Sets true if the values in an input array (data ) are in any of specified range (exclusive). More... | |
sakura_Status | sakura_SetTrueIfGreaterThanFloat (size_t num_data, float const data[], float threshold, bool result[]) |
Sets true if the values in the input array (data ) are greater than a threshold. More... | |
sakura_Status | sakura_SetTrueIfGreaterThanInt (size_t num_data, int const data[], int threshold, bool result[]) |
Sets true if the values in the input array (data ) are greater than a threshold. More... | |
sakura_Status | sakura_SetTrueIfGreaterThanOrEqualsFloat (size_t num_data, float const data[], float threshold, bool result[]) |
Sets true if the values in the input array (data ) are greater than or equal to a threshold. More... | |
sakura_Status | sakura_SetTrueIfGreaterThanOrEqualsInt (size_t num_data, int const data[], int threshold, bool result[]) |
Sets true if the values in the input array (data ) are greater than or equal to a threshold. More... | |
sakura_Status | sakura_SetTrueIfLessThanFloat (size_t num_data, float const data[], float threshold, bool result[]) |
Sets true if the values in the input array (data ) are less than a threshold. More... | |
sakura_Status | sakura_SetTrueIfLessThanInt (size_t num_data, int const data[], int threshold, bool result[]) |
Sets true if the values in the input array (data ) are less than a threshold. More... | |
sakura_Status | sakura_SetTrueIfLessThanOrEqualsFloat (size_t num_data, float const data[], float threshold, bool result[]) |
Sets true if the values in the input array (data ) are less than or equal to a threshold. More... | |
sakura_Status | sakura_SetTrueIfLessThanOrEqualsInt (size_t num_data, int const data[], int threshold, bool result[]) |
Sets true if the values in the input array (data ) are less than or equal to a threshold. More... | |
sakura_Status | sakura_SetFalseIfNanOrInfFloat (size_t num_data, float const data[], bool result[]) |
Sets false if the values in the input array (data ) are either NaN or infinity. More... | |
sakura_Status | sakura_Uint8ToBool (size_t num_data, uint8_t const data[], bool result[]) |
Convert an input array to a boolean array. More... | |
sakura_Status | sakura_Uint32ToBool (size_t num_data, uint32_t const data[], bool result[]) |
Convert an input array to a boolean array. More... | |
sakura_Status | sakura_InvertBool (size_t num_data, bool const data[], bool result[]) |
Invert a boolean array. More... | |
sakura_Status | sakura_OperateBitwiseAndUint8 (uint8_t bit_mask, size_t num_data, uint8_t const data[], bool const edit_mask[], uint8_t result[]) |
Invoke bit operation AND between a bit mask and an array. More... | |
sakura_Status | sakura_OperateBitwiseAndUint32 (uint32_t bit_mask, size_t num_data, uint32_t const data[], bool const edit_mask[], uint32_t result[]) |
Invoke bit operation AND between a bit mask and an array. More... | |
sakura_Status | sakura_OperateBitwiseConverseNonImplicationUint8 (uint8_t bit_mask, size_t num_data, uint8_t const data[], bool const edit_mask[], uint8_t result[]) |
Invoke bit operation, Converse Nonimplication, between a bit mask and an array. More... | |
sakura_Status | sakura_OperateBitwiseConverseNonImplicationUint32 (uint32_t bit_mask, size_t num_data, uint32_t const data[], bool const edit_mask[], uint32_t result[]) |
Invoke bit operation, Converse Nonimplication, between a bit mask and an array. More... | |
sakura_Status | sakura_OperateBitwiseImplicationUint8 (uint8_t bit_mask, size_t num_data, uint8_t const data[], bool const edit_mask[], uint8_t result[]) |
Invoke bit operation, Material Implication, between a bit mask and an array. More... | |
sakura_Status | sakura_OperateBitwiseImplicationUint32 (uint32_t bit_mask, size_t num_data, uint32_t const data[], bool const edit_mask[], uint32_t result[]) |
Invoke bit operation, Material Implication, between a bit mask and an array. More... | |
sakura_Status | sakura_OperateBitwiseNotUint8 (size_t num_data, uint8_t const data[], bool const edit_mask[], uint8_t result[]) |
Invoke bitwise NOT operation of an array. More... | |
sakura_Status | sakura_OperateBitwiseNotUint32 (size_t num_data, uint32_t const data[], bool const edit_mask[], uint32_t result[]) |
Invoke bitwise NOT operation of an array. More... | |
sakura_Status | sakura_OperateBitwiseOrUint8 (uint8_t bit_mask, size_t num_data, uint8_t const data[], bool const edit_mask[], uint8_t result[]) |
Invoke bit operation OR between a bit mask and an array. More... | |
sakura_Status | sakura_OperateBitwiseOrUint32 (uint32_t bit_mask, size_t num_data, uint32_t const data[], bool const edit_mask[], uint32_t result[]) |
Invoke bit operation OR between a bit mask and an array. More... | |
sakura_Status | sakura_OperateBitwiseXorUint8 (uint8_t bit_mask, size_t num_data, uint8_t const data[], bool const edit_mask[], uint8_t result[]) |
Invoke bit operation XOR between a bit mask and an array. More... | |
sakura_Status | sakura_OperateBitwiseXorUint32 (uint32_t bit_mask, size_t num_data, uint32_t const data[], bool const edit_mask[], uint32_t result[]) |
Invoke bit operation XOR between a bit mask and an array. More... | |
sakura_Status | sakura_InterpolateXAxisFloat (sakura_InterpolationMethod interpolation_method, uint8_t polynomial_order, size_t num_base, double const base_position[], size_t num_array, float const base_data[], bool const base_mask[], size_t num_interpolated, double const interpolated_position[], float interpolated_data[], bool interpolated_mask[]) |
Perform one-dimensional interpolation. More... | |
sakura_Status | sakura_InterpolateYAxisFloat (sakura_InterpolationMethod interpolation_method, uint8_t polynomial_order, size_t num_base, double const base_position[], size_t num_array, float const base_data[], bool const base_mask[], size_t num_interpolated, double const interpolated_position[], float interpolated_data[], bool interpolated_mask[]) |
Perform one-dimensional interpolation. More... | |
sakura_Status | sakura_CalibrateDataWithArrayScalingFloat (size_t num_data, float const scaling_factor[], float const data[], float const reference[], float result[]) |
Normalize data against reference value with scaling factor. More... | |
sakura_Status | sakura_CalibrateDataWithConstScalingFloat (float scaling_factor, size_t num_data, float const data[], float const reference[], float result[]) |
Normalize data against reference value with scaling factor. More... | |
sakura_Status | sakura_CreateGaussianKernelFloat (float peak_location, float kernel_width, size_t num_kernel, float kernel[]) |
Create Gaussian kernel. More... | |
sakura_Status | sakura_CreateConvolve1DContextFFTFloat (size_t num_kernel, float const kernel[], struct sakura_Convolve1DContextFloat **context) |
Create context for convolution using Fourier transformation. More... | |
sakura_Status | sakura_Convolve1DFloat (size_t num_kernel, float const kernel[], size_t num_data, float const input_data[], bool const input_mask[], float output_data[], float output_weight[]) |
Convolution is performed. More... | |
sakura_Status | sakura_Convolve1DFFTFloat (struct sakura_Convolve1DContextFloat const *context, size_t num_data, float const input_data[], float output_data[]) |
Convolution is performed using Fourier transformation. More... | |
sakura_Status | sakura_DestroyConvolve1DContextFloat (struct sakura_Convolve1DContextFloat *context) |
Destroy context for convolution. More... | |
sakura_Status | sakura_GetLSQCoefficientsDouble (size_t const num_data, float const data[], bool const mask[], size_t const num_model_bases, double const basis_data[], size_t const num_lsq_bases, size_t const use_bases_idx[], double lsq_matrix[], double lsq_vector[]) |
Compute coefficients of simultaneous equations used for Least-Square fitting. More... | |
sakura_Status | sakura_UpdateLSQCoefficientsDouble (size_t const num_data, float const data[], bool const mask[], size_t const num_exclude_indices, size_t const exclude_indices[], size_t const num_model_bases, double const basis_data[], size_t const num_lsq_bases, size_t const use_bases_idx[], double lsq_matrix[], double lsq_vector[]) |
Compute coefficients of simultaneous equations used for Least-Square fitting. More... | |
sakura_Status | sakura_SolveSimultaneousEquationsByLUDouble (size_t num_equations, double const in_matrix[], double const in_vector[], double out[]) |
Solve simultaneous equations via LU decomposition. More... | |
sakura_Status | sakura_LMFitGaussianFloat (size_t const num_data, float const data[], bool const mask[], size_t const num_peaks, double height[], double err_height[], double center[], double err_center[], double sigma[], double err_sigma[]) |
Fit Gaussian to the input data using Levenberg-Marquardt method. More... | |
sakura_Status | sakura_CreateLSQFitContextPolynomialFloat (sakura_LSQFitType const lsqfit_type, uint16_t order, size_t num_data, struct sakura_LSQFitContextFloat **context) |
Create an object containing model data for least-square fitting. More... | |
sakura_Status | sakura_CreateLSQFitContextCubicSplineFloat (uint16_t npiece, size_t num_data, struct sakura_LSQFitContextFloat **context) |
Create an object containing model data for least-square fitting. More... | |
sakura_Status | sakura_CreateLSQFitContextSinusoidFloat (uint16_t nwave, size_t num_data, struct sakura_LSQFitContextFloat **context) |
Create an object containing model data for least-square fitting. More... | |
sakura_Status | sakura_DestroyLSQFitContextFloat (struct sakura_LSQFitContextFloat *context) |
Destroy an object containing model data for least-square fitting. More... | |
sakura_Status | sakura_LSQFitPolynomialFloat (struct sakura_LSQFitContextFloat const *context, uint16_t order, size_t num_data, float const data[], bool const mask[], float clip_threshold_sigma, uint16_t num_fitting_max, size_t num_coeff, double coeff[], float best_fit[], float residual[], bool final_mask[], float *rms, sakura_LSQFitStatus *lsqfit_status) |
Fit a polynomial curve to input data. More... | |
sakura_Status | sakura_LSQFitCubicSplineFloat (struct sakura_LSQFitContextFloat const *context, size_t num_pieces, size_t num_data, float const data[], bool const mask[], float clip_threshold_sigma, uint16_t num_fitting_max, double coeff[][4], float best_fit[], float residual[], bool final_mask[], float *rms, size_t boundary[], sakura_LSQFitStatus *lsqfit_status) |
Fit a cubic spline curve to input data. More... | |
sakura_Status | sakura_LSQFitSinusoidFloat (struct sakura_LSQFitContextFloat const *context, size_t num_nwave, size_t const nwave[], size_t num_data, float const data[], bool const mask[], float clip_threshold_sigma, uint16_t num_fitting_max, size_t num_coeff, double coeff[], float best_fit[], float residual[], bool final_mask[], float *rms, sakura_LSQFitStatus *lsqfit_status) |
Fit a sinusoidal curve to input data. More... | |
sakura_Status | sakura_SubtractPolynomialFloat (struct sakura_LSQFitContextFloat const *context, size_t num_data, float const data[], size_t num_coeff, double const coeff[], float out[]) |
Subtract best-fit polynomial model from input data. More... | |
sakura_Status | sakura_SubtractCubicSplineFloat (struct sakura_LSQFitContextFloat const *context, size_t num_data, float const data[], size_t num_pieces, double const coeff[][4], size_t const boundary[], float out[]) |
Subtract best-fit cubic spline model from input data. More... | |
sakura_Status | sakura_SubtractSinusoidFloat (struct sakura_LSQFitContextFloat const *context, size_t num_data, float const data[], size_t num_nwave, size_t const nwave[], size_t num_coeff, double const coeff[], float out[]) |
Subtract best-fit sinusoidal model from input data. More... | |
sakura_Status | sakura_GetNumberOfCoefficientsFloat (struct sakura_LSQFitContextFloat const *context, uint16_t order, size_t *num_coeff) |
Return the number of basis functions used for least-square fitting. More... | |
sakura_Status | sakura_FlipArrayFloat (bool inner_most_untouched, size_t dims, size_t const elements[], float const src[], float dst[]) |
Copy elements in the src array into the dst array with flipping elements to reorder as some FFT library expects. More... | |
sakura_Status | sakura_UnflipArrayFloat (bool inner_most_untouched, size_t dims, size_t const elements[], float const src[], float dst[]) |
Copy elements in the src array into the dst array with unflipping elements to the original order. More... | |
sakura_Status | sakura_FlipArrayDouble (bool inner_most_untouched, size_t dims, size_t const elements[], double const src[], double dst[]) |
Copy elements in the src array into the dst array with flipping elements to reorder as some FFT library expects. More... | |
sakura_Status | sakura_UnflipArrayDouble (bool inner_most_untouched, size_t dims, size_t const elements[], double const src[], double dst[]) |
Copy elements in the src array into the dst array with unflipping elements to the original order. More... | |
sakura_Status | sakura_FlipArrayDouble2 (bool inner_most_untouched, size_t dims, size_t const elements[], double const src[][2], double dst[][2]) |
Same as sakura_FlipArrayFloat except the element type of the multidimensional arrays. More... | |
sakura_Status | sakura_UnflipArrayDouble2 (bool inner_most_untouched, size_t dims, size_t const elements[], double const src[][2], double dst[][2]) |
Same as sakura_UnflipArrayFloat except the element type of the multidimensional arrays. More... | |
sakura_Status | sakura_CreateMaskNearEdgeDouble (float fraction, double pixel_size, size_t num_data, double const x[], double const y[], double const *blc_x, double const *blc_y, double const *trc_x, double const *trc_y, bool mask[]) |
Create mask. More... | |
#define LIBSAKURA_NOEXCEPT /* noexcept */ |
#define LIBSAKURA_WARN_UNUSED_RESULT /* Don't ignore result value */ |
typedef void*(* sakura_UserAllocator) (size_t size) |
A type of the allocator function used by Sakura Library.
Implementation of the function of this type must be reentrant.
[in] | size | Size of required memory in bytes |
MT-safe
typedef void(* sakura_UserDeallocator) (void *pointer) |
A type of the deallocator function used by Sakura Library.
Implementation of the function of this type must be reentrant.
[in] | pointer | NULL or an address to be released which was allocated by the allocator of type sakura_UserAllocator . |
MT-safe
Enumerations to define interpolation types.
enum sakura_LSQFitStatus |
enum sakura_LSQFitType |
enum sakura_Status |
A result of function call.
void* sakura_AlignAny | ( | size_t | size_of_arena, |
void * | arena, | ||
size_t | size_required | ||
) |
Returns an aligned address close to arena by adding 0 or minimum offset.
It returns arena if arena is already aligned.
[in] | size_of_arena | Size of the memory region pointed by arena |
[in] | arena | Start address of a memory region |
[in] | size_required | Required size after alignment |
MT-safe
double* sakura_AlignDouble | ( | size_t | elements_in_arena, |
double * | arena, | ||
size_t | elements_required | ||
) |
Returns an aligned address close to arena by adding 0 or minimum offset.
It returns arena if arena is already aligned.
[in] | elements_in_arena | The number of elements in arena , not a size in bytes. |
[in] | arena | Start address of an array |
[in] | elements_required | Required number of elements after alignment |
MT-safe
float* sakura_AlignFloat | ( | size_t | elements_in_arena, |
float * | arena, | ||
size_t | elements_required | ||
) |
Returns an aligned address close to arena by adding 0 or minimum offset.
It returns arena if arena is already aligned.
[in] | elements_in_arena | The number of elements in arena , not a size in bytes. |
[in] | arena | Start address of an array |
[in] | elements_required | Required number of elements after alignment |
MT-safe
sakura_Status sakura_CalibrateDataWithArrayScalingFloat | ( | size_t | num_data, |
float const | scaling_factor[], | ||
float const | data[], | ||
float const | reference[], | ||
float | result[] | ||
) |
Normalize data against reference value with scaling factor.
Normalize the data against reference value with scaling factor. The function normalizes the data by the assumption that the value in reference is normalized to scaling_factor. Specifically, it will calculate,
result = scaling_factor * (data - reference) / reference
The function returns result status. For successful run, return value will be sakura_Status_kOK . On the other hand, it will return sakura_Status_kInvalidArgument if any invalid values are passed to arguments.
The function allows in-place calculation, i.e., result array can be either data or reference. In this case, data or reference will be overwritten.
Only difference from sakura_CalibrateDataWithConstScalingFloat is whether scaling_factor is array or scalar. Note that, due to this difference, order of the arguments is slightly different between the two.
[in] | num_data | Number of data |
[in] | scaling_factor | Scaling factor. This is an array and number of elements must be num_data. must-be-aligned |
[in] | data | Data to be normalized. Number of elements must be num_data. must-be-aligned |
[in] | reference | Reference data. Number of elements must be num_data. must-be-aligned |
[out] | result | Resulting normalized data. One can give same array with either data or reference for in-place calculation. Number of elements must be num_data. must-be-aligned |
MT-safe
sakura_Status sakura_CalibrateDataWithConstScalingFloat | ( | float | scaling_factor, |
size_t | num_data, | ||
float const | data[], | ||
float const | reference[], | ||
float | result[] | ||
) |
Normalize data against reference value with scaling factor.
Normalize the data against reference value with scaling factor. The function normalizes the data by the assumption that the value in reference is normalized to scaling_factor. Specifically, it will calculate,
result = scaling_factor * (data - reference) / reference
The function returns result status. For successful run, return value will be sakura_Status_kOK . On the other hand, it will return sakura_Status_kInvalidArgument if any invalid values are passed to arguments.
The function allows in-place calculation, i.e., result array can be either data or reference. In this case, data or reference will be overwritten.
See sakura_CalibrateDataWithArrayScalingFloat about the difference between those two similar functions.
[in] | scaling_factor | Scaling factor. This is a scalar. |
[in] | num_data | Number of data |
[in] | data | Data to be normalized. Number of elements must be num_data. must-be-aligned |
[in] | reference | Reference data. Number of elements must be num_data. must-be-aligned |
[out] | result | Resulting normalized data. One can give same array with either data or reference for in-place calculation. Number of elements must be num_data. must-be-aligned |
MT-safe
void sakura_CleanUp | ( | ) |
Cleans up Sakura Library.
When you call this function, no function of Sakura Library must be running.
MT-unsafe
sakura_Status sakura_ComputeAccurateStatisticsFloat | ( | size_t | num_data, |
float const | data[], | ||
bool const | is_valid[], | ||
sakura_StatisticsResultFloat * | result | ||
) |
Computes statistics.
Refer to sakura_StatisticsResultFloat to see what kind of statistics are computed.
[in] | num_data | The number of elements in data and is_valid . num_data <= INT32_MAX |
[in] | data | Data. If corresponding element in is_valid is true, the element in data must not be Inf nor NaN. must-be-aligned |
[in] | is_valid | Masks of data. If a value of element is false, the corresponding element in data is ignored. must-be-aligned |
[out] | result | An address where the result should be stored. Some fields may be set to NaN if it is impossible to figure out. If there is more than one occurrences of min or max value, it is undefined which index of the occurrences is selected for index_of_min or index_of_max. |
MT-safe
The result of this function is more accurate than that of sakura_ComputeStatisticsFloat if num_data is large. This function is slower than sakura_ComputeStatisticsFloat .
sakura_Status sakura_ComputeMedianAbsoluteDeviationFloat | ( | size_t | num_data, |
float const | data[], | ||
float | new_data[] | ||
) |
Computes median absolute deviation.
This function applies abs(data[i] - median) for each element in data and stores results to new_data. Then it sorts new_data in ascending order. data must be sorted in advance using e.g. sakura_SortValidValuesDenselyFloat . The median is a center value if num_data is odd. Otherwise, the median is an average of two center values. The caller is responsible to take a median value from new_data as a median absolute deviation as below.
mad = (new_data[num_data/2] + new_data[num_data/2 - num_data%2]) / 2 if num_data > 0
[in] | num_data | The number of elements in data and new_data . |
[in] | data | Each value in data must not be Inf nor NaN. data must be sorted. must-be-aligned |
[out] | new_data | An array where the results are stored. new_data may points data. must-be-aligned |
MT-safe
sakura_Status sakura_ComputeStatisticsFloat | ( | size_t | num_data, |
float const | data[], | ||
bool const | is_valid[], | ||
sakura_StatisticsResultFloat * | result | ||
) |
Computes statistics.
Refer to sakura_StatisticsResultFloat to see what kind of statistics are computed.
[in] | num_data | The number of elements in data and is_valid . num_data <= INT32_MAX |
[in] | data | Data. If corresponding element in is_valid is true, the element in data must not be Inf nor NaN. must-be-aligned |
[in] | is_valid | Masks of data. If a value of element is false, the corresponding element in data is ignored. must-be-aligned |
[out] | result | An address where the result should be stored. Some fields may be set to NaN if it is impossible to figure out. If there is more than one occurrences of min or max value, it is undefined which index of the occurrences is selected for index_of_min or index_of_max. |
MT-safe
sakura_Status sakura_Convolve1DFFTFloat | ( | struct sakura_Convolve1DContextFloat const * | context, |
size_t | num_data, | ||
float const | input_data[], | ||
float | output_data[] | ||
) |
Convolution is performed using Fourier transformation.
Convolution operation is performed using Fourier transformation of data and kernel arrays. The kernel is stored in a context in an internal format. The operation is carried out for all elements of data.
[in] | context | The context created by sakura_CreateConvolve1DContextFFTFloat. |
[in] | num_data | The number of elements in input_data and output_data. num_data must be equal to num_kernel in context . 0 < num_data <= INT_MAX |
[in] | input_data | Input data. All elements in input_data must not be Inf nor NaN. must-be-aligned |
[out] | output_data | Output data. The pointer of out is allowed to be equal to that of in (input_data == output_data), indicating in-place operation. must-be-aligned |
MT-safe
(But see sakura_CreateConvolve1DContextFFTFloat for detail about the thread-safety)
sakura_Status sakura_Convolve1DFloat | ( | size_t | num_kernel, |
float const | kernel[], | ||
size_t | num_data, | ||
float const | input_data[], | ||
bool const | input_mask[], | ||
float | output_data[], | ||
float | output_weight[] | ||
) |
Convolution is performed.
Convolution operation is performed by shifting a kernel over the input data. The operation is carried out only for elements of data whose mask values are true.
[in] | num_kernel | The number of elements in the kernel. num_kernel must be positive. 0 < num_kernel <= INT_MAX |
[in] | kernel | Convolution kernel. All elements in kernel must not be Inf nor NaN. must-be-aligned |
[in] | num_data | The number of elements in input_data, input_mask, output_data, and output_weight. 0 < num_data <= INT_MAX |
[in] | input_data | Input data. If corresponding element in input_mask is true, the element in input_data must not be Inf nor NaN. must-be-aligned |
[in] | input_mask | Mask of input data. input_data is used in operation if corresponding element of input_mask is true. If not, the corresponding elements in input_data is ignored. must-be-aligned |
[out] | output_data | Output data. must-be-aligned |
[out] | output_weight | Total weight of kernel summed up to corresponding elements of output_data in the convolution operation. must-be-aligned |
MT-safe
sakura_Status sakura_CreateConvolve1DContextFFTFloat | ( | size_t | num_kernel, |
float const | kernel[], | ||
struct sakura_Convolve1DContextFloat ** | context | ||
) |
Create context for convolution using Fourier transformation.
[in] | num_kernel | The number of elements in the kernel. num_kernel must be positive. 0 < num_kernel <= INT_MAX |
[in] | kernel | Convolution kernel. All elements in kernel must not be Inf nor NaN. must-be-aligned |
[out] | context | Context for convolution. The context can be shared between threads. It has to be destroyed by sakura_DestroyConvolve1DContextFloat after use by sakura_Convolve1DFFTFloat . Note also that null pointer will be set to *context in case this function fails. |
MT-unsafe
sakura_Status sakura_CreateGaussianKernelFloat | ( | float | peak_location, |
float | kernel_width, | ||
size_t | num_kernel, | ||
float | kernel[] | ||
) |
Create Gaussian kernel.
Create 1 dimensional Gaussian kernel according to a peak location and FWHM given by peak_location and kernel_width.
The value of the kernel for i-th element is calculated by integrating Gaussian from i-1/2 to i+1/2. Since definite integral of Gaussian doesn't have analytic form, this function is implemented using error function. Resulting kernel is normalized so that its sum to be 1.0. Actual formula for kernel is as follows.
s = kernel_width / sqrt(log(16)) l = (i - 1/2 - peak_location) / s r = (i + 1/2 - peak_location) / s kernel[i] = 1/2 * {erf(r) - erf(l)}
where sqrt and log are mathematical functions of square root and logarithm to natural base. The function erf is an error function. The variable s corresponds to sqrt(2) times standard deviation of the Gaussian. The above calculation is implemented using std::erf. See the link below for details about std::erf:
[in] | peak_location | the peak location of Gaussian |
[in] | kernel_width | FWHM (Full Width of Half Maximum) of Gaussian. kernel_width must be greater than 0. |
[in] | num_kernel | The number of elements in the kernel |
[out] | kernel | Output kernel array must-be-aligned |
sakura_Status sakura_CreateLSQFitContextCubicSplineFloat | ( | uint16_t | npiece, |
size_t | num_data, | ||
struct sakura_LSQFitContextFloat ** | context | ||
) |
Create an object containing model data for least-square fitting.
[in] | npiece | Number of spline pieces. It must be a positive value. |
[in] | num_data | Number of data to fit. It must be equal to or larger than the number of model bases ( npiece+3 ), thus the smallest allowed value of num_data is 4. Note also that ( npiece+3 ) * num_data must not exceed SIZE_MAX. |
[out] | context | Pointer to pointer to an object containing model data. When no longer used, the object must be destroyed by sakura_DestroyLSQFitContextFloat . Note also that null pointer will be set to *context in case this function fails. |
MT-safe
sakura_Status sakura_CreateLSQFitContextPolynomialFloat | ( | sakura_LSQFitType const | lsqfit_type, |
uint16_t | order, | ||
size_t | num_data, | ||
struct sakura_LSQFitContextFloat ** | context | ||
) |
Create an object containing model data for least-square fitting.
[in] | lsqfit_type | Type of basis function. It should be either of sakura_LSQFitType_kPolynomial or sakura_LSQFitType_kChebyshev . |
[in] | order | Polynomial order. It must be positive or zero. |
[in] | num_data | Number of data to fit. It must be equal to or larger than the number of model bases, which is ( order+1 ), thus the smallest allowed value of num_data is 1. Note also that ( order+1 ) * num_data must not exceed SIZE_MAX. |
[out] | context | Pointer to pointer to an object containing model data. When no longer used, the object must be destroyed by sakura_DestroyLSQFitContextFloat . Note also that null pointer will be set to *context in case this function fails. |
MT-safe
sakura_Status sakura_CreateLSQFitContextSinusoidFloat | ( | uint16_t | nwave, |
size_t | num_data, | ||
struct sakura_LSQFitContextFloat ** | context | ||
) |
Create an object containing model data for least-square fitting.
[in] | nwave | Maximum wave number of sinusoids. It must be positive or zero. |
[in] | num_data | Number of data to fit. It must be equal to or larger than the number of model bases plus one, which is ( nwave*2+2 ), thus the smallest allowed value of num_data is 2. Note also that ( nwave*2+2 ) * num_data must not exceed SIZE_MAX. |
[out] | context | Pointer to pointer to an object containing model data. When no longer used, the object must be destroyed by sakura_DestroyLSQFitContextFloat . Note also that null pointer will be set to *context in case this function fails. |
MT-safe
sakura_Status sakura_CreateMaskNearEdgeDouble | ( | float | fraction, |
double | pixel_size, | ||
size_t | num_data, | ||
double const | x[], | ||
double const | y[], | ||
double const * | blc_x, | ||
double const * | blc_y, | ||
double const * | trc_x, | ||
double const * | trc_y, | ||
bool | mask[] | ||
) |
Create mask.
Based on two dimensional position distribution given by x and y, CreateMaskNearEdgeDouble creates mask for each data point. The algorithm is as follows:
median_separation = sorted_separation[(num_data - 1) / 2];It is recommended to set pixel_size to zero. Please consider to set nonzero value only if zero pixel_size doesn't work on your data. Number of pixels along horizontal and vertical axes is determined by
(ceil(wx / pixel_size), ceil(wy / pixel_size))where wx and wy are evaluated from blc_x, blc_y, trc_x, and trc_y if they are all non-NULL. The formula is
wx = (*trc_x - *blc_x) * 1.1 wy = (*trc_y - *blc_y) * 1.1Otherwise, they are calculated by
wx = (max(x) - min(x)) * 1.1 wy = (max(y) - min(y)) * 1.1
[in] | fraction | Fraction of the data to be masked. Threshold for mask operation is evaluated by fraction times num_data. 0 <= fraction <= 1. |
[in] | pixel_size | Size of the pixel. If it is zero, pixel size is set to a half of median of the distance between neighboring two positions that are provided by x and y. pixel_size must be greater than or equal to 0. Setting 0 works most of the cases. |
[in] | num_data | Number of data. |
[in] | x | List of horizontal positions. The length must be num_data. must-be-aligned |
[in] | y | List of vertical positions. The length must be num_data. must-be-aligned |
[in] | blc_x | Horizontal limit for bottom-left corner of user-defined range. If non-NULL value is provided, the data having x value less than *blc_x will be ignored. NULL corresponds to no limit. |
[in] | blc_y | Vertical limit for bottom-left corner of user-defined range. If non-NULL value is provided, the data having y value less than *blc_y will be ignored. NULL corresponds to no limit. |
[in] | trc_x | Horizontal limit for top-right corner of user-defined range. If non-NULL value is provided, the data having x value greater than *trc_x will be ignored. NULL corresponds to no limit. |
[in] | trc_y | Vertical limit for top-right corner of user-defined range. If non-NULL value is provided, the data having y value greater than *trc_y will be ignored. NULL corresponds to no limit. |
[out] | mask | Output mask. The length must be num_data since mask is associating with each data points provided by x and y. The points near edge will be masked, i.e., mask value will be set to true for those points. must-be-aligned |
sakura_Status sakura_DestroyConvolve1DContextFloat | ( | struct sakura_Convolve1DContextFloat * | context | ) |
Destroy context for convolution.
[in] | context | The context created by sakura_CreateConvolve1DContextFFTFloat and to be destroyed. |
MT-unsafe
sakura_Status sakura_DestroyLSQFitContextFloat | ( | struct sakura_LSQFitContextFloat * | context | ) |
Destroy an object containing model data for least-square fitting.
[in] | context | A context created by sakura_CreateLSQFitContextPolynomialFloat or sakura_CreateLSQFitContextCubicSplineFloat or sakura_CreateLSQFitContextSinusoidFloat . |
MT-safe
sakura_Status sakura_FlipArrayDouble | ( | bool | inner_most_untouched, |
size_t | dims, | ||
size_t const | elements[], | ||
double const | src[], | ||
double | dst[] | ||
) |
Copy elements in the src array into the dst array with flipping elements to reorder as some FFT library expects.
When you provide inner_most_untouched = false, elements = {3, 4} and src = {
}, then you will get dst as below.
If inner_most_untouched = true, then you will get dst as below.
[in] | inner_most_untouched | If true, the order of the inner most dimension is untouched. |
[in] | dims | Dimensions of the array src and dst. In other words, the number of elements in elements. |
[in] | elements | Numbers of elements of each dimension of src and dst with the inner-to-outer order. For example, when you have float matrix[4][3];
{ 3, 4 }
|
[in] | src | Source multidimensional array. must-be-aligned |
[in] | dst | Destination multidimensional array. must-be-aligned |
MT-safe
sakura_Status sakura_FlipArrayDouble2 | ( | bool | inner_most_untouched, |
size_t | dims, | ||
size_t const | elements[], | ||
double const | src[][2], | ||
double | dst[][2] | ||
) |
Same as sakura_FlipArrayFloat except the element type of the multidimensional arrays.
sakura_Status sakura_FlipArrayFloat | ( | bool | inner_most_untouched, |
size_t | dims, | ||
size_t const | elements[], | ||
float const | src[], | ||
float | dst[] | ||
) |
Copy elements in the src array into the dst array with flipping elements to reorder as some FFT library expects.
When you provide inner_most_untouched = false, elements = {3, 4} and src = {
}, then you will get dst as below.
If inner_most_untouched = true, then you will get dst as below.
[in] | inner_most_untouched | If true, the order of the inner most dimension is untouched. |
[in] | dims | Dimensions of the array src and dst. In other words, the number of elements in elements. |
[in] | elements | Numbers of elements of each dimension of src and dst with the inner-to-outer order. For example, when you have float matrix[4][3];
{ 3, 4 }
|
[in] | src | Source multidimensional array. must-be-aligned |
[in] | dst | Destination multidimensional array. must-be-aligned |
MT-safe
size_t sakura_GetAlignment | ( | ) |
Returns an alignment that Sakura Library expects for arrays on which vector operations are performed.
MT-safe
sakura_Status sakura_GetLSQCoefficientsDouble | ( | size_t const | num_data, |
float const | data[], | ||
bool const | mask[], | ||
size_t const | num_model_bases, | ||
double const | basis_data[], | ||
size_t const | num_lsq_bases, | ||
size_t const | use_bases_idx[], | ||
double | lsq_matrix[], | ||
double | lsq_vector[] | ||
) |
Compute coefficients of simultaneous equations used for Least-Square fitting.
Suppose fitting ( num_data ) discrete data points yi with a linear combination of ( num_model_bases ) bases (ai, bi, ..., ni), which are also given as ( num_data ) discrete points. Assuming the best-fit model is given as (A * ai + B * bi + ... + N * ni), where (A, B, C, ...) are the coefficients to be solved, these values are connected via the following simultaneous equations known as normal equation:
Note that the summation means all the data points except masked ones are to be added. This function computes the coefficients of the above simultaneous equations, i.e., the elements of the matrix at the left side and of the vector at the right side.
[in] | num_data | The number of elements in the arrays data and mask, and also the number of elements in each model data (i.e., discrete values of basis function) consisting the entire model. It must be a positive number. |
[in] | data | Input data with length of num_data . must-be-aligned |
[in] | mask | The input mask data with length of num_data . The i th element of data is included in input data if the i th element of mask is true, while it is excluded from input data if the corresponding element of mask is false. must-be-aligned |
[in] | num_model_bases | Number of model basis functions. It must be in range 0 < num_model_bases <= num_data . |
[in] | basis_data | A 1D array containing values of all basis functions concatenated. Loop for basis index must be inside of that for data index, i.e., the n -th data of the m -th model should be stored at basis_data [ num_data * (n-1) + (m-1) ]. Its length must be equal to ( num_model_bases * num_data ). must-be-aligned |
[in] | num_lsq_bases | The number of model basis functions to be used in actual fitting. It must be a positive number and must not exceed num_model_bases . |
[in] | use_bases_idx | A 1D array containing indices of basis model that are to be used for fitting. As for fitting types other than sinusoidal, it should be always [0, 1, 2, ..., (num_lsq_bases-1)]. Element values must not be duplicated, and must be in ascending order. must-be-aligned |
[out] | lsq_matrix | A 1D array containing the values of a matrix at the left side of simultaneous equations for least-square fitting. Its length should therefore be equal to ( num_lsq_bases * num_lsq_bases ). Loop for columns comes inside that for rows, i.e., the value at the m -th row and n -th column is stored at lsq_matrix [ num_lsq_bases * ( m -1) + ( n -1)], though lsq_matrix is actually symmetric. must-be-aligned |
[out] | lsq_vector | The values of a vector at the right side of simultaneous equations for least-square fitting. Its length should be equal to num_model_bases . must-be-aligned |
MT-safe
sakura_Status sakura_GetNumberOfCoefficientsFloat | ( | struct sakura_LSQFitContextFloat const * | context, |
uint16_t | order, | ||
size_t * | num_coeff | ||
) |
Return the number of basis functions used for least-square fitting.
[in] | context | A context created by sakura_CreateLSQFitContextPolynomialFloat or sakura_CreateLSQFitContextCubicSplineFloat . |
[in] | order | Parameter for the specified function. It is the order (for polynomial and Chebyshev polynomial), or the number of pieces (for cubic spline). It must be positive for cubic spline, while other models accept zero value. The value should not exceed the order specified in creation of context . |
[out] | num_coeff | Number of basis functions to be used for least-square fitting. This value should be the actual number of normal equations. |
MT-safe
sakura_Status sakura_GridConvolvingFloat | ( | size_t | num_spectra, |
size_t | start_spectrum, | ||
size_t | end_spectrum, | ||
bool const | spectrum_mask[], | ||
double const | x[], | ||
double const | y[], | ||
size_t | support, | ||
size_t | sampling, | ||
size_t | num_polarizations, | ||
uint32_t const | polarization_map[], | ||
size_t | num_channels, | ||
uint32_t const | channel_map[], | ||
bool const | mask[], | ||
float const | value[], | ||
float const | weight[], | ||
bool | weight_only, | ||
size_t | num_convolution_table, | ||
float const | convolution_table[], | ||
size_t | num_polarizations_for_grid, | ||
size_t | num_channels_for_grid, | ||
size_t | width, | ||
size_t | height, | ||
double | weight_sum[], | ||
float | weight_of_grid[], | ||
float | grid[] | ||
) |
Grids data with convolution.
This function plot value from start_spectrum to end_spectrum whose location is at (x , y) onto grid as points which have width represented by convolution_table. The results of 2 * support pixels from the edges of the grids are not reliable due to implementation to gain speed. Polarizations and channels are mapped by polarization_map and channel_map when gridding. All floating-point values passed to this function must not be NaN nor +-Inf.
[in] | num_spectra | The number of spectra. It should be 0 <= start_spectrum <= end_spectrum <= num_spectra . |
[in] | start_spectrum | Starting index of spectrum to be processed. |
[in] | end_spectrum | An index next to the last index of spectrum to be processed. |
[in] | spectrum_mask | Masks to each spectrum. The number of elements must be num_spectra. If a value of element is false, the corresponding spectrum is ignored. must-be-aligned |
[in] | x | X position of the point projected to 2D plane of the grid. The number of elements must be num_spectra . Each element must be INT32_MIN < x[i] < INT32_MAX . must-be-aligned |
[in] | y | Y position of the point projected to 2D plane of the grid. The number of elements must be num_spectra . Each element must be INT32_MIN < y[i] < INT32_MAX . must-be-aligned |
[in] | support | A half width(the number of pixels from center to edge) of convolution kernel on the grid plane of size width x height. It must be 0 < support <= 256 . Note that this function may cause stack overflow because it allocates memory of a size proportional to support * sampling on stack when support * sampling is too large. |
[in] | sampling | A resolution of convolution kernel(samples/grid aka samples/pixel). It must be 0 < sampling <= INT32_MAX . |
[in] | num_polarizations | The number of polarizations. It must be 0 < num_polarizations <= INT32_MAX . |
[in] | polarization_map | The number of elements must be num_polarizations. Each element must be in range [0,num_polarizations_for_grid). must-be-aligned |
[in] | num_channels | The number of channels. It must be 0 < num_channels <= INT32_MAX . |
[in] | channel_map | The number of elements must be num_channels. Each element must be in range [0,num_channels_for_grid). must-be-aligned |
[in] | mask | Masks to each value . Its memory layout should be [num_spectra][num_polarizations][num_channels]. If a value of element is false, the corresponding combination of the spectrum, polarization and channel is ignored. must-be-aligned |
[in] | value | Values to be gridded. Its memory layout should be [num_spectra][num_polarizations][num_channels]. must-be-aligned |
[in] | weight | Weights for value. Its memory layout should be [num_spectra][num_channels]. must-be-aligned |
[in] | weight_only | True if you want to get a grid of weight itself rather than production of value and weight. Otherwise false. |
[in] | num_convolution_table | The number of elements of convolution_table. It should be ceil(sqrt(2.)*(support+1)*sampling) <= num_convolution_table <= INT32_MAX / 32 . |
[in] | convolution_table | An array which represents convolution kernel. The number of elements must be num_convolution_table. The first element corresponds to center of the point. must-be-aligned |
[in] | num_polarizations_for_grid | The number of polarizations on the grid. It should be 0 < num_polarizations_for_grid <= INT32_MAX . |
[in] | num_channels_for_grid | The number of channels on the grid. It should be 0 < num_channels_for_grid <= INT32_MAX . |
[in] | width | Width of the grid . It should be 0 < width <= INT32_MAX . |
[in] | height | Height of the grid . It should be 0 < height <= INT32_MAX . |
[out] | weight_sum | Sum of weights. Its memory layout should be [num_polarizations_for_grid][num_channels_for_grid]. must-be-aligned |
[out] | weight_of_grid | Weight for each grid. Its memory layout should be [height][width][num_polarizations_for_grid][num_channels_for_grid]. must-be-aligned |
[out] | grid | The resulting grid. Its memory layout should be [height][width][num_polarizations_for_grid][num_channels_for_grid]. must-be-aligned |
MT-safe
sakura_Status sakura_Initialize | ( | sakura_UserAllocator | allocator, |
sakura_UserDeallocator | deallocator | ||
) |
Initializes Sakura Library.
You must initialize libsakura by calling this function before calling any other function of Sakura Library.
Without calling sakura_CleanUp() , don't call this function again.
[in] | allocator | An allocator which is used when Sakura Library needs to allocate memory dynamically. posix_memalign(3) is used if NULL is provided. See sakura_UserAllocator . |
[in] | deallocator | A deallocator which is used when Sakura Library needs to free dynamically allocated memory. free(3) is used if NULL is provided. See sakura_UserDeallocator . |
MT-unsafe
sakura_Status sakura_InterpolateXAxisFloat | ( | sakura_InterpolationMethod | interpolation_method, |
uint8_t | polynomial_order, | ||
size_t | num_base, | ||
double const | base_position[], | ||
size_t | num_array, | ||
float const | base_data[], | ||
bool const | base_mask[], | ||
size_t | num_interpolated, | ||
double const | interpolated_position[], | ||
float | interpolated_data[], | ||
bool | interpolated_mask[] | ||
) |
Perform one-dimensional interpolation.
It performs one-dimensional interpolation based on two input arrays base_position and base_data. Size of input array is num_base for base_position and num_interpolated times num_array for base_data, where num_array is a number of arrays that are passed to the function so that interpolation on multiple arrays can be performed simultaneously. One can set boolean mask for base_data using base_mask, which has same array shape as base_data. Mask value is true if data is valid while the value is false if the data is invalid. Invalid data will be excluded from interpolation.
List of locations where interpolated value is evaluated have to be specified via interpolate_position whose length is num_interpolated. Interpolation result on each point specified by interpolate_position is stored to interpolated_data. No extrapolation will be performed. Instead, out of range points are filled by the value of nearest points. Output boolean mask is interpolated_mask. Data will be invalid if corresponding mask is false, i.e., the value is just a nominal one but a result of actual interpolation. The mask will be false when interpolation is skipped due to insufficient number of valid data elements.
The function returns result status. In the successful run, returned value is sakura_Status_kOK while appropriate error status will be returned for failure: sakura_Status_kInvalidArgument for invalid input arguments, sakura_Status_kNoMemory for memory allocation error for internal variables, and sakura_Status_kUnknownError for other unknown error. Possible reason for sakura_Status_kInvalidArgument is either of (1) invalid interpolation_method or (2) arrays are not aligned.
data0[0], data0[1], ..., data1[0], data1[1], ...On the other hand, the latter requires [num_base][num_array], i.e,
data0[0], data1[0], ..., data0[1], data1[1], ...Result array, interpolated_data, follows the memory layout required for base_data.
[in] | interpolation_method | Interpolation method. |
[in] | polynomial_order | Maximum polynomial order for polynomial interpolation. Actual order will be determined by a balance between polynomial_order and num_base. This parameter is effective only when interpolation_method is sakura_InterpolationMethod_kPolynomial . In other interpolation methods, it is ignored. |
[in] | num_base | Number of elements for data points. Its value must be greater than 0. |
[in] | base_position | Position of data points. Its length must be num_base. It must be sorted either ascending or descending. must-be-aligned |
[in] | num_array | Number of arrays given in base_data. |
[in] | base_data | Value of data points. Its length must be num_base times num_array. must-be-aligned |
[in] | base_mask | Boolean mask for data. Its length must be num_base times num_array. False points will be excluded from the interpolation must-be-aligned |
[in] | num_interpolated | Number of elements for points that wants to get interpolated value. |
[in] | interpolated_position | Location of points that wants to get interpolated value. Its length must be num_interpolated. must-be-aligned |
[out] | interpolated_data | Storage for interpolation result. Its length must be num_interpolated times num_array. must-be-aligned |
[out] | interpolated_mask | Boolean mask for interpolation result. Its length must be num_interpolated times num_array. must-be-aligned |
MT-safe
sakura_Status sakura_InterpolateYAxisFloat | ( | sakura_InterpolationMethod | interpolation_method, |
uint8_t | polynomial_order, | ||
size_t | num_base, | ||
double const | base_position[], | ||
size_t | num_array, | ||
float const | base_data[], | ||
bool const | base_mask[], | ||
size_t | num_interpolated, | ||
double const | interpolated_position[], | ||
float | interpolated_data[], | ||
bool | interpolated_mask[] | ||
) |
Perform one-dimensional interpolation.
It performs one-dimensional interpolation based on two input arrays base_position and base_data. Size of input array is num_base for base_position and num_interpolated times num_array for base_data, where num_array is a number of arrays that are passed to the function so that interpolation on multiple arrays can be performed simultaneously. One can set boolean mask for base_data using base_mask, which has same array shape as base_data. Mask value is true if data is valid while the value is false if the data is invalid. Invalid data will be excluded from interpolation.
List of locations where interpolated value is evaluated have to be specified via interpolate_position whose length is num_interpolated. Interpolation result on each point specified by interpolate_position is stored to interpolated_data. No extrapolation will be performed. Instead, out of range points are filled by the value of nearest points. Output boolean mask is interpolated_mask. Data will be invalid if corresponding mask is false, i.e., the value is just a nominal one but a result of actual interpolation. The mask will be false when interpolation is skipped due to insufficient number of valid data elements.
The function returns result status. In the successful run, returned value is sakura_Status_kOK while appropriate error status will be returned for failure: sakura_Status_kInvalidArgument for invalid input arguments, sakura_Status_kNoMemory for memory allocation error for internal variables, and sakura_Status_kUnknownError for other unknown error. Possible reason for sakura_Status_kInvalidArgument is either of (1) invalid interpolation_method or (2) arrays are not aligned.
data0[0], data0[1], ..., data1[0], data1[1], ...On the other hand, the latter requires [num_base][num_array], i.e,
data0[0], data1[0], ..., data0[1], data1[1], ...Result array, interpolated_data, follows the memory layout required for base_data.
[in] | interpolation_method | Interpolation method. |
[in] | polynomial_order | Maximum polynomial order for polynomial interpolation. Actual order will be determined by a balance between polynomial_order and num_base. This parameter is effective only when interpolation_method is sakura_InterpolationMethod_kPolynomial . In other interpolation methods, it is ignored. |
[in] | num_base | Number of elements for data points. Its value must be greater than 0. |
[in] | base_position | Position of data points. Its length must be num_base. It must be sorted either ascending or descending. must-be-aligned |
[in] | num_array | Number of arrays given in base_data. |
[in] | base_data | Value of data points. Its length must be num_base times num_array. must-be-aligned |
[in] | base_mask | Boolean mask for data. Its length must be num_base times num_array. False points will be excluded from the interpolation must-be-aligned |
[in] | num_interpolated | Number of elements for points that wants to get interpolated value. |
[in] | interpolated_position | Location of points that wants to get interpolated value. Its length must be num_interpolated. must-be-aligned |
[out] | interpolated_data | Storage for interpolation result. Its length must be num_interpolated times num_array. must-be-aligned |
[out] | interpolated_mask | Boolean mask for interpolation result. Its length must be num_interpolated times num_array. must-be-aligned |
MT-safe
sakura_Status sakura_InvertBool | ( | size_t | num_data, |
bool const | data[], | ||
bool | result[] | ||
) |
Invert a boolean array.
[in] | num_data | The number of elements in the arrays, data and result |
[in] | data | The input array of of size, num_data. must-be-aligned |
[out] | result | The output array of of size, num_data. The pointer of out is allowed to be equal to that of in (result == data), indicating in-place operation. must-be-aligned |
MT-safe
bool sakura_IsAligned | ( | void const * | ptr | ) |
Checks if ptr points the aligned address Sakura Library requires.
[in] | ptr | An address to be checked. NULL is allowed. |
MT-safe
sakura_Status sakura_LMFitGaussianFloat | ( | size_t const | num_data, |
float const | data[], | ||
bool const | mask[], | ||
size_t const | num_peaks, | ||
double | height[], | ||
double | err_height[], | ||
double | center[], | ||
double | err_center[], | ||
double | sigma[], | ||
double | err_sigma[] | ||
) |
Fit Gaussian to the input data using Levenberg-Marquardt method.
Fit Gaussian to the input data using Levenberg-Marquardt method. The input data can contain multiple Gaussian profiles, however, note that fitting multiple Gaussians at once is quite difficult and thus the results may not be correct.
[in] | num_data | Number of input data. It must be equal to or larger than the number of Gaussian parameters (3 * num_peaks ). |
[in] | data | The input data with length of num_data . |
[in] | mask | The input mask data with length of num_data . |
[in] | num_peaks | Number of Gaussians to be fitted. |
[in,out] | height | An array to store initial guess of Gaussian heights. Its length must be num_peaks and will be overwritten with the fitting results. |
[in,out] | center | An array to store initial guess of Gaussian centers. Its length must be num_peaks and will be overwritten with the fitting results. |
[in,out] | sigma | An array to store initial guess of Gaussian sigmas. Its length must be num_peaks and will be overwritten with the fitting results. |
[out] | err_height | An array to store errors of fitted Gaussian heights. Its length must be num_peaks . |
[out] | err_center | An array to store errors of fitted Gaussian centers. Its length must be num_peaks . |
[out] | err_sigma | An array to store errors of fitted Gaussian sigmas. Its length must be num_peaks . |
MT-safe
sakura_Status sakura_LSQFitCubicSplineFloat | ( | struct sakura_LSQFitContextFloat const * | context, |
size_t | num_pieces, | ||
size_t | num_data, | ||
float const | data[], | ||
bool const | mask[], | ||
float | clip_threshold_sigma, | ||
uint16_t | num_fitting_max, | ||
double | coeff[][4], | ||
float | best_fit[], | ||
float | residual[], | ||
bool | final_mask[], | ||
float * | rms, | ||
size_t | boundary[], | ||
sakura_LSQFitStatus * | lsqfit_status | ||
) |
Fit a cubic spline curve to input data.
A cubic spline curve is fitted to input data based on Least-Square method. As output, coefficients of bases of the best-fit curve, the best-fit curve value itself and the residuals (i.e., input - best-fit curve) are stored in coeff , best_fit and residual , respectively. If num_fitting_max greater than 1 is given, fitting is executed recursively. Once best-fit curve is subtracted from input data, mean and standard deviation of the residual is computed to update mask so that mask has false value at data points where residual value exceeds threshold defined by clip_threshold_sigma , then Least-Square fitting is again performed to the input data with updated mask to update coeff values. This procedure is repeatedly done for ( num_fitting_max-1 ) times or until the fitting result converges. Once fitting is done, the updated mask information is stored in final_mask .
[in] | context | A context created by sakura_CreateLSQFitContextCubicSplineFloat . |
[in] | num_pieces | Number of spline pieces. It must be positive and also must not exceed the number of spline pieces specified in creation of context . |
[in] | num_data | Number of elements in the arrays data , mask , best_fit , residual , and final_mask . It must be equal to num_data which was given to sakura_CreateLSQFitContextCubicSplineFloat to create context . |
[in] | data | Input data with length of num_data . must-be-aligned |
[in] | mask | Input mask data with length of num_data . The i th element of data is included in input spectrum if the i th element of mask is true, while it is excluded from input spectrum if the corresponding element of mask is false. must-be-aligned |
[in] | clip_threshold_sigma | Threshold of clipping in unit of sigma. It must be a positive value. |
[in] | num_fitting_max | Upper limit of how many times fitting is performed recursively. Before executing the second or later fitting, outlier in data is masked via clipping to be not used. If 1 is given, fitting is done just once and no clipping will be applied. If 0 is given, fitting is not executed and the values of residual should be identical with those of data . |
[out] | coeff | The coefficients of the best-fit cubic spline curve. It must be a 2D array with type of double[ num_pieces ][4]. coeff[i] is for the i th spline piece from the left side. The first element in each coeff[i] is of constant term, followed by those of first, second, and third orders. Null pointer can be given in case users do not need this value. must-be-aligned |
[out] | best_fit | The best-fit curve data, i.e., the result of least-square fitting itself. Its length must be num_data . data can be set to best_fit if users want to overwrite it in-place. Null pointer can be given in case users do not need this value. must-be-aligned |
[out] | residual | Residual (input - best-fit) data. Its length must be num_data . data can be set to residual if users want to overwrite it in-place. Null pointer can be given in case users do not need this value. must-be-aligned |
[out] | final_mask | The final status of mask data after recursive clipping procedure finish. Its length must be num_data . mask can be set to final_mask if users want to overwrite it in-place. must-be-aligned |
[out] | rms | The root-mean-square of residual . |
[out] | boundary | A 1D array containing the boundary positions of spline pieces, which are indices of data . Its length must be ( num_pieces +1). The element values will be stored in ascending order. The first element will always be zero, the left edge of the first (left-most) spline piece, while the last element will be num_data , which is the next of the right edge of the last (right-most) spline piece. must-be-aligned |
[out] | lsqfit_status | LSQFit-specific error code. |
MT-safe
sakura_Status sakura_LSQFitPolynomialFloat | ( | struct sakura_LSQFitContextFloat const * | context, |
uint16_t | order, | ||
size_t | num_data, | ||
float const | data[], | ||
bool const | mask[], | ||
float | clip_threshold_sigma, | ||
uint16_t | num_fitting_max, | ||
size_t | num_coeff, | ||
double | coeff[], | ||
float | best_fit[], | ||
float | residual[], | ||
bool | final_mask[], | ||
float * | rms, | ||
sakura_LSQFitStatus * | lsqfit_status | ||
) |
Fit a polynomial curve to input data.
A polynomial or Chebyshev polynomial is fitted to input data based on Least-Square method. As output, coefficients of bases of the best-fit curve, the best-fit curve value itself and the residuals (i.e., input - best-fit curve) are stored in coeff , best_fit and residual , respectively. If num_fitting_max greater than 1 is given, fitting is executed recursively. Once best-fit curve is subtracted from input data, mean and standard deviation of the residual is computed to update mask so that mask has false value at data points where residual value exceeds threshold defined by clip_threshold_sigma , then Least-Square fitting is again performed to the input data with updated mask to update coeff values. This procedure is repeatedly done for ( num_fitting_max-1 ) times or until the fitting result converges. Once fitting is done, the updated mask information is stored in final_mask .
[in] | context | A context created by sakura_CreateLSQFitContextPolynomialFloat . |
[in] | order | Polynomial order. It should not exceed the order specified in creation of context . |
[in] | num_data | Number of elements in the arrays data , mask , best_fit , residual , and final_mask . It must be equal to num_data which was given to sakura_CreateLSQFitContextPolynomialFloat to create context . |
[in] | data | Input data with length of num_data . must-be-aligned |
[in] | mask | Input mask data with length of num_data . The i th element of data is included in input spectrum if the i th element of mask is true, while it is excluded from input spectrum if the corresponding element of mask is false. must-be-aligned |
[in] | clip_threshold_sigma | Threshold of clipping in unit of sigma. It must be a positive value. |
[in] | num_fitting_max | Upper limit of how many times fitting is performed recursively. Before executing the second or later fitting, outlier in data is masked via clipping to be not used. If 1 is given, fitting is done just once and no clipping will be applied. If 0 is given, fitting is not executed and the values of residual should be identical with those of data . |
[in] | num_coeff | Number of elements in the array coeff . In case coeff is not null pointer, it must be equal to ( order +1), while the value is not checked when coeff is null pointer. |
[out] | coeff | Coefficients of the polynomial bases. Its length must be num_coeff . Coefficient values are stored in ascending order of polynomial order: coefficient of constant term comes first, followed by those of first order, second order, and so on. Null pointer can be given in case users do not need this value. must-be-aligned |
[out] | best_fit | The best-fit curve data, i.e., the result of least-square fitting itself. Its length must be num_data . data can be set to best_fit if users want to overwrite it in-place. Null pointer can be given in case users do not need this value. must-be-aligned |
[out] | residual | Residual (input - best-fit) data. Its length must be num_data . data can be set to residual if users want to overwrite it in-place. Null pointer can be given in case users do not need this value. must-be-aligned |
[out] | final_mask | The final status of mask data after recursive clipping procedure finish. Its length must be num_data . mask can be set to final_mask if users want to overwrite it in-place. must-be-aligned |
[out] | rms | The root-mean-square of residual . |
[out] | lsqfit_status | LSQFit-specific error code. |
MT-safe
sakura_Status sakura_LSQFitSinusoidFloat | ( | struct sakura_LSQFitContextFloat const * | context, |
size_t | num_nwave, | ||
size_t const | nwave[], | ||
size_t | num_data, | ||
float const | data[], | ||
bool const | mask[], | ||
float | clip_threshold_sigma, | ||
uint16_t | num_fitting_max, | ||
size_t | num_coeff, | ||
double | coeff[], | ||
float | best_fit[], | ||
float | residual[], | ||
bool | final_mask[], | ||
float * | rms, | ||
sakura_LSQFitStatus * | lsqfit_status | ||
) |
Fit a sinusoidal curve to input data.
A sinusoidal curve is fitted to input data based on Least-Square method. As output, coefficients of bases of the best-fit curve, the best-fit curve value itself and the residuals (i.e., input - best-fit curve) are stored in coeff , best_fit and residual , respectively. If num_fitting_max greater than 1 is given, fitting is executed recursively. Once best-fit curve is subtracted from input data, mean and standard deviation of the residual is computed to update mask so that mask has false value at data points where residual value exceeds threshold defined by clip_threshold_sigma , then Least-Square fitting is again performed to the input data with updated mask to update coeff values. This procedure is repeatedly done for ( num_fitting_max-1 ) times or until the fitting result converges. Once fitting is done, the updated mask information is stored in final_mask .
[in] | context | A context created by sakura_CreateLSQFitContextSinusoidFloat . |
[in] | num_nwave | The number of elements in the array nwave . |
[in] | nwave | Wave numbers within the index range of data to be used for sinusoidal fitting. The values must be positive or zero (for constant term), but not exceed the maximum wave number specified in creation of context . The values must be stored in ascending order and must not be duplicate. |
[in] | num_data | Number of elements in the arrays data , mask , best_fit , residual , and final_mask . It must be equal to num_data which was given to sakura_CreateLSQFitContextSinusoidFloat to create context . |
[in] | data | Input data with length of num_data . must-be-aligned |
[in] | mask | Input mask data with length of num_data . The i th element of data is included in input spectrum if the i th element of mask is true, while it is excluded from input spectrum if the corresponding element of mask is false. must-be-aligned |
[in] | clip_threshold_sigma | Threshold of clipping in unit of sigma. It must be a positive value. |
[in] | num_fitting_max | Upper limit of how many times fitting is performed recursively. Before executing the second or later fitting, outlier in data is masked via clipping to be not used. If 1 is given, fitting is done just once and no clipping will be applied. If 0 is given, fitting is not executed and the values of residual should be identical with those of data . |
[in] | num_coeff | The number of elements in the array coeff. If coeff is not null pointer, it must be ( num_nwave*2-1 ) or ( num_nwave*2 ) in cases nwave contains zero or not, respectively, and must not exceed num_data, while the value is not checked when coeff is null pointer. |
[out] | coeff | Coefficients of the sinusoidal fit. Its length must be num_coeff . If nwave contains zero, the first element is of the constant term. If not, the first and the second elements are for sine and cosine terms for the smallest wave number in nwave, respectively. Coefficients of sine and cosine terms for the second-smallest wave number come next, and so on. Null pointer can be given in case users do not need this value. must-be-aligned |
[out] | best_fit | The best-fit curve data, i.e., the result of least-square fitting itself. Its length must be num_data . data can be set to best_fit if users want to overwrite it in-place. Null pointer can be given in case users do not need this value. must-be-aligned |
[out] | residual | Residual (input - best-fit) data. Its length must be num_data . data can be set to residual if users want to overwrite it in-place. Null pointer can be given in case users do not need this value. must-be-aligned |
[out] | final_mask | The final status of mask data after recursive clipping procedure finish. Its length must be num_data . mask can be set to final_mask if users want to overwrite it in-place. must-be-aligned |
[out] | rms | The root-mean-square of residual . |
[out] | lsqfit_status | LSQFit-specific error code. |
MT-safe
sakura_Status sakura_OperateBitwiseAndUint32 | ( | uint32_t | bit_mask, |
size_t | num_data, | ||
uint32_t const | data[], | ||
bool const | edit_mask[], | ||
uint32_t | result[] | ||
) |
Invoke bit operation AND between a bit mask and an array.
Invokes the following bit operation to i- th element of result :
[in] | bit_mask | A bit mask. The bit operation is invoked between this value and the array, data. |
[in] | num_data | The number of elements in the arrays, data, edit_mask, and result. |
[in] | data | An input array of size, num_data. The bit operation is invoked between this array and bit_mask. must-be-aligned |
[in] | edit_mask | A boolean mask array of size, num_data. The bit operation is skipped for the elements with the value, false. must-be-aligned |
[out] | result | The output array of size, num_data. It stores the result of the bit operation between bit_mask and data. The bit operation is skipped and the value in array, data, is adopted for the elements where corresponding elements in edit_mask is false. The pointer of out is allowed to be equal to that of in (result == data), indicating in-place operation. must-be-aligned |
MT-safe
sakura_Status sakura_OperateBitwiseAndUint8 | ( | uint8_t | bit_mask, |
size_t | num_data, | ||
uint8_t const | data[], | ||
bool const | edit_mask[], | ||
uint8_t | result[] | ||
) |
Invoke bit operation AND between a bit mask and an array.
Invokes the following bit operation to i- th element of result :
[in] | bit_mask | A bit mask. The bit operation is invoked between this value and the array, data. |
[in] | num_data | The number of elements in the arrays, data, edit_mask, and result. |
[in] | data | An input array of size, num_data. The bit operation is invoked between this array and bit_mask. must-be-aligned |
[in] | edit_mask | A boolean mask array of size, num_data. The bit operation is skipped for the elements with the value, false. must-be-aligned |
[out] | result | The output array of size, num_data. It stores the result of the bit operation between bit_mask and data. The bit operation is skipped and the value in array, data, is adopted for the elements where corresponding elements in edit_mask is false. The pointer of out is allowed to be equal to that of in (result == data), indicating in-place operation. must-be-aligned |
MT-safe
sakura_Status sakura_OperateBitwiseConverseNonImplicationUint32 | ( | uint32_t | bit_mask, |
size_t | num_data, | ||
uint32_t const | data[], | ||
bool const | edit_mask[], | ||
uint32_t | result[] | ||
) |
Invoke bit operation, Converse Nonimplication, between a bit mask and an array.
Invokes the following bit operation to the i- th element of result :
[in] | bit_mask | A bit mask. The bit operation is invoked between this value and the array, data. |
[in] | num_data | The number of elements in the arrays, data, edit_mask, and result. |
[in] | data | An input array of size, num_data. The bit operation is invoked between this array and bit_mask. must-be-aligned |
[in] | edit_mask | A boolean mask array of size, num_data. The bit operation is skipped for the elements with the value, false. must-be-aligned |
[out] | result | The output array of size, num_data. It stores the result of the bit operation between bit_mask and data. The bit operation is skipped and the value in array, data, is adopted for the elements where corresponding elements in edit_mask is false. The pointer of out is allowed to be equal to that of in (result == data), indicating in-place operation. must-be-aligned |
MT-safe
sakura_Status sakura_OperateBitwiseConverseNonImplicationUint8 | ( | uint8_t | bit_mask, |
size_t | num_data, | ||
uint8_t const | data[], | ||
bool const | edit_mask[], | ||
uint8_t | result[] | ||
) |
Invoke bit operation, Converse Nonimplication, between a bit mask and an array.
Invokes the following bit operation to the i- th element of result :
[in] | bit_mask | A bit mask. The bit operation is invoked between this value and the array, data. |
[in] | num_data | The number of elements in the arrays, data, edit_mask, and result. |
[in] | data | An input array of size, num_data. The bit operation is invoked between this array and bit_mask. must-be-aligned |
[in] | edit_mask | A boolean mask array of size, num_data. The bit operation is skipped for the elements with the value, false. must-be-aligned |
[out] | result | The output array of size, num_data. It stores the result of the bit operation between bit_mask and data. The bit operation is skipped and the value in array, data, is adopted for the elements where corresponding elements in edit_mask is false. The pointer of out is allowed to be equal to that of in (result == data), indicating in-place operation. must-be-aligned |
MT-safe
sakura_Status sakura_OperateBitwiseImplicationUint32 | ( | uint32_t | bit_mask, |
size_t | num_data, | ||
uint32_t const | data[], | ||
bool const | edit_mask[], | ||
uint32_t | result[] | ||
) |
Invoke bit operation, Material Implication, between a bit mask and an array.
Invokes the following bit operation to the i- th element of result :
[in] | bit_mask | A bit mask. The bit operation is invoked between this value and the array, data. |
[in] | num_data | The number of elements in the arrays, data, edit_mask, and result. |
[in] | data | An input array of size, num_data. The bit operation is invoked between this array and bit_mask. must-be-aligned |
[in] | edit_mask | A boolean mask array of size, num_data. The bit operation is skipped for the elements with the value, false. must-be-aligned |
[out] | result | The output array of size, num_data. It stores the result of the bit operation between bit_mask and data. The bit operation is skipped and the value in array, data, is adopted for the elements where corresponding elements in edit_mask is false. The pointer of out is allowed to be equal to that of in (result == data), indicating in-place operation. must-be-aligned |
MT-safe
sakura_Status sakura_OperateBitwiseImplicationUint8 | ( | uint8_t | bit_mask, |
size_t | num_data, | ||
uint8_t const | data[], | ||
bool const | edit_mask[], | ||
uint8_t | result[] | ||
) |
Invoke bit operation, Material Implication, between a bit mask and an array.
Invokes the following bit operation to the i- th element of result :
[in] | bit_mask | A bit mask. The bit operation is invoked between this value and the array, data. |
[in] | num_data | The number of elements in the arrays, data, edit_mask, and result. |
[in] | data | An input array of size, num_data. The bit operation is invoked between this array and bit_mask. must-be-aligned |
[in] | edit_mask | A boolean mask array of size, num_data. The bit operation is skipped for the elements with the value, false. must-be-aligned |
[out] | result | The output array of size, num_data. It stores the result of the bit operation between bit_mask and data. The bit operation is skipped and the value in array, data, is adopted for the elements where corresponding elements in edit_mask is false. The pointer of out is allowed to be equal to that of in (result == data), indicating in-place operation. must-be-aligned |
MT-safe
sakura_Status sakura_OperateBitwiseNotUint32 | ( | size_t | num_data, |
uint32_t const | data[], | ||
bool const | edit_mask[], | ||
uint32_t | result[] | ||
) |
Invoke bitwise NOT operation of an array.
Invokes the following bit operation to i- th element of result :
[in] | num_data | The number of elements in the arrays, data, edit_mask, and result. |
[in] | data | An input array of size, num_data. The bit operation is invoked to this array. must-be-aligned |
[in] | edit_mask | A boolean mask array of size, num_data. The bit operation is skipped for the elements with the value, false. must-be-aligned |
[out] | result | The output array of size, num_data. It stores the result of the bit operation to data. The bit operation is skipped and the value in array, data, is adopted for the elements where corresponding elements in edit_mask is false. The pointer of out is allowed to be equal to that of in (result == data), indicating in-place operation. must-be-aligned |
MT-safe
sakura_Status sakura_OperateBitwiseNotUint8 | ( | size_t | num_data, |
uint8_t const | data[], | ||
bool const | edit_mask[], | ||
uint8_t | result[] | ||
) |
Invoke bitwise NOT operation of an array.
Invokes the following bit operation to i- th element of result :
[in] | num_data | The number of elements in the arrays, data, edit_mask, and result. |
[in] | data | An input array of size, num_data. The bit operation is invoked to this array. must-be-aligned |
[in] | edit_mask | A boolean mask array of size, num_data. The bit operation is skipped for the elements with the value, false. must-be-aligned |
[out] | result | The output array of size, num_data. It stores the result of the bit operation to data. The bit operation is skipped and the value in array, data, is adopted for the elements where corresponding elements in edit_mask is false. The pointer of out is allowed to be equal to that of in (result == data), indicating in-place operation. must-be-aligned |
MT-safe
sakura_Status sakura_OperateBitwiseOrUint32 | ( | uint32_t | bit_mask, |
size_t | num_data, | ||
uint32_t const | data[], | ||
bool const | edit_mask[], | ||
uint32_t | result[] | ||
) |
Invoke bit operation OR between a bit mask and an array.
Invokes the following bit operation to i- th element of result :
[in] | bit_mask | A bit mask. The bit operation is invoked between this value and the array, data. |
[in] | num_data | The number of elements in the arrays, data, edit_mask, and result. |
[in] | data | An input array of size, num_data. The bit operation is invoked between this array and bit_mask. must-be-aligned |
[in] | edit_mask | A boolean mask array of size, num_data. The bit operation is skipped for the elements with the value, false. must-be-aligned |
[out] | result | The output array of size, num_data. It stores the result of the bit operation between bit_mask and data. The bit operation is skipped and the value in array, data, is adopted for the elements where corresponding elements in edit_mask is false. The pointer of out is allowed to be equal to that of in (result == data), indicating in-place operation. must-be-aligned |
MT-safe
sakura_Status sakura_OperateBitwiseOrUint8 | ( | uint8_t | bit_mask, |
size_t | num_data, | ||
uint8_t const | data[], | ||
bool const | edit_mask[], | ||
uint8_t | result[] | ||
) |
Invoke bit operation OR between a bit mask and an array.
Invokes the following bit operation to i- th element of result :
[in] | bit_mask | A bit mask. The bit operation is invoked between this value and the array, data. |
[in] | num_data | The number of elements in the arrays, data, edit_mask, and result. |
[in] | data | An input array of size, num_data. The bit operation is invoked between this array and bit_mask. must-be-aligned |
[in] | edit_mask | A boolean mask array of size, num_data. The bit operation is skipped for the elements with the value, false. must-be-aligned |
[out] | result | The output array of size, num_data. It stores the result of the bit operation between bit_mask and data. The bit operation is skipped and the value in array, data, is adopted for the elements where corresponding elements in edit_mask is false. The pointer of out is allowed to be equal to that of in (result == data), indicating in-place operation. must-be-aligned |
MT-safe
sakura_Status sakura_OperateBitwiseXorUint32 | ( | uint32_t | bit_mask, |
size_t | num_data, | ||
uint32_t const | data[], | ||
bool const | edit_mask[], | ||
uint32_t | result[] | ||
) |
Invoke bit operation XOR between a bit mask and an array.
Invokes the following bit operation to i- th element of result :
[in] | bit_mask | A bit mask. The bit operation is invoked between this value and the array, data. |
[in] | num_data | The number of elements in the arrays, data, edit_mask, and result. |
[in] | data | An input array of size, num_data. The bit operation is invoked between this array and bit_mask. must-be-aligned |
[in] | edit_mask | A boolean mask array of size, num_data. The bit operation is skipped for the elements with the value, false. must-be-aligned |
[out] | result | The output array of size, num_data. It stores the result of the bit operation between bit_mask and data. The bit operation is skipped and the value in array, data, is adopted for the elements where corresponding elements in edit_mask is false. The pointer of out is allowed to be equal to that of in (result == data), indicating in-place operation. must-be-aligned |
MT-safe
sakura_Status sakura_OperateBitwiseXorUint8 | ( | uint8_t | bit_mask, |
size_t | num_data, | ||
uint8_t const | data[], | ||
bool const | edit_mask[], | ||
uint8_t | result[] | ||
) |
Invoke bit operation XOR between a bit mask and an array.
Invokes the following bit operation to i- th element of result :
[in] | bit_mask | A bit mask. The bit operation is invoked between this value and the array, data. |
[in] | num_data | The number of elements in the arrays, data, edit_mask, and result. |
[in] | data | An input array of size, num_data. The bit operation is invoked between this array and bit_mask. must-be-aligned |
[in] | edit_mask | A boolean mask array of size, num_data. The bit operation is skipped for the elements with the value, false. must-be-aligned |
[out] | result | The output array of size, num_data. It stores the result of the bit operation between bit_mask and data. The bit operation is skipped and the value in array, data, is adopted for the elements where corresponding elements in edit_mask is false. The pointer of out is allowed to be equal to that of in (result == data), indicating in-place operation. must-be-aligned |
MT-safe
sakura_Status sakura_SetFalseIfNanOrInfFloat | ( | size_t | num_data, |
float const | data[], | ||
bool | result[] | ||
) |
Sets false if the values in the input array (data ) are either NaN or infinity.
Elements of the output array are set to false if the corresponding element in the input array is not a number (NaN) or infinity. Otherwise, they are set to true.
[in] | num_data | The number of elements in the arrays, data and result |
[in] | data | The input array of of size, num_data. must-be-aligned |
[out] | result | The output array of of size, num_data. must-be-aligned |
MT-safe
sakura_Status sakura_SetTrueIfGreaterThanFloat | ( | size_t | num_data, |
float const | data[], | ||
float | threshold, | ||
bool | result[] | ||
) |
Sets true if the values in the input array (data ) are greater than a threshold.
Elements of the output array are set to true if the corresponding element in the input array is greater than a threshold,
otherwise they are set to false.
[in] | num_data | The number of elements in the arrays, data and result |
[in] | data | An input array of size, num_data. In case the array is floating-point type, the elements should not contain Inf nor NaN. must-be-aligned |
[in] | threshold | The threshold of evaluation. In case the parameter is floating-point type, the value should not be Inf nor NaN. |
[out] | result | The output array of size, num_data. must-be-aligned |
MT-safe
sakura_Status sakura_SetTrueIfGreaterThanInt | ( | size_t | num_data, |
int const | data[], | ||
int | threshold, | ||
bool | result[] | ||
) |
Sets true if the values in the input array (data ) are greater than a threshold.
Elements of the output array are set to true if the corresponding element in the input array is greater than a threshold,
otherwise they are set to false.
[in] | num_data | The number of elements in the arrays, data and result |
[in] | data | An input array of size, num_data. In case the array is floating-point type, the elements should not contain Inf nor NaN. must-be-aligned |
[in] | threshold | The threshold of evaluation. In case the parameter is floating-point type, the value should not be Inf nor NaN. |
[out] | result | The output array of size, num_data. must-be-aligned |
MT-safe
sakura_Status sakura_SetTrueIfGreaterThanOrEqualsFloat | ( | size_t | num_data, |
float const | data[], | ||
float | threshold, | ||
bool | result[] | ||
) |
Sets true if the values in the input array (data ) are greater than or equal to a threshold.
Elements of the output array are set to true if the corresponding element in the input array is greater than or equals to a threshold,
otherwise they are set to false.
[in] | num_data | The number of elements in the arrays, data and result |
[in] | data | An input array of size, num_data. In case the array is floating-point type, the elements should not contain Inf nor NaN. must-be-aligned |
[in] | threshold | The threshold of evaluation. In case the parameter is floating-point type, the value should not be Inf nor NaN. |
[out] | result | The output array of size, num_data. must-be-aligned |
MT-safe
sakura_Status sakura_SetTrueIfGreaterThanOrEqualsInt | ( | size_t | num_data, |
int const | data[], | ||
int | threshold, | ||
bool | result[] | ||
) |
Sets true if the values in the input array (data ) are greater than or equal to a threshold.
Elements of the output array are set to true if the corresponding element in the input array is greater than or equals to a threshold,
otherwise they are set to false.
[in] | num_data | The number of elements in the arrays, data and result |
[in] | data | An input array of size, num_data. In case the array is floating-point type, the elements should not contain Inf nor NaN. must-be-aligned |
[in] | threshold | The threshold of evaluation. In case the parameter is floating-point type, the value should not be Inf nor NaN. |
[out] | result | The output array of size, num_data. must-be-aligned |
MT-safe
sakura_Status sakura_SetTrueIfInRangesExclusiveFloat | ( | size_t | num_data, |
float const | data[], | ||
size_t | num_condition, | ||
float const | lower_bounds[], | ||
float const | upper_bounds[], | ||
bool | result[] | ||
) |
Sets true if the values in an input array (data ) are in any of specified range (exclusive).
Elements of the output array are set to true if the corresponding element in the input array is in range of upper and lower boundary pairs,
otherwise they are set to false.
The function takes more than one upper and lower boundary pairs as arrays, lower_bounds and upper_bounds.
[in] | num_data | The number of elements in the arrays, data and result |
[in] | data | An input array of size, num_data. In case the array is floating-point type, the elements should not contain Inf nor NaN. must-be-aligned |
[in] | num_condition | The number of elements in the arrays, lower_bounds and upper_bounds. |
[in] | lower_bounds | The input array of size, num_condition. In case the array is floating-point type, the elements should not contain Inf nor NaN. must-be-aligned |
[in] | upper_bounds | The input array of size, num_condition. In case the array is floating-point type, the elements should not contain Inf nor NaN. must-be-aligned |
[out] | result | The output array of size, num_data. must-be-aligned |
MT-safe
sakura_Status sakura_SetTrueIfInRangesExclusiveInt | ( | size_t | num_data, |
int const | data[], | ||
size_t | num_condition, | ||
int const | lower_bounds[], | ||
int const | upper_bounds[], | ||
bool | result[] | ||
) |
Sets true if the values in an input array (data ) are in any of specified range (exclusive).
Elements of the output array are set to true if the corresponding element in the input array is in range of upper and lower boundary pairs,
otherwise they are set to false.
The function takes more than one upper and lower boundary pairs as arrays, lower_bounds and upper_bounds.
[in] | num_data | The number of elements in the arrays, data and result |
[in] | data | An input array of size, num_data. In case the array is floating-point type, the elements should not contain Inf nor NaN. must-be-aligned |
[in] | num_condition | The number of elements in the arrays, lower_bounds and upper_bounds. |
[in] | lower_bounds | The input array of size, num_condition. In case the array is floating-point type, the elements should not contain Inf nor NaN. must-be-aligned |
[in] | upper_bounds | The input array of size, num_condition. In case the array is floating-point type, the elements should not contain Inf nor NaN. must-be-aligned |
[out] | result | The output array of size, num_data. must-be-aligned |
MT-safe
sakura_Status sakura_SetTrueIfInRangesInclusiveFloat | ( | size_t | num_data, |
float const | data[], | ||
size_t | num_condition, | ||
float const | lower_bounds[], | ||
float const | upper_bounds[], | ||
bool | result[] | ||
) |
Sets true if the values in an input array (data ) are in any of specified range (inclusive).
Elements of the output array are set to true if the corresponding element in the input array is in range of upper and lower boundary pairs,
otherwise they are set to false.
The function takes more than one upper and lower boundary pairs as arrays, lower_bounds and upper_bounds.
[in] | num_data | The number of elements in the arrays, data and result . |
[in] | data | An input array of size, num_data. In case the array is floating-point type, the elements should not contain Inf nor NaN. must-be-aligned |
[in] | num_condition | The number of elements in the arrays, lower_bounds and upper_bounds. |
[in] | lower_bounds | The input array of size, num_condition. In case the array is floating-point type, the elements should not contain Inf nor NaN. must-be-aligned |
[in] | upper_bounds | The input array of size, num_condition. In case the array is floating-point type, the elements should not contain Inf nor NaN. must-be-aligned |
[out] | result | The output array of size, num_data. must-be-aligned |
MT-safe
sakura_Status sakura_SetTrueIfInRangesInclusiveInt | ( | size_t | num_data, |
int const | data[], | ||
size_t | num_condition, | ||
int const | lower_bounds[], | ||
int const | upper_bounds[], | ||
bool | result[] | ||
) |
Sets true if the values in an input array (data ) are in any of specified range (inclusive).
Elements of the output array are set to true if the corresponding element in the input array is in range of upper and lower boundary pairs,
otherwise they are set to false.
The function takes more than one upper and lower boundary pairs as arrays, lower_bounds and upper_bounds.
[in] | num_data | The number of elements in the arrays, data and result . |
[in] | data | An input array of size, num_data. In case the array is floating-point type, the elements should not contain Inf nor NaN. must-be-aligned |
[in] | num_condition | The number of elements in the arrays, lower_bounds and upper_bounds. |
[in] | lower_bounds | The input array of size, num_condition. In case the array is floating-point type, the elements should not contain Inf nor NaN. must-be-aligned |
[in] | upper_bounds | The input array of size, num_condition. In case the array is floating-point type, the elements should not contain Inf nor NaN. must-be-aligned |
[out] | result | The output array of size, num_data. must-be-aligned |
MT-safe
sakura_Status sakura_SetTrueIfLessThanFloat | ( | size_t | num_data, |
float const | data[], | ||
float | threshold, | ||
bool | result[] | ||
) |
Sets true if the values in the input array (data ) are less than a threshold.
Elements of the output array are set to true if the corresponding element in the input array is less than a threshold,
otherwise they are set to false.
[in] | num_data | The number of elements in the arrays, data and result |
[in] | data | An input array of size, num_data. In case the array is floating-point type, the elements should not contain Inf nor NaN. must-be-aligned |
[in] | threshold | The threshold of evaluation. In case the parameter is floating-point type, the value should not be Inf nor NaN. |
[out] | result | The output array of size, num_data. must-be-aligned |
MT-safe
sakura_Status sakura_SetTrueIfLessThanInt | ( | size_t | num_data, |
int const | data[], | ||
int | threshold, | ||
bool | result[] | ||
) |
Sets true if the values in the input array (data ) are less than a threshold.
Elements of the output array are set to true if the corresponding element in the input array is less than a threshold,
otherwise they are set to false.
[in] | num_data | The number of elements in the arrays, data and result |
[in] | data | An input array of size, num_data. In case the array is floating-point type, the elements should not contain Inf nor NaN. must-be-aligned |
[in] | threshold | The threshold of evaluation. In case the parameter is floating-point type, the value should not be Inf nor NaN. |
[out] | result | The output array of size, num_data. must-be-aligned |
MT-safe
sakura_Status sakura_SetTrueIfLessThanOrEqualsFloat | ( | size_t | num_data, |
float const | data[], | ||
float | threshold, | ||
bool | result[] | ||
) |
Sets true if the values in the input array (data ) are less than or equal to a threshold.
Elements of the output array are set to true if the corresponding element in the input array is less than or equals to a threshold,
otherwise they are set to false.
[in] | num_data | The number of elements in the arrays, data and result |
[in] | data | An input array of size, num_data. In case the array is floating-point type, the elements should not contain Inf nor NaN. must-be-aligned |
[in] | threshold | The threshold of evaluation. In case the parameter is floating-point type, the value should not be Inf nor NaN. |
[out] | result | The output array of size, num_data. must-be-aligned |
MT-safe
sakura_Status sakura_SetTrueIfLessThanOrEqualsInt | ( | size_t | num_data, |
int const | data[], | ||
int | threshold, | ||
bool | result[] | ||
) |
Sets true if the values in the input array (data ) are less than or equal to a threshold.
Elements of the output array are set to true if the corresponding element in the input array is less than or equals to a threshold,
otherwise they are set to false.
[in] | num_data | The number of elements in the arrays, data and result |
[in] | data | An input array of size, num_data. In case the array is floating-point type, the elements should not contain Inf nor NaN. must-be-aligned |
[in] | threshold | The threshold of evaluation. In case the parameter is floating-point type, the value should not be Inf nor NaN. |
[out] | result | The output array of size, num_data. must-be-aligned |
MT-safe
sakura_Status sakura_SolveSimultaneousEquationsByLUDouble | ( | size_t | num_equations, |
double const | in_matrix[], | ||
double const | in_vector[], | ||
double | out[] | ||
) |
Solve simultaneous equations via LU decomposition.
Suppose solving simultaneous equations A x = y to derive x, where A is a square matrix of num_equations rows and columns, and x and y are vectors with length of num_equations . Given A and y values, this function computes x values using LU decomposition of A.
[in] | num_equations | Number of equations. |
[in] | in_matrix | A 1D array containing values of the matrix A in the left side of the above simultaneous equations. Loop for columns comes inside that for rows, i.e., the value at the m -th row and n -th column is stored at in_matrix [ num_equations * ( m -1) + ( n -1)]. Its length must be (num_equations * num_equations). must-be-aligned |
[in] | in_vector | A 1D array containing values of the vector y in the right side of the above simultaneous equations. Its length must be num_equations . must-be-aligned |
[out] | out | The solution (x in the above equations). Its length must be num_equations . must-be-aligned |
MT-safe
sakura_Status sakura_SortValidValuesDenselyFloat | ( | size_t | num_data, |
bool const | is_valid[], | ||
float | data[], | ||
size_t * | new_num_data | ||
) |
Sorts only valid data in ascending order.
[in] | num_data | The number of elements in data and is_valid . |
[in] | is_valid | Masks of data. If a value of element is false, the corresponding element in data is ignored. |
[in,out] | data | Data to be sorted. Since data is sorted in place, contents of this array are not preserved. If corresponding element in is_valid is true, the element in data must not be Inf nor NaN. |
[out] | new_num_data | The number of sorted elements that don't include invalid data( <= num_data ) is stored here. |
MT-safe
sakura_Status sakura_SubtractCubicSplineFloat | ( | struct sakura_LSQFitContextFloat const * | context, |
size_t | num_data, | ||
float const | data[], | ||
size_t | num_pieces, | ||
double const | coeff[][4], | ||
size_t const | boundary[], | ||
float | out[] | ||
) |
Subtract best-fit cubic spline model from input data.
[in] | context | A context created by sakura_CreateLSQFitContextCubicSplineFloat . |
[in] | num_data | The number of elements in data and out . It must be equal to num_data which was given to sakura_CreateLSQFitContextCubicSplineFloat to create context . |
[in] | data | The input data with length of num_data . must-be-aligned |
[in] | num_pieces | The number of spline pieces. If zero is given, no subtraction executed. |
[in] | coeff | Coefficients of cubic spline curve. It must be a 2D array with type of double[ num_pieces ][4]. coeff[i] is for the i th spline piece from the left side. The first element in each coeff[i] must be of constant term, followed by the ones of first, second, and third orders. must-be-aligned |
[in] | boundary | A 1D array containing the boundary positions, which are indices of parameters data and out , of spline pieces. Its length must be ( num_pieces +1). The element values should be stored in ascending order. The first element must always be zero, the left edge of the first (left-most) spline piece, while the last element must be num_data , which is the next of the right edge of the last (right-most) spline piece. must-be-aligned |
[out] | out | The output data. Its length must be num_data . must-be-aligned |
MT-safe
sakura_Status sakura_SubtractPolynomialFloat | ( | struct sakura_LSQFitContextFloat const * | context, |
size_t | num_data, | ||
float const | data[], | ||
size_t | num_coeff, | ||
double const | coeff[], | ||
float | out[] | ||
) |
Subtract best-fit polynomial model from input data.
[in] | context | A context created by sakura_CreateLSQFitContextPolynomialFloat . |
[in] | num_data | The number of elements in data and out. It must be equal to num_data which was given to sakura_CreateLSQFitContextPolynomialFloat to create context . |
[in] | data | The input data with length of num_data . must-be-aligned |
[in] | num_coeff | The number of elements in coeff . It must be in range (0 < num_coeff <= order+1 ), where order is the polynomial or Chebyshev polynomial order specified in creation of context . |
[in] | coeff | Coefficients of model data. Its length must be num_coeff. Its first element must be of constant term, followed by those of first order, second order, and so on. must-be-aligned |
[out] | out | The output data. Its length must be num_data . must-be-aligned |
MT-safe
sakura_Status sakura_SubtractSinusoidFloat | ( | struct sakura_LSQFitContextFloat const * | context, |
size_t | num_data, | ||
float const | data[], | ||
size_t | num_nwave, | ||
size_t const | nwave[], | ||
size_t | num_coeff, | ||
double const | coeff[], | ||
float | out[] | ||
) |
Subtract best-fit sinusoidal model from input data.
[in] | context | A context created by sakura_CreateLSQFitContextSinusoidFloat . |
[in] | num_data | The number of elements in data and out. It must be equal to num_data which was given to sakura_CreateLSQFitContextSinusoidFloat to create context . |
[in] | data | The input data with length of num_data . must-be-aligned |
[in] | num_nwave | The number of elements in the array nwave . |
[in] | nwave | Wave numbers within the index range of data to be used for sinusoidal fitting. The values must be positive or zero (for constant term), but not exceed the maximum wave number specified in creation of context . The values must be stored in ascending order and must not be duplicate. |
[in] | num_coeff | The number of elements in coeff . It must be in range (0 < num_coeff <= num_model_bases ), where num_model_bases is ( num_nwave*2-1 ) or ( num_nwave*2 ) in cases nwave contains zero or not, respectively. Also it must not exceed num_data . |
[in] | coeff | Coefficients of model data. Its length must be num_coeff. If nwave contains zero, the first element must be of the constant term. If not, the first and the second elements are for sine and cosine terms for the smallest wave number in nwave, respectively. Coefficients of sine and cosine terms for the second-smallest wave number come next, and so on. must-be-aligned |
[out] | out | The output data. Its length must be num_data . must-be-aligned |
MT-safe
sakura_Status sakura_Uint32ToBool | ( | size_t | num_data, |
uint32_t const | data[], | ||
bool | result[] | ||
) |
Convert an input array to a boolean array.
Returns true if the corresponding element in input array != 0.
[in] | num_data | The number of elements in the arrays, data and result |
[in] | data | The input array of of size, num_data. must-be-aligned |
[out] | result | The output array of of size, num_data. must-be-aligned |
MT-safe
sakura_Status sakura_Uint8ToBool | ( | size_t | num_data, |
uint8_t const | data[], | ||
bool | result[] | ||
) |
Convert an input array to a boolean array.
Returns true if the corresponding element in input array != 0.
[in] | num_data | The number of elements in the arrays, data and result |
[in] | data | The input array of of size, num_data. must-be-aligned |
[out] | result | The output array of of size, num_data. must-be-aligned |
MT-safe
sakura_Status sakura_UnflipArrayDouble | ( | bool | inner_most_untouched, |
size_t | dims, | ||
size_t const | elements[], | ||
double const | src[], | ||
double | dst[] | ||
) |
Copy elements in the src array into the dst array with unflipping elements to the original order.
When you provide inner_most_untouched = false, elements = {3, 4} and src = {
}, then you will get dst as below.
[in] | inner_most_untouched | If true, the order of the inner most dimension is untouched. |
[in] | dims | Dimensions of the array src and dst. In other words, the number of elements in elements. |
[in] | elements | Numbers of elements of each dimension of src and dst with the inner-to-outer order. For example, when you have float matrix[4][3];
{ 3, 4 }
|
[in] | src | Source multidimensional array. must-be-aligned |
[in] | dst | Destination multidimensional array. must-be-aligned |
MT-safe
sakura_Status sakura_UnflipArrayDouble2 | ( | bool | inner_most_untouched, |
size_t | dims, | ||
size_t const | elements[], | ||
double const | src[][2], | ||
double | dst[][2] | ||
) |
Same as sakura_UnflipArrayFloat except the element type of the multidimensional arrays.
sakura_Status sakura_UnflipArrayFloat | ( | bool | inner_most_untouched, |
size_t | dims, | ||
size_t const | elements[], | ||
float const | src[], | ||
float | dst[] | ||
) |
Copy elements in the src array into the dst array with unflipping elements to the original order.
When you provide inner_most_untouched = false, elements = {3, 4} and src = {
}, then you will get dst as below.
[in] | inner_most_untouched | If true, the order of the inner most dimension is untouched. |
[in] | dims | Dimensions of the array src and dst. In other words, the number of elements in elements. |
[in] | elements | Numbers of elements of each dimension of src and dst with the inner-to-outer order. For example, when you have float matrix[4][3];
{ 3, 4 }
|
[in] | src | Source multidimensional array. must-be-aligned |
[in] | dst | Destination multidimensional array. must-be-aligned |
MT-safe
sakura_Status sakura_UpdateLSQCoefficientsDouble | ( | size_t const | num_data, |
float const | data[], | ||
bool const | mask[], | ||
size_t const | num_exclude_indices, | ||
size_t const | exclude_indices[], | ||
size_t const | num_model_bases, | ||
double const | basis_data[], | ||
size_t const | num_lsq_bases, | ||
size_t const | use_bases_idx[], | ||
double | lsq_matrix[], | ||
double | lsq_vector[] | ||
) |
Compute coefficients of simultaneous equations used for Least-Square fitting.
This function updates the coefficients of normal equation for LSQ fitting created by sakura_GetLSQCoefficientsDouble , by subtracting values corresponding to data points which have been used in the previous calculation but not this time. This is faster than newly calculating coefficients if the number of points to be excluded this time is less than half of those previously used.
[in] | num_data | The number of elements in the array data and the number of elements in each model data (i.e., discrete values of basis function) consisting the entire model. It must be a positive number. |
[in] | data | Input data with length of num_data . must-be-aligned |
[in] | mask | The input mask data with length of num_data . The i th element of data is included in input data if the i th element of mask is true, while it is excluded from input data if the corresponding element of mask is false. must-be-aligned |
[in] | num_exclude_indices | The number of data points to be excluded this time. The range of allowed value is between 0 and num_data . |
[in] | exclude_indices | An array containing indices of data points (the row index of basis_data ) to be excluded this time. The indices must be stored in the first num_exclude_indices elements. Its length should be num_exclude_indices . must-be-aligned |
[in] | num_model_bases | Number of model basis functions. It must be in range 0 < num_model_bases <= num_data . |
[in] | basis_data | A 1D array containing values of all basis functions concatenated. Loop for basis index must be inside of that for data index, i.e., the n -th data of the m -th model should be stored at basis_data [ num_data * (n-1) + (m-1) ]. Its length must be equal to ( num_model_bases * num_data ). must-be-aligned |
[in] | num_lsq_bases | The number of model basis functions to be used in actual fitting. It must be in range 0 < num_lsq_bases <= num_model_bases. |
[in] | use_bases_idx | A 1D array containing indices of basis model that are to be used for fitting. As for fitting types other than sinusoidal, it should be always [0, 1, 2, ..., (num_lsq_bases-1)]. Element values must be in ascending order. must-be-aligned |
[in,out] | lsq_matrix | A 1D array containing the values of a matrix at the left side of simultaneous equations for least-square fitting. Its length should therefore be equal to ( num_lsq_bases * num_lsq_bases ). Loop for columns comes inside that for rows, i.e., the value at the m -th row and n -th column is stored at lsq_matrix [ num_lsq_bases * ( m -1) + ( n -1)], though lsq_matrix is actually symmetric. |
[in,out] | lsq_vector | The values of a vector at the right side of simultaneous equations for least-square fitting. Its length should be equal to num_lsq_bases . must-be-aligned |
MT-safe