File indexing completed on 2024-04-21 03:48:14
0001 /* 0002 File : NSLSmoothTest.cpp 0003 Project : LabPlot 0004 Description : NSL Tests for smoothing 0005 -------------------------------------------------------------------- 0006 SPDX-FileCopyrightText: 2019 Stefan Gerlach <stefan.gerlach@uni.kn> 0007 0008 SPDX-License-Identifier: GPL-2.0-or-later 0009 */ 0010 0011 #include "NSLSmoothTest.h" 0012 0013 extern "C" { 0014 #include "backend/nsl/nsl_smooth.h" 0015 } 0016 0017 // ############################################################################## 0018 // ################# moving average tests 0019 // ############################################################################## 0020 0021 const int N = 9; 0022 const int mapoints = 5; 0023 const nsl_smooth_weight_type weight = nsl_smooth_weight_uniform; 0024 0025 void NSLSmoothTest::testMA_padnone() { 0026 double data[] = {2, 2, 5, 2, 1, 0, 1, 4, 9}; 0027 const double result[] = {2, 3, 2.4, 2, 1.8, 1.6, 3, 14 / 3., 9}; 0028 0029 int status = nsl_smooth_moving_average(data, N, mapoints, weight, nsl_smooth_pad_none); 0030 QCOMPARE(status, 0); 0031 for (int i = 0; i < N; i++) 0032 QCOMPARE(data[i], result[i]); 0033 } 0034 0035 void NSLSmoothTest::testMA_padmirror() { 0036 double data[] = {2, 2, 5, 2, 1, 0, 1, 4, 9}; 0037 const double result[] = {3.2, 2.6, 2.4, 2, 1.8, 1.6, 3, 3.6, 3.8}; 0038 0039 int status = nsl_smooth_moving_average(data, N, mapoints, weight, nsl_smooth_pad_mirror); 0040 QCOMPARE(status, 0); 0041 for (int i = 0; i < N; i++) 0042 QCOMPARE(data[i], result[i]); 0043 } 0044 0045 void NSLSmoothTest::testMA_padnearest() { 0046 double data[] = {2, 2, 5, 2, 1, 0, 1, 4, 9}; 0047 const double result[] = {2.6, 2.6, 2.4, 2, 1.8, 1.6, 3, 4.6, 6.4}; 0048 0049 int status = nsl_smooth_moving_average(data, N, mapoints, weight, nsl_smooth_pad_nearest); 0050 QCOMPARE(status, 0); 0051 for (int i = 0; i < N; i++) 0052 QCOMPARE(data[i], result[i]); 0053 } 0054 0055 void NSLSmoothTest::testMA_padconstant() { 0056 double data[] = {2, 2, 5, 2, 1, 0, 1, 4, 9}; 0057 const double result[] = {1.8, 2.2, 2.4, 2, 1.8, 1.6, 3, 2.8, 2.8}; 0058 0059 int status = nsl_smooth_moving_average(data, N, mapoints, weight, nsl_smooth_pad_constant); 0060 QCOMPARE(status, 0); 0061 for (int i = 0; i < N; i++) 0062 QCOMPARE(data[i], result[i]); 0063 } 0064 0065 void NSLSmoothTest::testMA_padperiodic() { 0066 double data[] = {2, 2, 5, 2, 1, 0, 1, 4, 9}; 0067 const double result[] = {4.4, 4, 2.4, 2, 1.8, 1.6, 3, 3.2, 3.6}; 0068 0069 int status = nsl_smooth_moving_average(data, N, mapoints, weight, nsl_smooth_pad_periodic); 0070 QCOMPARE(status, 0); 0071 for (int i = 0; i < N; i++) 0072 QCOMPARE(data[i], result[i]); 0073 } 0074 0075 // ############################################################################## 0076 // ################# lagged moving average tests 0077 // ############################################################################## 0078 0079 void NSLSmoothTest::testMAL_padnone() { 0080 double data[] = {2, 2, 5, 2, 1, 0, 1, 4, 9}; 0081 const double result[] = {2, 2, 3, 2.75, 2.4, 2, 1.8, 1.6, 3}; 0082 0083 int status = nsl_smooth_moving_average_lagged(data, N, mapoints, weight, nsl_smooth_pad_none); 0084 QCOMPARE(status, 0); 0085 for (int i = 0; i < N; i++) 0086 QCOMPARE(data[i], result[i]); 0087 } 0088 0089 void NSLSmoothTest::testMAL_padmirror() { 0090 double data[] = {2, 2, 5, 2, 1, 0, 1, 4, 9}; 0091 const double result[] = {2.4, 2.6, 3.2, 2.6, 2.4, 2, 1.8, 1.6, 3}; 0092 0093 int status = nsl_smooth_moving_average_lagged(data, N, mapoints, weight, nsl_smooth_pad_mirror); 0094 QCOMPARE(status, 0); 0095 for (int i = 0; i < N; i++) 0096 QCOMPARE(data[i], result[i]); 0097 } 0098 0099 void NSLSmoothTest::testMAL_padnearest() { 0100 double data[] = {2, 2, 5, 2, 1, 0, 1, 4, 9}; 0101 const double result[] = {2., 2., 2.6, 2.6, 2.4, 2, 1.8, 1.6, 3}; 0102 0103 int status = nsl_smooth_moving_average_lagged(data, N, mapoints, weight, nsl_smooth_pad_nearest); 0104 QCOMPARE(status, 0); 0105 for (int i = 0; i < N; i++) 0106 QCOMPARE(data[i], result[i]); 0107 } 0108 0109 void NSLSmoothTest::testMAL_padconstant() { 0110 double data[] = {2, 2, 5, 2, 1, 0, 1, 4, 9}; 0111 const double result[] = {0.4, 0.8, 1.8, 2.2, 2.4, 2, 1.8, 1.6, 3}; 0112 0113 int status = nsl_smooth_moving_average_lagged(data, N, mapoints, weight, nsl_smooth_pad_constant); 0114 QCOMPARE(status, 0); 0115 for (int i = 0; i < N; i++) 0116 QCOMPARE(data[i], result[i]); 0117 } 0118 0119 void NSLSmoothTest::testMAL_padperiodic() { 0120 double data[] = {2, 2, 5, 2, 1, 0, 1, 4, 9}; 0121 const double result[] = {3.2, 3.6, 4.4, 4, 2.4, 2, 1.8, 1.6, 3}; 0122 0123 int status = nsl_smooth_moving_average_lagged(data, N, mapoints, weight, nsl_smooth_pad_periodic); 0124 QCOMPARE(status, 0); 0125 for (int i = 0; i < N; i++) 0126 QCOMPARE(data[i], result[i]); 0127 } 0128 0129 // ############################################################################## 0130 // ################# percentile tests 0131 // ############################################################################## 0132 0133 const double percentile = 0.5; 0134 0135 void NSLSmoothTest::testPercentile_padnone() { 0136 double data[] = {2, 2, 5, 2, 1, 0, 1, 4, 9}; 0137 const double result[] = {2, 2, 2, 2, 1, 1, 1, 4, 9}; 0138 0139 int status = nsl_smooth_percentile(data, N, mapoints, percentile, nsl_smooth_pad_none); 0140 QCOMPARE(status, 0); 0141 for (int i = 0; i < N; i++) 0142 QCOMPARE(data[i], result[i]); 0143 } 0144 0145 void NSLSmoothTest::testPercentile_padmirror() { 0146 double data[] = {2, 2, 5, 2, 1, 0, 1, 4, 9}; 0147 const double result[] = {2, 2, 2, 2, 1, 1, 1, 4, 4}; 0148 0149 int status = nsl_smooth_percentile(data, N, mapoints, percentile, nsl_smooth_pad_mirror); 0150 QCOMPARE(status, 0); 0151 for (int i = 0; i < N; i++) 0152 QCOMPARE(data[i], result[i]); 0153 } 0154 0155 void NSLSmoothTest::testPercentile_padnearest() { 0156 double data[] = {2, 2, 5, 2, 1, 0, 1, 4, 9}; 0157 const double result[] = {2, 2, 2, 2, 1, 1, 1, 4, 9}; 0158 0159 int status = nsl_smooth_percentile(data, N, mapoints, percentile, nsl_smooth_pad_nearest); 0160 QCOMPARE(status, 0); 0161 for (int i = 0; i < N; i++) 0162 QCOMPARE(data[i], result[i]); 0163 } 0164 0165 void NSLSmoothTest::testPercentile_padconstant() { 0166 double data[] = {2, 2, 5, 2, 1, 0, 1, 4, 9}; 0167 const double result[] = {2, 2, 2, 2, 1, 1, 1, 1, 1}; 0168 0169 int status = nsl_smooth_percentile(data, N, mapoints, percentile, nsl_smooth_pad_constant); 0170 QCOMPARE(status, 0); 0171 for (int i = 0; i < N; i++) 0172 QCOMPARE(data[i], result[i]); 0173 } 0174 0175 void NSLSmoothTest::testPercentile_padperiodic() { 0176 double data[] = {2, 2, 5, 2, 1, 0, 1, 4, 9}; 0177 const double result[] = {4, 2, 2, 2, 1, 1, 1, 2, 2}; 0178 0179 int status = nsl_smooth_percentile(data, N, mapoints, percentile, nsl_smooth_pad_periodic); 0180 QCOMPARE(status, 0); 0181 for (int i = 0; i < N; i++) 0182 QCOMPARE(data[i], result[i]); 0183 } 0184 0185 // ############################################################################## 0186 // ################# Savitzky-Golay coeff tests 0187 // ############################################################################## 0188 0189 void NSLSmoothTest::testSG_coeff31() { 0190 int points = 3, order = 1; 0191 gsl_matrix* h = gsl_matrix_alloc(points, points); 0192 const double result[] = {1, 1, 1}; 0193 0194 int status = nsl_smooth_savgol_coeff(points, order, h); 0195 QCOMPARE(status, 0); 0196 for (int i = 0; i < points; i++) 0197 QCOMPARE(points * gsl_matrix_get(h, (points - 1) / 2, i), result[i]); 0198 0199 gsl_matrix_free(h); 0200 } 0201 0202 void NSLSmoothTest::testSG_coeff51() { 0203 int points = 5, order = 1; 0204 gsl_matrix* h = gsl_matrix_alloc(points, points); 0205 const double result[] = {1, 1, 1, 1, 1}; 0206 0207 int status = nsl_smooth_savgol_coeff(points, order, h); 0208 QCOMPARE(status, 0); 0209 for (int i = 0; i < points; i++) 0210 QCOMPARE(points * gsl_matrix_get(h, (points - 1) / 2, i), result[i]); 0211 0212 gsl_matrix_free(h); 0213 } 0214 0215 void NSLSmoothTest::testSG_coeff53() { 0216 int points = 5, order = 3; 0217 gsl_matrix* h = gsl_matrix_alloc(points, points); 0218 const double result[] = {-3, 12, 17, 12, -3}; 0219 0220 int status = nsl_smooth_savgol_coeff(points, order, h); 0221 QCOMPARE(status, 0); 0222 for (int i = 0; i < points; i++) 0223 QCOMPARE(35 * gsl_matrix_get(h, (points - 1) / 2, i), result[i]); 0224 0225 gsl_matrix_free(h); 0226 } 0227 0228 void NSLSmoothTest::testSG_coeff73() { 0229 int points = 7, order = 3; 0230 gsl_matrix* h = gsl_matrix_alloc(points, points); 0231 const double result[] = {-2, 3, 6, 7, 6, 3, -2}; 0232 0233 int status = nsl_smooth_savgol_coeff(points, order, h); 0234 QCOMPARE(status, 0); 0235 for (int i = 0; i < points; i++) 0236 QCOMPARE(21 * gsl_matrix_get(h, (points - 1) / 2, i), result[i]); 0237 0238 gsl_matrix_free(h); 0239 } 0240 0241 void NSLSmoothTest::testSG_coeff74() { 0242 int points = 7, order = 4; 0243 gsl_matrix* h = gsl_matrix_alloc(points, points); 0244 const double result[] = {5, -30, 75, 131, 75, -30, 5}; 0245 0246 int status = nsl_smooth_savgol_coeff(points, order, h); 0247 QCOMPARE(status, 0); 0248 0249 for (int i = 0; i < points; i++) { 0250 FuzzyCompare(result[i], 231 * gsl_matrix_get(h, (points - 1) / 2, i), 1.e-10); 0251 // QCOMPARE(231 * gsl_matrix_get(h,(points-1)/2, i), result[i]); 0252 } 0253 0254 gsl_matrix_free(h); 0255 } 0256 0257 void NSLSmoothTest::testSG_coeff92() { 0258 int points = 9, order = 2; 0259 gsl_matrix* h = gsl_matrix_alloc(points, points); 0260 const double result[] = {-21, 14, 39, 54, 59, 54, 39, 14, -21}; 0261 0262 int status = nsl_smooth_savgol_coeff(points, order, h); 0263 QCOMPARE(status, 0); 0264 0265 for (int i = 0; i < points; i++) { 0266 // FuzzyCompare(result[i], 231 * gsl_matrix_get(h,(points-1)/2, i), 1.e-10); 0267 QCOMPARE(231 * gsl_matrix_get(h, (points - 1) / 2, i), result[i]); 0268 } 0269 0270 gsl_matrix_free(h); 0271 } 0272 0273 void NSLSmoothTest::testSG_coeff94() { 0274 int points = 9, order = 4; 0275 gsl_matrix* h = gsl_matrix_alloc(points, points); 0276 const double result[] = {15, -55, 30, 135, 179, 135, 30, -55, 15}; 0277 0278 int status = nsl_smooth_savgol_coeff(points, order, h); 0279 QCOMPARE(status, 0); 0280 0281 for (int i = 0; i < points; i++) { 0282 FuzzyCompare(result[i], 429 * gsl_matrix_get(h, (points - 1) / 2, i), 1.e-10); 0283 // QCOMPARE(429 * gsl_matrix_get(h,(points-1)/2, i), result[i]); 0284 } 0285 0286 gsl_matrix_free(h); 0287 } 0288 0289 // ############################################################################## 0290 // ################# Savitzky-Golay modes 0291 // ############################################################################## 0292 0293 const int n = 9, m = 5, sgorder = 2; 0294 0295 void NSLSmoothTest::testSG_mode_interp() { 0296 double data[] = {2, 2, 5, 2, 1, 0, 1, 4, 9}; 0297 const double result[] = {1.65714285714286, 3.17142857142857, 3.54285714285714, 2.85714285714285, 0.65714285714287, 0.17142857142858, 1., 4., 9.}; 0298 0299 int status = nsl_smooth_savgol(data, n, m, sgorder, nsl_smooth_pad_interp); 0300 QCOMPARE(status, 0); 0301 0302 for (int i = 0; i < n; i++) 0303 QCOMPARE(data[i], result[i]); 0304 // FuzzyCompare(data[i], result[i], 1.e-11); 0305 } 0306 0307 void NSLSmoothTest::testSG_mode_mirror() { 0308 double data[] = {2, 2, 5, 2, 1, 0, 1, 4, 9}; 0309 const double result[] = 0310 {1.48571428571430, 3.02857142857143, 3.542857142857, 2.857142857143, 0.657142857143, 0.171428571429, 1., 5.02857142857142, 6.94285714285713}; 0311 0312 int status = nsl_smooth_savgol(data, n, m, sgorder, nsl_smooth_pad_mirror); 0313 QCOMPARE(status, 0); 0314 0315 for (int i = 0; i < n; i++) 0316 FuzzyCompare(data[i], result[i], 1.e-11); 0317 // QCOMPARE(data[i], result[i]); 0318 } 0319 0320 void NSLSmoothTest::testSG_mode_nearest() { 0321 double data[] = {2, 2, 5, 2, 1, 0, 1, 4, 9}; 0322 const double result[] = {1.74285714285715, 3.02857142857143, 3.542857142857, 2.857142857143, 0.657142857143, 0.171428571429, 1., 4.6, 7.97142857142856}; 0323 0324 int status = nsl_smooth_savgol(data, n, m, sgorder, nsl_smooth_pad_nearest); 0325 QCOMPARE(status, 0); 0326 0327 for (int i = 0; i < n; i++) 0328 FuzzyCompare(data[i], result[i], 1.e-11); 0329 // QCOMPARE(data[i], result[i]); 0330 } 0331 0332 void NSLSmoothTest::testSG_mode_constant() { 0333 double data[] = {2, 2, 5, 2, 1, 0, 1, 4, 9}; 0334 const double result[] = {1.22857142857143, 3.2, 3.542857142857, 2.857142857143, 0.657142857143, 0.171428571429, 1., 5.37142857142855, 5.65714285714284}; 0335 0336 int status = nsl_smooth_savgol(data, n, m, sgorder, nsl_smooth_pad_constant); 0337 QCOMPARE(status, 0); 0338 0339 for (int i = 0; i < n; i++) 0340 FuzzyCompare(data[i], result[i], 1.e-11); 0341 // QCOMPARE(data[i], result[i]); 0342 } 0343 0344 void NSLSmoothTest::testSG_mode_periodic() { 0345 double data[] = {2, 2, 5, 2, 1, 0, 1, 4, 9}; 0346 const double result[] = {3.97142857142858, 2.42857142857144, 3.542857142857, 2.857142857143, 0.657142857143, 0.171428571429, 1., 5.2, 6.17142857142856}; 0347 0348 int status = nsl_smooth_savgol(data, n, m, sgorder, nsl_smooth_pad_periodic); 0349 QCOMPARE(status, 0); 0350 0351 for (int i = 0; i < n; i++) 0352 FuzzyCompare(data[i], result[i], 1.e-11); 0353 // QCOMPARE(data[i], result[i]); 0354 } 0355 0356 // ############################################################################## 0357 // ################# performance 0358 // ############################################################################## 0359 0360 const int nn = 1e6; 0361 0362 void NSLSmoothTest::testPerformance_interp() { 0363 QScopedArrayPointer<double> data(new double[nn]); 0364 0365 QBENCHMARK { 0366 for (int i = 0; i < nn; i++) 0367 data[i] = i; 0368 int status = nsl_smooth_savgol(data.data(), nn, m, sgorder, nsl_smooth_pad_interp); 0369 QCOMPARE(status, 0); 0370 } 0371 } 0372 0373 void NSLSmoothTest::testPerformance_mirror() { 0374 QScopedArrayPointer<double> data(new double[nn]); 0375 0376 QBENCHMARK { 0377 for (int i = 0; i < nn; i++) 0378 data[i] = i; 0379 int status = nsl_smooth_savgol(data.data(), nn, m, sgorder, nsl_smooth_pad_mirror); 0380 QCOMPARE(status, 0); 0381 } 0382 } 0383 0384 void NSLSmoothTest::testPerformance_nearest() { 0385 QScopedArrayPointer<double> data(new double[nn]); 0386 0387 QBENCHMARK { 0388 for (int i = 0; i < nn; i++) 0389 data[i] = i; 0390 int status = nsl_smooth_savgol(data.data(), nn, m, sgorder, nsl_smooth_pad_nearest); 0391 QCOMPARE(status, 0); 0392 } 0393 } 0394 0395 void NSLSmoothTest::testPerformance_constant() { 0396 QScopedArrayPointer<double> data(new double[nn]); 0397 0398 QBENCHMARK { 0399 for (int i = 0; i < nn; i++) 0400 data[i] = i; 0401 int status = nsl_smooth_savgol(data.data(), nn, m, sgorder, nsl_smooth_pad_constant); 0402 QCOMPARE(status, 0); 0403 } 0404 } 0405 0406 void NSLSmoothTest::testPerformance_periodic() { 0407 QScopedArrayPointer<double> data(new double[nn]); 0408 0409 QBENCHMARK { 0410 for (int i = 0; i < nn; i++) 0411 data[i] = i; 0412 int status = nsl_smooth_savgol(data.data(), nn, m, sgorder, nsl_smooth_pad_periodic); 0413 QCOMPARE(status, 0); 0414 } 0415 } 0416 0417 QTEST_MAIN(NSLSmoothTest)