1namespace Eigen { 2 3/** \eigenManualPage InplaceDecomposition Inplace matrix decompositions 4 5Starting from %Eigen 3.3, the LU, Cholesky, and QR decompositions can operate \em inplace, that is, directly within the given input matrix. 6This feature is especially useful when dealing with huge matrices, and or when the available memory is very limited (embedded systems). 7 8To this end, the respective decomposition class must be instantiated with a Ref<> matrix type, and the decomposition object must be constructed with the input matrix as argument. As an example, let us consider an inplace LU decomposition with partial pivoting. 9 10Let's start with the basic inclusions, and declaration of a 2x2 matrix \c A: 11 12<table class="example"> 13<tr><th>code</th><th>output</th></tr> 14<tr> 15 <td>\snippet TutorialInplaceLU.cpp init 16 </td> 17 <td>\snippet TutorialInplaceLU.out init 18 </td> 19</tr> 20</table> 21 22No surprise here! Then, let's declare our inplace LU object \c lu, and check the content of the matrix \c A: 23 24<table class="example"> 25<tr> 26 <td>\snippet TutorialInplaceLU.cpp declaration 27 </td> 28 <td>\snippet TutorialInplaceLU.out declaration 29 </td> 30</tr> 31</table> 32 33Here, the \c lu object computes and stores the \c L and \c U factors within the memory held by the matrix \c A. 34The coefficients of \c A have thus been destroyed during the factorization, and replaced by the L and U factors as one can verify: 35 36<table class="example"> 37<tr> 38 <td>\snippet TutorialInplaceLU.cpp matrixLU 39 </td> 40 <td>\snippet TutorialInplaceLU.out matrixLU 41 </td> 42</tr> 43</table> 44 45Then, one can use the \c lu object as usual, for instance to solve the Ax=b problem: 46<table class="example"> 47<tr> 48 <td>\snippet TutorialInplaceLU.cpp solve 49 </td> 50 <td>\snippet TutorialInplaceLU.out solve 51 </td> 52</tr> 53</table> 54 55Here, since the content of the original matrix \c A has been lost, we had to declared a new matrix \c A0 to verify the result. 56 57Since the memory is shared between \c A and \c lu, modifying the matrix \c A will make \c lu invalid. 58This can easily be verified by modifying the content of \c A and trying to solve the initial problem again: 59 60<table class="example"> 61<tr> 62 <td>\snippet TutorialInplaceLU.cpp modifyA 63 </td> 64 <td>\snippet TutorialInplaceLU.out modifyA 65 </td> 66</tr> 67</table> 68 69Note that there is no shared pointer under the hood, it is the \b responsibility \b of \b the \b user to keep the input matrix \c A in life as long as \c lu is living. 70 71If one wants to update the factorization with the modified A, one has to call the compute method as usual: 72<table class="example"> 73<tr> 74 <td>\snippet TutorialInplaceLU.cpp recompute 75 </td> 76 <td>\snippet TutorialInplaceLU.out recompute 77 </td> 78</tr> 79</table> 80 81Note that calling compute does not change the memory which is referenced by the \c lu object. Therefore, if the compute method is called with another matrix \c A1 different than \c A, then the content of \c A1 won't be modified. This is still the content of \c A that will be used to store the L and U factors of the matrix \c A1. 82This can easily be verified as follows: 83<table class="example"> 84<tr> 85 <td>\snippet TutorialInplaceLU.cpp recompute_bis0 86 </td> 87 <td>\snippet TutorialInplaceLU.out recompute_bis0 88 </td> 89</tr> 90</table> 91The matrix \c A1 is unchanged, and one can thus solve A1*x=b, and directly check the residual without any copy of \c A1: 92<table class="example"> 93<tr> 94 <td>\snippet TutorialInplaceLU.cpp recompute_bis1 95 </td> 96 <td>\snippet TutorialInplaceLU.out recompute_bis1 97 </td> 98</tr> 99</table> 100 101 102Here is the list of matrix decompositions supporting this inplace mechanism: 103 104- class LLT 105- class LDLT 106- class PartialPivLU 107- class FullPivLU 108- class HouseholderQR 109- class ColPivHouseholderQR 110- class FullPivHouseholderQR 111- class CompleteOrthogonalDecomposition 112 113*/ 114 115}