1 #include <Eigen/Sparse>
2 #include <vector>
3 #include <QImage>
4 
5 typedef Eigen::SparseMatrix<double> SpMat; // declares a column-major sparse matrix type of double
6 typedef Eigen::Triplet<double> T;
7 
insertCoefficient(int id,int i,int j,double w,std::vector<T> & coeffs,Eigen::VectorXd & b,const Eigen::VectorXd & boundary)8 void insertCoefficient(int id, int i, int j, double w, std::vector<T>& coeffs,
9                        Eigen::VectorXd& b, const Eigen::VectorXd& boundary)
10 {
11   int n = boundary.size();
12   int id1 = i+j*n;
13 
14         if(i==-1 || i==n) b(id) -= w * boundary(j); // constrained coefficient
15   else  if(j==-1 || j==n) b(id) -= w * boundary(i); // constrained coefficient
16   else  coeffs.push_back(T(id,id1,w));              // unknown coefficient
17 }
18 
buildProblem(std::vector<T> & coefficients,Eigen::VectorXd & b,int n)19 void buildProblem(std::vector<T>& coefficients, Eigen::VectorXd& b, int n)
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 }
36 
saveAsBitmap(const Eigen::VectorXd & x,int n,const char * filename)37 void saveAsBitmap(const Eigen::VectorXd& x, int n, const char* filename)
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 }
45