Skip to content

Documentation for Parameters Estimation

Least Squares Methods for parameter estimation

Estimators

Ordinary Least Squares for linear parameter estimation

Source code in sysidentpy\parameter_estimation\estimators.py
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
class Estimators:
    """Ordinary Least Squares for linear parameter estimation"""

    def __init__(
        self,
        max_lag=1,
        lam=0.98,
        delta=0.01,
        offset_covariance=0.2,
        mu=0.01,
        eps=np.finfo(np.float64).eps,
        ridge_param=np.finfo(np.float64).eps,  # for regularized ridge regression
        gama=0.2,
        weight=0.02,
        basis_function=None,
    ):
        self.eps = eps
        self.ridge_param = ridge_param  # for regularized ridge regression
        self.mu = mu
        self.offset_covariance = offset_covariance
        self.max_lag = max_lag
        self.lam = lam
        self.delta = delta
        self.gama = gama
        self.weight = weight  # <0  e <1
        self.xi = None
        self.theta_evolution = None
        self.basis_function = basis_function
        self._validate_params()

    def _validate_params(self):
        """Validate input params."""
        attributes = {
            "max_lag": self.max_lag,
            "lam": self.lam,
            "delta": self.delta,
            "offset_covariance": self.offset_covariance,
            "mu": self.mu,
            "eps": self.eps,
            "ridge_param": self.ridge_param,
            "gama": self.gama,
            "weight": self.weight,
        }
        for attribute, value in attributes.items():
            if not isinstance(value, (np.integer, int, float)):
                raise ValueError(
                    f"{attribute} must be int or float (positive).Got {type(attribute)}"
                )

            if attribute in ["lam", "weight", "offset_covariance"]:
                if value > 1 or value < 0:
                    raise ValueError(
                        f"{attribute} must lies on [0 1] range. Got {value}"
                    )

            if value < 0:
                raise ValueError(
                    f"{attribute} must be positive. Got {value}"
                    "Check the documentation for allowed values"
                )

    def _check_linear_dependence_rows(self, psi):
        if np.linalg.matrix_rank(psi) != psi.shape[1]:
            warnings.warn(
                "Psi matrix might have linearly dependent rows."
                "Be careful and check your data",
                stacklevel=2,
            )

    def least_squares(self, psi, y):
        """Estimate the model parameters using Least Squares method.

        Parameters
        ----------
        psi : ndarray of floats
            The information matrix of the model.
        y : array-like of shape = y_training
            The data used to training the model.

        Returns
        -------
        theta : array-like of shape = number_of_model_elements
            The estimated parameters of the model.

        References
        ----------
        - Manuscript: Sorenson, H. W. (1970). Least-squares estimation:
           from Gauss to Kalman. IEEE spectrum, 7(7), 63-68.
           http://pzs.dstu.dp.ua/DataMining/mls/bibl/Gauss2Kalman.pdf
        - Book (Portuguese): Aguirre, L. A. (2007). Introdução identificação
           de sistemas: técnicas lineares e não-lineares aplicadas a sistemas
           reais. Editora da UFMG. 3a edição.
        - Manuscript: Markovsky, I., & Van Huffel, S. (2007).
           Overview of total least-squares methods.
           Signal processing, 87(10), 2283-2302.
           https://eprints.soton.ac.uk/263855/1/tls_overview.pdf
        - Wikipedia entry on Least Squares
           https://en.wikipedia.org/wiki/Least_squares

        """
        self._check_linear_dependence_rows(psi)

        y = y[self.max_lag :, 0].reshape(-1, 1)
        theta = np.linalg.lstsq(psi, y, rcond=None)[0]
        return theta

    def ridge_regression(self, psi, y):
        """Estimate the model parameters using the regularized least squares method
           known as ridge regression.  Based on the least_squares module and uses
           the same data format but you need to pass ridge_param in the call to
           FROLS

        Parameters
        ----------
        psi : ndarray of floats
            The information matrix of the model.
        y : array-like of shape = y_training
            The data used to training the model.
        ridge_param : ridge regression parameter that regularizes the algorithm
            to prevent over fitting.  If the input is a noisy signal, the ridge
            parameter is likely to be set close to the noise level, at least
            as a starting point.  Entered through the self data structure.

        Returns
        -------
        theta : array-like of shape = number_of_model_elements
            The estimated parameters of the model.

        References
        ----------
        - Wikipedia entry on ridge regression
          https://en.wikipedia.org/wiki/Ridge_regression

        ridge_parm multiplied by the identity matrix (np.eye) favors models (theta) that
        have small size using an L2 norm.  This prevents over fitting of the model.
        For applications where preventing overfitting is important, see, for example,
        D. J. Gauthier, E. Bollt, A. Griffith, W. A. S. Barbosa, ‘Next generation
        reservoir computing,’ Nat. Commun. 12, 5564 (2021).
        https://www.nature.com/articles/s41467-021-25801-2

        """
        self._check_linear_dependence_rows(psi)

        y = y[self.max_lag :, 0].reshape(-1, 1)
        theta = (
            np.linalg.pinv(psi.T @ psi + self.ridge_param * np.eye(psi.shape[1]))
            @ psi.T
            @ y
        )
        return theta

    def _unbiased_estimator(self, psi, y, theta, elag, max_lag, estimator):
        """Estimate the model parameters using Extended Least Squares method.

        Parameters
        ----------
        psi : ndarray of floats
            The information matrix of the model.
        X : ndarray of floats
            The input data to be used in the training process.
        y : array-like of shape = y_training
            The data used to training the model.
        biased_theta : array-like of shape = number_of_model_elements
            The estimated biased parameters of the model.

        Returns
        -------
        theta : array-like of shape = number_of_model_elements
            The estimated unbiased parameters of the model.

        References
        ----------
        - Manuscript: Sorenson, H. W. (1970). Least-squares estimation:
           from Gauss to Kalman. IEEE spectrum, 7(7), 63-68.
           http://pzs.dstu.dp.ua/DataMining/mls/bibl/Gauss2Kalman.pdf
        - Book (Portuguese): Aguirre, L. A. (2007). Introdução a identificação
           de sistemas: técnicas lineares e não-lineares aplicadas a sistemas
           reais. Editora da UFMG. 3a edição.
        - Manuscript: Markovsky, I., & Van Huffel, S. (2007).
           Overview of total least-squares methods.
           Signal processing, 87(10), 2283-2302.
            https://eprints.soton.ac.uk/263855/1/tls_overview.pdf
        - Wikipedia entry on Least Squares
           https://en.wikipedia.org/wiki/Least_squares

        """
        e = y[max_lag:, 0].reshape(-1, 1) - np.dot(psi, theta)
        im = InformationMatrix(ylag=elag)
        for _ in range(30):
            e = np.concatenate([np.zeros([max_lag, 1]), e], axis=0)

            lagged_data = im.build_output_matrix(None, e)

            e_regressors = self.basis_function.fit(
                lagged_data, max_lag, predefined_regressors=None
            )

            psi_extended = np.concatenate([psi, e_regressors], axis=1)
            unbiased_theta = getattr(self, estimator)(psi_extended, y)
            e = y[max_lag:, 0].reshape(-1, 1) - np.dot(
                psi_extended, unbiased_theta.reshape(-1, 1)
            )

        return unbiased_theta[0 : theta.shape[0], 0].reshape(-1, 1)

    def total_least_squares(self, psi, y):
        """Estimate the model parameters using Total Least Squares method.

        Parameters
        ----------
        psi : ndarray of floats
            The information matrix of the model.
        y : array-like of shape = y_training
            The data used to training the model.

        Returns
        -------
        theta : array-like of shape = number_of_model_elements
            The estimated parameters of the model.

        References
        ----------
        - Manuscript: Golub, G. H., & Van Loan, C. F. (1980).
           An analysis of the total least squares problem.
           SIAM journal on numerical analysis, 17(6), 883-893.
        - Manuscript: Markovsky, I., & Van Huffel, S. (2007).
           Overview of total least-squares methods.
           Signal processing, 87(10), 2283-2302.
           https://eprints.soton.ac.uk/263855/1/tls_overview.pdf
        - Wikipedia entry on Total Least Squares
           https://en.wikipedia.org/wiki/Total_least_squares

        """
        y = y[self.max_lag :, 0].reshape(-1, 1)
        full = np.hstack((psi, y))
        n = psi.shape[1]
        _, _, v = np.linalg.svd(full, full_matrices=True)
        theta = -v.T[:n, n:] / v.T[n:, n:]
        return theta.reshape(-1, 1)

    def _initial_values(self, y, psi):
        y = y[self.max_lag :, 0].reshape(-1, 1)
        n_theta = psi.shape[1]
        n = len(psi)
        theta = np.zeros([n_theta, n])
        xi = np.zeros([n, 1])
        return y, n_theta, n, theta, xi

    def recursive_least_squares(self, psi, y):
        """Estimate the model parameters using the Recursive Least Squares method.

        The implementation consider the forgetting factor.
        Parameters
        ----------
        psi : ndarray of floats
            The information matrix of the model.
        y : array-like of shape = y_training
            The data used to training the model.

        Returns
        -------
        theta : array-like of shape = number_of_model_elements
            The estimated parameters of the model.

        Notes
        -----
        A more in-depth documentation of all methods for parameters estimation
        will be available soon. For now, please refer to the mentioned
        references.

        References
        ----------
        - Book (Portuguese): Aguirre, L. A. (2007). Introdução identificação
           de sistemas: técnicas lineares e não-lineares aplicadas a sistemas
           reais. Editora da UFMG. 3a edição.

        """
        y, n_theta, n, theta, self.xi = self._initial_values(y, psi)

        p = np.eye(n_theta) / self.delta

        for i in range(2, n):
            psi_tmp = psi[i, :].reshape(-1, 1)
            k_numerator = self.lam ** (-1) * p.dot(psi_tmp)
            k_denominator = 1 + self.lam ** (-1) * psi_tmp.T.dot(p).dot(psi_tmp)
            k = np.divide(k_numerator, k_denominator)
            self.xi[i, 0] = y[i, 0] - np.dot(psi_tmp.T, theta[:, i - 1])
            tmp_list = theta[:, i - 1].reshape(-1, 1) + k.dot(self.xi[i, 0])
            theta[:, i] = tmp_list.flatten()

            p1 = p.dot(psi[i, :].reshape(-1, 1)).dot(psi[i, :].reshape(-1, 1).T).dot(p)
            p2 = (
                psi[i, :].reshape(-1, 1).T.dot(p).dot(psi[i, :].reshape(-1, 1))
                + self.lam
            )

            p_numerator = p - np.divide(p1, p2)
            p = np.divide(p_numerator, self.lam)

        self.theta_evolution = theta.copy()
        return theta[:, -1].reshape(-1, 1)

    def affine_least_mean_squares(self, psi, y):
        """Estimate the model parameters using the Affine Least Mean Squares.

        Parameters
        ----------
        psi : ndarray of floats
            The information matrix of the model.
        y : array-like of shape = y_training
            The data used to training the model.

        Returns
        -------
        theta : array-like of shape = number_of_model_elements
            The estimated parameters of the model.

        Notes
        -----
        A more in-depth documentation of all methods for parameters estimation
        will be available soon. For now, please refer to the mentioned
        references.

        References
        ----------
        - Book: Poularikas, A. D. (2017). Adaptive filtering: Fundamentals
           of least mean squares with MATLAB®. CRC Press.

        """
        y, n_theta, n, theta, self.xi = self._initial_values(y, psi)

        for i in range(n_theta, n):
            self.xi = y - psi.dot(theta[:, i - 1].reshape(-1, 1))
            aux = (
                self.mu
                * psi
                @ np.linalg.pinv(psi.T @ psi + self.offset_covariance * np.eye(n_theta))
            )
            tmp_list = theta[:, i - 1].reshape(-1, 1) + aux.T.dot(self.xi)
            theta[:, i] = tmp_list.flatten()

        return theta[:, -1].reshape(-1, 1)

    def least_mean_squares(self, psi, y):
        """Estimate the model parameters using the Least Mean Squares filter.

        Parameters
        ----------
        psi : ndarray of floats
            The information matrix of the model.
        y : array-like of shape = y_training
            The data used to training the model.

        Returns
        -------
        theta : array-like of shape = number_of_model_elements
            The estimated parameters of the model.

        Notes
        -----
        A more in-depth documentation of all methods for parameters estimation
        will be available soon. For now, please refer to the mentioned
        references.

        References
        ----------
        - Book: Haykin, S., & Widrow, B. (Eds.). (2003). Least-mean-square
           adaptive filters (Vol. 31). John Wiley & Sons.
        - Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação,
           análise estatística e novas estratégias de algoritmos LMS de passo
           variável.
        - Wikipedia entry on Least Mean Squares
           https://en.wikipedia.org/wiki/Least_mean_squares_filter

        """
        y, n_theta, n, theta, self.xi = self._initial_values(y, psi)

        for i in range(n_theta, n):
            psi_tmp = psi[i, :].reshape(-1, 1)
            self.xi[i, 0] = y[i, 0] - np.dot(psi_tmp.T, theta[:, i - 1])
            tmp_list = (
                theta[:, i - 1].reshape(-1, 1) + 2 * self.mu * self.xi[i, 0] * psi_tmp
            )
            theta[:, i] = tmp_list.flatten()

        return theta[:, -1].reshape(-1, 1)

    def least_mean_squares_sign_error(self, psi, y):
        """Parameter estimation using the Sign-Error  Least Mean Squares filter.

        The sign-error LMS algorithm uses the sign of the error vector
        to change the filter coefficients.

        Parameters
        ----------
        psi : ndarray of floats
            The information matrix of the model.
        y : array-like of shape = y_training
            The data used to training the model.

        Returns
        -------
        theta : array-like of shape = number_of_model_elements
            The estimated parameters of the model.

        Notes
        -----
        A more in-depth documentation of all methods for parameters estimation
        will be available soon. For now, please refer to the mentioned
        references.

        References
        ----------
        - Book: Hayes, M. H. (2009). Statistical digital signal processing
           and modeling. John Wiley & Sons.
        - Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação,
           análise estatística e novas estratégias de algoritmos LMS de passo
           variável.
        - Wikipedia entry on Least Mean Squares
           https://en.wikipedia.org/wiki/Least_mean_squares_filter

        """
        y, n_theta, n, theta, self.xi = self._initial_values(y, psi)

        for i in range(n_theta, n):
            psi_tmp = psi[i, :].reshape(-1, 1)
            self.xi[i, 0] = y[i, 0] - np.dot(psi_tmp.T, theta[:, i - 1])
            tmp_list = (
                theta[:, i - 1].reshape(-1, 1)
                + self.mu * np.sign(self.xi[i, 0]) * psi_tmp
            )
            theta[:, i] = tmp_list.flatten()

        return theta[:, -1].reshape(-1, 1)

    def normalized_least_mean_squares(self, psi, y):
        """Parameter estimation using the Normalized Least Mean Squares filter.

        The normalization is used to avoid numerical instability when updating
        the estimated parameters.

        Parameters
        ----------
        psi : ndarray of floats
            The information matrix of the model.
        y : array-like of shape = y_training
            The data used to training the model.

        Returns
        -------
        theta : array-like of shape = number_of_model_elements
            The estimated parameters of the model.

        Notes
        -----
        A more in-depth documentation of all methods for parameters estimation
        will be available soon. For now, please refer to the mentioned
        references.

        References
        ----------
        - Book: Hayes, M. H. (2009). Statistical digital signal processing
           and modeling. John Wiley & Sons.
        - Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação,
           análise estatística e novas estratégias de algoritmos LMS de passo
           variável.
        - Wikipedia entry on Least Mean Squares
           https://en.wikipedia.org/wiki/Least_mean_squares_filter

        """
        y, n_theta, n, theta, self.xi = self._initial_values(y, psi)

        for i in range(n_theta, n):
            psi_tmp = psi[i, :].reshape(-1, 1)
            self.xi[i, 0] = y[i, 0] - np.dot(psi_tmp.T, theta[:, i - 1])
            tmp_list = theta[:, i - 1].reshape(-1, 1) + 2 * self.mu * self.xi[i, 0] * (
                psi_tmp / (self.eps + np.dot(psi_tmp.T, psi_tmp))
            )
            theta[:, i] = tmp_list.flatten()

        return theta[:, -1].reshape(-1, 1)

    def least_mean_squares_normalized_sign_error(self, psi, y):
        """Parameter estimation using the Normalized Sign-Error LMS filter.

        The normalization is used to avoid numerical instability when updating
        the estimated parameters and the sign of the error vector is used to
        to change the filter coefficients.

        Parameters
        ----------
        psi : ndarray of floats
            The information matrix of the model.
        y : array-like of shape = y_training
            The data used to training the model.

        Returns
        -------
        theta : array-like of shape = number_of_model_elements
            The estimated parameters of the model.

        Notes
        -----
        A more in-depth documentation of all methods for parameters estimation
        will be available soon. For now, please refer to the mentioned
        references.

        References
        ----------
        - Book: Hayes, M. H. (2009). Statistical digital signal processing
           and modeling. John Wiley & Sons.
        - Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação,
           análise estatística e novas estratégias de algoritmos LMS de passo
           variável.
        - Wikipedia entry on Least Mean Squares
           https://en.wikipedia.org/wiki/Least_mean_squares_filter

        """
        y, n_theta, n, theta, self.xi = self._initial_values(y, psi)

        for i in range(n_theta, n):
            psi_tmp = psi[i, :].reshape(-1, 1)
            self.xi[i, 0] = y[i, 0] - np.dot(psi_tmp.T, theta[:, i - 1])
            tmp_list = theta[:, i - 1].reshape(-1, 1) + 2 * self.mu * np.sign(
                self.xi[i, 0]
            ) * (psi_tmp / (self.eps + np.dot(psi_tmp.T, psi_tmp)))
            theta[:, i] = tmp_list.flatten()

        return theta[:, -1].reshape(-1, 1)

    def least_mean_squares_sign_regressor(self, psi, y):
        """Parameter estimation using the  Sign-Regressor LMS filter.

        The sign-regressor LMS algorithm uses the sign of the matrix
        information to change the filter coefficients.

        Parameters
        ----------
        psi : ndarray of floats
            The information matrix of the model.
        y : array-like of shape = y_training
            The data used to training the model.

        Returns
        -------
        theta : array-like of shape = number_of_model_elements
            The estimated parameters of the model.

        Notes
        -----
        A more in-depth documentation of all methods for parameters estimation
        will be available soon. For now, please refer to the mentioned
        references.

        References
        ----------
        - Book: Hayes, M. H. (2009). Statistical digital signal processing
           and modeling. John Wiley & Sons.
        - Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação,
           análise estatística e novas estratégias de algoritmos LMS de passo
           variável.
        - Wikipedia entry on Least Mean Squares
           https://en.wikipedia.org/wiki/Least_mean_squares_filter

        """
        y, n_theta, n, theta, self.xi = self._initial_values(y, psi)

        for i in range(n_theta, n):
            psi_tmp = psi[i, :].reshape(-1, 1)
            self.xi[i, 0] = y[i, 0] - np.dot(psi_tmp.T, theta[:, i - 1])
            tmp_list = theta[:, i - 1].reshape(-1, 1) + self.mu * self.xi[
                i, 0
            ] * np.sign(psi_tmp)
            theta[:, i] = tmp_list.flatten()

        return theta[:, -1].reshape(-1, 1)

    def least_mean_squares_normalized_sign_regressor(self, psi, y):
        """Parameter estimation using the Normalized Sign-Regressor LMS filter.

        The normalization is used to avoid numerical instability when updating
        the estimated parameters and the sign of the information matrix is
        used to change the filter coefficients.

        Parameters
        ----------
        psi : ndarray of floats
            The information matrix of the model.
        y : array-like of shape = y_training
            The data used to training the model.

        Returns
        -------
        theta : array-like of shape = number_of_model_elements
            The estimated parameters of the model.

        Notes
        -----
        A more in-depth documentation of all methods for parameters estimation
        will be available soon. For now, please refer to the mentioned
        references.

        References
        ----------
        .. [1] Book: Hayes, M. H. (2009). Statistical digital signal processing
           and modeling. John Wiley & Sons.
        .. [2] Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação,
           análise estatística e novas estratégias de algoritmos LMS de passo
           variável.
        .. [3] Wikipedia entry on Least Mean Squares
           https://en.wikipedia.org/wiki/Least_mean_squares_filter

        """
        y, n_theta, n, theta, self.xi = self._initial_values(y, psi)

        for i in range(n_theta, n):
            psi_tmp = psi[i, :].reshape(-1, 1)
            self.xi[i, 0] = y[i, 0] - np.dot(psi_tmp.T, theta[:, i - 1])
            tmp_list = theta[:, i - 1].reshape(-1, 1) + self.mu * self.xi[i, 0] * (
                np.sign(psi_tmp) / (self.eps + np.dot(psi_tmp.T, psi_tmp))
            )
            theta[:, i] = tmp_list.flatten()

        return theta[:, -1].reshape(-1, 1)

    def least_mean_squares_sign_sign(self, psi, y):
        """Parameter estimation using the  Sign-Sign LMS filter.

        The sign-regressor LMS algorithm uses both the sign of the matrix
        information and the sign of the error vector to change the filter
        coefficients.

        Parameters
        ----------
        psi : ndarray of floats
            The information matrix of the model.
        y : array-like of shape = y_training
            The data used to training the model.

        Returns
        -------
        theta : array-like of shape = number_of_model_elements
            The estimated parameters of the model.

        Notes
        -----
        A more in-depth documentation of all methods for parameters estimation
        will be available soon. For now, please refer to the mentioned
        references.

        References
        ----------
        - Book: Hayes, M. H. (2009). Statistical digital signal processing
           and modeling. John Wiley & Sons.
        - Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação,
           análise estatística e novas estratégias de algoritmos LMS de passo
           variável.
        - Wikipedia entry on Least Mean Squares
           https://en.wikipedia.org/wiki/Least_mean_squares_filter

        """
        y, n_theta, n, theta, self.xi = self._initial_values(y, psi)

        for i in range(n_theta, n):
            psi_tmp = psi[i, :].reshape(-1, 1)
            self.xi[i, 0] = y[i, 0] - np.dot(psi_tmp.T, theta[:, i - 1])
            tmp_list = theta[:, i - 1].reshape(-1, 1) + 2 * self.mu * np.sign(
                self.xi[i, 0]
            ) * np.sign(psi_tmp)
            theta[:, i] = tmp_list.flatten()

        return theta[:, -1].reshape(-1, 1)

    def least_mean_squares_normalized_sign_sign(self, psi, y):
        """Parameter estimation using the Normalized Sign-Sign LMS filter.

        The normalization is used to avoid numerical instability when updating
        the estimated parameters and both the sign of the information matrix
        and the sign of the error vector are used to change the filter
        coefficients.

        Parameters
        ----------
        psi : ndarray of floats
            The information matrix of the model.
        y : array-like of shape = y_training
            The data used to training the model.

        Returns
        -------
        theta : array-like of shape = number_of_model_elements
            The estimated parameters of the model.

        Notes
        -----
        A more in-depth documentation of all methods for parameters estimation
        will be available soon. For now, please refer to the mentioned
        references.

        References
        ----------
        - Book: Hayes, M. H. (2009). Statistical digital signal processing
           and modeling. John Wiley & Sons.
        - Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação,
           análise estatística e novas estratégias de algoritmos LMS de passo
           variável.
        - Wikipedia entry on Least Mean Squares
           https://en.wikipedia.org/wiki/Least_mean_squares_filter

        """
        y, n_theta, n, theta, self.xi = self._initial_values(y, psi)

        for i in range(n_theta, n):
            psi_tmp = psi[i, :].reshape(-1, 1)
            self.xi[i, 0] = y[i, 0] - np.dot(psi_tmp.T, theta[:, i - 1])
            tmp_list = theta[:, i - 1].reshape(-1, 1) + 2 * self.mu * np.sign(
                self.xi[i, 0]
            ) * (np.sign(psi_tmp) / (self.eps + np.dot(psi_tmp.T, psi_tmp)))
            theta[:, i] = tmp_list.flatten()

        return theta[:, -1].reshape(-1, 1)

    def least_mean_squares_normalized_leaky(self, psi, y):
        """Parameter estimation using the  Normalized Leaky LMS filter.

        When the leakage factor, gama, is set to 0 then there is no
        leakage in the estimation process.

        Parameters
        ----------
        psi : ndarray of floats
            The information matrix of the model.
        y : array-like of shape = y_training
            The data used to training the model.

        Returns
        -------
        theta : array-like of shape = number_of_model_elements
            The estimated parameters of the model.

        Notes
        -----
        A more in-depth documentation of all methods for parameters estimation
        will be available soon. For now, please refer to the mentioned
        references.

        References
        ----------
        - Book: Hayes, M. H. (2009). Statistical digital signal processing
           and modeling. John Wiley & Sons.
        - Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação,
           análise estatística e novas estratégias de algoritmos LMS de passo
           variável.
        - Wikipedia entry on Least Mean Squares
           https://en.wikipedia.org/wiki/Least_mean_squares_filter

        """
        y, n_theta, n, theta, self.xi = self._initial_values(y, psi)

        for i in range(n_theta, n):
            psi_tmp = psi[i, :].reshape(-1, 1)
            self.xi[i, 0] = y[i, 0] - np.dot(psi_tmp.T, theta[:, i - 1])
            tmp_list = theta[:, i - 1].reshape(-1, 1) * (
                1 - self.mu * self.gama
            ) + self.mu * self.xi[i, 0] * psi_tmp / (
                self.eps + np.dot(psi_tmp.T, psi_tmp)
            )
            theta[:, i] = tmp_list.flatten()

        return theta[:, -1].reshape(-1, 1)

    def least_mean_squares_leaky(self, psi, y):
        """Parameter estimation using the  Leaky LMS filter.

        When the leakage factor, gama, is set to 0 then there is no
        leakage in the estimation process.

        Parameters
        ----------
        psi : ndarray of floats
            The information matrix of the model.
        y : array-like of shape = y_training
            The data used to training the model.

        Returns
        -------
        theta : array-like of shape = number_of_model_elements
            The estimated parameters of the model.

        Notes
        -----
        A more in-depth documentation of all methods for parameters estimation
        will be available soon. For now, please refer to the mentioned
        references.

        References
        ----------
        - Book: Hayes, M. H. (2009). Statistical digital signal processing
           and modeling. John Wiley & Sons.
        - Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação,
           análise estatística e novas estratégias de algoritmos LMS de passo
           variável.
        - Wikipedia entry on Least Mean Squares
           https://en.wikipedia.org/wiki/Least_mean_squares_filter

        """
        y, n_theta, n, theta, self.xi = self._initial_values(y, psi)

        for i in range(n_theta, n):
            psi_tmp = psi[i, :].reshape(-1, 1)
            self.xi[i, 0] = y[i, 0] - np.dot(psi_tmp.T, theta[:, i - 1])
            tmp_list = (
                theta[:, i - 1].reshape(-1, 1) * (1 - self.mu * self.gama)
                + self.mu * self.xi[i, 0] * psi_tmp
            )
            theta[:, i] = tmp_list.flatten()

        return theta[:, -1].reshape(-1, 1)

    def least_mean_squares_fourth(self, psi, y):
        """Parameter estimation using the  LMS Fourth filter.

        When the leakage factor, gama, is set to 0 then there is no
        leakage in the estimation process.

        Parameters
        ----------
        psi : ndarray of floats
            The information matrix of the model.
        y : array-like of shape = y_training
            The data used to training the model.

        Returns
        -------
        theta : array-like of shape = number_of_model_elements
            The estimated parameters of the model.

        Notes
        -----
        A more in-depth documentation of all methods for parameters estimation
        will be available soon. For now, please refer to the mentioned
        references.

        References
        ----------
        - Book: Hayes, M. H. (2009). Statistical digital signal processing
           and modeling. John Wiley & Sons.
        - Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação,
           análise estatística e novas estratégias de algoritmos LMS de passo
           variável.
        - Manuscript:Gui, G., Mehbodniya, A., & Adachi, F. (2013).
           Least mean square/fourth algorithm with application to sparse
           channel estimation. arXiv preprint arXiv:1304.3911.
           https://arxiv.org/pdf/1304.3911.pdf
        - Manuscript: Nascimento, V. H., & Bermudez, J. C. M. (2005, March).
           When is the least-mean fourth algorithm mean-square stable?
           In Proceedings.(ICASSP'05). IEEE International Conference on
           Acoustics, Speech, and Signal Processing, 2005.
           (Vol. 4, pp. iv-341). IEEE.
           http://www.lps.usp.br/vitor/artigos/icassp05.pdf
        - Wikipedia entry on Least Mean Squares
           https://en.wikipedia.org/wiki/Least_mean_squares_filter

        """
        y, n_theta, n, theta, self.xi = self._initial_values(y, psi)

        for i in range(n_theta, n):
            psi_tmp = psi[i, :].reshape(-1, 1)
            self.xi[i, 0] = y[i, 0] - np.dot(psi_tmp.T, theta[:, i - 1])
            tmp_list = (
                theta[:, i - 1].reshape(-1, 1) + self.mu * psi_tmp * self.xi[i, 0] ** 3
            )
            theta[:, i] = tmp_list.flatten()

        return theta[:, -1].reshape(-1, 1)

    def least_mean_squares_mixed_norm(self, psi, y):
        """Parameter estimation using the Mixed-norm LMS filter.

        The weight factor controls the proportions of the error norms
        and offers an extra degree of freedom within the adaptation.

        Parameters
        ----------
        psi : ndarray of floats
            The information matrix of the model.
        y : array-like of shape = y_training
            The data used to training the model.

        Returns
        -------
        theta : array-like of shape = number_of_model_elements
            The estimated parameters of the model.

        Notes
        -----
        A more in-depth documentation of all methods for parameters estimation
        will be available soon. For now, please refer to the mentioned
        references.

        References
        ----------
        - Chambers, J. A., Tanrikulu, O., & Constantinides, A. G. (1994).
           Least mean mixed-norm adaptive filtering.
           Electronics letters, 30(19), 1574-1575.
           https://ieeexplore.ieee.org/document/326382
        - Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação,
           análise estatística e novas estratégias de algoritmos LMS de passo
           variável.
        - Wikipedia entry on Least Mean Squares
           https://en.wikipedia.org/wiki/Least_mean_squares_filter

        """
        y, n_theta, n, theta, self.xi = self._initial_values(y, psi)

        for i in range(n_theta, n):
            psi_tmp = psi[i, :].reshape(-1, 1)
            self.xi[i, 0] = y[i, 0] - np.dot(psi_tmp.T, theta[:, i - 1])
            tmp_list = theta[:, i - 1].reshape(-1, 1) + self.mu * psi_tmp * self.xi[
                i, 0
            ] * (self.weight + (1 - self.weight) * self.xi[i, 0] ** 2)
            theta[:, i] = tmp_list.flatten()

        return theta[:, -1].reshape(-1, 1)

affine_least_mean_squares(psi, y)

Estimate the model parameters using the Affine Least Mean Squares.

Parameters:

Name Type Description Default
psi ndarray of floats

The information matrix of the model.

required
y array-like of shape = y_training

The data used to training the model.

required

Returns:

Name Type Description
theta array-like of shape = number_of_model_elements

The estimated parameters of the model.

Notes

A more in-depth documentation of all methods for parameters estimation will be available soon. For now, please refer to the mentioned references.

References
  • Book: Poularikas, A. D. (2017). Adaptive filtering: Fundamentals of least mean squares with MATLAB®. CRC Press.
Source code in sysidentpy\parameter_estimation\estimators.py
def affine_least_mean_squares(self, psi, y):
    """Estimate the model parameters using the Affine Least Mean Squares.

    Parameters
    ----------
    psi : ndarray of floats
        The information matrix of the model.
    y : array-like of shape = y_training
        The data used to training the model.

    Returns
    -------
    theta : array-like of shape = number_of_model_elements
        The estimated parameters of the model.

    Notes
    -----
    A more in-depth documentation of all methods for parameters estimation
    will be available soon. For now, please refer to the mentioned
    references.

    References
    ----------
    - Book: Poularikas, A. D. (2017). Adaptive filtering: Fundamentals
       of least mean squares with MATLAB®. CRC Press.

    """
    y, n_theta, n, theta, self.xi = self._initial_values(y, psi)

    for i in range(n_theta, n):
        self.xi = y - psi.dot(theta[:, i - 1].reshape(-1, 1))
        aux = (
            self.mu
            * psi
            @ np.linalg.pinv(psi.T @ psi + self.offset_covariance * np.eye(n_theta))
        )
        tmp_list = theta[:, i - 1].reshape(-1, 1) + aux.T.dot(self.xi)
        theta[:, i] = tmp_list.flatten()

    return theta[:, -1].reshape(-1, 1)

least_mean_squares(psi, y)

Estimate the model parameters using the Least Mean Squares filter.

Parameters:

Name Type Description Default
psi ndarray of floats

The information matrix of the model.

required
y array-like of shape = y_training

The data used to training the model.

required

Returns:

Name Type Description
theta array-like of shape = number_of_model_elements

The estimated parameters of the model.

Notes

A more in-depth documentation of all methods for parameters estimation will be available soon. For now, please refer to the mentioned references.

References
  • Book: Haykin, S., & Widrow, B. (Eds.). (2003). Least-mean-square adaptive filters (Vol. 31). John Wiley & Sons.
  • Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação, análise estatística e novas estratégias de algoritmos LMS de passo variável.
  • Wikipedia entry on Least Mean Squares https://en.wikipedia.org/wiki/Least_mean_squares_filter
Source code in sysidentpy\parameter_estimation\estimators.py
def least_mean_squares(self, psi, y):
    """Estimate the model parameters using the Least Mean Squares filter.

    Parameters
    ----------
    psi : ndarray of floats
        The information matrix of the model.
    y : array-like of shape = y_training
        The data used to training the model.

    Returns
    -------
    theta : array-like of shape = number_of_model_elements
        The estimated parameters of the model.

    Notes
    -----
    A more in-depth documentation of all methods for parameters estimation
    will be available soon. For now, please refer to the mentioned
    references.

    References
    ----------
    - Book: Haykin, S., & Widrow, B. (Eds.). (2003). Least-mean-square
       adaptive filters (Vol. 31). John Wiley & Sons.
    - Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação,
       análise estatística e novas estratégias de algoritmos LMS de passo
       variável.
    - Wikipedia entry on Least Mean Squares
       https://en.wikipedia.org/wiki/Least_mean_squares_filter

    """
    y, n_theta, n, theta, self.xi = self._initial_values(y, psi)

    for i in range(n_theta, n):
        psi_tmp = psi[i, :].reshape(-1, 1)
        self.xi[i, 0] = y[i, 0] - np.dot(psi_tmp.T, theta[:, i - 1])
        tmp_list = (
            theta[:, i - 1].reshape(-1, 1) + 2 * self.mu * self.xi[i, 0] * psi_tmp
        )
        theta[:, i] = tmp_list.flatten()

    return theta[:, -1].reshape(-1, 1)

least_mean_squares_fourth(psi, y)

Parameter estimation using the LMS Fourth filter.

When the leakage factor, gama, is set to 0 then there is no leakage in the estimation process.

Parameters:

Name Type Description Default
psi ndarray of floats

The information matrix of the model.

required
y array-like of shape = y_training

The data used to training the model.

required

Returns:

Name Type Description
theta array-like of shape = number_of_model_elements

The estimated parameters of the model.

Notes

A more in-depth documentation of all methods for parameters estimation will be available soon. For now, please refer to the mentioned references.

References
  • Book: Hayes, M. H. (2009). Statistical digital signal processing and modeling. John Wiley & Sons.
  • Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação, análise estatística e novas estratégias de algoritmos LMS de passo variável.
  • Manuscript:Gui, G., Mehbodniya, A., & Adachi, F. (2013). Least mean square/fourth algorithm with application to sparse channel estimation. arXiv preprint arXiv:1304.3911. https://arxiv.org/pdf/1304.3911.pdf
  • Manuscript: Nascimento, V. H., & Bermudez, J. C. M. (2005, March). When is the least-mean fourth algorithm mean-square stable? In Proceedings.(ICASSP'05). IEEE International Conference on Acoustics, Speech, and Signal Processing, 2005. (Vol. 4, pp. iv-341). IEEE. http://www.lps.usp.br/vitor/artigos/icassp05.pdf
  • Wikipedia entry on Least Mean Squares https://en.wikipedia.org/wiki/Least_mean_squares_filter
Source code in sysidentpy\parameter_estimation\estimators.py
def least_mean_squares_fourth(self, psi, y):
    """Parameter estimation using the  LMS Fourth filter.

    When the leakage factor, gama, is set to 0 then there is no
    leakage in the estimation process.

    Parameters
    ----------
    psi : ndarray of floats
        The information matrix of the model.
    y : array-like of shape = y_training
        The data used to training the model.

    Returns
    -------
    theta : array-like of shape = number_of_model_elements
        The estimated parameters of the model.

    Notes
    -----
    A more in-depth documentation of all methods for parameters estimation
    will be available soon. For now, please refer to the mentioned
    references.

    References
    ----------
    - Book: Hayes, M. H. (2009). Statistical digital signal processing
       and modeling. John Wiley & Sons.
    - Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação,
       análise estatística e novas estratégias de algoritmos LMS de passo
       variável.
    - Manuscript:Gui, G., Mehbodniya, A., & Adachi, F. (2013).
       Least mean square/fourth algorithm with application to sparse
       channel estimation. arXiv preprint arXiv:1304.3911.
       https://arxiv.org/pdf/1304.3911.pdf
    - Manuscript: Nascimento, V. H., & Bermudez, J. C. M. (2005, March).
       When is the least-mean fourth algorithm mean-square stable?
       In Proceedings.(ICASSP'05). IEEE International Conference on
       Acoustics, Speech, and Signal Processing, 2005.
       (Vol. 4, pp. iv-341). IEEE.
       http://www.lps.usp.br/vitor/artigos/icassp05.pdf
    - Wikipedia entry on Least Mean Squares
       https://en.wikipedia.org/wiki/Least_mean_squares_filter

    """
    y, n_theta, n, theta, self.xi = self._initial_values(y, psi)

    for i in range(n_theta, n):
        psi_tmp = psi[i, :].reshape(-1, 1)
        self.xi[i, 0] = y[i, 0] - np.dot(psi_tmp.T, theta[:, i - 1])
        tmp_list = (
            theta[:, i - 1].reshape(-1, 1) + self.mu * psi_tmp * self.xi[i, 0] ** 3
        )
        theta[:, i] = tmp_list.flatten()

    return theta[:, -1].reshape(-1, 1)

least_mean_squares_leaky(psi, y)

Parameter estimation using the Leaky LMS filter.

When the leakage factor, gama, is set to 0 then there is no leakage in the estimation process.

Parameters:

Name Type Description Default
psi ndarray of floats

The information matrix of the model.

required
y array-like of shape = y_training

The data used to training the model.

required

Returns:

Name Type Description
theta array-like of shape = number_of_model_elements

The estimated parameters of the model.

Notes

A more in-depth documentation of all methods for parameters estimation will be available soon. For now, please refer to the mentioned references.

References
  • Book: Hayes, M. H. (2009). Statistical digital signal processing and modeling. John Wiley & Sons.
  • Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação, análise estatística e novas estratégias de algoritmos LMS de passo variável.
  • Wikipedia entry on Least Mean Squares https://en.wikipedia.org/wiki/Least_mean_squares_filter
Source code in sysidentpy\parameter_estimation\estimators.py
def least_mean_squares_leaky(self, psi, y):
    """Parameter estimation using the  Leaky LMS filter.

    When the leakage factor, gama, is set to 0 then there is no
    leakage in the estimation process.

    Parameters
    ----------
    psi : ndarray of floats
        The information matrix of the model.
    y : array-like of shape = y_training
        The data used to training the model.

    Returns
    -------
    theta : array-like of shape = number_of_model_elements
        The estimated parameters of the model.

    Notes
    -----
    A more in-depth documentation of all methods for parameters estimation
    will be available soon. For now, please refer to the mentioned
    references.

    References
    ----------
    - Book: Hayes, M. H. (2009). Statistical digital signal processing
       and modeling. John Wiley & Sons.
    - Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação,
       análise estatística e novas estratégias de algoritmos LMS de passo
       variável.
    - Wikipedia entry on Least Mean Squares
       https://en.wikipedia.org/wiki/Least_mean_squares_filter

    """
    y, n_theta, n, theta, self.xi = self._initial_values(y, psi)

    for i in range(n_theta, n):
        psi_tmp = psi[i, :].reshape(-1, 1)
        self.xi[i, 0] = y[i, 0] - np.dot(psi_tmp.T, theta[:, i - 1])
        tmp_list = (
            theta[:, i - 1].reshape(-1, 1) * (1 - self.mu * self.gama)
            + self.mu * self.xi[i, 0] * psi_tmp
        )
        theta[:, i] = tmp_list.flatten()

    return theta[:, -1].reshape(-1, 1)

least_mean_squares_mixed_norm(psi, y)

Parameter estimation using the Mixed-norm LMS filter.

The weight factor controls the proportions of the error norms and offers an extra degree of freedom within the adaptation.

Parameters:

Name Type Description Default
psi ndarray of floats

The information matrix of the model.

required
y array-like of shape = y_training

The data used to training the model.

required

Returns:

Name Type Description
theta array-like of shape = number_of_model_elements

The estimated parameters of the model.

Notes

A more in-depth documentation of all methods for parameters estimation will be available soon. For now, please refer to the mentioned references.

References
Source code in sysidentpy\parameter_estimation\estimators.py
def least_mean_squares_mixed_norm(self, psi, y):
    """Parameter estimation using the Mixed-norm LMS filter.

    The weight factor controls the proportions of the error norms
    and offers an extra degree of freedom within the adaptation.

    Parameters
    ----------
    psi : ndarray of floats
        The information matrix of the model.
    y : array-like of shape = y_training
        The data used to training the model.

    Returns
    -------
    theta : array-like of shape = number_of_model_elements
        The estimated parameters of the model.

    Notes
    -----
    A more in-depth documentation of all methods for parameters estimation
    will be available soon. For now, please refer to the mentioned
    references.

    References
    ----------
    - Chambers, J. A., Tanrikulu, O., & Constantinides, A. G. (1994).
       Least mean mixed-norm adaptive filtering.
       Electronics letters, 30(19), 1574-1575.
       https://ieeexplore.ieee.org/document/326382
    - Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação,
       análise estatística e novas estratégias de algoritmos LMS de passo
       variável.
    - Wikipedia entry on Least Mean Squares
       https://en.wikipedia.org/wiki/Least_mean_squares_filter

    """
    y, n_theta, n, theta, self.xi = self._initial_values(y, psi)

    for i in range(n_theta, n):
        psi_tmp = psi[i, :].reshape(-1, 1)
        self.xi[i, 0] = y[i, 0] - np.dot(psi_tmp.T, theta[:, i - 1])
        tmp_list = theta[:, i - 1].reshape(-1, 1) + self.mu * psi_tmp * self.xi[
            i, 0
        ] * (self.weight + (1 - self.weight) * self.xi[i, 0] ** 2)
        theta[:, i] = tmp_list.flatten()

    return theta[:, -1].reshape(-1, 1)

least_mean_squares_normalized_leaky(psi, y)

Parameter estimation using the Normalized Leaky LMS filter.

When the leakage factor, gama, is set to 0 then there is no leakage in the estimation process.

Parameters:

Name Type Description Default
psi ndarray of floats

The information matrix of the model.

required
y array-like of shape = y_training

The data used to training the model.

required

Returns:

Name Type Description
theta array-like of shape = number_of_model_elements

The estimated parameters of the model.

Notes

A more in-depth documentation of all methods for parameters estimation will be available soon. For now, please refer to the mentioned references.

References
  • Book: Hayes, M. H. (2009). Statistical digital signal processing and modeling. John Wiley & Sons.
  • Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação, análise estatística e novas estratégias de algoritmos LMS de passo variável.
  • Wikipedia entry on Least Mean Squares https://en.wikipedia.org/wiki/Least_mean_squares_filter
Source code in sysidentpy\parameter_estimation\estimators.py
def least_mean_squares_normalized_leaky(self, psi, y):
    """Parameter estimation using the  Normalized Leaky LMS filter.

    When the leakage factor, gama, is set to 0 then there is no
    leakage in the estimation process.

    Parameters
    ----------
    psi : ndarray of floats
        The information matrix of the model.
    y : array-like of shape = y_training
        The data used to training the model.

    Returns
    -------
    theta : array-like of shape = number_of_model_elements
        The estimated parameters of the model.

    Notes
    -----
    A more in-depth documentation of all methods for parameters estimation
    will be available soon. For now, please refer to the mentioned
    references.

    References
    ----------
    - Book: Hayes, M. H. (2009). Statistical digital signal processing
       and modeling. John Wiley & Sons.
    - Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação,
       análise estatística e novas estratégias de algoritmos LMS de passo
       variável.
    - Wikipedia entry on Least Mean Squares
       https://en.wikipedia.org/wiki/Least_mean_squares_filter

    """
    y, n_theta, n, theta, self.xi = self._initial_values(y, psi)

    for i in range(n_theta, n):
        psi_tmp = psi[i, :].reshape(-1, 1)
        self.xi[i, 0] = y[i, 0] - np.dot(psi_tmp.T, theta[:, i - 1])
        tmp_list = theta[:, i - 1].reshape(-1, 1) * (
            1 - self.mu * self.gama
        ) + self.mu * self.xi[i, 0] * psi_tmp / (
            self.eps + np.dot(psi_tmp.T, psi_tmp)
        )
        theta[:, i] = tmp_list.flatten()

    return theta[:, -1].reshape(-1, 1)

least_mean_squares_normalized_sign_error(psi, y)

Parameter estimation using the Normalized Sign-Error LMS filter.

The normalization is used to avoid numerical instability when updating the estimated parameters and the sign of the error vector is used to to change the filter coefficients.

Parameters:

Name Type Description Default
psi ndarray of floats

The information matrix of the model.

required
y array-like of shape = y_training

The data used to training the model.

required

Returns:

Name Type Description
theta array-like of shape = number_of_model_elements

The estimated parameters of the model.

Notes

A more in-depth documentation of all methods for parameters estimation will be available soon. For now, please refer to the mentioned references.

References
  • Book: Hayes, M. H. (2009). Statistical digital signal processing and modeling. John Wiley & Sons.
  • Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação, análise estatística e novas estratégias de algoritmos LMS de passo variável.
  • Wikipedia entry on Least Mean Squares https://en.wikipedia.org/wiki/Least_mean_squares_filter
Source code in sysidentpy\parameter_estimation\estimators.py
def least_mean_squares_normalized_sign_error(self, psi, y):
    """Parameter estimation using the Normalized Sign-Error LMS filter.

    The normalization is used to avoid numerical instability when updating
    the estimated parameters and the sign of the error vector is used to
    to change the filter coefficients.

    Parameters
    ----------
    psi : ndarray of floats
        The information matrix of the model.
    y : array-like of shape = y_training
        The data used to training the model.

    Returns
    -------
    theta : array-like of shape = number_of_model_elements
        The estimated parameters of the model.

    Notes
    -----
    A more in-depth documentation of all methods for parameters estimation
    will be available soon. For now, please refer to the mentioned
    references.

    References
    ----------
    - Book: Hayes, M. H. (2009). Statistical digital signal processing
       and modeling. John Wiley & Sons.
    - Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação,
       análise estatística e novas estratégias de algoritmos LMS de passo
       variável.
    - Wikipedia entry on Least Mean Squares
       https://en.wikipedia.org/wiki/Least_mean_squares_filter

    """
    y, n_theta, n, theta, self.xi = self._initial_values(y, psi)

    for i in range(n_theta, n):
        psi_tmp = psi[i, :].reshape(-1, 1)
        self.xi[i, 0] = y[i, 0] - np.dot(psi_tmp.T, theta[:, i - 1])
        tmp_list = theta[:, i - 1].reshape(-1, 1) + 2 * self.mu * np.sign(
            self.xi[i, 0]
        ) * (psi_tmp / (self.eps + np.dot(psi_tmp.T, psi_tmp)))
        theta[:, i] = tmp_list.flatten()

    return theta[:, -1].reshape(-1, 1)

least_mean_squares_normalized_sign_regressor(psi, y)

Parameter estimation using the Normalized Sign-Regressor LMS filter.

The normalization is used to avoid numerical instability when updating the estimated parameters and the sign of the information matrix is used to change the filter coefficients.

Parameters:

Name Type Description Default
psi ndarray of floats

The information matrix of the model.

required
y array-like of shape = y_training

The data used to training the model.

required

Returns:

Name Type Description
theta array-like of shape = number_of_model_elements

The estimated parameters of the model.

Notes

A more in-depth documentation of all methods for parameters estimation will be available soon. For now, please refer to the mentioned references.

References

.. [1] Book: Hayes, M. H. (2009). Statistical digital signal processing and modeling. John Wiley & Sons. .. [2] Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação, análise estatística e novas estratégias de algoritmos LMS de passo variável. .. [3] Wikipedia entry on Least Mean Squares https://en.wikipedia.org/wiki/Least_mean_squares_filter

Source code in sysidentpy\parameter_estimation\estimators.py
def least_mean_squares_normalized_sign_regressor(self, psi, y):
    """Parameter estimation using the Normalized Sign-Regressor LMS filter.

    The normalization is used to avoid numerical instability when updating
    the estimated parameters and the sign of the information matrix is
    used to change the filter coefficients.

    Parameters
    ----------
    psi : ndarray of floats
        The information matrix of the model.
    y : array-like of shape = y_training
        The data used to training the model.

    Returns
    -------
    theta : array-like of shape = number_of_model_elements
        The estimated parameters of the model.

    Notes
    -----
    A more in-depth documentation of all methods for parameters estimation
    will be available soon. For now, please refer to the mentioned
    references.

    References
    ----------
    .. [1] Book: Hayes, M. H. (2009). Statistical digital signal processing
       and modeling. John Wiley & Sons.
    .. [2] Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação,
       análise estatística e novas estratégias de algoritmos LMS de passo
       variável.
    .. [3] Wikipedia entry on Least Mean Squares
       https://en.wikipedia.org/wiki/Least_mean_squares_filter

    """
    y, n_theta, n, theta, self.xi = self._initial_values(y, psi)

    for i in range(n_theta, n):
        psi_tmp = psi[i, :].reshape(-1, 1)
        self.xi[i, 0] = y[i, 0] - np.dot(psi_tmp.T, theta[:, i - 1])
        tmp_list = theta[:, i - 1].reshape(-1, 1) + self.mu * self.xi[i, 0] * (
            np.sign(psi_tmp) / (self.eps + np.dot(psi_tmp.T, psi_tmp))
        )
        theta[:, i] = tmp_list.flatten()

    return theta[:, -1].reshape(-1, 1)

least_mean_squares_normalized_sign_sign(psi, y)

Parameter estimation using the Normalized Sign-Sign LMS filter.

The normalization is used to avoid numerical instability when updating the estimated parameters and both the sign of the information matrix and the sign of the error vector are used to change the filter coefficients.

Parameters:

Name Type Description Default
psi ndarray of floats

The information matrix of the model.

required
y array-like of shape = y_training

The data used to training the model.

required

Returns:

Name Type Description
theta array-like of shape = number_of_model_elements

The estimated parameters of the model.

Notes

A more in-depth documentation of all methods for parameters estimation will be available soon. For now, please refer to the mentioned references.

References
  • Book: Hayes, M. H. (2009). Statistical digital signal processing and modeling. John Wiley & Sons.
  • Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação, análise estatística e novas estratégias de algoritmos LMS de passo variável.
  • Wikipedia entry on Least Mean Squares https://en.wikipedia.org/wiki/Least_mean_squares_filter
Source code in sysidentpy\parameter_estimation\estimators.py
def least_mean_squares_normalized_sign_sign(self, psi, y):
    """Parameter estimation using the Normalized Sign-Sign LMS filter.

    The normalization is used to avoid numerical instability when updating
    the estimated parameters and both the sign of the information matrix
    and the sign of the error vector are used to change the filter
    coefficients.

    Parameters
    ----------
    psi : ndarray of floats
        The information matrix of the model.
    y : array-like of shape = y_training
        The data used to training the model.

    Returns
    -------
    theta : array-like of shape = number_of_model_elements
        The estimated parameters of the model.

    Notes
    -----
    A more in-depth documentation of all methods for parameters estimation
    will be available soon. For now, please refer to the mentioned
    references.

    References
    ----------
    - Book: Hayes, M. H. (2009). Statistical digital signal processing
       and modeling. John Wiley & Sons.
    - Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação,
       análise estatística e novas estratégias de algoritmos LMS de passo
       variável.
    - Wikipedia entry on Least Mean Squares
       https://en.wikipedia.org/wiki/Least_mean_squares_filter

    """
    y, n_theta, n, theta, self.xi = self._initial_values(y, psi)

    for i in range(n_theta, n):
        psi_tmp = psi[i, :].reshape(-1, 1)
        self.xi[i, 0] = y[i, 0] - np.dot(psi_tmp.T, theta[:, i - 1])
        tmp_list = theta[:, i - 1].reshape(-1, 1) + 2 * self.mu * np.sign(
            self.xi[i, 0]
        ) * (np.sign(psi_tmp) / (self.eps + np.dot(psi_tmp.T, psi_tmp)))
        theta[:, i] = tmp_list.flatten()

    return theta[:, -1].reshape(-1, 1)

least_mean_squares_sign_error(psi, y)

Parameter estimation using the Sign-Error Least Mean Squares filter.

The sign-error LMS algorithm uses the sign of the error vector to change the filter coefficients.

Parameters:

Name Type Description Default
psi ndarray of floats

The information matrix of the model.

required
y array-like of shape = y_training

The data used to training the model.

required

Returns:

Name Type Description
theta array-like of shape = number_of_model_elements

The estimated parameters of the model.

Notes

A more in-depth documentation of all methods for parameters estimation will be available soon. For now, please refer to the mentioned references.

References
  • Book: Hayes, M. H. (2009). Statistical digital signal processing and modeling. John Wiley & Sons.
  • Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação, análise estatística e novas estratégias de algoritmos LMS de passo variável.
  • Wikipedia entry on Least Mean Squares https://en.wikipedia.org/wiki/Least_mean_squares_filter
Source code in sysidentpy\parameter_estimation\estimators.py
def least_mean_squares_sign_error(self, psi, y):
    """Parameter estimation using the Sign-Error  Least Mean Squares filter.

    The sign-error LMS algorithm uses the sign of the error vector
    to change the filter coefficients.

    Parameters
    ----------
    psi : ndarray of floats
        The information matrix of the model.
    y : array-like of shape = y_training
        The data used to training the model.

    Returns
    -------
    theta : array-like of shape = number_of_model_elements
        The estimated parameters of the model.

    Notes
    -----
    A more in-depth documentation of all methods for parameters estimation
    will be available soon. For now, please refer to the mentioned
    references.

    References
    ----------
    - Book: Hayes, M. H. (2009). Statistical digital signal processing
       and modeling. John Wiley & Sons.
    - Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação,
       análise estatística e novas estratégias de algoritmos LMS de passo
       variável.
    - Wikipedia entry on Least Mean Squares
       https://en.wikipedia.org/wiki/Least_mean_squares_filter

    """
    y, n_theta, n, theta, self.xi = self._initial_values(y, psi)

    for i in range(n_theta, n):
        psi_tmp = psi[i, :].reshape(-1, 1)
        self.xi[i, 0] = y[i, 0] - np.dot(psi_tmp.T, theta[:, i - 1])
        tmp_list = (
            theta[:, i - 1].reshape(-1, 1)
            + self.mu * np.sign(self.xi[i, 0]) * psi_tmp
        )
        theta[:, i] = tmp_list.flatten()

    return theta[:, -1].reshape(-1, 1)

least_mean_squares_sign_regressor(psi, y)

Parameter estimation using the Sign-Regressor LMS filter.

The sign-regressor LMS algorithm uses the sign of the matrix information to change the filter coefficients.

Parameters:

Name Type Description Default
psi ndarray of floats

The information matrix of the model.

required
y array-like of shape = y_training

The data used to training the model.

required

Returns:

Name Type Description
theta array-like of shape = number_of_model_elements

The estimated parameters of the model.

Notes

A more in-depth documentation of all methods for parameters estimation will be available soon. For now, please refer to the mentioned references.

References
  • Book: Hayes, M. H. (2009). Statistical digital signal processing and modeling. John Wiley & Sons.
  • Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação, análise estatística e novas estratégias de algoritmos LMS de passo variável.
  • Wikipedia entry on Least Mean Squares https://en.wikipedia.org/wiki/Least_mean_squares_filter
Source code in sysidentpy\parameter_estimation\estimators.py
def least_mean_squares_sign_regressor(self, psi, y):
    """Parameter estimation using the  Sign-Regressor LMS filter.

    The sign-regressor LMS algorithm uses the sign of the matrix
    information to change the filter coefficients.

    Parameters
    ----------
    psi : ndarray of floats
        The information matrix of the model.
    y : array-like of shape = y_training
        The data used to training the model.

    Returns
    -------
    theta : array-like of shape = number_of_model_elements
        The estimated parameters of the model.

    Notes
    -----
    A more in-depth documentation of all methods for parameters estimation
    will be available soon. For now, please refer to the mentioned
    references.

    References
    ----------
    - Book: Hayes, M. H. (2009). Statistical digital signal processing
       and modeling. John Wiley & Sons.
    - Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação,
       análise estatística e novas estratégias de algoritmos LMS de passo
       variável.
    - Wikipedia entry on Least Mean Squares
       https://en.wikipedia.org/wiki/Least_mean_squares_filter

    """
    y, n_theta, n, theta, self.xi = self._initial_values(y, psi)

    for i in range(n_theta, n):
        psi_tmp = psi[i, :].reshape(-1, 1)
        self.xi[i, 0] = y[i, 0] - np.dot(psi_tmp.T, theta[:, i - 1])
        tmp_list = theta[:, i - 1].reshape(-1, 1) + self.mu * self.xi[
            i, 0
        ] * np.sign(psi_tmp)
        theta[:, i] = tmp_list.flatten()

    return theta[:, -1].reshape(-1, 1)

least_mean_squares_sign_sign(psi, y)

Parameter estimation using the Sign-Sign LMS filter.

The sign-regressor LMS algorithm uses both the sign of the matrix information and the sign of the error vector to change the filter coefficients.

Parameters:

Name Type Description Default
psi ndarray of floats

The information matrix of the model.

required
y array-like of shape = y_training

The data used to training the model.

required

Returns:

Name Type Description
theta array-like of shape = number_of_model_elements

The estimated parameters of the model.

Notes

A more in-depth documentation of all methods for parameters estimation will be available soon. For now, please refer to the mentioned references.

References
  • Book: Hayes, M. H. (2009). Statistical digital signal processing and modeling. John Wiley & Sons.
  • Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação, análise estatística e novas estratégias de algoritmos LMS de passo variável.
  • Wikipedia entry on Least Mean Squares https://en.wikipedia.org/wiki/Least_mean_squares_filter
Source code in sysidentpy\parameter_estimation\estimators.py
def least_mean_squares_sign_sign(self, psi, y):
    """Parameter estimation using the  Sign-Sign LMS filter.

    The sign-regressor LMS algorithm uses both the sign of the matrix
    information and the sign of the error vector to change the filter
    coefficients.

    Parameters
    ----------
    psi : ndarray of floats
        The information matrix of the model.
    y : array-like of shape = y_training
        The data used to training the model.

    Returns
    -------
    theta : array-like of shape = number_of_model_elements
        The estimated parameters of the model.

    Notes
    -----
    A more in-depth documentation of all methods for parameters estimation
    will be available soon. For now, please refer to the mentioned
    references.

    References
    ----------
    - Book: Hayes, M. H. (2009). Statistical digital signal processing
       and modeling. John Wiley & Sons.
    - Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação,
       análise estatística e novas estratégias de algoritmos LMS de passo
       variável.
    - Wikipedia entry on Least Mean Squares
       https://en.wikipedia.org/wiki/Least_mean_squares_filter

    """
    y, n_theta, n, theta, self.xi = self._initial_values(y, psi)

    for i in range(n_theta, n):
        psi_tmp = psi[i, :].reshape(-1, 1)
        self.xi[i, 0] = y[i, 0] - np.dot(psi_tmp.T, theta[:, i - 1])
        tmp_list = theta[:, i - 1].reshape(-1, 1) + 2 * self.mu * np.sign(
            self.xi[i, 0]
        ) * np.sign(psi_tmp)
        theta[:, i] = tmp_list.flatten()

    return theta[:, -1].reshape(-1, 1)

least_squares(psi, y)

Estimate the model parameters using Least Squares method.

Parameters:

Name Type Description Default
psi ndarray of floats

The information matrix of the model.

required
y array-like of shape = y_training

The data used to training the model.

required

Returns:

Name Type Description
theta array-like of shape = number_of_model_elements

The estimated parameters of the model.

References
Source code in sysidentpy\parameter_estimation\estimators.py
def least_squares(self, psi, y):
    """Estimate the model parameters using Least Squares method.

    Parameters
    ----------
    psi : ndarray of floats
        The information matrix of the model.
    y : array-like of shape = y_training
        The data used to training the model.

    Returns
    -------
    theta : array-like of shape = number_of_model_elements
        The estimated parameters of the model.

    References
    ----------
    - Manuscript: Sorenson, H. W. (1970). Least-squares estimation:
       from Gauss to Kalman. IEEE spectrum, 7(7), 63-68.
       http://pzs.dstu.dp.ua/DataMining/mls/bibl/Gauss2Kalman.pdf
    - Book (Portuguese): Aguirre, L. A. (2007). Introdução identificação
       de sistemas: técnicas lineares e não-lineares aplicadas a sistemas
       reais. Editora da UFMG. 3a edição.
    - Manuscript: Markovsky, I., & Van Huffel, S. (2007).
       Overview of total least-squares methods.
       Signal processing, 87(10), 2283-2302.
       https://eprints.soton.ac.uk/263855/1/tls_overview.pdf
    - Wikipedia entry on Least Squares
       https://en.wikipedia.org/wiki/Least_squares

    """
    self._check_linear_dependence_rows(psi)

    y = y[self.max_lag :, 0].reshape(-1, 1)
    theta = np.linalg.lstsq(psi, y, rcond=None)[0]
    return theta

normalized_least_mean_squares(psi, y)

Parameter estimation using the Normalized Least Mean Squares filter.

The normalization is used to avoid numerical instability when updating the estimated parameters.

Parameters:

Name Type Description Default
psi ndarray of floats

The information matrix of the model.

required
y array-like of shape = y_training

The data used to training the model.

required

Returns:

Name Type Description
theta array-like of shape = number_of_model_elements

The estimated parameters of the model.

Notes

A more in-depth documentation of all methods for parameters estimation will be available soon. For now, please refer to the mentioned references.

References
  • Book: Hayes, M. H. (2009). Statistical digital signal processing and modeling. John Wiley & Sons.
  • Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação, análise estatística e novas estratégias de algoritmos LMS de passo variável.
  • Wikipedia entry on Least Mean Squares https://en.wikipedia.org/wiki/Least_mean_squares_filter
Source code in sysidentpy\parameter_estimation\estimators.py
def normalized_least_mean_squares(self, psi, y):
    """Parameter estimation using the Normalized Least Mean Squares filter.

    The normalization is used to avoid numerical instability when updating
    the estimated parameters.

    Parameters
    ----------
    psi : ndarray of floats
        The information matrix of the model.
    y : array-like of shape = y_training
        The data used to training the model.

    Returns
    -------
    theta : array-like of shape = number_of_model_elements
        The estimated parameters of the model.

    Notes
    -----
    A more in-depth documentation of all methods for parameters estimation
    will be available soon. For now, please refer to the mentioned
    references.

    References
    ----------
    - Book: Hayes, M. H. (2009). Statistical digital signal processing
       and modeling. John Wiley & Sons.
    - Dissertation (Portuguese): Zipf, J. G. F. (2011). Classificação,
       análise estatística e novas estratégias de algoritmos LMS de passo
       variável.
    - Wikipedia entry on Least Mean Squares
       https://en.wikipedia.org/wiki/Least_mean_squares_filter

    """
    y, n_theta, n, theta, self.xi = self._initial_values(y, psi)

    for i in range(n_theta, n):
        psi_tmp = psi[i, :].reshape(-1, 1)
        self.xi[i, 0] = y[i, 0] - np.dot(psi_tmp.T, theta[:, i - 1])
        tmp_list = theta[:, i - 1].reshape(-1, 1) + 2 * self.mu * self.xi[i, 0] * (
            psi_tmp / (self.eps + np.dot(psi_tmp.T, psi_tmp))
        )
        theta[:, i] = tmp_list.flatten()

    return theta[:, -1].reshape(-1, 1)

recursive_least_squares(psi, y)

Estimate the model parameters using the Recursive Least Squares method.

The implementation consider the forgetting factor.

Parameters:

Name Type Description Default
psi ndarray of floats

The information matrix of the model.

required
y array-like of shape = y_training

The data used to training the model.

required

Returns:

Name Type Description
theta array-like of shape = number_of_model_elements

The estimated parameters of the model.

Notes

A more in-depth documentation of all methods for parameters estimation will be available soon. For now, please refer to the mentioned references.

References
  • Book (Portuguese): Aguirre, L. A. (2007). Introdução identificação de sistemas: técnicas lineares e não-lineares aplicadas a sistemas reais. Editora da UFMG. 3a edição.
Source code in sysidentpy\parameter_estimation\estimators.py
def recursive_least_squares(self, psi, y):
    """Estimate the model parameters using the Recursive Least Squares method.

    The implementation consider the forgetting factor.
    Parameters
    ----------
    psi : ndarray of floats
        The information matrix of the model.
    y : array-like of shape = y_training
        The data used to training the model.

    Returns
    -------
    theta : array-like of shape = number_of_model_elements
        The estimated parameters of the model.

    Notes
    -----
    A more in-depth documentation of all methods for parameters estimation
    will be available soon. For now, please refer to the mentioned
    references.

    References
    ----------
    - Book (Portuguese): Aguirre, L. A. (2007). Introdução identificação
       de sistemas: técnicas lineares e não-lineares aplicadas a sistemas
       reais. Editora da UFMG. 3a edição.

    """
    y, n_theta, n, theta, self.xi = self._initial_values(y, psi)

    p = np.eye(n_theta) / self.delta

    for i in range(2, n):
        psi_tmp = psi[i, :].reshape(-1, 1)
        k_numerator = self.lam ** (-1) * p.dot(psi_tmp)
        k_denominator = 1 + self.lam ** (-1) * psi_tmp.T.dot(p).dot(psi_tmp)
        k = np.divide(k_numerator, k_denominator)
        self.xi[i, 0] = y[i, 0] - np.dot(psi_tmp.T, theta[:, i - 1])
        tmp_list = theta[:, i - 1].reshape(-1, 1) + k.dot(self.xi[i, 0])
        theta[:, i] = tmp_list.flatten()

        p1 = p.dot(psi[i, :].reshape(-1, 1)).dot(psi[i, :].reshape(-1, 1).T).dot(p)
        p2 = (
            psi[i, :].reshape(-1, 1).T.dot(p).dot(psi[i, :].reshape(-1, 1))
            + self.lam
        )

        p_numerator = p - np.divide(p1, p2)
        p = np.divide(p_numerator, self.lam)

    self.theta_evolution = theta.copy()
    return theta[:, -1].reshape(-1, 1)

ridge_regression(psi, y)

Estimate the model parameters using the regularized least squares method known as ridge regression. Based on the least_squares module and uses the same data format but you need to pass ridge_param in the call to FROLS

Parameters:

Name Type Description Default
psi ndarray of floats

The information matrix of the model.

required
y array-like of shape = y_training

The data used to training the model.

required
ridge_param ridge regression parameter that regularizes the algorithm

to prevent over fitting. If the input is a noisy signal, the ridge parameter is likely to be set close to the noise level, at least as a starting point. Entered through the self data structure.

required

Returns:

Name Type Description
theta array-like of shape = number_of_model_elements

The estimated parameters of the model.

References

ridge_parm multiplied by the identity matrix (np.eye) favors models (theta) that have small size using an L2 norm. This prevents over fitting of the model. For applications where preventing overfitting is important, see, for example, D. J. Gauthier, E. Bollt, A. Griffith, W. A. S. Barbosa, ‘Next generation reservoir computing,’ Nat. Commun. 12, 5564 (2021). https://www.nature.com/articles/s41467-021-25801-2

Source code in sysidentpy\parameter_estimation\estimators.py
def ridge_regression(self, psi, y):
    """Estimate the model parameters using the regularized least squares method
       known as ridge regression.  Based on the least_squares module and uses
       the same data format but you need to pass ridge_param in the call to
       FROLS

    Parameters
    ----------
    psi : ndarray of floats
        The information matrix of the model.
    y : array-like of shape = y_training
        The data used to training the model.
    ridge_param : ridge regression parameter that regularizes the algorithm
        to prevent over fitting.  If the input is a noisy signal, the ridge
        parameter is likely to be set close to the noise level, at least
        as a starting point.  Entered through the self data structure.

    Returns
    -------
    theta : array-like of shape = number_of_model_elements
        The estimated parameters of the model.

    References
    ----------
    - Wikipedia entry on ridge regression
      https://en.wikipedia.org/wiki/Ridge_regression

    ridge_parm multiplied by the identity matrix (np.eye) favors models (theta) that
    have small size using an L2 norm.  This prevents over fitting of the model.
    For applications where preventing overfitting is important, see, for example,
    D. J. Gauthier, E. Bollt, A. Griffith, W. A. S. Barbosa, ‘Next generation
    reservoir computing,’ Nat. Commun. 12, 5564 (2021).
    https://www.nature.com/articles/s41467-021-25801-2

    """
    self._check_linear_dependence_rows(psi)

    y = y[self.max_lag :, 0].reshape(-1, 1)
    theta = (
        np.linalg.pinv(psi.T @ psi + self.ridge_param * np.eye(psi.shape[1]))
        @ psi.T
        @ y
    )
    return theta

total_least_squares(psi, y)

Estimate the model parameters using Total Least Squares method.

Parameters:

Name Type Description Default
psi ndarray of floats

The information matrix of the model.

required
y array-like of shape = y_training

The data used to training the model.

required

Returns:

Name Type Description
theta array-like of shape = number_of_model_elements

The estimated parameters of the model.

References
Source code in sysidentpy\parameter_estimation\estimators.py
def total_least_squares(self, psi, y):
    """Estimate the model parameters using Total Least Squares method.

    Parameters
    ----------
    psi : ndarray of floats
        The information matrix of the model.
    y : array-like of shape = y_training
        The data used to training the model.

    Returns
    -------
    theta : array-like of shape = number_of_model_elements
        The estimated parameters of the model.

    References
    ----------
    - Manuscript: Golub, G. H., & Van Loan, C. F. (1980).
       An analysis of the total least squares problem.
       SIAM journal on numerical analysis, 17(6), 883-893.
    - Manuscript: Markovsky, I., & Van Huffel, S. (2007).
       Overview of total least-squares methods.
       Signal processing, 87(10), 2283-2302.
       https://eprints.soton.ac.uk/263855/1/tls_overview.pdf
    - Wikipedia entry on Total Least Squares
       https://en.wikipedia.org/wiki/Total_least_squares

    """
    y = y[self.max_lag :, 0].reshape(-1, 1)
    full = np.hstack((psi, y))
    n = psi.shape[1]
    _, _, v = np.linalg.svd(full, full_matrices=True)
    theta = -v.T[:n, n:] / v.T[n:, n:]
    return theta.reshape(-1, 1)