sablib
Loading...
Searching...
No Matches
airpls.cpp
Go to the documentation of this file.
1
6
7#include <cmath>
8
9#include "airpls.h"
10
11namespace sablib {
12
13//
14// Implementation of BaselineAirPLS() function
15//
16const std::vector<double> BaselineAirPLS(
17 std::vector<double> & y, const double lambda, const unsigned int s,
18 const unsigned int loop, const double eps
19)
20{
21 if(y.size() == 0) {
22 throw std::invalid_argument("BaselineAirPLS(): the length of y is zero.");
23 }
24
25 if(lambda <= 0) {
26 throw std::invalid_argument("BaselineAirPLS(): non-positive lambda value is given.");
27 }
28
29 if(s == 0 || s > 3) {
30 throw std::invalid_argument("BaselineAirPLS(): s must be 1, 2 or 3.");
31 }
32
33 if(loop == 0) {
34 throw std::invalid_argument("BaselineAirPLS(): loop is zero.");
35 }
36
37 if(eps <= 0) {
38 throw std::invalid_argument("BaselineAirPLS(): non-positive eps value is given.");
39 }
40
41 size_t m = y.size();
42 Eigen::VectorXd yy, w, z, d;
43 Eigen::SparseMatrix<double> I, D, lambdaDTD;
44 double y_abs_sum, d_sum_abs;
45
46 yy = Eigen::VectorXd::Map(y.data(), m);
47
48 w.setOnes(m);
49 y_abs_sum = yy.array().abs().matrix().sum();
50
51 I.resize(m, m);
52 I.setIdentity();
53 D = Diff(I, s);
54 lambdaDTD = lambda * (D.transpose() * D);
55
56 for(unsigned int i = 0; i < loop; i++) {
57 z = Whittaker(yy, w, lambdaDTD);
58
59 d = (yy.array() >= z.array()).select(0, yy - z);
60 d_sum_abs = std::fabs(d.sum());
61
62 if (d_sum_abs < eps * y_abs_sum) {
63 break;
64 }
65
66 w = (yy.array() >= z.array()).select(0, ((loop * d.array().abs()) / d_sum_abs).exp());
67 w(0) = w(w.size() - 1) = std::exp((loop * d.maxCoeff() / d_sum_abs));
68 }
69
70 std::vector<double> result(z.size());
71
72 Eigen::VectorXd::Map(result.data(), result.size()) = z;
73
74 return result;
75}
76
77}; // namespace sablib
const std::vector< double > BaselineAirPLS(std::vector< double > &y, const double lambda, const unsigned int s, const unsigned int loop, const double eps)
Performs baseline estimation using adaptive iteratively reweighted Penalized Least Squares(airPLS).
Definition airpls.cpp:16
Baseline estimation using adaptive iteratively reweighted Penalized Least Squares(airPLS).
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