sablib
Loading...
Searching...
No Matches
arpls.cpp
Go to the documentation of this file.
1
6
7#include <cmath>
8
10
11#include "arpls.h"
12
13namespace sablib {
14
15//
16// Implementation of BaselineArPLS() function
17//
18const std::vector<double> BaselineArPLS(
19 std::vector<double> & y, const double lambda, const unsigned int s,
20 const unsigned int loop, const double eps
21)
22{
23 if(y.size() == 0) {
24 throw std::invalid_argument("BaselineArPLS(): the length of y is zero.");
25 }
26
27 if(lambda <= 0) {
28 throw std::invalid_argument("BaselineArPLS(): non-positive lambda value is given.");
29 }
30
31 if(s == 0 || s > 3) {
32 throw std::invalid_argument("BaselineArPLS(): s must be 1, 2 or 3.");
33 }
34
35 if(loop == 0) {
36 throw std::invalid_argument("BaselineArPLS(): loop is zero.");
37 }
38
39 if(eps <= 0) {
40 throw std::invalid_argument("BaselineArPLS(): non-positive eps value is given.");
41 }
42
43 size_t m = y.size();
44 Eigen::VectorXd yy, w, wt, z, d, tmp;
45 Eigen::SparseMatrix<double> I, D, lambdaDTD;
46 double mean, sd;
47 Eigen::Index ct;
48
49 auto logistic = [](const double xx) {
50 if(xx < 0) {
51 double t = std::exp(xx);
52 return t / (t + 1);
53 }
54 else {
55 return 1 / (1 + std::exp(-xx));
56 }
57 };
58
59 yy = Eigen::VectorXd::Map(y.data(), m);
60
61 w.setOnes(m);
62
63 I.resize(m, m);
64 I.setIdentity();
65 D = Diff(I, s);
66 lambdaDTD = lambda * (D.transpose() * D);
67
68 for(unsigned int i = 0; i < loop; i++) {
69 z = Whittaker(yy, w, lambdaDTD);
70
71 d = yy - z;
72
73 ct = (d.array() < 0).count();
74 mean = (d.array() < 0).select(d, 0).sum() / ct;
75 sd = std::sqrt((d.array() < 0).select((d.array() - mean).matrix(), 0).squaredNorm() / (ct - 1));
76
77 tmp = (-2 * (d.array() + mean - 2 * sd) / sd).matrix().unaryExpr(logistic);
78 wt = (d.array() >= 0).select(tmp, 1);
79
80 if ((w - wt).norm() / w.norm() < eps) {
81 break;
82 }
83
84 w = wt;
85 }
86
87 std::vector<double> result(z.size());
88 Eigen::VectorXd::Map(result.data(), result.size()) = z;
89
90 return result;
91}
92
93}; // namespace sablib
const std::vector< double > BaselineArPLS(std::vector< double > &y, const double lambda, const unsigned int s, const unsigned int loop, const double eps)
Performs baseline estimation using asymmetrically reweighted Penalized Least Squares(arPLS).
Definition arpls.cpp:18
Baseline estimation using asymmetrically reweighted Penalized Least Squares(arPLS).
const Derived::PlainObject Diff(const Eigen::MatrixBase< Derived > &m0, const int n=1, const Dir dir=Dir::RowWise)
Calculates the n-th discrete difference along the given axis.
Definition diff.h:32
const std::vector< double > Whittaker(const std::vector< double > &y, const std::vector< double > &w, const double lambda, const unsigned int s)
Performs Whittaker smoothing (std::vector<double> version, with weights).
Definition whittaker.cpp:14
Smoothing using Whittaker smoother.