Public Member Functions | List of all members
WireCell::LassoModel Class Reference

#include <LassoModel.h>

Inheritance diagram for WireCell::LassoModel:
WireCell::ElasticNetModel WireCell::LinearModel

Public Member Functions

 LassoModel (double lambda=1., int max_iter=100000, double TOL=1e-3, bool non_negtive=true)
 
 ~LassoModel ()
 
void Fit ()
 
void Set_init_values (std::vector< double > values)
 
double chi2_l1 ()
 
- Public Member Functions inherited from WireCell::ElasticNetModel
 ElasticNetModel (double lambda=1., double alpha=1., int max_iter=100000, double TOL=1e-3, bool non_negtive=true)
 
 ~ElasticNetModel ()
 
void SetLambdaWeight (Eigen::VectorXd w)
 
void SetLambdaWeight (int i, double weight)
 
void SetX (Eigen::MatrixXd X)
 
- Public Member Functions inherited from WireCell::LinearModel
 LinearModel ()
 
virtual ~LinearModel ()
 
Eigen::VectorXd & Gety ()
 
Eigen::MatrixXd & GetX ()
 
Eigen::VectorXd & Getbeta ()
 
virtual void SetData (Eigen::MatrixXd X, Eigen::VectorXd y)
 
virtual void Sety (Eigen::VectorXd y)
 
virtual void Setbeta (Eigen::VectorXd beta)
 
Eigen::VectorXd Predict ()
 
double chi2_base ()
 
double MeanResidual ()
 

Additional Inherited Members

- Public Attributes inherited from WireCell::ElasticNetModel
double lambda
 
double alpha
 
int max_iter
 
double TOL
 
bool non_negtive
 
Eigen::VectorXd lambda_weight
 
- Public Attributes inherited from WireCell::LinearModel
std::string name
 
- Protected Member Functions inherited from WireCell::ElasticNetModel
double _soft_thresholding (double x, double lambda_)
 
- Protected Attributes inherited from WireCell::ElasticNetModel
std::vector< bool_active_beta
 
- Protected Attributes inherited from WireCell::LinearModel
Eigen::VectorXd _y
 
Eigen::MatrixXd _X
 
Eigen::VectorXd _beta
 

Detailed Description

Definition at line 8 of file LassoModel.h.

Constructor & Destructor Documentation

WireCell::LassoModel::LassoModel ( double  lambda = 1.,
int  max_iter = 100000,
double  TOL = 1e-3,
bool  non_negtive = true 
)

Definition at line 14 of file LassoModel.cxx.

16 {
17  name = "Lasso";
18 }
ElasticNetModel(double lambda=1., double alpha=1., int max_iter=100000, double TOL=1e-3, bool non_negtive=true)
WireCell::LassoModel::~LassoModel ( )

Definition at line 20 of file LassoModel.cxx.

21 {}

Member Function Documentation

double WireCell::LassoModel::chi2_l1 ( )

Definition at line 125 of file LassoModel.cxx.

126 {
127  return 2 * lambda * Getbeta().lpNorm<1>() * Gety().size() ;
128 }
Eigen::VectorXd & Gety()
Definition: LinearModel.h:14
Eigen::VectorXd & Getbeta()
Definition: LinearModel.h:16
void WireCell::LassoModel::Fit ( )
virtual

Reimplemented from WireCell::ElasticNetModel.

Definition at line 32 of file LassoModel.cxx.

33 {
34  // initialize solution to zero unless user set beta already
35  Eigen::VectorXd beta = _beta;
36  if (0 == beta.size()) {
37  beta = VectorXd::Zero(_X.cols());
38  }
39 
40  // initialize active_beta to true
41  int nbeta = beta.size();
42  _active_beta = vector<bool>(nbeta, true);
43 
44  Eigen::VectorXd y = Gety();
45  Eigen::MatrixXd X = GetX();
46 
47  // cooridate decsent
48  //int N = y.size();
49  VectorXd norm(nbeta);
50  for (int j=0; j<nbeta; j++) {
51  norm(j) = X.col(j).squaredNorm();
52  if (norm(j) < 1e-6) {
53  cerr << "warning: the " << j << "th variable is not used, please consider removing it." << endl;
54  norm(j) = 1;
55  }
56  }
57  double tol2 = TOL*TOL*nbeta;
58 
59  // calculate the inner product
60  Eigen::VectorXd ydX(nbeta);
61  Eigen::SparseMatrix<double> XdX(nbeta,nbeta);
62  for (int i=0;i!=nbeta;i++){
63  ydX(i) = y.dot(X.col(i));
64  //beta(i) = ydX(i) / norm(i); // first time result saved here ...
65  for (int j=0;j!=nbeta;j++){
66  double value = X.col(i).dot(X.col(j));
67  if (value!=0)
68  XdX.insert(i,j) = value;
69  //std::cout << i << " " << j << " " << value << " " << XdX.coeffRef(i,j) << std::endl;
70  }
71  //std::cout << ydX(i) << " " << norm(i) << std::endl;
72  }
73  // std::cout << "Xin " << nbeta << std::endl;
74 
75 
76  // start interation ...
77  int double_check = 0;
78  for (int i =0; i< max_iter; i++){
79  VectorXd betalast = beta;
80 
81  //loop through sparse matrix ...
82  for (int j=0;j!=nbeta;j++){
83  if (!_active_beta[j]) {continue;}
84  beta(j) = ydX(j);
85  for (SparseMatrix<double>::InnerIterator it(XdX,j); it; ++it){
86  //if (it.row()!=j && beta(it.row())!=0)
87  if (it.row()!=j)
88  beta(j) -= it.value() * beta(it.row());
89  }
90  beta(j) = _soft_thresholding( beta(j)/norm(j), lambda * lambda_weight(j));
91 
92  if(fabs(beta(j)) < 1e-6) { _active_beta[j] = false; }
93  // VectorXd X_j = X.col(j);
94  // VectorXd beta_tmp = betalast;
95  // beta_tmp(j) = 0;
96  // VectorXd r_j = (y - X * beta_tmp);
97  // double delta_j = X_j.dot(r_j);
98  //std::cout << i << " " << j << " " << beta(j) << " " << betalast(j) << std::endl;
99 
100  }
101  double_check++;
102  // cout << endl;
103  VectorXd diff = beta - betalast;
104 
105  // std::cout << i << " " << diff.squaredNorm() << " " << tol2 << std::endl;
106 
107  if (diff.squaredNorm()<tol2) {
108  if (double_check!=1) {
109  double_check = 0;
110  for (int k=0; k<nbeta; k++) {
111  _active_beta[k] = true;
112  }
113  }else {
114  //cout << "found minimum at iteration: " << i << " " << flag_initial_values << endl;
115  break;
116  }
117  }
118  }
119 
120 
121  // save results in the model
122  Setbeta(beta);
123 }
Eigen::MatrixXd _X
Definition: LinearModel.h:33
std::vector< bool > _active_beta
Eigen::MatrixXd & GetX()
Definition: LinearModel.h:15
virtual void Setbeta(Eigen::VectorXd beta)
Definition: LinearModel.h:21
Eigen::VectorXd _beta
Definition: LinearModel.h:34
double y
const double e
Eigen::VectorXd & Gety()
Definition: LinearModel.h:14
double _soft_thresholding(double x, double lambda_)
auto norm(Vector const &v)
Return norm of the specified vector.
const GenericPointer< typename T::ValueType > T2 value
Definition: pointer.h:1225
Eigen::VectorXd lambda_weight
QTextStream & endl(QTextStream &s)
void WireCell::LassoModel::Set_init_values ( std::vector< double >  values)

Definition at line 23 of file LassoModel.cxx.

24 {
25  const size_t nvals = values.size();
26  _beta = VectorXd::Zero(nvals);
27  for (size_t i=0; i != nvals; ++i) {
28  _beta(i) = values.at(i);
29  }
30 }
Eigen::VectorXd _beta
Definition: LinearModel.h:34

The documentation for this class was generated from the following files: