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