Classes | Typedefs | Enumerations | Functions
WireCell::Ress Namespace Reference

Classes

struct  Params
 

Typedefs

typedef Eigen::VectorXd vector_t
 
typedef Eigen::MatrixXd matrix_t
 

Enumerations

enum  Model { unknown =0, lasso, elnet }
 

Functions

vector_t solve (matrix_t response, vector_t measured, const Params &params=Params(), vector_t source=Eigen::VectorXd(), vector_t weights=Eigen::VectorXd())
 
vector_t predict (matrix_t response, vector_t source)
 
double chi2 (vector_t measured, vector_t predicted)
 
double mean_residual (vector_t measured, vector_t predicted)
 
double chi2_l1 (vector_t measured, vector_t solved, double lambda=1.0)
 

Typedef Documentation

typedef Eigen::MatrixXd WireCell::Ress::matrix_t

Definition at line 13 of file Ress.h.

typedef Eigen::VectorXd WireCell::Ress::vector_t

Definition at line 12 of file Ress.h.

Enumeration Type Documentation

Enumerator
unknown 
lasso 
elnet 

Definition at line 15 of file Ress.h.

15  {
16  unknown=0,
17  lasso, // Lasso model
18  elnet // elastic net
19  };

Function Documentation

double WireCell::Ress::chi2 ( vector_t  measured,
vector_t  predicted 
)

Definition at line 57 of file Ress.cxx.

58 {
59  return ( measured- predicted ).squaredNorm();
60 }
double WireCell::Ress::chi2_l1 ( vector_t  measured,
vector_t  solved,
double  lambda = 1.0 
)

Definition at line 47 of file Ress.cxx.

48 {
49  return 2 * lambda * solved.lpNorm<1>() * measured.size();
50 }
double WireCell::Ress::mean_residual ( vector_t  measured,
vector_t  predicted 
)

Definition at line 62 of file Ress.cxx.

63 {
64  return ( measured - predicted ).norm() / measured.size();
65 }
Ress::vector_t WireCell::Ress::predict ( matrix_t  response,
vector_t  source 
)

Definition at line 52 of file Ress.cxx.

53 {
54  return response * source;
55 }
const CharType(& source)[N]
Definition: pointer.h:1147
Ress::vector_t WireCell::Ress::solve ( Ress::matrix_t  matrix,
Ress::vector_t  measured,
const Params params = Params(),
Ress::vector_t  initial = Eigen::VectorXd(),
Ress::vector_t  weights = Eigen::VectorXd() 
)

Definition at line 7 of file Ress.cxx.

12 {
13  // Provide a uniform interface to RESS solving models. RESS
14  // *almost* already provides this.
15 
16  if (params.model == Ress::lasso) {
17  WireCell::LassoModel model(params.lambda,
18  params.max_iter, params.tolerance, params.non_negative);
19  if (initial.size()) {
20  model.Setbeta(initial);
21  }
22  if (weights.size()) {
23  model.SetLambdaWeight(weights);
24  }
25  model.SetData(matrix, measured);
26  model.Fit();
27  return model.Getbeta();
28  }
29 
30  if (params.model == Ress::elnet) {
31  WireCell::ElasticNetModel model(params.lambda, params.alpha,
32  params.max_iter, params.tolerance, params.non_negative);
33  if (initial.size()) {
34  model.Setbeta(initial);
35  }
36  if (weights.size()) {
37  model.SetLambdaWeight(weights);
38  }
39  model.SetData(matrix, measured);
40  model.Fit();
41  return model.Getbeta();
42  }
43 
44  return Ress::vector_t();
45 }
Eigen::VectorXd vector_t
Definition: Ress.h:12