Libpolycomp  1.0
A compression/decompression library that implements the polynomial compression and other simple compression schemes
libpolycomp.h
Go to the documentation of this file.
1 /* libpolycomp.h - interface to "polycomp", a compression library
2  * aimed to smooth, noise-free data
3  *
4  * Copyright (c) 2015 Maurizio Tomasi
5  *
6  * Permission is hereby granted, free of charge, to any person
7  * obtaining a copy of this software and associated documentation
8  * files (the "Software"), to deal in the Software without
9  * restriction, including without limitation the rights to use, copy,
10  * modify, merge, publish, distribute, sublicense, and/or sell copies
11  * of the Software, and to permit persons to whom the Software is
12  * furnished to do so, subject to the following conditions:
13  *
14  * The above copyright notice and this permission notice shall be
15  * included in all copies or substantial portions of the Software.
16  *
17  * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND,
18  * EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF
19  * MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND
20  * NONINFRINGEMENT. IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS
21  * BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN
22  * ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN
23  * CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE
24  * SOFTWARE.
25  */
26 
46 #ifndef LIBPOLYCOMP_H_GUARD
47 #define LIBPOLYCOMP_H_GUARD
48 
49 #include <stddef.h>
50 #include <stdint.h>
51 
52 /***********************************************************************
53  * Error codes
54  */
55 #define PCOMP_STAT_SUCCESS 0
56 #define PCOMP_STAT_INVALID_ENCODING \
57  1
59 #define PCOMP_STAT_INVALID_BUFFER \
60  2
61 #define PCOMP_STAT_INVALID_FIT 3
63 /***********************************************************************
64  * Version information
65  */
66 
67 /* Return the major/minor version of the library */
68 void pcomp_version(int* major, int* minor);
69 
70 /***********************************************************************
71  * RLE compression functions
72  *
73  * All of the following function perform a direct/inverse run-length
74  * encoding of the input data. The buffer "output_buf" must be
75  * preallocated with a number of elements which is at least equal to
76  * the value returned by the function "pcomp_rle_bufsize".
77  *
78  * The inverse encoding (decompression) routines require the user to
79  * know how large "output_buf" it will be. This can be achieved by
80  * saving this information somewhere during compression.
81  */
82 
83 /* Return the minimum number of elements necessary for the buffer that
84  * will contain the compressed data, if the input data has
85  * "input_size" elements. */
86 size_t pcomp_rle_bufsize(size_t input_size);
87 
88 /* Compression of integer sequences */
89 int pcomp_compress_rle_int8(int8_t* output_buf, size_t* output_size,
90  const int8_t* input_buf, size_t input_size);
91 int pcomp_compress_rle_int16(int16_t* output_buf, size_t* output_size,
92  const int16_t* input_buf,
93  size_t input_size);
94 int pcomp_compress_rle_int32(int32_t* output_buf, size_t* output_size,
95  const int32_t* input_buf,
96  size_t input_size);
97 int pcomp_compress_rle_int64(int64_t* output_buf, size_t* output_size,
98  const int64_t* input_buf,
99  size_t input_size);
100 
101 /* Compression of unsigned integer sequences */
102 int pcomp_compress_rle_uint8(uint8_t* output_buf, size_t* output_size,
103  const uint8_t* input_buf,
104  size_t input_size);
105 int pcomp_compress_rle_uint16(uint16_t* output_buf, size_t* output_size,
106  const uint16_t* input_buf,
107  size_t input_size);
108 int pcomp_compress_rle_uint32(uint32_t* output_buf, size_t* output_size,
109  const uint32_t* input_buf,
110  size_t input_size);
111 int pcomp_compress_rle_uint64(uint64_t* output_buf, size_t* output_size,
112  const uint64_t* input_buf,
113  size_t input_size);
114 
115 /* Decompression of integer sequences */
116 int pcomp_decompress_rle_int8(int8_t* output_buf, size_t output_size,
117  const int8_t* input_buf,
118  size_t input_size);
119 int pcomp_decompress_rle_int16(int16_t* output_buf, size_t output_size,
120  const int16_t* input_buf,
121  size_t input_size);
122 int pcomp_decompress_rle_int32(int32_t* output_buf, size_t output_size,
123  const int32_t* input_buf,
124  size_t input_size);
125 int pcomp_decompress_rle_int64(int64_t* output_buf, size_t output_size,
126  const int64_t* input_buf,
127  size_t input_size);
128 
129 /* Decompression of unsigned integer sequences */
130 int pcomp_decompress_rle_uint8(uint8_t* output_buf, size_t output_size,
131  const uint8_t* input_buf,
132  size_t input_size);
133 int pcomp_decompress_rle_uint16(uint16_t* output_buf,
134  size_t output_size,
135  const uint16_t* input_buf,
136  size_t input_size);
137 int pcomp_decompress_rle_uint32(uint32_t* output_buf,
138  size_t output_size,
139  const uint32_t* input_buf,
140  size_t input_size);
141 int pcomp_decompress_rle_uint64(uint64_t* output_buf,
142  size_t output_size,
143  const uint64_t* input_buf,
144  size_t input_size);
145 
146 /***********************************************************************
147  * Differenced RLE compression functions
148  *
149  * All of the following function perform a direct/inverse run-length
150  * encoding of consecutive differences in the samples of the input
151  * data. The buffer "output_buf" must be preallocated with a number of
152  * elements which is at least equal to the value returned by the
153  * function "pcomp_diffrle_bufsize".
154  *
155  * The inverse encoding (decompression) routines require the user to
156  * know how large "output_buf" it will be. This can be achieved by
157  * saving this information somewhere during compression.
158  */
159 
160 /* Return the minimum number of elements necessary for the buffer that
161  * will contain the compressed data. */
162 size_t pcomp_diffrle_bufsize(size_t input_size);
163 
164 /* Compression of integer sequences */
165 int pcomp_compress_diffrle_int8(int8_t* output_buf, size_t* output_size,
166  const int8_t* input_buf,
167  size_t input_size);
168 int pcomp_compress_diffrle_int16(int16_t* output_buf,
169  size_t* output_size,
170  const int16_t* input_buf,
171  size_t input_size);
172 int pcomp_compress_diffrle_int32(int32_t* output_buf,
173  size_t* output_size,
174  const int32_t* input_buf,
175  size_t input_size);
176 int pcomp_compress_diffrle_int64(int64_t* output_buf,
177  size_t* output_size,
178  const int64_t* input_buf,
179  size_t input_size);
180 
181 /* Compression of unsigned integer sequences */
182 int pcomp_compress_diffrle_uint8(uint8_t* output_buf,
183  size_t* output_size,
184  const uint8_t* input_buf,
185  size_t input_size);
186 int pcomp_compress_diffrle_uint16(uint16_t* output_buf,
187  size_t* output_size,
188  const uint16_t* input_buf,
189  size_t input_size);
190 int pcomp_compress_diffrle_uint32(uint32_t* output_buf,
191  size_t* output_size,
192  const uint32_t* input_buf,
193  size_t input_size);
194 int pcomp_compress_diffrle_uint64(uint64_t* output_buf,
195  size_t* output_size,
196  const uint64_t* input_buf,
197  size_t input_size);
198 
199 /* Decompression of integer sequences */
200 int pcomp_decompress_diffrle_int8(int8_t* output_buf,
201  size_t output_size,
202  const int8_t* input_buf,
203  size_t input_size);
204 int pcomp_decompress_diffrle_int16(int16_t* output_buf,
205  size_t output_size,
206  const int16_t* input_buf,
207  size_t input_size);
208 int pcomp_decompress_diffrle_int32(int32_t* output_buf,
209  size_t output_size,
210  const int32_t* input_buf,
211  size_t input_size);
212 int pcomp_decompress_diffrle_int64(int64_t* output_buf,
213  size_t output_size,
214  const int64_t* input_buf,
215  size_t input_size);
216 
217 /* Decompression of unsigned integer sequences */
218 int pcomp_decompress_diffrle_uint8(uint8_t* output_buf,
219  size_t output_size,
220  const uint8_t* input_buf,
221  size_t input_size);
222 int pcomp_decompress_diffrle_uint16(uint16_t* output_buf,
223  size_t output_size,
224  const uint16_t* input_buf,
225  size_t input_size);
226 int pcomp_decompress_diffrle_uint32(uint32_t* output_buf,
227  size_t output_size,
228  const uint32_t* input_buf,
229  size_t input_size);
230 int pcomp_decompress_diffrle_uint64(uint64_t* output_buf,
231  size_t output_size,
232  const uint64_t* input_buf,
233  size_t input_size);
234 
235 /***********************************************************************
236  * Quantization
237  *
238  * The following functions perform a quantization/dequantization of
239  * the input data. The buffer "output_buf" must have been preallocated
240  * with a number of bytes which is at least equal to the value of the
241  * function.
242  *
243  * The inverse encoding (decompression) routines require the user to
244  * know how large "output_buf" it will be. This can be achieved by
245  * saving this information somewhere during compression.
246  */
247 
248 /* Parameters used in the quantization. */
251 
252 pcomp_quant_params_t* pcomp_init_quant_params(size_t element_size,
253  size_t bits_per_sample);
255 
256 size_t pcomp_quant_element_size(const pcomp_quant_params_t* params);
258 double pcomp_quant_normalization(const pcomp_quant_params_t* params);
259 double pcomp_quant_offset(const pcomp_quant_params_t* params);
260 
262  double normalization, double offset);
263 
264 /* Return the minimum number of bytes necessary for the buffer that
265  * will contain the compressed data, if the input data has
266  * "input_size" elements, each requiring "element_size" bytes, and if
267  * quantization will use a number of bits per sample equal to
268  * "bits_per_sample". Return zero if the input is invalid. */
269 size_t pcomp_quant_bufsize(size_t input_size,
270  const pcomp_quant_params_t* params);
271 
272 int pcomp_compress_quant_float(void* output_buf, size_t* output_size,
273  const float* input_buf,
274  size_t input_size,
275  pcomp_quant_params_t* params);
276 int pcomp_compress_quant_double(void* output_buf, size_t* output_size,
277  const double* input_buf,
278  size_t input_size,
279  pcomp_quant_params_t* params);
280 
281 int pcomp_decompress_quant_float(float* output_buf, size_t output_size,
282  const void* input_buf,
283  size_t input_size,
284  const pcomp_quant_params_t* params);
285 int pcomp_decompress_quant_double(double* output_buf,
286  size_t output_size,
287  const void* input_buf,
288  size_t input_size,
289  const pcomp_quant_params_t* params);
290 
291 /***********************************************************************
292  * Polynomial fitting routines
293  */
294 
297 
298 pcomp_poly_fit_data_t* pcomp_init_poly_fit(size_t num_of_samples,
299  size_t num_of_coeffs);
301 size_t
303 size_t
305 
306 int pcomp_run_poly_fit(pcomp_poly_fit_data_t* poly_fit, double* coeffs,
307  const double* points);
308 
309 /***********************************************************************
310  * Chebyshev transform routines
311  */
312 
313 struct __pcomp_chebyshev_t;
315 
320 typedef enum {
328 
330 pcomp_init_chebyshev(size_t num_of_samples,
337  pcomp_transform_direction_t dir, double* output,
338  const double* input);
339 const double* pcomp_chebyshev_input(const pcomp_chebyshev_t* plan);
340 const double* pcomp_chebyshev_output(const pcomp_chebyshev_t* plan);
341 
342 /***********************************************************************
343  * Polynomial compression (low-level functions)
344  */
345 
346 typedef uint8_t pcomp_poly_size_t;
347 typedef uint16_t pcomp_chunk_size_t;
348 
355 typedef enum {
363 
364 struct __pcomp_polycomp_t;
365 typedef struct __pcomp_polycomp_t pcomp_polycomp_t;
366 
369 
371 pcomp_init_chunk(pcomp_chunk_size_t num_of_samples);
373 pcomp_init_uncompressed_chunk(pcomp_chunk_size_t num_of_samples,
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);
380 
382 
383 pcomp_chunk_size_t
387 const double*
389 
390 pcomp_poly_size_t
392 const double*
394 
395 pcomp_chunk_size_t
397 const double*
399 size_t pcomp_chunk_cheby_mask_size(pcomp_chunk_size_t chunk_size);
400 const uint8_t*
402 
403 void pcomp_straighten(double* output, const double* input,
404  size_t num_of_samples, double period);
405 
407 pcomp_init_polycomp(pcomp_chunk_size_t samples_per_chunk,
408  pcomp_poly_size_t num_of_coeffs,
409  double max_allowable_error,
410  pcomp_polycomp_algorithm_t algorithm);
412 
414  const double* input,
415  pcomp_chunk_size_t num_of_samples,
416  pcomp_polycomp_chunk_t* chunk,
417  double* max_error);
418 
419 int pcomp_decompress_polycomp_chunk(double* output,
420  const pcomp_polycomp_chunk_t* chunk,
421  pcomp_chebyshev_t* inv_chebyshev);
422 
423 /***********************************************************************
424  * Polynomial compression (high-level functions)
425  *
426  * The following functions implement the "polynomial compression",
427  * which is based on a mixture of polynomial least-square fitting and
428  * Chebyshev transform techniques. Unlike the other functions defined
429  * above, this group handles memory allocation autonomously. This
430  * means that "output_buf" must not be preallocated before calling one
431  * of the *_compress_* functions, and the pre-existing value of
432  * "output_size" is ignored in the call (it is only set on exit).
433  */
434 
435 pcomp_chunk_size_t
437 pcomp_poly_size_t
439 double pcomp_polycomp_max_error(const pcomp_polycomp_t* params);
446 double pcomp_polycomp_period(const pcomp_polycomp_t* params);
447 
448 void pcomp_polycomp_set_period(pcomp_polycomp_t* params, double period);
449 
451  double* coeffs, double* cheby_residuals,
452  const double* input,
453  double* max_residual);
454 
455 int pcomp_mask_get_bit(const uint8_t* mask, size_t pos);
456 void pcomp_mask_set_bit(uint8_t* mask, size_t pos, int value);
458  pcomp_chebyshev_t* inv_chebyshev,
459  double max_allowable_error,
460  uint8_t* mask, double* max_error);
461 
463  size_t* num_of_chunks,
464  const double* input_buf, size_t input_size,
465  const pcomp_polycomp_t* params);
466 
467 size_t
469  size_t num_of_chunks);
470 
471 /* The buffer "output_buf" must have room for a number of "double"
472  * values which is at least the number returned by
473  * "pcomp_total_num_of_samples". */
475  double* output_buf, pcomp_polycomp_chunk_t* const chunk_array[],
476  size_t num_of_chunks);
477 
478 /* Free the heap memory allocated for an array of chunks */
479 void pcomp_free_chunks(pcomp_polycomp_chunk_t* chunk_array[],
480  size_t num_of_chunks);
481 
483  size_t num_of_chunks);
484 
485 int pcomp_encode_chunks(void* buf, size_t* buf_size,
486  pcomp_polycomp_chunk_t* const chunk_array[],
487  size_t num_of_chunks);
488 
489 int pcomp_decode_chunks(pcomp_polycomp_chunk_t** chunk_array[],
490  size_t* num_of_chunks, const void* buf);
491 
492 #endif /* LIBPOLYCOMP_H_GUARD */
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...
Definition: poly.c:1266
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.
Definition: poly.c:365
int pcomp_run_poly_fit(pcomp_poly_fit_data_t *poly_fit, double *coeffs, const double *points)
Calculates a polynomial least-squares fit.
Definition: poly.c:301
double pcomp_quant_normalization(const pcomp_quant_params_t *params)
Return the normalization constant used for converting floating-point numbers into quantized integers...
Definition: quant.c:170
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...
Definition: poly.c:773
pcomp_transform_direction_t pcomp_chebyshev_direction(const pcomp_chebyshev_t *plan)
Return the direction of a Chebyshev transform.
Definition: poly.c:435
int pcomp_decompress_polycomp_chunk(double *output, const pcomp_polycomp_chunk_t *chunk, pcomp_chebyshev_t *inv_chebyshev)
Decompress the data in a chunk.
Definition: poly.c:1833
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...
Definition: poly.c:677
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.
Definition: poly.c:2193
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. ...
Definition: libpolycomp.h:326
pcomp_polycomp_algorithm_t
Kind of algorithm used for the polynomial compression.
Definition: libpolycomp.h:355
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.
Definition: poly.c:2287
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)
Definition: libpolycomp.h:323
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...
Definition: poly.c:1230
size_t pcomp_quant_element_size(const pcomp_quant_params_t *params)
Return the size (in bytes) of the elements to be quantized.
Definition: quant.c:142
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...
Definition: poly.c:1033
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.
Definition: poly.c:695
void pcomp_free_polycomp(pcomp_polycomp_t *params)
Free the memory allocated by pcomp_init_polycomp for a pcomp_polycomp_t structure.
Definition: poly.c:628
pcomp_polycomp_chunk_t * pcomp_init_chunk(pcomp_chunk_size_t num_of_samples)
Allocate memory for a pcomp_polycomp_chunk_t object.
Definition: poly.c:936
size_t pcomp_rle_bufsize(size_t input_size)
Calculate an upper limit for the size of a buffer holding RLE-encoded streams.
Definition: rle.c:95
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.
Definition: poly.c:1194
pcomp_chunk_size_t pcomp_chunk_num_of_samples(const pcomp_polycomp_chunk_t *chunk)
Return the number of samples in a chunk.
Definition: poly.c:1124
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...
Definition: poly.c:733
pcomp_chunk_size_t pcomp_polycomp_samples_per_chunk(const pcomp_polycomp_t *params)
Return the number of samples per chunk.
Definition: poly.c:655
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...
Definition: poly.c:1212
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...
Definition: poly.c:1522
When needed, apply the Chebyshev transform to the residuals of the polynomial fit.
Definition: libpolycomp.h:358
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.
Definition: poly.c:1450
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...
Definition: diff_rle.c:78
void pcomp_free_chebyshev(pcomp_chebyshev_t *plan)
Free the memory allocated by a previous call to pcomp_init_chebyshev.
Definition: poly.c:390
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.
Definition: poly.c:416
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.
Definition: poly.c:1312
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.
Definition: poly.c:1357
const double * pcomp_chebyshev_output(const pcomp_chebyshev_t *plan)
Return the output (Chebyshev transform) of the last call to pcomp_run_chebyshev.
Definition: poly.c:558
void pcomp_straighten(double *output, const double *input, size_t num_of_samples, double period)
Remove sudden jumps from input.
Definition: poly.c:864
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.
Definition: poly.c:2022
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...
Definition: poly.c:975
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.
Definition: quant.c:155
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.
Definition: poly.c:1143
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.
Definition: poly.c:253
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 .
Definition: poly.c:1466
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.
Definition: poly.c:1666
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.
Definition: poly.c:1295
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.
Definition: poly.c:173
pcomp_quant_params_t * pcomp_init_quant_params(size_t element_size, size_t bits_per_sample)
Initialize a pcomp_quant_params_t structure.
Definition: quant.c:104
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.
Definition: poly.c:799
double pcomp_quant_offset(const pcomp_quant_params_t *params)
Return the additive constant used for converting floating-point numbers into quantized integers...
Definition: quant.c:185
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...
Definition: poly.c:210
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.
Definition: poly.c:1095
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...
Definition: poly.c:752
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.
Definition: poly.c:598
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.
Definition: poly.c:2059
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 ...
Definition: quant.c:244
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...
Definition: quant.c:212
void pcomp_free_quant_params(pcomp_quant_params_t *params)
Free a pcomp_quant_params_t structure.
Definition: quant.c:128
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.
Definition: libpolycomp.h:320
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.
Definition: poly.c:714
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.
Definition: poly.c:1951
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.
Definition: poly.c:2138
int pcomp_chunk_is_compressed(const pcomp_polycomp_chunk_t *chunk)
Return nonzero if the chunk holds data in uncompressed form.
Definition: poly.c:1180
If the absolute value of the residuals of a polynomial fit are too large, store the data in uncompres...
Definition: libpolycomp.h:361
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.
Definition: poly.c:234
void pcomp_free_chunks(pcomp_polycomp_chunk_t *chunk_array[], size_t num_of_chunks)
Free an array of chunks.
Definition: poly.c:2106
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.
Definition: poly.c:491
const double * pcomp_chebyshev_input(const pcomp_chebyshev_t *plan)
Return the input data used in the last call to pcomp_run_chebyshev.
Definition: poly.c:537
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...
Definition: poly.c:1248
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.