TutorialInplaceLU.cpp File Reference
#include <iostream>
#include <Eigen/Dense>

Classes

struct  init
 

Functions

int main ()
 

Variables

init init_obj
 

Function Documentation

◆ main()

int main ( )
9  {
10  Eigen::MatrixXd A(2, 2);
11  A << 2, -1, 1, 3;
12  std::cout << "Here is the input matrix A before decomposition:\n" << A << "\n";
13  std::cout << "[init]\n";
14 
15  std::cout << "[declaration]\n";
17  std::cout << "Here is the input matrix A after decomposition:\n" << A << "\n";
18  std::cout << "[declaration]\n";
19 
20  std::cout << "[matrixLU]\n";
21  std::cout << "Here is the matrix storing the L and U factors:\n" << lu.matrixLU() << "\n";
22  std::cout << "[matrixLU]\n";
23 
24  std::cout << "[solve]\n";
25  Eigen::MatrixXd A0(2, 2);
26  A0 << 2, -1, 1, 3;
27  Eigen::VectorXd b(2);
28  b << 1, 2;
29  Eigen::VectorXd x = lu.solve(b);
30  std::cout << "Residual: " << (A0 * x - b).norm() << "\n";
31  std::cout << "[solve]\n";
32 
33  std::cout << "[modifyA]\n";
34  A << 3, 4, -2, 1;
35  x = lu.solve(b);
36  std::cout << "Residual: " << (A0 * x - b).norm() << "\n";
37  std::cout << "[modifyA]\n";
38 
39  std::cout << "[recompute]\n";
40  A0 = A; // save A
41  lu.compute(A);
42  x = lu.solve(b);
43  std::cout << "Residual: " << (A0 * x - b).norm() << "\n";
44  std::cout << "[recompute]\n";
45 
46  std::cout << "[recompute_bis0]\n";
47  Eigen::MatrixXd A1(2, 2);
48  A1 << 5, -2, 3, 4;
49  lu.compute(A1);
50  std::cout << "Here is the input matrix A1 after decomposition:\n" << A1 << "\n";
51  std::cout << "[recompute_bis0]\n";
52 
53  std::cout << "[recompute_bis1]\n";
54  x = lu.solve(b);
55  std::cout << "Residual: " << (A1 * x - b).norm() << "\n";
56  std::cout << "[recompute_bis1]\n";
57 }
Scalar * b
Definition: benchVecAdd.cpp:17
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition: bench_gemm.cpp:47
LU decomposition of a matrix with partial pivoting, and related features.
Definition: PartialPivLU.h:77
cout<< "Here is the matrix m:"<< endl<< m<< endl;Eigen::FullPivLU< Matrix5x3 > lu(m)
list x
Definition: plotDoE.py:28

References b, lu(), and plotDoE::x.

Variable Documentation

◆ init_obj

init init_obj