File indexing completed on 2024-05-12 15:27:10
0001 /*************************************************************************** 0002 File : nsl_smooth.h 0003 Project : LabPlot 0004 Description : NSL smooth functions 0005 -------------------------------------------------------------------- 0006 Copyright : (C) 2016 by Stefan Gerlach (stefan.gerlach@uni.kn) 0007 0008 ***************************************************************************/ 0009 0010 /*************************************************************************** 0011 * * 0012 * This program is free software; you can redistribute it and/or modify * 0013 * it under the terms of the GNU General Public License as published by * 0014 * the Free Software Foundation; either version 2 of the License, or * 0015 * (at your option) any later version. * 0016 * * 0017 * This program is distributed in the hope that it will be useful, * 0018 * but WITHOUT ANY WARRANTY; without even the implied warranty of * 0019 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the * 0020 * GNU General Public License for more details. * 0021 * * 0022 * You should have received a copy of the GNU General Public License * 0023 * along with this program; if not, write to the Free Software * 0024 * Foundation, Inc., 51 Franklin Street, Fifth Floor, * 0025 * Boston, MA 02110-1301 USA * 0026 * * 0027 ***************************************************************************/ 0028 0029 #ifndef NSL_SMOOTH_H 0030 #define NSL_SMOOTH_H 0031 0032 #include <gsl/gsl_matrix.h> 0033 0034 #define NSL_SMOOTH_TYPE_COUNT 4 0035 typedef enum {nsl_smooth_type_moving_average, nsl_smooth_type_moving_average_lagged, nsl_smooth_type_percentile, 0036 nsl_smooth_type_savitzky_golay} nsl_smooth_type; 0037 extern const char* nsl_smooth_type_name[]; 0038 /* TODO: LOWESS/etc., Bezier, B-Spline, (FFT Filter) */ 0039 0040 /* mode of extension for padding signal 0041 * none: reduce points at edges 0042 * interp: polynomial interpolation 0043 * mirror: 3 2 | 1 2 3 4 5 | 4 3 (reflect) 0044 * nearest: 1 1 | 1 2 3 4 5 | 5 5 (repeat) 0045 * constant: L L | 1 2 3 4 5 | R R 0046 * periodic: 4 5 | 1 2 3 4 5 | 1 2 (wrap) 0047 */ 0048 #define NSL_SMOOTH_PAD_MODE_COUNT 6 0049 typedef enum {nsl_smooth_pad_none, nsl_smooth_pad_interp, nsl_smooth_pad_mirror, nsl_smooth_pad_nearest, 0050 nsl_smooth_pad_constant, nsl_smooth_pad_periodic} nsl_smooth_pad_mode; 0051 extern const char* nsl_smooth_pad_mode_name[]; 0052 0053 #define NSL_SMOOTH_WEIGHT_TYPE_COUNT 8 0054 typedef enum {nsl_smooth_weight_uniform, nsl_smooth_weight_triangular, nsl_smooth_weight_binomial, 0055 nsl_smooth_weight_parabolic, nsl_smooth_weight_quartic, nsl_smooth_weight_triweight, 0056 nsl_smooth_weight_tricube, nsl_smooth_weight_cosine } nsl_smooth_weight_type; 0057 extern const char* nsl_smooth_weight_type_name[]; 0058 /*TODO: IIR: exponential, Gaussian, see nsl_sf_kernel */ 0059 0060 /* values used for constant padding */ 0061 extern double nsl_smooth_pad_constant_lvalue, nsl_smooth_pad_constant_rvalue; 0062 0063 /********* Smoothing algorithms **********/ 0064 0065 /* Moving average */ 0066 int nsl_smooth_moving_average(double *data, size_t n, size_t points, nsl_smooth_weight_type weight, nsl_smooth_pad_mode mode); 0067 0068 /* Lagged moving average */ 0069 int nsl_smooth_moving_average_lagged(double *data, size_t n, size_t points, nsl_smooth_weight_type weight, nsl_smooth_pad_mode mode); 0070 0071 /* Percentile filter */ 0072 int nsl_smooth_percentile(double *data, size_t n, size_t points, double percentile, nsl_smooth_pad_mode mode); 0073 0074 /* Savitzky-Golay coefficients */ 0075 /** 0076 * \brief Compute Savitzky-Golay coefficients and store them into #h. 0077 * 0078 * This function follows GSL conventions in that it writes its result into a matrix allocated by 0079 * the caller and returns a non-zero result on error. 0080 * 0081 * The coefficient matrix is defined as the matrix H mapping a set of input values to the values 0082 * of the polynomial of order #polynom_order which minimizes squared deviations from the input 0083 * values. It is computed using the formula \$H=V(V^TV)^(-1)V^T\$, where \$V\$ is the Vandermonde 0084 * matrix of the point indices. 0085 * 0086 * For a short description of the mathematical background, see 0087 * http://www.statistics4u.info/fundstat_eng/cc_filter_savgol_math.html 0088 */ 0089 int nsl_smooth_savgol_coeff(size_t points, int order, gsl_matrix *h); 0090 0091 /* set values for constant padding */ 0092 void nsl_smooth_pad_constant_set(double lvalue, double rvalue); 0093 0094 /* Savitzky-Golay smoothing */ 0095 /** 0096 * \brief Savitzky-Golay smoothing of (uniformly distributed) data. 0097 * 0098 * When the data is not uniformly distributed, Savitzky-Golay looses its interesting conservation 0099 * properties. On the other hand, a central point of the algorithm is that for uniform data, the 0100 * operation can be implemented as a convolution. This is considerably more efficient than a more 0101 * generic method able to handle non-uniform input data. 0102 */ 0103 int nsl_smooth_savgol(double *data, size_t n, size_t points, int order, nsl_smooth_pad_mode mode); 0104 0105 /* Savitzky-Golay default smoothing (interp) */ 0106 int nsl_smooth_savgol_default(double *data, size_t n, size_t points, int order); 0107 0108 /* TODO SmoothFilter::smoothModifiedSavGol(double *x_in, double *y_inout) 0109 see SmoothFilter.cpp of libscidavis 0110 */ 0111 0112 #endif /* NSL_SMOOTH_H */ 0113 0114