46 #ifndef LIBPOLYCOMP_H_GUARD
47 #define LIBPOLYCOMP_H_GUARD
55 #define PCOMP_STAT_SUCCESS 0
56 #define PCOMP_STAT_INVALID_ENCODING \
59 #define PCOMP_STAT_INVALID_BUFFER \
61 #define PCOMP_STAT_INVALID_FIT 3
68 void pcomp_version(
int* major,
int* minor);
90 const int8_t* input_buf,
size_t input_size);
92 const int16_t* input_buf,
95 const int32_t* input_buf,
98 const int64_t* input_buf,
103 const uint8_t* input_buf,
106 const uint16_t* input_buf,
109 const uint32_t* input_buf,
112 const uint64_t* input_buf,
116 int pcomp_decompress_rle_int8(int8_t* output_buf,
size_t output_size,
117 const int8_t* input_buf,
119 int pcomp_decompress_rle_int16(int16_t* output_buf,
size_t output_size,
120 const int16_t* input_buf,
122 int pcomp_decompress_rle_int32(int32_t* output_buf,
size_t output_size,
123 const int32_t* input_buf,
125 int pcomp_decompress_rle_int64(int64_t* output_buf,
size_t output_size,
126 const int64_t* input_buf,
130 int pcomp_decompress_rle_uint8(uint8_t* output_buf,
size_t output_size,
131 const uint8_t* input_buf,
133 int pcomp_decompress_rle_uint16(uint16_t* output_buf,
135 const uint16_t* input_buf,
137 int pcomp_decompress_rle_uint32(uint32_t* output_buf,
139 const uint32_t* input_buf,
141 int pcomp_decompress_rle_uint64(uint64_t* output_buf,
143 const uint64_t* input_buf,
166 const int8_t* input_buf,
170 const int16_t* input_buf,
174 const int32_t* input_buf,
178 const int64_t* input_buf,
184 const uint8_t* input_buf,
188 const uint16_t* input_buf,
192 const uint32_t* input_buf,
196 const uint64_t* input_buf,
202 const int8_t* input_buf,
206 const int16_t* input_buf,
210 const int32_t* input_buf,
214 const int64_t* input_buf,
220 const uint8_t* input_buf,
224 const uint16_t* input_buf,
228 const uint32_t* input_buf,
232 const uint64_t* input_buf,
253 size_t bits_per_sample);
262 double normalization,
double offset);
273 const float* input_buf,
277 const double* input_buf,
281 int pcomp_decompress_quant_float(
float* output_buf,
size_t output_size,
282 const void* input_buf,
285 int pcomp_decompress_quant_double(
double* output_buf,
287 const void* input_buf,
299 size_t num_of_coeffs);
307 const double* points);
338 const double* input);
346 typedef uint8_t pcomp_poly_size_t;
347 typedef uint16_t pcomp_chunk_size_t;
374 const double* samples);
376 pcomp_chunk_size_t num_of_samples,
377 pcomp_poly_size_t num_of_poly_coeffs,
const double* poly_coeffs,
378 pcomp_chunk_size_t num_of_cheby_coeffs,
379 const uint8_t* chebyshev_mask,
const double* cheby_coeffs);
404 size_t num_of_samples,
double period);
408 pcomp_poly_size_t num_of_coeffs,
409 double max_allowable_error,
415 pcomp_chunk_size_t num_of_samples,
451 double* coeffs,
double* cheby_residuals,
453 double* max_residual);
459 double max_allowable_error,
460 uint8_t* mask,
double* max_error);
463 size_t* num_of_chunks,
464 const double* input_buf,
size_t input_size,
469 size_t num_of_chunks);
476 size_t num_of_chunks);
480 size_t num_of_chunks);
483 size_t num_of_chunks);
487 size_t num_of_chunks);
490 size_t* num_of_chunks,
const void* buf);
int pcomp_decompress_diffrle_uint32(uint32_t *output_buf, size_t output_size, const uint32_t *input_buf, size_t input_size)
Compress an array of uint32_t values using the diffRLE compression.
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 transf...
int pcomp_compress_rle_int16(int16_t *output_buf, size_t *output_size, const int16_t *input_buf, size_t input_size)
Compress an array of int16_t values using the RLE compression.
pcomp_chebyshev_t * pcomp_init_chebyshev(size_t num_of_samples, pcomp_transform_direction_t dir)
Allocate a new instance of the pcomp_chebyshev_t structure on the heap.
int pcomp_run_poly_fit(pcomp_poly_fit_data_t *poly_fit, double *coeffs, const double *points)
Calculates a polynomial least-squares fit.
double pcomp_quant_normalization(const pcomp_quant_params_t *params)
Return the normalization constant used for converting floating-point numbers into quantized integers...
double pcomp_polycomp_period(const pcomp_polycomp_t *params)
Return the period of the input data, or a number less than or equal to 0 if the data have no periodic...
pcomp_transform_direction_t pcomp_chebyshev_direction(const pcomp_chebyshev_t *plan)
Return the direction of a Chebyshev transform.
int pcomp_decompress_polycomp_chunk(double *output, const pcomp_polycomp_chunk_t *chunk, pcomp_chebyshev_t *inv_chebyshev)
Decompress the data in a chunk.
pcomp_poly_size_t pcomp_polycomp_num_of_poly_coeffs(const pcomp_polycomp_t *params)
Return the number of coefficients for the fitting polynomial used in the polynomial compression...
int pcomp_encode_chunks(void *buf, size_t *buf_size, pcomp_polycomp_chunk_t *const chunk_array[], size_t num_of_chunks)
Encode a list of chunks into a sequence of raw bytes.
int pcomp_compress_diffrle_uint8(uint8_t *output_buf, size_t *output_size, const uint8_t *input_buf, size_t input_size)
Compress an array of uint8_t values using the diffRLE compression.
int pcomp_decompress_diffrle_int16(int16_t *output_buf, size_t output_size, const int16_t *input_buf, size_t input_size)
Compress an array of int16_t values using the diffRLE compression.
Compute a backward Chebyshev transform, with a normalization factor equal to one. ...
pcomp_polycomp_algorithm_t
Kind of algorithm used for the polynomial compression.
int pcomp_decode_chunks(pcomp_polycomp_chunk_t **chunk_array[], size_t *num_of_chunks, const void *buf)
Decode a byte sequence created by pcomp_encode_chunks into an array of chunks.
int pcomp_compress_rle_int32(int32_t *output_buf, size_t *output_size, const int32_t *input_buf, size_t input_size)
Compress an array of int32_t values using the RLE compression.
Compute a forward Chebyshev transform, with a normalization factor 1/(N + 1)
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 coeffic...
size_t pcomp_quant_element_size(const pcomp_quant_params_t *params)
Return the size (in bytes) of the elements to be quantized.
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 *chebyshev_mask, const double *cheby_coeffs)
Allocate memory for a pcomp_polycomp_chunk_t object and fill it with data compressed using the polyno...
int pcomp_compress_diffrle_int32(int32_t *output_buf, size_t *output_size, const int32_t *input_buf, size_t input_size)
Compress an array of int32_t values using the diffRLE compression.
double pcomp_polycomp_max_error(const pcomp_polycomp_t *params)
Return the upper bound on the error of the polynomial compression.
void pcomp_free_polycomp(pcomp_polycomp_t *params)
Free the memory allocated by pcomp_init_polycomp for a pcomp_polycomp_t structure.
pcomp_polycomp_chunk_t * pcomp_init_chunk(pcomp_chunk_size_t num_of_samples)
Allocate memory for a pcomp_polycomp_chunk_t object.
size_t pcomp_rle_bufsize(size_t input_size)
Calculate an upper limit for the size of a buffer holding RLE-encoded streams.
int pcomp_decompress_diffrle_uint8(uint8_t *output_buf, size_t output_size, const uint8_t *input_buf, size_t input_size)
Compress an array of uint8_t values using the diffRLE compression.
int pcomp_compress_rle_uint64(uint64_t *output_buf, size_t *output_size, const uint64_t *input_buf, size_t input_size)
Compress an array of uint64_t values using the RLE compression.
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_chunk_size_t pcomp_chunk_num_of_samples(const pcomp_polycomp_chunk_t *chunk)
Return the number of samples in a chunk.
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...
pcomp_chunk_size_t pcomp_polycomp_samples_per_chunk(const pcomp_polycomp_t *params)
Return the number of samples per chunk.
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 comp...
int pcomp_compress_diffrle_uint16(uint16_t *output_buf, size_t *output_size, const uint16_t *input_buf, size_t input_size)
Compress an array of uint16_t values using the diffRLE compression.
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...
When needed, apply the Chebyshev transform to the residuals of the polynomial fit.
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.
size_t pcomp_diffrle_bufsize(size_t input_size)
Calculate an upper limit for the size of a buffer holding streams encoded using differenced RLE...
void pcomp_free_chebyshev(pcomp_chebyshev_t *plan)
Free the memory allocated by a previous call to pcomp_init_chebyshev.
int pcomp_decompress_diffrle_int32(int32_t *output_buf, size_t output_size, const int32_t *input_buf, size_t input_size)
Compress an array of int32_t values using the diffRLE compression.
size_t pcomp_chebyshev_num_of_samples(const pcomp_chebyshev_t *plan)
Return the number of samples in a Chebyshev transform.
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.
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.
const double * pcomp_chebyshev_output(const pcomp_chebyshev_t *plan)
Return the output (Chebyshev transform) of the last call to pcomp_run_chebyshev.
void pcomp_straighten(double *output, const double *input, size_t num_of_samples, double period)
Remove sudden jumps from input.
size_t pcomp_total_num_of_samples(pcomp_polycomp_chunk_t *const chunk_array[], size_t num_of_chunks)
Compute the sum of the number of samples encoded in chunk_array.
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...
size_t pcomp_quant_bits_per_sample(const pcomp_quant_params_t *params)
Return the number of bits that must be used for each quantized sample.
int pcomp_decompress_diffrle_uint16(uint16_t *output_buf, size_t output_size, const uint16_t *input_buf, size_t input_size)
Compress an array of uint16_t values using the diffRLE compression.
int pcomp_compress_quant_float(void *output_buf, size_t *output_size, const float *input_buf, size_t input_size, pcomp_quant_params_t *params)
Quantize a stream of 32-bit floating point numbers.
int pcomp_compress_quant_double(void *output_buf, size_t *output_size, const double *input_buf, size_t input_size, pcomp_quant_params_t *params)
Quantize a stream of 64-bit floating point numbers.
int pcomp_compress_rle_uint16(uint16_t *output_buf, size_t *output_size, const uint16_t *input_buf, size_t input_size)
Compress an array of uint16_t values using the RLE compression.
size_t pcomp_chunk_num_of_bytes(const pcomp_polycomp_chunk_t *chunk)
Return the number of bytes necessary to encode a chunk.
size_t pcomp_poly_fit_num_of_coeffs(const pcomp_poly_fit_data_t *poly_fit)
Return the number of coefficients of the least-squares fitting polynomial.
int pcomp_compress_diffrle_uint32(uint32_t *output_buf, size_t *output_size, const uint32_t *input_buf, size_t input_size)
Compress an array of uint32_t values using the diffRLE compression.
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 .
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.
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.
pcomp_poly_fit_data_t * pcomp_init_poly_fit(size_t num_of_samples, size_t num_of_coeffs)
Allocate a new instance of the pcomp_poly_fit_data_t structure on the heap.
pcomp_quant_params_t * pcomp_init_quant_params(size_t element_size, size_t bits_per_sample)
Initialize a pcomp_quant_params_t structure.
int pcomp_decompress_diffrle_int8(int8_t *output_buf, size_t output_size, const int8_t *input_buf, size_t input_size)
Decompress an array of int8_t values encoded using the diffRLE compression.
void pcomp_polycomp_set_period(pcomp_polycomp_t *params, double period)
Set the periodicity of the data to be compressed.
double pcomp_quant_offset(const pcomp_quant_params_t *params)
Return the additive constant used for converting floating-point numbers into quantized integers...
void pcomp_free_poly_fit(pcomp_poly_fit_data_t *poly_fit)
Free an instance of the pcomp_poly_fit_data_t that has been allocated via a call to pcomp_init_poly_f...
int pcomp_compress_diffrle_int16(int16_t *output_buf, size_t *output_size, const int16_t *input_buf, size_t input_size)
Compress an array of int16_t values using the diffRLE compression.
void pcomp_free_chunk(pcomp_polycomp_chunk_t *chunk)
Free memory associated with a pcomp_poly_chunk_t.
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...
int pcomp_compress_rle_int64(int64_t *output_buf, size_t *output_size, const int64_t *input_buf, size_t input_size)
Compress an array of int64_t values using the RLE compression.
pcomp_polycomp_t * pcomp_init_polycomp(pcomp_chunk_size_t samples_per_chunk, pcomp_poly_size_t num_of_coeffs, double max_allowable_error, pcomp_polycomp_algorithm_t algorithm)
Allocate space for a pcomp_polycomp_t structure.
int pcomp_compress_rle_int8(int8_t *output_buf, size_t *output_size, const int8_t *input_buf, size_t input_size)
Compress an array of int8_t values using the RLE compression.
int pcomp_decompress_polycomp(double *output_buf, pcomp_polycomp_chunk_t *const chunk_array[], size_t num_of_chunks)
Decompress a sequence of chunks.
int pcomp_decompress_diffrle_uint64(uint64_t *output_buf, size_t output_size, const uint64_t *input_buf, size_t input_size)
Compress an array of uint64_t values using the diffRLE compression.
int pcomp_compress_rle_uint32(uint32_t *output_buf, size_t *output_size, const uint32_t *input_buf, size_t input_size)
Compress an array of uint32_t values using the RLE compression.
size_t pcomp_quant_bufsize(size_t input_size, const pcomp_quant_params_t *params)
Return the size (in bytes) of the buffer that will contain a quantized stream of input_size floating ...
int pcomp_decompress_diffrle_int64(int64_t *output_buf, size_t output_size, const int64_t *input_buf, size_t input_size)
Compress an array of int64_t values using the diffRLE compression.
void pcomp_quant_set_normalization(pcomp_quant_params_t *params, double normalization, double offset)
Set the normalization constants (multiplicative and additive) used to quantize floating-point numbers...
void pcomp_free_quant_params(pcomp_quant_params_t *params)
Free a pcomp_quant_params_t structure.
int pcomp_compress_diffrle_int8(int8_t *output_buf, size_t *output_size, const int8_t *input_buf, size_t input_size)
Compress an array of int8_t values using the diffRLE compression.
int pcomp_compress_rle_uint8(uint8_t *output_buf, size_t *output_size, const uint8_t *input_buf, size_t input_size)
Compress an array of uint8_t values using the RLE compression.
pcomp_transform_direction_t
Direction of a Chebyshev transform.
int pcomp_compress_diffrle_int64(int64_t *output_buf, size_t *output_size, const int64_t *input_buf, size_t input_size)
Compress an array of int64_t values using the diffRLE compression.
pcomp_polycomp_algorithm_t pcomp_polycomp_algorithm(const pcomp_polycomp_t *params)
Return the kind of algorithm used for a polynomial compression.
int pcomp_compress_polycomp(pcomp_polycomp_chunk_t **chunk_array[], size_t *num_of_chunks, const double *input_buf, size_t input_size, const pcomp_polycomp_t *params)
Compress the array input_buf using polynomial compression.
size_t pcomp_chunks_num_of_bytes(pcomp_polycomp_chunk_t *const chunks[], size_t num_of_chunks)
Number of bytes required by pcomp_encode_chunks.
int pcomp_chunk_is_compressed(const pcomp_polycomp_chunk_t *chunk)
Return nonzero if the chunk holds data in uncompressed form.
If the absolute value of the residuals of a polynomial fit are too large, store the data in uncompres...
size_t pcomp_poly_fit_num_of_samples(const pcomp_poly_fit_data_t *poly_fit)
Return the number of samples to be used in a polynomial fit.
void pcomp_free_chunks(pcomp_polycomp_chunk_t *chunk_array[], size_t num_of_chunks)
Free an array of chunks.
int pcomp_run_chebyshev(pcomp_chebyshev_t *plan, pcomp_transform_direction_t dir, double *output, const double *input)
Compute a forward/backward Chebyshev discrete transform.
const double * pcomp_chebyshev_input(const pcomp_chebyshev_t *plan)
Return the input data used in the last call to pcomp_run_chebyshev.
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 t...
int pcomp_compress_diffrle_uint64(uint64_t *output_buf, size_t *output_size, const uint64_t *input_buf, size_t input_size)
Compress an array of uint64_t values using the diffRLE compression.