Libpolycomp  1.0
A compression/decompression library that implements the polynomial compression and other simple compression schemes
test_polycomp_low_level.c
1 /* test_polycomp_low_level.c - Tests for low-level polynomial
2  * compression functions
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 
27 #include <libpolycomp.h>
28 #include <assert.h>
29 #include <math.h>
30 #include <stdlib.h>
31 
32 #include <stdio.h>
33 
34 #define EPSILON 1.0e-7
35 #define MAX_ALLOWABLE_ERROR 0.3
36 
37 int test_chunk_creation(void)
38 {
39  const double samples[] = { 1.0, 2.0, 3.0, 4.0 };
40  const size_t num_of_samples = sizeof(samples) / sizeof(samples[0]);
41 
42  const double poly[] = { 3.0, 2.0 };
43  const size_t num_of_poly = sizeof(poly) / sizeof(poly[0]);
44 
45  const uint8_t cheby_mask[] = { 0xFF };
46 
47  const double cheby[] = { -1.0, -2.0, -3.0, -4.0 };
48  const size_t num_of_cheby = sizeof(cheby) / sizeof(cheby[0]);
49 
50  const double* values;
51  size_t idx;
52 
53  pcomp_polycomp_chunk_t* chunk = NULL;
54 
55  chunk = pcomp_init_uncompressed_chunk(num_of_samples, &samples[0]);
56  assert(chunk != NULL);
57  assert(!pcomp_chunk_is_compressed(chunk));
58  assert(pcomp_chunk_num_of_samples(chunk) == num_of_samples);
59  pcomp_free_chunk(chunk);
60 
61  assert(pcomp_chunk_cheby_mask_size(num_of_samples) == 1);
62  chunk = pcomp_init_compressed_chunk(num_of_samples, num_of_poly,
63  &poly[0], num_of_cheby,
64  &cheby_mask[0], &cheby[0]);
65  assert(chunk != NULL);
66  assert(pcomp_chunk_is_compressed(chunk));
67 
68  assert(pcomp_chunk_num_of_poly_coeffs(chunk) == num_of_poly);
69  values = pcomp_chunk_poly_coeffs(chunk);
70  for (idx = 0; idx < pcomp_chunk_num_of_poly_coeffs(chunk); ++idx) {
71  assert(values[idx] == poly[idx]);
72  }
73 
74  assert(pcomp_chunk_num_of_cheby_coeffs(chunk) == num_of_cheby);
75  values = pcomp_chunk_cheby_coeffs(chunk);
76  for (idx = 0; idx < pcomp_chunk_num_of_cheby_coeffs(chunk); ++idx) {
77  assert(values[idx] == cheby[idx]);
78  }
79 
80  pcomp_free_chunk(chunk);
81 
82  return 0;
83 }
84 
85 int test_no_compression(void)
86 {
87  /* It is impossible to compress these data using the polynomial
88  * compression algorithm, as they contain too many jumps */
89  double input[]
90  = { 1.0, -2.0, 3.0, -4.0, 5.0, 6.0, 7.0, -8.0, 9.0, -10.0 };
91  size_t input_size = sizeof(input) / sizeof(input[0]);
92  double* decompr = malloc(sizeof(double) * input_size);
93  pcomp_polycomp_chunk_t* chunk = pcomp_init_chunk(input_size);
94  double max_error;
95  pcomp_chebyshev_t* inv_chebyshev
98  input_size, 2, MAX_ALLOWABLE_ERROR, PCOMP_ALG_USE_CHEBYSHEV);
99  size_t idx;
100 
101  pcomp_run_polycomp_on_chunk(polycomp, input, input_size, chunk,
102  &max_error);
103  assert(pcomp_chunk_num_of_samples(chunk) == 10);
104  assert(!pcomp_chunk_is_compressed(chunk));
105  assert(pcomp_chunk_poly_coeffs(chunk) == NULL);
106  assert(pcomp_chunk_cheby_coeffs(chunk) == NULL);
107 
108  pcomp_decompress_polycomp_chunk(decompr, chunk, inv_chebyshev);
109  for (idx = 0; idx < input_size; ++idx) {
110  assert(input[idx] == decompr[idx]);
111  }
112 
113  free(decompr);
114  pcomp_free_chunk(chunk);
115  pcomp_free_chebyshev(inv_chebyshev);
116  pcomp_free_polycomp(polycomp);
117 
118  return 0;
119 }
120 
121 int test_no_chebyshev(void)
122 {
123  /* These data are fitted perfectly by a straight line */
124  double input[]
125  = { 1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0, 8.0, 9.0, 10.0 };
126  size_t input_size = sizeof(input) / sizeof(input[0]);
127  double* decompr = malloc(sizeof(double) * input_size);
128  pcomp_polycomp_chunk_t* chunk = pcomp_init_chunk(input_size);
129  double max_error;
130  pcomp_chebyshev_t* inv_chebyshev
131  = pcomp_init_chebyshev(input_size, PCOMP_TD_INVERSE);
133  input_size, 2, MAX_ALLOWABLE_ERROR, PCOMP_ALG_USE_CHEBYSHEV);
134  size_t idx;
135 
136  pcomp_run_polycomp_on_chunk(polycomp, input, input_size, chunk,
137  &max_error);
138  assert(pcomp_chunk_num_of_samples(chunk) == 10);
139  assert(pcomp_chunk_is_compressed(chunk));
140  assert(pcomp_chunk_num_of_poly_coeffs(chunk) == 2);
141  assert(pcomp_chunk_poly_coeffs(chunk) != NULL);
142  assert(fabs(pcomp_chunk_poly_coeffs(chunk)[0] - 0.0) < EPSILON);
143  assert(fabs(pcomp_chunk_poly_coeffs(chunk)[1] - 1.0) < EPSILON);
144  assert(pcomp_chunk_cheby_coeffs(chunk) == NULL);
145 
146  pcomp_decompress_polycomp_chunk(decompr, chunk, inv_chebyshev);
147  for (idx = 0; idx < input_size; ++idx) {
148  assert(fabs(input[idx] - decompr[idx]) < EPSILON);
149  }
150 
151  free(decompr);
152  pcomp_free_chunk(chunk);
153  pcomp_free_chebyshev(inv_chebyshev);
154  pcomp_free_polycomp(polycomp);
155 
156  return 0;
157 }
158 
159 int test_complete_compression_and_decompression(void)
160 {
161  double input[]
162  = { 1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0, 8.0, 9.0, 11.0 };
163  size_t input_size = sizeof(input) / sizeof(input[0]);
164  double* decompr = malloc(sizeof(double) * input_size);
165  pcomp_polycomp_chunk_t* chunk = pcomp_init_chunk(input_size);
166  double max_error;
167  pcomp_chebyshev_t* inv_chebyshev
168  = pcomp_init_chebyshev(input_size, PCOMP_TD_INVERSE);
170  input_size, 2, MAX_ALLOWABLE_ERROR, PCOMP_ALG_USE_CHEBYSHEV);
171  size_t idx;
172 
173  pcomp_run_polycomp_on_chunk(polycomp, input, input_size, chunk,
174  &max_error);
175  assert(pcomp_chunk_num_of_samples(chunk) == 10);
176  assert(pcomp_chunk_is_compressed(chunk));
177  assert(pcomp_chunk_num_of_poly_coeffs(chunk) == 2);
178  assert(pcomp_chunk_poly_coeffs(chunk) != NULL);
179  assert(pcomp_chunk_num_of_cheby_coeffs(chunk) < 10);
180 
181  pcomp_decompress_polycomp_chunk(decompr, chunk, inv_chebyshev);
182  for (idx = 0; idx < input_size; ++idx) {
183  assert(fabs(input[idx] - decompr[idx]) < MAX_ALLOWABLE_ERROR);
184  }
185 
186  free(decompr);
187  pcomp_free_chunk(chunk);
188  pcomp_free_chebyshev(inv_chebyshev);
189  pcomp_free_polycomp(polycomp);
190 
191  return 0;
192 }
193 
194 int main(void)
195 {
196  int result;
197 
198  result = test_chunk_creation();
199  if (result != 0)
200  return result;
201 
202  result = test_no_compression();
203  if (result != 0)
204  return result;
205 
206  result = test_no_chebyshev();
207  if (result != 0)
208  return result;
209 
210  result = test_complete_compression_and_decompression();
211  if (result != 0)
212  return result;
213 
214  return 0;
215 }
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
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_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
Compute a backward Chebyshev transform, with a normalization factor equal to one. ...
Definition: libpolycomp.h:326
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
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
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
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_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
When needed, apply the Chebyshev transform to the residuals of the polynomial fit.
Definition: libpolycomp.h:358
void pcomp_free_chebyshev(pcomp_chebyshev_t *plan)
Free the memory allocated by a previous call to pcomp_init_chebyshev.
Definition: poly.c:390
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
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
void pcomp_free_chunk(pcomp_polycomp_chunk_t *chunk)
Free memory associated with a pcomp_poly_chunk_t.
Definition: poly.c:1095
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_chunk_is_compressed(const pcomp_polycomp_chunk_t *chunk)
Return nonzero if the chunk holds data in uncompressed form.
Definition: poly.c:1180
Header file for Libpolycomp.
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