Tutorial_sparse_example.cpp File Reference

Go to the source code of this file.

Typedefs

typedef Eigen::SparseMatrix< double > SpMat
 
typedef Eigen::Triplet< double > T
 

Functions

void buildProblem (std::vector< T > &coefficients, Eigen::VectorXd &b, int n)
 
int main (int argc, char **argv)
 
void saveAsBitmap (const Eigen::VectorXd &x, int n, const char *filename)
 

Typedef Documentation

◆ SpMat

typedef Eigen::SparseMatrix<double> SpMat

Definition at line 5 of file Tutorial_sparse_example.cpp.

◆ T

typedef Eigen::Triplet<double> T

Definition at line 6 of file Tutorial_sparse_example.cpp.

Function Documentation

◆ buildProblem()

void buildProblem ( std::vector< T > &  coefficients,
Eigen::VectorXd b,
int  n 
)

Definition at line 19 of file Tutorial_sparse_example_details.cpp.

20 {
21  b.setZero();
22  Eigen::ArrayXd boundary = Eigen::ArrayXd::LinSpaced(n, 0,M_PI).sin().pow(2);
23  for(int j=0; j<n; ++j)
24  {
25  for(int i=0; i<n; ++i)
26  {
27  int id = i+j*n;
28  insertCoefficient(id, i-1,j, -1, coefficients, b, boundary);
29  insertCoefficient(id, i+1,j, -1, coefficients, b, boundary);
30  insertCoefficient(id, i,j-1, -1, coefficients, b, boundary);
31  insertCoefficient(id, i,j+1, -1, coefficients, b, boundary);
32  insertCoefficient(id, i,j, 4, coefficients, b, boundary);
33  }
34  }
35 }
Array< int, 3, 1 > b
int n
void insertCoefficient(int id, int i, int j, double w, std::vector< T > &coeffs, Eigen::VectorXd &b, const Eigen::VectorXd &boundary)
General-purpose arrays with easy API for coefficient-wise operations.
Definition: Array.h:49
static EIGEN_DEPRECATED const RandomAccessLinSpacedReturnType LinSpaced(Sequential_t, Index size, const Scalar &low, const Scalar &high)
std::ptrdiff_t j

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 11 of file Tutorial_sparse_example.cpp.

12 {
13  if(argc!=2) {
14  std::cerr << "Error: expected one and only one argument.\n";
15  return -1;
16  }
17 
18  int n = 300; // size of the image
19  int m = n*n; // number of unknowns (=number of pixels)
20 
21  // Assembly:
22  std::vector<T> coefficients; // list of non-zeros coefficients
23  Eigen::VectorXd b(m); // the right hand side-vector resulting from the constraints
24  buildProblem(coefficients, b, n);
25 
26  SpMat A(m,m);
27  A.setFromTriplets(coefficients.begin(), coefficients.end());
28 
29  // Solving:
30  Eigen::SimplicialCholesky<SpMat> chol(A); // performs a Cholesky factorization of A
31  Eigen::VectorXd x = chol.solve(b); // use the factorization to solve for the given right hand side
32 
33  // Export the result to a file:
34  saveAsBitmap(x, n, argv[1]);
35 
36  return 0;
37 }
Matrix3f m
MatrixXcf A
void buildProblem(std::vector< T > &coefficients, Eigen::VectorXd &b, int n)
void saveAsBitmap(const Eigen::VectorXd &x, int n, const char *filename)
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:182

◆ saveAsBitmap()

void saveAsBitmap ( const Eigen::VectorXd x,
int  n,
const char *  filename 
)

Definition at line 37 of file Tutorial_sparse_example_details.cpp.

38 {
39  Eigen::Array<unsigned char,Eigen::Dynamic,Eigen::Dynamic> bits = (x*255).cast<unsigned char>();
40  QImage img(bits.data(), n,n,QImage::Format_Indexed8);
41  img.setColorCount(256);
42  for(int i=0;i<256;i++) img.setColor(i,qRgb(i,i,i));
43  img.save(filename);
44 }
const Scalar * data() const