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