File indexing completed on 2024-04-28 15:14:12

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