Libpolycomp
1.0
A compression/decompression library that implements the polynomial compression and other simple compression schemes
|
Functions | |
pcomp_chebyshev_t * | pcomp_polycomp_forward_cheby (const pcomp_polycomp_t *params) |
Return a pointer to a pcomp_chebyshev_t structure representing the forward Chebyshev transform. More... | |
pcomp_chebyshev_t * | pcomp_polycomp_backward_cheby (const pcomp_polycomp_t *params) |
Return a pointer to a pcomp_chebyshev_t structure representing the forward Chebyshev transform. More... | |
void | pcomp_straighten (double *output, const double *input, size_t num_of_samples, double period) |
Remove sudden jumps from input. More... | |
pcomp_polycomp_chunk_t * | pcomp_init_chunk (pcomp_chunk_size_t num_of_samples) |
Allocate memory for a pcomp_polycomp_chunk_t object. More... | |
pcomp_polycomp_chunk_t * | pcomp_init_uncompressed_chunk (pcomp_chunk_size_t num_of_samples, const double *samples) |
Allocate memory for a pcomp_polycomp_chunk_t object and fill it with data in uncompressed form. More... | |
pcomp_polycomp_chunk_t * | pcomp_init_compressed_chunk (pcomp_chunk_size_t num_of_samples, pcomp_poly_size_t num_of_poly_coeffs, const double *poly_coeffs, pcomp_chunk_size_t num_of_cheby_coeffs, const uint8_t *cheby_mask, const double *cheby_coeffs) |
Allocate memory for a pcomp_polycomp_chunk_t object and fill it with data compressed using the polynomial compression algorithm. More... | |
void | pcomp_free_chunk (pcomp_polycomp_chunk_t *chunk) |
Free memory associated with a pcomp_poly_chunk_t. More... | |
pcomp_chunk_size_t | pcomp_chunk_num_of_samples (const pcomp_polycomp_chunk_t *chunk) |
Return the number of samples in a chunk. More... | |
size_t | pcomp_chunk_num_of_bytes (const pcomp_polycomp_chunk_t *chunk) |
Return the number of bytes necessary to encode a chunk. More... | |
int | pcomp_chunk_is_compressed (const pcomp_polycomp_chunk_t *chunk) |
Return nonzero if the chunk holds data in uncompressed form. | |
const double * | pcomp_chunk_uncompressed_data (const pcomp_polycomp_chunk_t *chunk) |
If the chunks contain uncompressed data, returns a pointer to the first element. Otherwise, return NULL . | |
pcomp_poly_size_t | pcomp_chunk_num_of_poly_coeffs (const pcomp_polycomp_chunk_t *chunk) |
If the chunks contain compressed data, returns the number of polynomial coefficients used in the compression. Otherwise, return zero. | |
const double * | pcomp_chunk_poly_coeffs (const pcomp_polycomp_chunk_t *chunk) |
If the chunks contain compressed data, returns a pointer to the first element of the array of coefficients of the interpolating polynomial. Otherwise, return NULL . | |
pcomp_chunk_size_t | pcomp_chunk_num_of_cheby_coeffs (const pcomp_polycomp_chunk_t *chunk) |
If the chunks contain compressed data, returns the number of nonzero Chebyshev coefficients held in the chunk. Otherwise, return zero. | |
const double * | pcomp_chunk_cheby_coeffs (const pcomp_polycomp_chunk_t *chunk) |
If the chunks contain compressed data, returns a pointer to the first element of the Chebyshev transform of the fit residuals. Otherwise, return NULL . | |
size_t | pcomp_chunk_cheby_mask_size (pcomp_chunk_size_t chunk_size) |
Return the number of bytes required for the bitmask of nonzero Chebyshev coefficients. More... | |
const uint8_t * | pcomp_chunk_cheby_mask (const pcomp_polycomp_chunk_t *chunk) |
Return a pointer to the bitmask of nonzero Chebyshev coefficients for a chunk. More... | |
int | pcomp_polyfit_and_chebyshev (pcomp_polycomp_t *params, double *coeffs, double *cheby_residuals, const double *input, double *max_residual) |
Compute a polynomial fit of the data in input and a Chebyshev transform of the residuals. More... | |
int | pcomp_mask_get_bit (const uint8_t *mask, size_t pos) |
Return the value of the bit at the position pos in the bitmask mask. More... | |
void | pcomp_mask_set_bit (uint8_t *mask, size_t pos, int value) |
Set the value of the bit at position pos in the bitmask . More... | |
size_t | pcomp_find_chebyshev_mask (pcomp_chebyshev_t *chebyshev, pcomp_chebyshev_t *inv_chebyshev, double max_allowable_error, uint8_t *mask, double *max_error) |
Find the smallest subset of Chebyshev coefficients that can approximate a Chebyshev transform with an error less than max_allowable_error. More... | |
int | pcomp_run_polycomp_on_chunk (pcomp_polycomp_t *params, const double *input, pcomp_chunk_size_t num_of_samples, pcomp_polycomp_chunk_t *chunk, double *max_error) |
Compress the first num_of_samples elements in input and store them in chunk. More... | |
int | pcomp_decompress_polycomp_chunk (double *output, const pcomp_polycomp_chunk_t *chunk, pcomp_chebyshev_t *inv_chebyshev) |
Decompress the data in a chunk. More... | |
const uint8_t* pcomp_chunk_cheby_mask | ( | const pcomp_polycomp_chunk_t * | chunk | ) |
Return a pointer to the bitmask of nonzero Chebyshev coefficients for a chunk.
[in] | chunk | Pointer to the chunk |
size_t pcomp_chunk_cheby_mask_size | ( | pcomp_chunk_size_t | chunk_size | ) |
Return the number of bytes required for the bitmask of nonzero Chebyshev coefficients.
The polynomial compression compresses Chebyshev transforms by saving only those coefficients that are significantly different from zero. In order to keep track of the position of such coefficients in the full array, a bit mask is used. This function determines how many bytes are required for such mask, which is internally represented by Libpolycomp as an array of uint8_t
values.
[in] | chunk_size | Number of samples in the chunk |
uint8_t
values) required for the mask. size_t pcomp_chunk_num_of_bytes | ( | const pcomp_polycomp_chunk_t * | chunk | ) |
Return the number of bytes necessary to encode a chunk.
Refer to pcomp_encode_chunks and pcomp_decode_chunks for further details.
[in] | chunk | Pointer to the chunk data |
pcomp_chunk_size_t pcomp_chunk_num_of_samples | ( | const pcomp_polycomp_chunk_t * | chunk | ) |
int pcomp_decompress_polycomp_chunk | ( | double * | output, |
const pcomp_polycomp_chunk_t * | chunk, | ||
pcomp_chebyshev_t * | inv_chebyshev | ||
) |
Decompress the data in a chunk.
This function performs the decompression of a chunk, and it is the counterpart of pcomp_run_polycomp_on_chunk. Here is an example:
[out] | output | Pointer to the array that will contain the uncompressed data |
[in] | chunk | The chunk to decompress |
[in] | inv_chebyshev | Pointer to a pcomp_chebyshev_t object that performs the inverse Chebyshev transform. The function does not allocate an object of this kind because in this way such objects can be reused on subsequent calls to pcomp_decompress_polycomp_chunk. |
size_t pcomp_find_chebyshev_mask | ( | pcomp_chebyshev_t * | chebyshev, |
pcomp_chebyshev_t * | inv_chebyshev, | ||
double | max_allowable_error, | ||
uint8_t * | mask, | ||
double * | max_error | ||
) |
Find the smallest subset of Chebyshev coefficients that can approximate a Chebyshev transform with an error less than max_allowable_error.
On exit, the bits in bitmask will be set to 1 in correspondence of every Chebyshev coefficient that must be retained. The function returns the number of Chebyshev coefficients to retain (i.e., the number of bits in mask that have been set to 1).
[in] | chebyshev | Pointer to a pcomp_chebyshev_t structure used to compute the forward Chebyshev transform (from the space of the fit residuals to the Chebyshev space) |
[in] | inv_chebyshev | Pointer to a pcomp_chebyshev_t structure representing the inverse transform of chebyshev. |
[in] | max_allowable_error | The maximum allowed discrepancy for the chopped Chebyshev, as measured in the space of the fit residuals. |
[out] | mask | Bitmask that will contain the position of the unchopped Chebyshev terms of the transform of the fit residuals. Use pcomp_mask_get_bit to access each element. |
[out] | max_error | If not NULL , it will be set to the maximum error due to the chopping of the Chebyshev transform represented by mask. |
void pcomp_free_chunk | ( | pcomp_polycomp_chunk_t * | chunk | ) |
Free memory associated with a pcomp_poly_chunk_t.
This function releases the memory allocated by one of the following functions:
[in] | chunk | Pointer to the object to be freed. |
pcomp_polycomp_chunk_t* pcomp_init_chunk | ( | pcomp_chunk_size_t | num_of_samples | ) |
Allocate memory for a pcomp_polycomp_chunk_t object.
[in] | num_of_samples | Number of samples that the chunk will be capable to hold. |
pcomp_polycomp_chunk_t* pcomp_init_compressed_chunk | ( | pcomp_chunk_size_t | num_of_samples, |
pcomp_poly_size_t | num_of_poly_coeffs, | ||
const double * | poly_coeffs, | ||
pcomp_chunk_size_t | num_of_cheby_coeffs, | ||
const uint8_t * | cheby_mask, | ||
const double * | cheby_coeffs | ||
) |
Allocate memory for a pcomp_polycomp_chunk_t object and fill it with data compressed using the polynomial compression algorithm.
[in] | num_of_samples | Number of samples that the chunk will be capable to hold. |
[in] | num_of_poly_coeffs | Number of coefficients of the interpolating polynomial. |
[in] | poly_coeffs | Pointer to the coefficients of the interpolating polynomial. Their number must be equal to the parameter num_of_poly_coeffs. |
[in] | num_of_cheby_coeffs | Number of nonzero Chebyshev coefficients associated with the polynomial fit. This number is always less than num_of_samples. Zero is allowed. |
[in] | cheby_mask | Bitmask representing the position of the nonzero coefficients in cheby_coeffs within the full sequence. (Use pcomp_mask_get_bit and pcomp_mask_set_bit to read/write bits in the sequence.) |
[in] | cheby_coeffs | Array of nonzero Chebyshev coefficients. Their number must be equal to num_of_cheby_coeffs. |
pcomp_polycomp_chunk_t* pcomp_init_uncompressed_chunk | ( | pcomp_chunk_size_t | num_of_samples, |
const double * | samples | ||
) |
Allocate memory for a pcomp_polycomp_chunk_t object and fill it with data in uncompressed form.
[in] | num_of_samples | Number of samples that the chunk will be capable to hold. |
[in] | samples | The (uncompressed) samples to copy into the chunk. After the call, input is no longer needed and can be freed without invalidating the pointer returned by the function. |
int pcomp_mask_get_bit | ( | const uint8_t * | mask, |
size_t | pos | ||
) |
void pcomp_mask_set_bit | ( | uint8_t * | mask, |
size_t | pos, | ||
int | value | ||
) |
pcomp_chebyshev_t* pcomp_polycomp_backward_cheby | ( | const pcomp_polycomp_t * | params | ) |
pcomp_chebyshev_t* pcomp_polycomp_forward_cheby | ( | const pcomp_polycomp_t * | params | ) |
int pcomp_polyfit_and_chebyshev | ( | pcomp_polycomp_t * | params, |
double * | coeffs, | ||
double * | cheby_residuals, | ||
const double * | input, | ||
double * | max_residual | ||
) |
Compute a polynomial fit of the data in input and a Chebyshev transform of the residuals.
Note that this function always computes the Chebyshev transform of the data, even if there is a perfect fit between the polynomial and the input data.
[in] | params | Pointer to a pcomp_polycomp_t structure initialized by pcomp_init_polycomp. |
[out] | coeffs | Pointer to the array that on exit will hold the coefficients of the best-fit polynomial. It must have enough room for a number of elements equal to the return value of pcomp_polycomp_num_of_poly_coeffs. |
[out] | cheby_residuals | Pointer to an array that on exit will contain the Chebyshev transform of the residuals of the fit. It can be NULL ; in any case, these numbers can be obtained by the use of a call to pcomp_polycomp_forward_cheby and pcomp_chebyshev_output. |
[in] | input | Pointer to the array of values to be transformed. The number of values used is equal to the return value of the function pcomp_polycomp_samples_per_chunk. |
[out] | max_residual | Pointer to a variable that will hold the maximum absolute value of the discrepancy between each sample in input and the polynomial fit. It can be NULL . |
int pcomp_run_polycomp_on_chunk | ( | pcomp_polycomp_t * | params, |
const double * | input, | ||
pcomp_chunk_size_t | num_of_samples, | ||
pcomp_polycomp_chunk_t * | chunk, | ||
double * | max_error | ||
) |
Compress the first num_of_samples elements in input and store them in chunk.
The function determines if the data in input can be efficiently compressed using polynomial compression with the parameters described by params. If it is so, it stores the compressed data in chunk. If the compression ratio is small or equal to one, or if the compression error is too large, the function copies the data in input into chunk in uncompressed format.
The following example shows how to use this function together with pcomp_init_chunk:
[in] | params | Pointer to a pcomp_polycomp_t structure (created using pcomp_init_polycomp) which provides the parameters of the compression. |
[in] | input | The sequence of numbers to compress. Their number is equal to the parameter num_of_samples |
[in] | num_of_samples | Number of values in input to compress |
[in,out] | chunk | The chunk that will contain the data, either in compressed or uncompressed format. It must have already been initialized via a call to pcomp_init_chunk. |
[out] | max_error | On exit, the function writes the compression error here. It can be NULL . |
void pcomp_straighten | ( | double * | output, |
const double * | input, | ||
size_t | num_of_samples, | ||
double | period | ||
) |
Remove sudden jumps from input.
Assuming that the data in the array input have a periodicity equal to period, the function copies them to output while applying a positive/negative offset equal to a multiple of period.
It is ok for input and output to point to the same memory location.
[out] | output | Pointer to the array that will contain the result. It must have room for at least num_of_samples values. |
[in] | input | Array of num_of_samples values to process. |
[in] | num_of_samples | Number of samples to process in input |
[in] | period | Periodicity of the data. If less or equal to zero, input is copied verbatim to output. |