<!DOCTYPE html PUBLIC "-//W3C//DTD XHTML 1.0 Transitional//EN" "http://www.w3.org/TR/xhtml1/DTD/xhtml1-transitional.dtd"> <html xmlns="http://www.w3.org/1999/xhtml"> <head> <meta http-equiv="Content-Type" content="text/xhtml;charset=UTF-8"/> <title>linbox: LQUPMatrix< Field > Class Template Reference</title> <link href="tabs.css" rel="stylesheet" type="text/css"/> <link href="doxygen.css" rel="stylesheet" type="text/css"/> </head> <body> <!-- Generated by Doxygen 1.7.4 --> <script type="text/javascript"> function hasClass(ele,cls) { return ele.className.match(new RegExp('(\\s|^)'+cls+'(\\s|$)')); } function addClass(ele,cls) { if (!this.hasClass(ele,cls)) ele.className += " "+cls; } function removeClass(ele,cls) { if (hasClass(ele,cls)) { var reg = new RegExp('(\\s|^)'+cls+'(\\s|$)'); ele.className=ele.className.replace(reg,' '); } } function toggleVisibility(linkObj) { var base = linkObj.getAttribute('id'); var summary = document.getElementById(base + '-summary'); var content = document.getElementById(base + '-content'); var trigger = document.getElementById(base + '-trigger'); if ( hasClass(linkObj,'closed') ) { summary.style.display = 'none'; content.style.display = 'block'; trigger.src = 'open.png'; removeClass(linkObj,'closed'); addClass(linkObj,'opened'); } else if ( hasClass(linkObj,'opened') ) { summary.style.display = 'block'; content.style.display = 'none'; trigger.src = 'closed.png'; removeClass(linkObj,'opened'); addClass(linkObj,'closed'); } return false; } </script> <div id="top"> <div id="titlearea"> <table cellspacing="0" cellpadding="0"> <tbody> <tr style="height: 56px;"> <td style="padding-left: 0.5em;"> <div id="projectname">linbox</div> </td> </tr> </tbody> </table> </div> <div id="navrow1" class="tabs"> <ul class="tablist"> <li><a href="index.html"><span>Main Page</span></a></li> <li><a href="pages.html"><span>Related Pages</span></a></li> <li><a href="modules.html"><span>Modules</span></a></li> <li><a href="namespaces.html"><span>Namespaces</span></a></li> <li class="current"><a href="annotated.html"><span>Data Structures</span></a></li> <li><a href="files.html"><span>Files</span></a></li> <li><a href="dirs.html"><span>Directories</span></a></li> <li><a href="examples.html"><span>Examples</span></a></li> </ul> </div> <div id="navrow2" class="tabs2"> <ul class="tablist"> <li><a href="annotated.html"><span>Data Structures</span></a></li> <li><a href="hierarchy.html"><span>Class Hierarchy</span></a></li> <li><a href="functions.html"><span>Data Fields</span></a></li> </ul> </div> <div id="nav-path" class="navpath"> <ul> <li class="navelem"><a class="el" href="namespace_lin_box.html">LinBox</a> </li> <li class="navelem"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html">LQUPMatrix</a> </li> </ul> </div> </div> <div class="header"> <div class="summary"> <a href="#pub-methods">Public Member Functions</a> </div> <div class="headertitle"> <div class="title">LQUPMatrix< Field > Class Template Reference</div> </div> </div> <div class="contents"> <!-- doxytag: class="LinBox::LQUPMatrix" --> <p>LQUP factorisation. <a href="class_lin_box_1_1_l_q_u_p_matrix.html#details">More...</a></p> <p><code>#include <factorized-matrix.h></code></p> <table class="memberdecls"> <tr><td colspan="2"><h2><a name="pub-methods"></a> Public Member Functions</h2></td></tr> <tr><td class="memItemLeft" align="right" valign="top"><a class="anchor" id="ace8099aa462c72ca4c77e3055ef9d871"></a><!-- doxytag: member="LinBox::LQUPMatrix::LQUPMatrix" ref="ace8099aa462c72ca4c77e3055ef9d871" args="(const Field &F, const BlasMatrix< Element > &A)" -->  </td><td class="memItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#ace8099aa462c72ca4c77e3055ef9d871">LQUPMatrix</a> (const <a class="el" href="class_lin_box_1_1_modular_3_01uint32__t_01_4.html">Field</a> &F, const <a class="el" href="class_lin_box_1_1_blas_matrix.html">BlasMatrix</a>< Element > &A)</td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">Contruction of LQUP factorization of A (making a copy of A) <br/></td></tr> <tr><td class="memItemLeft" align="right" valign="top"><a class="anchor" id="a05f33626a0cd74d38cfacbe3d8685981"></a><!-- doxytag: member="LinBox::LQUPMatrix::LQUPMatrix" ref="a05f33626a0cd74d38cfacbe3d8685981" args="(const Field &F, BlasMatrix< Element > &A)" -->  </td><td class="memItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#a05f33626a0cd74d38cfacbe3d8685981">LQUPMatrix</a> (const <a class="el" href="class_lin_box_1_1_modular_3_01uint32__t_01_4.html">Field</a> &F, <a class="el" href="class_lin_box_1_1_blas_matrix.html">BlasMatrix</a>< Element > &A)</td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">Contruction of LQUP factorization of A (in-place in A) <br/></td></tr> <tr><td class="memItemLeft" align="right" valign="top"><a class="anchor" id="a0749c939e0aa50828847bd400e5ab42f"></a><!-- doxytag: member="LinBox::LQUPMatrix::LQUPMatrix" ref="a0749c939e0aa50828847bd400e5ab42f" args="(const BlasBlackbox< Field > &A)" -->  </td><td class="memItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#a0749c939e0aa50828847bd400e5ab42f">LQUPMatrix</a> (const <a class="el" href="class_lin_box_1_1_blas_blackbox.html">BlasBlackbox</a>< <a class="el" href="class_lin_box_1_1_modular_3_01uint32__t_01_4.html">Field</a> > &A)</td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">Contruction of LQUP factorization of A (making a copy of A) <br/></td></tr> <tr><td class="memItemLeft" align="right" valign="top"><a class="anchor" id="aacdece45f9e769fc6deb52b0263a5b4f"></a><!-- doxytag: member="LinBox::LQUPMatrix::LQUPMatrix" ref="aacdece45f9e769fc6deb52b0263a5b4f" args="(BlasBlackbox< Field > &A)" -->  </td><td class="memItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#aacdece45f9e769fc6deb52b0263a5b4f">LQUPMatrix</a> (<a class="el" href="class_lin_box_1_1_blas_blackbox.html">BlasBlackbox</a>< <a class="el" href="class_lin_box_1_1_modular_3_01uint32__t_01_4.html">Field</a> > &A)</td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">Contruction of LQUP factorization of A (in-place in A) <br/></td></tr> <tr><td class="memItemLeft" align="right" valign="top"> </td><td class="memItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#a199ad8081b66b01d769c5047f92f2a1d">LQUPMatrix</a> (const <a class="el" href="class_lin_box_1_1_modular_3_01uint32__t_01_4.html">Field</a> &F, const <a class="el" href="class_lin_box_1_1_blas_matrix.html">BlasMatrix</a>< Element > &A, <a class="el" href="class_lin_box_1_1_blas_permutation.html">BlasPermutation</a>< size_t > &P, <a class="el" href="class_lin_box_1_1_blas_permutation.html">BlasPermutation</a>< size_t > &Q)</td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">Contruction of LQUP factorization of A (making a copy of A). <a href="#a199ad8081b66b01d769c5047f92f2a1d"></a><br/></td></tr> <tr><td class="memItemLeft" align="right" valign="top"> </td><td class="memItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#a7d3f799e47f858367ac448272de524a7">LQUPMatrix</a> (const <a class="el" href="class_lin_box_1_1_modular_3_01uint32__t_01_4.html">Field</a> &F, <a class="el" href="class_lin_box_1_1_blas_matrix.html">BlasMatrix</a>< Element > &A, <a class="el" href="class_lin_box_1_1_blas_permutation.html">BlasPermutation</a>< size_t > &P, <a class="el" href="class_lin_box_1_1_blas_permutation.html">BlasPermutation</a>< size_t > &Q)</td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">Contruction of LQUP factorization of A (in-place in A). <a href="#a7d3f799e47f858367ac448272de524a7"></a><br/></td></tr> <tr><td class="memItemLeft" align="right" valign="top"> </td><td class="memItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#ab58544c11a196ff59a72d59f3f186347">LQUPMatrix</a> (const <a class="el" href="class_lin_box_1_1_blas_blackbox.html">BlasBlackbox</a>< <a class="el" href="class_lin_box_1_1_modular_3_01uint32__t_01_4.html">Field</a> > &A, <a class="el" href="class_lin_box_1_1_blas_permutation.html">BlasPermutation</a>< size_t > &P, <a class="el" href="class_lin_box_1_1_blas_permutation.html">BlasPermutation</a>< size_t > &Q)</td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">Contruction of LQUP factorization of A (making a copy of A). <a href="#ab58544c11a196ff59a72d59f3f186347"></a><br/></td></tr> <tr><td class="memItemLeft" align="right" valign="top"> </td><td class="memItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#aa244b808d07cee4eba16953e42e0583f">LQUPMatrix</a> (<a class="el" href="class_lin_box_1_1_blas_blackbox.html">BlasBlackbox</a>< <a class="el" href="class_lin_box_1_1_modular_3_01uint32__t_01_4.html">Field</a> > &A, <a class="el" href="class_lin_box_1_1_blas_permutation.html">BlasPermutation</a>< size_t > &P, <a class="el" href="class_lin_box_1_1_blas_permutation.html">BlasPermutation</a>< size_t > &Q)</td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">Contruction of LQUP factorization of A (in-place in A). <a href="#aa244b808d07cee4eba16953e42e0583f"></a><br/></td></tr> <tr><td class="memItemLeft" align="right" valign="top"><a class="anchor" id="af8008f5fa004411c9ed8cafaa4262472"></a><!-- doxytag: member="LinBox::LQUPMatrix::~LQUPMatrix" ref="af8008f5fa004411c9ed8cafaa4262472" args="()" -->  </td><td class="memItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#af8008f5fa004411c9ed8cafaa4262472">~LQUPMatrix</a> ()</td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">destructor. <br/></td></tr> <tr><td class="memItemLeft" align="right" valign="top"><a class="anchor" id="a05eca9e10246b4737dd67c2a6cb58592"></a><!-- doxytag: member="LinBox::LQUPMatrix::field" ref="a05eca9e10246b4737dd67c2a6cb58592" args="()" --> <a class="el" href="class_lin_box_1_1_modular_3_01uint32__t_01_4.html">Field</a> & </td><td class="memItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#a05eca9e10246b4737dd67c2a6cb58592">field</a> ()</td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">get the field on which the factorization is done <br/></td></tr> <tr><td class="memItemLeft" align="right" valign="top"><a class="anchor" id="a0fee90e6b8ef7fcd8d76a7f3e1f268e6"></a><!-- doxytag: member="LinBox::LQUPMatrix::rowdim" ref="a0fee90e6b8ef7fcd8d76a7f3e1f268e6" args="() const " --> size_t </td><td class="memItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#a0fee90e6b8ef7fcd8d76a7f3e1f268e6">rowdim</a> () const </td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">get row dimension <br/></td></tr> <tr><td class="memItemLeft" align="right" valign="top"><a class="anchor" id="a32edb490d3597f5553152d14b102e227"></a><!-- doxytag: member="LinBox::LQUPMatrix::coldim" ref="a32edb490d3597f5553152d14b102e227" args="() const " --> size_t </td><td class="memItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#a32edb490d3597f5553152d14b102e227">coldim</a> () const </td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">get column dimension <br/></td></tr> <tr><td class="memItemLeft" align="right" valign="top"><a class="anchor" id="a8161e33322b767f76d94da6b4db30034"></a><!-- doxytag: member="LinBox::LQUPMatrix::getrank" ref="a8161e33322b767f76d94da6b4db30034" args="() const " --> size_t </td><td class="memItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#a8161e33322b767f76d94da6b4db30034">getrank</a> () const </td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">get the rank of matrix <br/></td></tr> <tr><td class="memItemLeft" align="right" valign="top">const <a class="el" href="class_lin_box_1_1_blas_permutation.html">BlasPermutation</a>< size_t > & </td><td class="memItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#a0161f9158ac37f2553e084ef2aa7ca6c">getP</a> () const </td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">get the permutation P. <a href="#a0161f9158ac37f2553e084ef2aa7ca6c"></a><br/></td></tr> <tr><td class="memItemLeft" align="right" valign="top"><a class="el" href="class_lin_box_1_1_blas_permutation.html">BlasPermutation</a>< size_t > & </td><td class="memItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#ac9ce865f44bee3b064fd08112e8c5870">getP</a> (<a class="el" href="class_lin_box_1_1_blas_permutation.html">BlasPermutation</a>< size_t > &P) const </td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">get the permutation P. <a href="#ac9ce865f44bee3b064fd08112e8c5870"></a><br/></td></tr> <tr><td class="memItemLeft" align="right" valign="top">const <a class="el" href="class_lin_box_1_1_blas_permutation.html">BlasPermutation</a>< size_t > & </td><td class="memItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#a5002c0c098bdc55d1ad298fed01c7d76">getQ</a> () const </td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">Get the <em>transpose</em> of the permutation <code>Q</code>. <a href="#a5002c0c098bdc55d1ad298fed01c7d76"></a><br/></td></tr> <tr><td class="memItemLeft" align="right" valign="top"><a class="el" href="class_lin_box_1_1_blas_permutation.html">BlasPermutation</a>< size_t > & </td><td class="memItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#aba559e3bebe2ecc86d3f1602cf374323">getQ</a> (<a class="el" href="class_lin_box_1_1_blas_permutation.html">BlasPermutation</a>< size_t > &Qt) const </td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">get the permutation Qt. <a href="#aba559e3bebe2ecc86d3f1602cf374323"></a><br/></td></tr> <tr><td class="memItemLeft" align="right" valign="top"><a class="el" href="class_lin_box_1_1_triangular_blas_matrix.html">TriangularBlasMatrix</a>< Element > & </td><td class="memItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#abe14dc4d0857287c11c06ec1ff10f8e8">getL</a> (<a class="el" href="class_lin_box_1_1_triangular_blas_matrix.html">TriangularBlasMatrix</a>< Element > &L, bool _QLUP=false) const </td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">get the Matrix <code>L</code>. <a href="#abe14dc4d0857287c11c06ec1ff10f8e8"></a><br/></td></tr> <tr><td class="memItemLeft" align="right" valign="top"><a class="el" href="class_lin_box_1_1_triangular_blas_matrix.html">TriangularBlasMatrix</a>< Element > & </td><td class="memItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#a9285cdfbe78257552e71eed2c319ec3c">getU</a> (<a class="el" href="class_lin_box_1_1_triangular_blas_matrix.html">TriangularBlasMatrix</a>< Element > &U) const </td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">get the matrix <code>U</code>. <a href="#a9285cdfbe78257552e71eed2c319ec3c"></a><br/></td></tr> <tr><td class="memItemLeft" align="right" valign="top"><a class="el" href="class_lin_box_1_1_blas_matrix.html">BlasMatrix</a>< Element > & </td><td class="memItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#a9242c61339a9675bb14e3f5e35436199">getS</a> (<a class="el" href="class_lin_box_1_1_blas_matrix.html">BlasMatrix</a>< Element > &S) const </td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">get the matrix S. <a href="#a9242c61339a9675bb14e3f5e35436199"></a><br/></td></tr> <tr><td colspan="2"><div class="groupHeader"></div></td></tr> <tr><td class="memTemplParams" colspan="2"><a class="anchor" id="ad6a16839616999ffdc3e5bc153e6caaf"></a><!-- doxytag: member="LinBox::LQUPMatrix::left_solve" ref="ad6a16839616999ffdc3e5bc153e6caaf" args="(Operand &X, const Operand &B) const " --> template<class Operand > </td></tr> <tr><td class="memTemplItemLeft" align="right" valign="top">Operand & </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#ad6a16839616999ffdc3e5bc153e6caaf">left_solve</a> (Operand &X, const Operand &B) const </td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">Solvers with matrices or vectors Operand can be a <a class="el" href="class_lin_box_1_1_blas_matrix.html">BlasMatrix<Element></a> or a std::vector<Element> <br/></td></tr> <tr><td class="memTemplParams" colspan="2"><a class="anchor" id="a31f31fdce8a0bcb3b192973cfa525562"></a><!-- doxytag: member="LinBox::LQUPMatrix::left_solve" ref="a31f31fdce8a0bcb3b192973cfa525562" args="(Operand &B) const " --> template<class Operand > </td></tr> <tr><td class="memTemplItemLeft" align="right" valign="top">Operand & </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#a31f31fdce8a0bcb3b192973cfa525562">left_solve</a> (Operand &B) const </td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">Solvers with matrices or vectors Operand can be a <a class="el" href="class_lin_box_1_1_blas_matrix.html">BlasMatrix<Element></a> or a std::vector<Element> <br/></td></tr> <tr><td class="memTemplParams" colspan="2"><a class="anchor" id="ab8fb1f3b7812db60492a59b6ef90d9a2"></a><!-- doxytag: member="LinBox::LQUPMatrix::right_solve" ref="ab8fb1f3b7812db60492a59b6ef90d9a2" args="(Operand &X, const Operand &B) const " --> template<class Operand > </td></tr> <tr><td class="memTemplItemLeft" align="right" valign="top">Operand & </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#ab8fb1f3b7812db60492a59b6ef90d9a2">right_solve</a> (Operand &X, const Operand &B) const </td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">Solvers with matrices or vectors Operand can be a <a class="el" href="class_lin_box_1_1_blas_matrix.html">BlasMatrix<Element></a> or a std::vector<Element> <br/></td></tr> <tr><td class="memTemplParams" colspan="2"><a class="anchor" id="a9b6a1cb332d7fe9d3866a593ebba1ba8"></a><!-- doxytag: member="LinBox::LQUPMatrix::right_solve" ref="a9b6a1cb332d7fe9d3866a593ebba1ba8" args="(Operand &B) const " --> template<class Operand > </td></tr> <tr><td class="memTemplItemLeft" align="right" valign="top">Operand & </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#a9b6a1cb332d7fe9d3866a593ebba1ba8">right_solve</a> (Operand &B) const </td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">Solvers with matrices or vectors Operand can be a <a class="el" href="class_lin_box_1_1_blas_matrix.html">BlasMatrix<Element></a> or a std::vector<Element> <br/></td></tr> <tr><td class="memTemplParams" colspan="2"><a class="anchor" id="a32c830804b303f098737c1330d7f0fc1"></a><!-- doxytag: member="LinBox::LQUPMatrix::left_Lsolve" ref="a32c830804b303f098737c1330d7f0fc1" args="(Operand &X, const Operand &B) const " --> template<class Operand > </td></tr> <tr><td class="memTemplItemLeft" align="right" valign="top">Operand & </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#a32c830804b303f098737c1330d7f0fc1">left_Lsolve</a> (Operand &X, const Operand &B) const </td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">Solvers with matrices or vectors Operand can be a <a class="el" href="class_lin_box_1_1_blas_matrix.html">BlasMatrix<Element></a> or a std::vector<Element> <br/></td></tr> <tr><td class="memTemplParams" colspan="2"><a class="anchor" id="a9bb58b0d7925b2bde6e67a16a7d73ddc"></a><!-- doxytag: member="LinBox::LQUPMatrix::left_Lsolve" ref="a9bb58b0d7925b2bde6e67a16a7d73ddc" args="(Operand &B) const " --> template<class Operand > </td></tr> <tr><td class="memTemplItemLeft" align="right" valign="top">Operand & </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#a9bb58b0d7925b2bde6e67a16a7d73ddc">left_Lsolve</a> (Operand &B) const </td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">Solvers with matrices or vectors Operand can be a <a class="el" href="class_lin_box_1_1_blas_matrix.html">BlasMatrix<Element></a> or a std::vector<Element> <br/></td></tr> <tr><td class="memTemplParams" colspan="2"><a class="anchor" id="a8f88e31c93e13e5038c778faefd0a11d"></a><!-- doxytag: member="LinBox::LQUPMatrix::right_Lsolve" ref="a8f88e31c93e13e5038c778faefd0a11d" args="(Operand &X, const Operand &B) const " --> template<class Operand > </td></tr> <tr><td class="memTemplItemLeft" align="right" valign="top">Operand & </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#a8f88e31c93e13e5038c778faefd0a11d">right_Lsolve</a> (Operand &X, const Operand &B) const </td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">Solvers with matrices or vectors Operand can be a <a class="el" href="class_lin_box_1_1_blas_matrix.html">BlasMatrix<Element></a> or a std::vector<Element> <br/></td></tr> <tr><td class="memTemplParams" colspan="2"><a class="anchor" id="a71c2fe6d0ab31798437e2cafae110dea"></a><!-- doxytag: member="LinBox::LQUPMatrix::right_Lsolve" ref="a71c2fe6d0ab31798437e2cafae110dea" args="(Operand &B) const " --> template<class Operand > </td></tr> <tr><td class="memTemplItemLeft" align="right" valign="top">Operand & </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#a71c2fe6d0ab31798437e2cafae110dea">right_Lsolve</a> (Operand &B) const </td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">Solvers with matrices or vectors Operand can be a <a class="el" href="class_lin_box_1_1_blas_matrix.html">BlasMatrix<Element></a> or a std::vector<Element> <br/></td></tr> <tr><td class="memTemplParams" colspan="2"><a class="anchor" id="a1df1bf8335ae23f18f98d028fde7f2f7"></a><!-- doxytag: member="LinBox::LQUPMatrix::left_Usolve" ref="a1df1bf8335ae23f18f98d028fde7f2f7" args="(Operand &X, const Operand &B) const " --> template<class Operand > </td></tr> <tr><td class="memTemplItemLeft" align="right" valign="top">Operand & </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#a1df1bf8335ae23f18f98d028fde7f2f7">left_Usolve</a> (Operand &X, const Operand &B) const </td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">Solvers with matrices or vectors Operand can be a <a class="el" href="class_lin_box_1_1_blas_matrix.html">BlasMatrix<Element></a> or a std::vector<Element> <br/></td></tr> <tr><td class="memTemplParams" colspan="2"><a class="anchor" id="aa8e267c9082b4e64e0b30376b3cfd34a"></a><!-- doxytag: member="LinBox::LQUPMatrix::rleft_Usolve" ref="aa8e267c9082b4e64e0b30376b3cfd34a" args="(Operand &B) const " --> template<class Operand > </td></tr> <tr><td class="memTemplItemLeft" align="right" valign="top">Operand & </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#aa8e267c9082b4e64e0b30376b3cfd34a">rleft_Usolve</a> (Operand &B) const </td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">Solvers with matrices or vectors Operand can be a <a class="el" href="class_lin_box_1_1_blas_matrix.html">BlasMatrix<Element></a> or a std::vector<Element> <br/></td></tr> <tr><td class="memTemplParams" colspan="2"><a class="anchor" id="aa1c481aaade4db3eebd2d061fc30982b"></a><!-- doxytag: member="LinBox::LQUPMatrix::right_Usolve" ref="aa1c481aaade4db3eebd2d061fc30982b" args="(Operand &X, const Operand &B) const " --> template<class Operand > </td></tr> <tr><td class="memTemplItemLeft" align="right" valign="top">Operand & </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#aa1c481aaade4db3eebd2d061fc30982b">right_Usolve</a> (Operand &X, const Operand &B) const </td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">Solvers with matrices or vectors Operand can be a <a class="el" href="class_lin_box_1_1_blas_matrix.html">BlasMatrix<Element></a> or a std::vector<Element> <br/></td></tr> <tr><td class="memTemplParams" colspan="2"><a class="anchor" id="a3cc40de05940887f9934490ebb7a2634"></a><!-- doxytag: member="LinBox::LQUPMatrix::right_Usolve" ref="a3cc40de05940887f9934490ebb7a2634" args="(Operand &B) const " --> template<class Operand > </td></tr> <tr><td class="memTemplItemLeft" align="right" valign="top">Operand & </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#a3cc40de05940887f9934490ebb7a2634">right_Usolve</a> (Operand &B) const </td></tr> <tr><td class="mdescLeft"> </td><td class="mdescRight">Solvers with matrices or vectors Operand can be a <a class="el" href="class_lin_box_1_1_blas_matrix.html">BlasMatrix<Element></a> or a std::vector<Element> <br/></td></tr> </table> <hr/><a name="details" id="details"></a><h2>Detailed Description</h2> <div class="textblock"><h3>template<class Field><br/> class LinBox::LQUPMatrix< Field ></h3> <p>LQUP factorisation. </p> <p>This is a class to ease the use LU factorisation (see LUdivine in FFPACK)</p> <p>The factorisation is <img class="formulaInl" alt="$ A = L Q U P $" src="form_56.png"/> with <code>L</code> lower unit triangular, <code>U</code> upper non-unit triangular, <code>P</code> and <code>Q</code> permutations.</p> <p>There are two kind of contructors (with and without permutations) and they build a LQUP factorisation of a BlasMatrix/BlasBlackbox on a finite field. There are methods for retrieving L,Q,U and P matrices and methods for solving systems. </p> </div><hr/><h2>Constructor & Destructor Documentation</h2> <a class="anchor" id="a199ad8081b66b01d769c5047f92f2a1d"></a><!-- doxytag: member="LinBox::LQUPMatrix::LQUPMatrix" ref="a199ad8081b66b01d769c5047f92f2a1d" args="(const Field &F, const BlasMatrix< Element > &A, BlasPermutation< size_t > &P, BlasPermutation< size_t > &Q)" --> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html">LQUPMatrix</a> </td> <td>(</td> <td class="paramtype">const <a class="el" href="class_lin_box_1_1_modular_3_01uint32__t_01_4.html">Field</a> & </td> <td class="paramname"><em>F</em>, </td> </tr> <tr> <td class="paramkey"></td> <td></td> <td class="paramtype">const <a class="el" href="class_lin_box_1_1_blas_matrix.html">BlasMatrix</a>< Element > & </td> <td class="paramname"><em>A</em>, </td> </tr> <tr> <td class="paramkey"></td> <td></td> <td class="paramtype"><a class="el" href="class_lin_box_1_1_blas_permutation.html">BlasPermutation</a>< size_t > & </td> <td class="paramname"><em>P</em>, </td> </tr> <tr> <td class="paramkey"></td> <td></td> <td class="paramtype"><a class="el" href="class_lin_box_1_1_blas_permutation.html">BlasPermutation</a>< size_t > & </td> <td class="paramname"><em>Q</em> </td> </tr> <tr> <td></td> <td>)</td> <td></td><td><code> [inline]</code></td> </tr> </table> </div> <div class="memdoc"> <p>Contruction of LQUP factorization of A (making a copy of A). </p> <p>P and Q are arguments ! </p> </div> </div> <a class="anchor" id="a7d3f799e47f858367ac448272de524a7"></a><!-- doxytag: member="LinBox::LQUPMatrix::LQUPMatrix" ref="a7d3f799e47f858367ac448272de524a7" args="(const Field &F, BlasMatrix< Element > &A, BlasPermutation< size_t > &P, BlasPermutation< size_t > &Q)" --> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html">LQUPMatrix</a> </td> <td>(</td> <td class="paramtype">const <a class="el" href="class_lin_box_1_1_modular_3_01uint32__t_01_4.html">Field</a> & </td> <td class="paramname"><em>F</em>, </td> </tr> <tr> <td class="paramkey"></td> <td></td> <td class="paramtype"><a class="el" href="class_lin_box_1_1_blas_matrix.html">BlasMatrix</a>< Element > & </td> <td class="paramname"><em>A</em>, </td> </tr> <tr> <td class="paramkey"></td> <td></td> <td class="paramtype"><a class="el" href="class_lin_box_1_1_blas_permutation.html">BlasPermutation</a>< size_t > & </td> <td class="paramname"><em>P</em>, </td> </tr> <tr> <td class="paramkey"></td> <td></td> <td class="paramtype"><a class="el" href="class_lin_box_1_1_blas_permutation.html">BlasPermutation</a>< size_t > & </td> <td class="paramname"><em>Q</em> </td> </tr> <tr> <td></td> <td>)</td> <td></td><td><code> [inline]</code></td> </tr> </table> </div> <div class="memdoc"> <p>Contruction of LQUP factorization of A (in-place in A). </p> <p>P and Q are arguments ! </p> </div> </div> <a class="anchor" id="ab58544c11a196ff59a72d59f3f186347"></a><!-- doxytag: member="LinBox::LQUPMatrix::LQUPMatrix" ref="ab58544c11a196ff59a72d59f3f186347" args="(const BlasBlackbox< Field > &A, BlasPermutation< size_t > &P, BlasPermutation< size_t > &Q)" --> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html">LQUPMatrix</a> </td> <td>(</td> <td class="paramtype">const <a class="el" href="class_lin_box_1_1_blas_blackbox.html">BlasBlackbox</a>< <a class="el" href="class_lin_box_1_1_modular_3_01uint32__t_01_4.html">Field</a> > & </td> <td class="paramname"><em>A</em>, </td> </tr> <tr> <td class="paramkey"></td> <td></td> <td class="paramtype"><a class="el" href="class_lin_box_1_1_blas_permutation.html">BlasPermutation</a>< size_t > & </td> <td class="paramname"><em>P</em>, </td> </tr> <tr> <td class="paramkey"></td> <td></td> <td class="paramtype"><a class="el" href="class_lin_box_1_1_blas_permutation.html">BlasPermutation</a>< size_t > & </td> <td class="paramname"><em>Q</em> </td> </tr> <tr> <td></td> <td>)</td> <td></td><td><code> [inline]</code></td> </tr> </table> </div> <div class="memdoc"> <p>Contruction of LQUP factorization of A (making a copy of A). </p> <p>P and Q are arguments ! </p> </div> </div> <a class="anchor" id="aa244b808d07cee4eba16953e42e0583f"></a><!-- doxytag: member="LinBox::LQUPMatrix::LQUPMatrix" ref="aa244b808d07cee4eba16953e42e0583f" args="(BlasBlackbox< Field > &A, BlasPermutation< size_t > &P, BlasPermutation< size_t > &Q)" --> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname"><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html">LQUPMatrix</a> </td> <td>(</td> <td class="paramtype"><a class="el" href="class_lin_box_1_1_blas_blackbox.html">BlasBlackbox</a>< <a class="el" href="class_lin_box_1_1_modular_3_01uint32__t_01_4.html">Field</a> > & </td> <td class="paramname"><em>A</em>, </td> </tr> <tr> <td class="paramkey"></td> <td></td> <td class="paramtype"><a class="el" href="class_lin_box_1_1_blas_permutation.html">BlasPermutation</a>< size_t > & </td> <td class="paramname"><em>P</em>, </td> </tr> <tr> <td class="paramkey"></td> <td></td> <td class="paramtype"><a class="el" href="class_lin_box_1_1_blas_permutation.html">BlasPermutation</a>< size_t > & </td> <td class="paramname"><em>Q</em> </td> </tr> <tr> <td></td> <td>)</td> <td></td><td><code> [inline]</code></td> </tr> </table> </div> <div class="memdoc"> <p>Contruction of LQUP factorization of A (in-place in A). </p> <p>P and Q are arguments ! </p> </div> </div> <hr/><h2>Member Function Documentation</h2> <a class="anchor" id="a0161f9158ac37f2553e084ef2aa7ca6c"></a><!-- doxytag: member="LinBox::LQUPMatrix::getP" ref="a0161f9158ac37f2553e084ef2aa7ca6c" args="() const " --> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">const <a class="el" href="class_lin_box_1_1_blas_permutation.html">BlasPermutation</a><size_t>& getP </td> <td>(</td> <td class="paramname"></td><td>)</td> <td> const<code> [inline]</code></td> </tr> </table> </div> <div class="memdoc"> <p>get the permutation P. </p> <p>(no copy) </p> </div> </div> <a class="anchor" id="ac9ce865f44bee3b064fd08112e8c5870"></a><!-- doxytag: member="LinBox::LQUPMatrix::getP" ref="ac9ce865f44bee3b064fd08112e8c5870" args="(BlasPermutation< size_t > &P) const " --> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname"><a class="el" href="class_lin_box_1_1_blas_permutation.html">BlasPermutation</a><size_t>& getP </td> <td>(</td> <td class="paramtype"><a class="el" href="class_lin_box_1_1_blas_permutation.html">BlasPermutation</a>< size_t > & </td> <td class="paramname"><em>P</em></td><td>)</td> <td> const<code> [inline]</code></td> </tr> </table> </div> <div class="memdoc"> <p>get the permutation P. </p> <p>(copy) </p> </div> </div> <a class="anchor" id="a5002c0c098bdc55d1ad298fed01c7d76"></a><!-- doxytag: member="LinBox::LQUPMatrix::getQ" ref="a5002c0c098bdc55d1ad298fed01c7d76" args="() const " --> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">const <a class="el" href="class_lin_box_1_1_blas_permutation.html">BlasPermutation</a><size_t>& getQ </td> <td>(</td> <td class="paramname"></td><td>)</td> <td> const<code> [inline]</code></td> </tr> </table> </div> <div class="memdoc"> <p>Get the <em>transpose</em> of the permutation <code>Q</code>. </p> <dl class="warning"><dt><b>Warning:</b></dt><dd>This does not return <code>Q</code> itself! (because it is more difficult to compute) If needed, <code>Q</code> can be obtained as a <code><a class="el" href="class_lin_box_1_1_transposed_blas_matrix.html" title="TransposedBlasMatrix.">TransposedBlasMatrix</a></code> from the return value One reason this confusion exists is that left-multiplying by a permuation matrix corresponds to a row permuation <img class="formulaInl" alt="$\pi \in S_n$" src="form_57.png"/>, while right-multiplying by the same matrix corresponds to the inverse column permutation <img class="formulaInl" alt="$\pi^{-1}$" src="form_58.png"/>! Usually this is handled intelligently (eg by <code>applyP</code>) but you must be careful with <code><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#a5002c0c098bdc55d1ad298fed01c7d76" title="Get the transpose of the permutation Q.">getQ()</a></code>. </dd></dl> </div> </div> <a class="anchor" id="aba559e3bebe2ecc86d3f1602cf374323"></a><!-- doxytag: member="LinBox::LQUPMatrix::getQ" ref="aba559e3bebe2ecc86d3f1602cf374323" args="(BlasPermutation< size_t > &Qt) const " --> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname"><a class="el" href="class_lin_box_1_1_blas_permutation.html">BlasPermutation</a><size_t>& getQ </td> <td>(</td> <td class="paramtype"><a class="el" href="class_lin_box_1_1_blas_permutation.html">BlasPermutation</a>< size_t > & </td> <td class="paramname"><em>Qt</em></td><td>)</td> <td> const<code> [inline]</code></td> </tr> </table> </div> <div class="memdoc"> <p>get the permutation Qt. </p> <p>(copy) </p> <dl class="warning"><dt><b>Warning:</b></dt><dd>see <code><a class="el" href="class_lin_box_1_1_l_q_u_p_matrix.html#a5002c0c098bdc55d1ad298fed01c7d76" title="Get the transpose of the permutation Q.">LQUPMatrix::getQ()</a></code> </dd></dl> </div> </div> <a class="anchor" id="abe14dc4d0857287c11c06ec1ff10f8e8"></a><!-- doxytag: member="LinBox::LQUPMatrix::getL" ref="abe14dc4d0857287c11c06ec1ff10f8e8" args="(TriangularBlasMatrix< Element > &L, bool _QLUP=false) const " --> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname"><a class="el" href="class_lin_box_1_1_triangular_blas_matrix.html">TriangularBlasMatrix</a>< typename Field::Element > & getL </td> <td>(</td> <td class="paramtype"><a class="el" href="class_lin_box_1_1_triangular_blas_matrix.html">TriangularBlasMatrix</a>< Element > & </td> <td class="paramname"><em>L</em>, </td> </tr> <tr> <td class="paramkey"></td> <td></td> <td class="paramtype">bool </td> <td class="paramname"><em>_QLUP</em> = <code>false</code> </td> </tr> <tr> <td></td> <td>)</td> <td></td><td> const<code> [inline]</code></td> </tr> </table> </div> <div class="memdoc"> <p>get the Matrix <code>L</code>. </p> <dl><dt><b>Parameters:</b></dt><dd> <table class="params"> <tr><td class="paramdir">[out]</td><td class="paramname">L</td><td></td></tr> <tr><td class="paramdir"></td><td class="paramname">_QLUP</td><td>if true then <code>L</code> form <code>QLUP</code> decomposition, else <code>L</code> is form <code>LQUP</code> decomposition. </td></tr> </table> </dd> </dl> <dl class="pre"><dt><b>Precondition:</b></dt><dd><code>L</code> has unit diagonal </dd></dl> </div> </div> <a class="anchor" id="a9285cdfbe78257552e71eed2c319ec3c"></a><!-- doxytag: member="LinBox::LQUPMatrix::getU" ref="a9285cdfbe78257552e71eed2c319ec3c" args="(TriangularBlasMatrix< Element > &U) const " --> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname"><a class="el" href="class_lin_box_1_1_triangular_blas_matrix.html">TriangularBlasMatrix</a>< typename Field::Element > & getU </td> <td>(</td> <td class="paramtype"><a class="el" href="class_lin_box_1_1_triangular_blas_matrix.html">TriangularBlasMatrix</a>< Element > & </td> <td class="paramname"><em>U</em></td><td>)</td> <td> const<code> [inline]</code></td> </tr> </table> </div> <div class="memdoc"> <p>get the matrix <code>U</code>. </p> <dl class="pre"><dt><b>Precondition:</b></dt><dd><code>U</code> has non-unit diagonal </dd></dl> </div> </div> <a class="anchor" id="a9242c61339a9675bb14e3f5e35436199"></a><!-- doxytag: member="LinBox::LQUPMatrix::getS" ref="a9242c61339a9675bb14e3f5e35436199" args="(BlasMatrix< Element > &S) const " --> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname"><a class="el" href="class_lin_box_1_1_blas_matrix.html">BlasMatrix</a>< typename Field::Element > & getS </td> <td>(</td> <td class="paramtype"><a class="el" href="class_lin_box_1_1_blas_matrix.html">BlasMatrix</a>< Element > & </td> <td class="paramname"><em>S</em></td><td>)</td> <td> const<code> [inline]</code></td> </tr> </table> </div> <div class="memdoc"> <p>get the matrix S. </p> <p>from the LSP factorization of A deduced from LQUP) </p> </div> </div> <hr/>The documentation for this class was generated from the following files:<ul> <li>factorized-matrix.h</li> <li>factorized-matrix.inl</li> </ul> </div> <hr class="footer"/><address class="footer"><small>Generated on Tue Aug 30 2011 for linbox by  <a href="http://www.doxygen.org/index.html"> <img class="footer" src="doxygen.png" alt="doxygen"/></a> 1.7.4 </small></address> </body> </html>