<!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"/> <meta http-equiv="X-UA-Compatible" content="IE=9"/> <meta name="generator" content="Doxygen 1.8.5"/> <title>Eigen: Reductions, visitors and broadcasting</title> <link href="tabs.css" rel="stylesheet" type="text/css"/> <script type="text/javascript" src="jquery.js"></script> <script type="text/javascript" src="dynsections.js"></script> <link href="navtree.css" rel="stylesheet" type="text/css"/> <script type="text/javascript" src="resize.js"></script> <script type="text/javascript" src="navtree.js"></script> <script type="text/javascript"> $(document).ready(initResizable); $(window).load(resizeHeight); </script> <link href="search/search.css" rel="stylesheet" type="text/css"/> <script type="text/javascript" src="search/search.js"></script> <script type="text/javascript"> $(document).ready(function() { searchBox.OnSelectItem(0); }); </script> <link href="doxygen.css" rel="stylesheet" type="text/css" /> <link href="eigendoxy.css" rel="stylesheet" type="text/css"> <!-- --> <script type="text/javascript" src="eigen_navtree_hacks.js"></script> <!-- <script type="text/javascript"> --> <!-- </script> --> </head> <body> <div id="top"><!-- do not remove this div, it is closed by doxygen! --> <!-- <a name="top"></a> --> <div id="titlearea"> <table cellspacing="0" cellpadding="0"> <tbody> <tr style="height: 56px;"> <td id="projectlogo"><img alt="Logo" src="Eigen_Silly_Professor_64x64.png"/></td> <td style="padding-left: 0.5em;"> <div id="projectname"><a href="http://eigen.tuxfamily.org">Eigen</a>  <span id="projectnumber">3.2.0</span> </div> </td> <td> <div id="MSearchBox" class="MSearchBoxInactive"> <span class="left"> <img id="MSearchSelect" src="search/mag_sel.png" onmouseover="return searchBox.OnSearchSelectShow()" onmouseout="return searchBox.OnSearchSelectHide()" alt=""/> <input type="text" id="MSearchField" value="Search" accesskey="S" onfocus="searchBox.OnSearchFieldFocus(true)" onblur="searchBox.OnSearchFieldFocus(false)" onkeyup="searchBox.OnSearchFieldChange(event)"/> </span><span class="right"> <a id="MSearchClose" href="javascript:searchBox.CloseResultsWindow()"><img id="MSearchCloseImg" border="0" src="search/close.png" alt=""/></a> </span> </div> </td> </tr> </tbody> </table> </div> <!-- end header part --> <!-- Generated by Doxygen 1.8.5 --> <script type="text/javascript"> var searchBox = new SearchBox("searchBox", "search",false,'Search'); </script> </div><!-- top --> <div id="side-nav" class="ui-resizable side-nav-resizable"> <div id="nav-tree"> <div id="nav-tree-contents"> <div id="nav-sync" class="sync"></div> </div> </div> <div id="splitbar" style="-moz-user-select:none;" class="ui-resizable-handle"> </div> </div> <script type="text/javascript"> $(document).ready(function(){initNavTree('group__TutorialReductionsVisitorsBroadcasting.html','');}); </script> <div id="doc-content"> <!-- window showing the filter options --> <div id="MSearchSelectWindow" onmouseover="return searchBox.OnSearchSelectShow()" onmouseout="return searchBox.OnSearchSelectHide()" onkeydown="return searchBox.OnSearchSelectKey(event)"> <a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(0)"><span class="SelectionMark"> </span>All</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(1)"><span class="SelectionMark"> </span>Classes</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(2)"><span class="SelectionMark"> </span>Namespaces</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(3)"><span class="SelectionMark"> </span>Functions</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(4)"><span class="SelectionMark"> </span>Variables</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(5)"><span class="SelectionMark"> </span>Typedefs</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(6)"><span class="SelectionMark"> </span>Enumerations</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(7)"><span class="SelectionMark"> </span>Enumerator</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(8)"><span class="SelectionMark"> </span>Friends</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(9)"><span class="SelectionMark"> </span>Groups</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(10)"><span class="SelectionMark"> </span>Pages</a></div> <!-- iframe showing the search results (closed by default) --> <div id="MSearchResultsWindow"> <iframe src="javascript:void(0)" frameborder="0" name="MSearchResults" id="MSearchResults"> </iframe> </div> <div class="header"> <div class="headertitle"> <div class="title">Reductions, visitors and broadcasting<div class="ingroups"><a class="el" href="group__DenseMatrixManipulation__chapter.html">Dense matrix and array manipulation</a></div></div> </div> </div><!--header--> <div class="contents"> <p>This page explains <a class="el" href="namespaceEigen.html" title="Namespace containing all symbols from the Eigen library. ">Eigen</a>'s reductions, visitors and broadcasting and how they are used with <a class="el" href="classEigen_1_1MatrixBase.html">matrices </a> and <a class="el" href="classEigen_1_1ArrayBase.html">arrays </a>.</p> <h1><a class="anchor" id="TutorialReductionsVisitorsBroadcastingReductions"></a> Reductions</h1> <p>In <a class="el" href="namespaceEigen.html" title="Namespace containing all symbols from the Eigen library. ">Eigen</a>, a reduction is a function taking a matrix or array, and returning a single scalar value. One of the most used reductions is <a class="el" href="classEigen_1_1DenseBase.html#a3a3b3fb530d3364ecef0bf9c9daf0983">.sum() </a>, returning the sum of all the coefficients inside a given matrix or array.</p> <table class="example"> <tr> <th>Example:</th><th>Output: </th></tr> <tr> <td><div class="fragment"><div class="line"><span class="preprocessor">#include <iostream></span></div> <div class="line"><span class="preprocessor">#include <Eigen/Dense></span></div> <div class="line"></div> <div class="line"><span class="keyword">using namespace </span>std;</div> <div class="line"><span class="keywordtype">int</span> main()</div> <div class="line">{</div> <div class="line"> <a class="code" href="classEigen_1_1Matrix.html">Eigen::Matrix2d</a> mat;</div> <div class="line"> mat << 1, 2,</div> <div class="line"> 3, 4;</div> <div class="line"> cout << <span class="stringliteral">"Here is mat.sum(): "</span> << mat.<a class="code" href="classEigen_1_1DenseBase.html#a3a3b3fb530d3364ecef0bf9c9daf0983">sum</a>() << endl;</div> <div class="line"> cout << <span class="stringliteral">"Here is mat.prod(): "</span> << mat.<a class="code" href="classEigen_1_1DenseBase.html#a6bdcbfa7e3b07d3246ad80de7170b0f5">prod</a>() << endl;</div> <div class="line"> cout << <span class="stringliteral">"Here is mat.mean(): "</span> << mat.<a class="code" href="classEigen_1_1DenseBase.html#a0af2b3991862a079e3efaef3e4d17d96">mean</a>() << endl;</div> <div class="line"> cout << <span class="stringliteral">"Here is mat.minCoeff(): "</span> << mat.<a class="code" href="classEigen_1_1DenseBase.html#add6cb2d85282829eb9adc9565ce784d6">minCoeff</a>() << endl;</div> <div class="line"> cout << <span class="stringliteral">"Here is mat.maxCoeff(): "</span> << mat.<a class="code" href="classEigen_1_1DenseBase.html#a878f0dae18b28d8158c5f1c232edced2">maxCoeff</a>() << endl;</div> <div class="line"> cout << <span class="stringliteral">"Here is mat.trace(): "</span> << mat.<a class="code" href="classEigen_1_1MatrixBase.html#a71696dd0adbf4731561fd60e55c3a96e">trace</a>() << endl;</div> <div class="line">}</div> </div><!-- fragment --> </td><td><pre class="fragment">Here is mat.sum(): 10 Here is mat.prod(): 24 Here is mat.mean(): 2.5 Here is mat.minCoeff(): 1 Here is mat.maxCoeff(): 4 Here is mat.trace(): 5 </pre> </td></tr> </table> <p>The <em>trace</em> of a matrix, as returned by the function <code>trace()</code>, is the sum of the diagonal coefficients and can equivalently be computed <code>a.diagonal().sum()</code>.</p> <h2><a class="anchor" id="TutorialReductionsVisitorsBroadcastingReductionsNorm"></a> Norm computations</h2> <p>The (Euclidean a.k.a. <img class="formulaInl" alt="$\ell^2$" src="form_190.png"/>) squared norm of a vector can be obtained <a class="el" href="classEigen_1_1MatrixBase.html#a229bd5cc6237359a1d85401743476ede">squaredNorm() </a>. It is equal to the dot product of the vector by itself, and equivalently to the sum of squared absolute values of its coefficients.</p> <p><a class="el" href="namespaceEigen.html" title="Namespace containing all symbols from the Eigen library. ">Eigen</a> also provides the <a class="el" href="classEigen_1_1MatrixBase.html#a0be1b433c65ce9d92c81a4718daf54e5">norm() </a> method, which returns the square root of <a class="el" href="classEigen_1_1MatrixBase.html#a229bd5cc6237359a1d85401743476ede">squaredNorm() </a>.</p> <p>These operations can also operate on matrices; in that case, a n-by-p matrix is seen as a vector of size (n*p), so for example the <a class="el" href="classEigen_1_1MatrixBase.html#a0be1b433c65ce9d92c81a4718daf54e5">norm() </a> method returns the "Frobenius" or "Hilbert-Schmidt" norm. We refrain from speaking of the <img class="formulaInl" alt="$\ell^2$" src="form_190.png"/> norm of a matrix because that can mean different things.</p> <p>If you want other <img class="formulaInl" alt="$\ell^p$" src="form_191.png"/> norms, use the <a class="el" href="">lpNnorm<p>() </a> method. The template parameter <em>p</em> can take the special value <em>Infinity</em> if you want the <img class="formulaInl" alt="$\ell^\infty$" src="form_192.png"/> norm, which is the maximum of the absolute values of the coefficients.</p> <p>The following example demonstrates these methods.</p> <table class="example"> <tr> <th>Example:</th><th>Output: </th></tr> <tr> <td><div class="fragment"><div class="line"><span class="preprocessor">#include <Eigen/Dense></span></div> <div class="line"><span class="preprocessor">#include <iostream></span></div> <div class="line"></div> <div class="line"><span class="keyword">using namespace </span>std;</div> <div class="line"><span class="keyword">using namespace </span>Eigen;</div> <div class="line"></div> <div class="line"><span class="keywordtype">int</span> main()</div> <div class="line">{</div> <div class="line"> VectorXf v(2);</div> <div class="line"> <a class="code" href="classEigen_1_1Matrix.html">MatrixXf</a> m(2,2), n(2,2);</div> <div class="line"> </div> <div class="line"> v << -1,</div> <div class="line"> 2;</div> <div class="line"> </div> <div class="line"> m << 1,-2,</div> <div class="line"> -3,4;</div> <div class="line"></div> <div class="line"> cout << <span class="stringliteral">"v.squaredNorm() = "</span> << v.<a class="code" href="classEigen_1_1MatrixBase.html#a229bd5cc6237359a1d85401743476ede">squaredNorm</a>() << endl;</div> <div class="line"> cout << <span class="stringliteral">"v.norm() = "</span> << v.norm() << endl;</div> <div class="line"> cout << <span class="stringliteral">"v.lpNorm<1>() = "</span> << v.lpNorm<1>() << endl;</div> <div class="line"> cout << <span class="stringliteral">"v.lpNorm<Infinity>() = "</span> << v.lpNorm<<a class="code" href="namespaceEigen.html#a6cff718ef6bce746aa5e100036226625">Infinity</a>>() << endl;</div> <div class="line"></div> <div class="line"> cout << endl;</div> <div class="line"> cout << <span class="stringliteral">"m.squaredNorm() = "</span> << m.<a class="code" href="classEigen_1_1MatrixBase.html#a229bd5cc6237359a1d85401743476ede">squaredNorm</a>() << endl;</div> <div class="line"> cout << <span class="stringliteral">"m.norm() = "</span> << m.<a class="code" href="classEigen_1_1MatrixBase.html#a0be1b433c65ce9d92c81a4718daf54e5">norm</a>() << endl;</div> <div class="line"> cout << <span class="stringliteral">"m.lpNorm<1>() = "</span> << m.lpNorm<1>() << endl;</div> <div class="line"> cout << <span class="stringliteral">"m.lpNorm<Infinity>() = "</span> << m.lpNorm<<a class="code" href="namespaceEigen.html#a6cff718ef6bce746aa5e100036226625">Infinity</a>>() << endl;</div> <div class="line">}</div> </div><!-- fragment --> </td><td><pre class="fragment">v.squaredNorm() = 5 v.norm() = 2.23607 v.lpNorm<1>() = 3 v.lpNorm<Infinity>() = 2 m.squaredNorm() = 30 m.norm() = 5.47723 m.lpNorm<1>() = 10 m.lpNorm<Infinity>() = 4 </pre> </td></tr> </table> <h2><a class="anchor" id="TutorialReductionsVisitorsBroadcastingReductionsBool"></a> Boolean reductions</h2> <p>The following reductions operate on boolean values:</p> <ul> <li><a class="el" href="classEigen_1_1DenseBase.html#aea914316af61df197f21629e14e7870a">all() </a> returns <b>true</b> if all of the coefficients in a given <a class="el" href="classEigen_1_1Matrix.html" title="The matrix class, also used for vectors and row-vectors. ">Matrix</a> or <a class="el" href="classEigen_1_1Array.html" title="General-purpose arrays with easy API for coefficient-wise operations. ">Array</a> evaluate to <b>true</b> .</li> <li><a class="el" href="classEigen_1_1DenseBase.html#a42571e028736ca9103bac8b50f269824">any() </a> returns <b>true</b> if at least one of the coefficients in a given <a class="el" href="classEigen_1_1Matrix.html" title="The matrix class, also used for vectors and row-vectors. ">Matrix</a> or <a class="el" href="classEigen_1_1Array.html" title="General-purpose arrays with easy API for coefficient-wise operations. ">Array</a> evaluates to <b>true</b> .</li> <li><a class="el" href="classEigen_1_1DenseBase.html#aa671b5ea336ba21a7644d3fa6577ee00">count() </a> returns the number of coefficients in a given <a class="el" href="classEigen_1_1Matrix.html" title="The matrix class, also used for vectors and row-vectors. ">Matrix</a> or <a class="el" href="classEigen_1_1Array.html" title="General-purpose arrays with easy API for coefficient-wise operations. ">Array</a> that evaluate to <b>true</b>.</li> </ul> <p>These are typically used in conjunction with the coefficient-wise comparison and equality operators provided by <a class="el" href="classEigen_1_1Array.html" title="General-purpose arrays with easy API for coefficient-wise operations. ">Array</a>. For instance, <code>array > 0</code> is an Array of the same size as <code>array</code> , with <b>true</b> at those positions where the corresponding coefficient of <code>array</code> is positive. Thus, <code>(array > 0).all()</code> tests whether all coefficients of <code>array</code> are positive. This can be seen in the following example:</p> <table class="example"> <tr> <th>Example:</th><th>Output: </th></tr> <tr> <td><div class="fragment"><div class="line"><span class="preprocessor">#include <Eigen/Dense></span></div> <div class="line"><span class="preprocessor">#include <iostream></span></div> <div class="line"></div> <div class="line"><span class="keyword">using namespace </span>std;</div> <div class="line"><span class="keyword">using namespace </span>Eigen;</div> <div class="line"></div> <div class="line"><span class="keywordtype">int</span> main()</div> <div class="line">{</div> <div class="line"> ArrayXXf a(2,2);</div> <div class="line"> </div> <div class="line"> a << 1,2,</div> <div class="line"> 3,4;</div> <div class="line"></div> <div class="line"> cout << <span class="stringliteral">"(a > 0).all() = "</span> << (a > 0).all() << endl;</div> <div class="line"> cout << <span class="stringliteral">"(a > 0).any() = "</span> << (a > 0).any() << endl;</div> <div class="line"> cout << <span class="stringliteral">"(a > 0).count() = "</span> << (a > 0).count() << endl;</div> <div class="line"> cout << endl;</div> <div class="line"> cout << <span class="stringliteral">"(a > 2).all() = "</span> << (a > 2).all() << endl;</div> <div class="line"> cout << <span class="stringliteral">"(a > 2).any() = "</span> << (a > 2).any() << endl;</div> <div class="line"> cout << <span class="stringliteral">"(a > 2).count() = "</span> << (a > 2).count() << endl;</div> <div class="line">}</div> </div><!-- fragment --> </td><td><pre class="fragment">(a > 0).all() = 1 (a > 0).any() = 1 (a > 0).count() = 4 (a > 2).all() = 0 (a > 2).any() = 1 (a > 2).count() = 2 </pre> </td></tr> </table> <h2><a class="anchor" id="TutorialReductionsVisitorsBroadcastingReductionsUserdefined"></a> User defined reductions</h2> <p>TODO</p> <p>In the meantime you can have a look at the DenseBase::redux() function.</p> <h1><a class="anchor" id="TutorialReductionsVisitorsBroadcastingVisitors"></a> Visitors</h1> <p>Visitors are useful when one wants to obtain the location of a coefficient inside a <a class="el" href="classEigen_1_1Matrix.html" title="The matrix class, also used for vectors and row-vectors. ">Matrix</a> or <a class="el" href="classEigen_1_1Array.html" title="General-purpose arrays with easy API for coefficient-wise operations. ">Array</a>. The simplest examples are <a class="el" href="classEigen_1_1DenseBase.html#a878f0dae18b28d8158c5f1c232edced2">maxCoeff(&x,&y) </a> and <a class="el" href="classEigen_1_1DenseBase.html#add6cb2d85282829eb9adc9565ce784d6">minCoeff(&x,&y)</a>, which can be used to find the location of the greatest or smallest coefficient in a <a class="el" href="classEigen_1_1Matrix.html" title="The matrix class, also used for vectors and row-vectors. ">Matrix</a> or <a class="el" href="classEigen_1_1Array.html" title="General-purpose arrays with easy API for coefficient-wise operations. ">Array</a>.</p> <p>The arguments passed to a visitor are pointers to the variables where the row and column position are to be stored. These variables should be of type <a class="el" href="classEigen_1_1DenseBase.html#a4d4873e91be950c079f067fa97fd5c40">Index </a>, as shown below:</p> <table class="example"> <tr> <th>Example:</th><th>Output: </th></tr> <tr> <td><div class="fragment"><div class="line"><span class="preprocessor">#include <iostream></span></div> <div class="line"><span class="preprocessor">#include <Eigen/Dense></span></div> <div class="line"></div> <div class="line"><span class="keyword">using namespace </span>std;</div> <div class="line"><span class="keyword">using namespace </span>Eigen;</div> <div class="line"></div> <div class="line"><span class="keywordtype">int</span> main()</div> <div class="line">{</div> <div class="line"> <a class="code" href="classEigen_1_1Matrix.html">Eigen::MatrixXf</a> m(2,2);</div> <div class="line"> </div> <div class="line"> m << 1, 2,</div> <div class="line"> 3, 4;</div> <div class="line"></div> <div class="line"> <span class="comment">//get location of maximum</span></div> <div class="line"> <a class="code" href="classEigen_1_1DenseBase.html#a4d4873e91be950c079f067fa97fd5c40">MatrixXf::Index</a> maxRow, maxCol;</div> <div class="line"> <span class="keywordtype">float</span> max = m.<a class="code" href="classEigen_1_1DenseBase.html#a878f0dae18b28d8158c5f1c232edced2">maxCoeff</a>(&maxRow, &maxCol);</div> <div class="line"></div> <div class="line"> <span class="comment">//get location of minimum</span></div> <div class="line"> <a class="code" href="classEigen_1_1DenseBase.html#a4d4873e91be950c079f067fa97fd5c40">MatrixXf::Index</a> minRow, minCol;</div> <div class="line"> <span class="keywordtype">float</span> min = m.<a class="code" href="classEigen_1_1DenseBase.html#add6cb2d85282829eb9adc9565ce784d6">minCoeff</a>(&minRow, &minCol);</div> <div class="line"></div> <div class="line"> cout << <span class="stringliteral">"Max: "</span> << max << <span class="stringliteral">", at: "</span> <<</div> <div class="line"> maxRow << <span class="stringliteral">","</span> << maxCol << endl;</div> <div class="line"> cout << <span class="stringliteral">"Min: "</span> << min << <span class="stringliteral">", at: "</span> <<</div> <div class="line"> minRow << <span class="stringliteral">","</span> << minCol << endl;</div> <div class="line">}</div> </div><!-- fragment --> </td><td><pre class="fragment">Max: 4, at: 1,1 Min: 1, at: 0,0 </pre> </td></tr> </table> <p>Note that both functions also return the value of the minimum or maximum coefficient if needed, as if it was a typical reduction operation.</p> <h1><a class="anchor" id="TutorialReductionsVisitorsBroadcastingPartialReductions"></a> Partial reductions</h1> <p>Partial reductions are reductions that can operate column- or row-wise on a <a class="el" href="classEigen_1_1Matrix.html" title="The matrix class, also used for vectors and row-vectors. ">Matrix</a> or <a class="el" href="classEigen_1_1Array.html" title="General-purpose arrays with easy API for coefficient-wise operations. ">Array</a>, applying the reduction operation on each column or row and returning a column or row-vector with the corresponding values. Partial reductions are applied with <a class="el" href="classEigen_1_1DenseBase.html#abe7ae69362c464b6721adbb47c655874">colwise() </a> or <a class="el" href="classEigen_1_1DenseBase.html#a8a7fd1e8004d4bd93a7ea36957aa8e99">rowwise() </a>.</p> <p>A simple example is obtaining the maximum of the elements in each column in a given matrix, storing the result in a row-vector:</p> <table class="example"> <tr> <th>Example:</th><th>Output: </th></tr> <tr> <td><div class="fragment"><div class="line"><span class="preprocessor">#include <iostream></span></div> <div class="line"><span class="preprocessor">#include <Eigen/Dense></span></div> <div class="line"></div> <div class="line"><span class="keyword">using namespace </span>std;</div> <div class="line"><span class="keywordtype">int</span> main()</div> <div class="line">{</div> <div class="line"> <a class="code" href="classEigen_1_1Matrix.html">Eigen::MatrixXf</a> mat(2,4);</div> <div class="line"> mat << 1, 2, 6, 9,</div> <div class="line"> 3, 1, 7, 2;</div> <div class="line"> </div> <div class="line"> std::cout << <span class="stringliteral">"Column's maximum: "</span> << std::endl</div> <div class="line"> << mat.<a class="code" href="classEigen_1_1DenseBase.html#a49a617f24129ca31a27fe8a67ec20370">colwise</a>().<a class="code" href="classEigen_1_1VectorwiseOp.html#ae45bc1f13fd426ca7bea74ae4b8f99d6">maxCoeff</a>() << std::endl;</div> <div class="line">}</div> </div><!-- fragment --> </td><td><pre class="fragment">Column's maximum: 3 2 7 9 </pre> </td></tr> </table> <p>The same operation can be performed row-wise:</p> <table class="example"> <tr> <th>Example:</th><th>Output: </th></tr> <tr> <td><div class="fragment"><div class="line"><span class="preprocessor">#include <iostream></span></div> <div class="line"><span class="preprocessor">#include <Eigen/Dense></span></div> <div class="line"></div> <div class="line"><span class="keyword">using namespace </span>std;</div> <div class="line"><span class="keywordtype">int</span> main()</div> <div class="line">{</div> <div class="line"> <a class="code" href="classEigen_1_1Matrix.html">Eigen::MatrixXf</a> mat(2,4);</div> <div class="line"> mat << 1, 2, 6, 9,</div> <div class="line"> 3, 1, 7, 2;</div> <div class="line"> </div> <div class="line"> std::cout << <span class="stringliteral">"Row's maximum: "</span> << std::endl</div> <div class="line"> << mat.<a class="code" href="classEigen_1_1DenseBase.html#a3af2f03b1d2affcec24e0748edf892cd">rowwise</a>().<a class="code" href="classEigen_1_1VectorwiseOp.html#ae45bc1f13fd426ca7bea74ae4b8f99d6">maxCoeff</a>() << std::endl;</div> <div class="line">}</div> </div><!-- fragment --> </td><td><pre class="fragment">Row's maximum: 9 7 </pre> </td></tr> </table> <p><b>Note that column-wise operations return a 'row-vector' while row-wise operations return a 'column-vector'</b></p> <h2><a class="anchor" id="TutorialReductionsVisitorsBroadcastingPartialReductionsCombined"></a> Combining partial reductions with other operations</h2> <p>It is also possible to use the result of a partial reduction to do further processing. Here is another example that finds the column whose sum of elements is the maximum within a matrix. With column-wise partial reductions this can be coded as:</p> <table class="example"> <tr> <th>Example:</th><th>Output: </th></tr> <tr> <td><div class="fragment"><div class="line"><span class="preprocessor">#include <iostream></span></div> <div class="line"><span class="preprocessor">#include <Eigen/Dense></span></div> <div class="line"></div> <div class="line"><span class="keyword">using namespace </span>std;</div> <div class="line"><span class="keyword">using namespace </span>Eigen;</div> <div class="line"><span class="keywordtype">int</span> main()</div> <div class="line">{</div> <div class="line"> <a class="code" href="classEigen_1_1Matrix.html">MatrixXf</a> mat(2,4);</div> <div class="line"> mat << 1, 2, 6, 9,</div> <div class="line"> 3, 1, 7, 2;</div> <div class="line"> </div> <div class="line"> <a class="code" href="classEigen_1_1DenseBase.html#a4d4873e91be950c079f067fa97fd5c40">MatrixXf::Index</a> maxIndex;</div> <div class="line"> <span class="keywordtype">float</span> maxNorm = mat.<a class="code" href="classEigen_1_1DenseBase.html#a49a617f24129ca31a27fe8a67ec20370">colwise</a>().<a class="code" href="classEigen_1_1VectorwiseOp.html#af27cb8c9043ed6408e64f94645439362">sum</a>().maxCoeff(&maxIndex);</div> <div class="line"> </div> <div class="line"> std::cout << <span class="stringliteral">"Maximum sum at position "</span> << maxIndex << std::endl;</div> <div class="line"></div> <div class="line"> std::cout << <span class="stringliteral">"The corresponding vector is: "</span> << std::endl;</div> <div class="line"> std::cout << mat.<a class="code" href="classEigen_1_1DenseBase.html#a58c77695de3b33405f01f2fdf3dc389d">col</a>( maxIndex ) << std::endl;</div> <div class="line"> std::cout << <span class="stringliteral">"And its sum is is: "</span> << maxNorm << std::endl;</div> <div class="line">}</div> </div><!-- fragment --> </td><td><pre class="fragment">Maximum sum at position 2 The corresponding vector is: 6 7 And its sum is is: 13 </pre> </td></tr> </table> <p>The previous example applies the <a class="el" href="classEigen_1_1DenseBase.html#a3a3b3fb530d3364ecef0bf9c9daf0983">sum() </a> reduction on each column though the <a class="el" href="classEigen_1_1DenseBase.html#abe7ae69362c464b6721adbb47c655874">colwise() </a> visitor, obtaining a new matrix whose size is 1x4.</p> <p>Therefore, if </p> <p class="formulaDsp"> <img class="formulaDsp" alt="\[ \mbox{m} = \begin{bmatrix} 1 & 2 & 6 & 9 \\ 3 & 1 & 7 & 2 \end{bmatrix} \]" src="form_193.png"/> </p> <p>then</p> <p class="formulaDsp"> <img class="formulaDsp" alt="\[ \mbox{m.colwise().sum()} = \begin{bmatrix} 4 & 3 & 13 & 11 \end{bmatrix} \]" src="form_194.png"/> </p> <p>The <a class="el" href="classEigen_1_1DenseBase.html#a878f0dae18b28d8158c5f1c232edced2">maxCoeff() </a> reduction is finally applied to obtain the column index where the maximum sum is found, which is the column index 2 (third column) in this case.</p> <h1><a class="anchor" id="TutorialReductionsVisitorsBroadcastingBroadcasting"></a> Broadcasting</h1> <p>The concept behind broadcasting is similar to partial reductions, with the difference that broadcasting constructs an expression where a vector (column or row) is interpreted as a matrix by replicating it in one direction.</p> <p>A simple example is to add a certain column-vector to each column in a matrix. This can be accomplished with:</p> <table class="example"> <tr> <th>Example:</th><th>Output: </th></tr> <tr> <td><div class="fragment"><div class="line"><span class="preprocessor">#include <iostream></span></div> <div class="line"><span class="preprocessor">#include <Eigen/Dense></span></div> <div class="line"></div> <div class="line"><span class="keyword">using namespace </span>std;</div> <div class="line"><span class="keywordtype">int</span> main()</div> <div class="line">{</div> <div class="line"> <a class="code" href="classEigen_1_1Matrix.html">Eigen::MatrixXf</a> mat(2,4);</div> <div class="line"> <a class="code" href="classEigen_1_1Matrix.html">Eigen::VectorXf</a> v(2);</div> <div class="line"> </div> <div class="line"> mat << 1, 2, 6, 9,</div> <div class="line"> 3, 1, 7, 2;</div> <div class="line"> </div> <div class="line"> v << 0,</div> <div class="line"> 1;</div> <div class="line"> </div> <div class="line"> <span class="comment">//add v to each column of m</span></div> <div class="line"> mat.<a class="code" href="classEigen_1_1DenseBase.html#a49a617f24129ca31a27fe8a67ec20370">colwise</a>() += v;</div> <div class="line"> </div> <div class="line"> std::cout << <span class="stringliteral">"Broadcasting result: "</span> << std::endl;</div> <div class="line"> std::cout << mat << std::endl;</div> <div class="line">}</div> </div><!-- fragment --> </td><td><pre class="fragment">Broadcasting result: 1 2 6 9 4 2 8 3 </pre> </td></tr> </table> <p>We can interpret the instruction <code>mat.colwise() += v</code> in two equivalent ways. It adds the vector <code>v</code> to every column of the matrix. Alternatively, it can be interpreted as repeating the vector <code>v</code> four times to form a four-by-two matrix which is then added to <code>mat:</code> </p> <p class="formulaDsp"> <img class="formulaDsp" alt="\[ \begin{bmatrix} 1 & 2 & 6 & 9 \\ 3 & 1 & 7 & 2 \end{bmatrix} + \begin{bmatrix} 0 & 0 & 0 & 0 \\ 1 & 1 & 1 & 1 \end{bmatrix} = \begin{bmatrix} 1 & 2 & 6 & 9 \\ 4 & 2 & 8 & 3 \end{bmatrix}. \]" src="form_195.png"/> </p> <p> The operators <code>-=</code>, <code>+</code> and <code>-</code> can also be used column-wise and row-wise. On arrays, we can also use the operators <code>*=</code>, <code>/=</code>, <code>*</code> and <code>/</code> to perform coefficient-wise multiplication and division column-wise or row-wise. These operators are not available on matrices because it is not clear what they would do. If you want multiply column 0 of a matrix <code>mat</code> with <code>v(0)</code>, column 1 with <code>v(1)</code>, and so on, then use <code>mat = mat * v.asDiagonal()</code>.</p> <p>It is important to point out that the vector to be added column-wise or row-wise must be of type Vector, and cannot be a <a class="el" href="classEigen_1_1Matrix.html" title="The matrix class, also used for vectors and row-vectors. ">Matrix</a>. If this is not met then you will get compile-time error. This also means that broadcasting operations can only be applied with an object of type Vector, when operating with <a class="el" href="classEigen_1_1Matrix.html" title="The matrix class, also used for vectors and row-vectors. ">Matrix</a>. The same applies for the <a class="el" href="classEigen_1_1Array.html" title="General-purpose arrays with easy API for coefficient-wise operations. ">Array</a> class, where the equivalent for VectorXf is ArrayXf. As always, you should not mix arrays and matrices in the same expression.</p> <p>To perform the same operation row-wise we can do:</p> <table class="example"> <tr> <th>Example:</th><th>Output: </th></tr> <tr> <td><div class="fragment"><div class="line"><span class="preprocessor">#include <iostream></span></div> <div class="line"><span class="preprocessor">#include <Eigen/Dense></span></div> <div class="line"></div> <div class="line"><span class="keyword">using namespace </span>std;</div> <div class="line"><span class="keywordtype">int</span> main()</div> <div class="line">{</div> <div class="line"> <a class="code" href="classEigen_1_1Matrix.html">Eigen::MatrixXf</a> mat(2,4);</div> <div class="line"> <a class="code" href="classEigen_1_1Matrix.html">Eigen::VectorXf</a> v(4);</div> <div class="line"> </div> <div class="line"> mat << 1, 2, 6, 9,</div> <div class="line"> 3, 1, 7, 2;</div> <div class="line"> </div> <div class="line"> v << 0,1,2,3;</div> <div class="line"> </div> <div class="line"> <span class="comment">//add v to each row of m</span></div> <div class="line"> mat.<a class="code" href="classEigen_1_1DenseBase.html#a3af2f03b1d2affcec24e0748edf892cd">rowwise</a>() += v.transpose();</div> <div class="line"> </div> <div class="line"> std::cout << <span class="stringliteral">"Broadcasting result: "</span> << std::endl;</div> <div class="line"> std::cout << mat << std::endl;</div> <div class="line">}</div> </div><!-- fragment --> </td><td><pre class="fragment">Broadcasting result: 1 3 8 12 3 2 9 5 </pre> </td></tr> </table> <h2><a class="anchor" id="TutorialReductionsVisitorsBroadcastingBroadcastingCombined"></a> Combining broadcasting with other operations</h2> <p>Broadcasting can also be combined with other operations, such as <a class="el" href="classEigen_1_1Matrix.html" title="The matrix class, also used for vectors and row-vectors. ">Matrix</a> or <a class="el" href="classEigen_1_1Array.html" title="General-purpose arrays with easy API for coefficient-wise operations. ">Array</a> operations, reductions and partial reductions.</p> <p>Now that broadcasting, reductions and partial reductions have been introduced, we can dive into a more advanced example that finds the nearest neighbour of a vector <code>v</code> within the columns of matrix <code>m</code>. The Euclidean distance will be used in this example, computing the squared Euclidean distance with the partial reduction named <a class="el" href="classEigen_1_1MatrixBase.html#a229bd5cc6237359a1d85401743476ede">squaredNorm() </a>:</p> <table class="example"> <tr> <th>Example:</th><th>Output: </th></tr> <tr> <td><div class="fragment"><div class="line"><span class="preprocessor">#include <iostream></span></div> <div class="line"><span class="preprocessor">#include <Eigen/Dense></span></div> <div class="line"></div> <div class="line"><span class="keyword">using namespace </span>std;</div> <div class="line"><span class="keyword">using namespace </span>Eigen;</div> <div class="line"></div> <div class="line"><span class="keywordtype">int</span> main()</div> <div class="line">{</div> <div class="line"> <a class="code" href="classEigen_1_1Matrix.html">Eigen::MatrixXf</a> m(2,4);</div> <div class="line"> <a class="code" href="classEigen_1_1Matrix.html">Eigen::VectorXf</a> v(2);</div> <div class="line"> </div> <div class="line"> m << 1, 23, 6, 9,</div> <div class="line"> 3, 11, 7, 2;</div> <div class="line"> </div> <div class="line"> v << 2,</div> <div class="line"> 3;</div> <div class="line"></div> <div class="line"> <a class="code" href="classEigen_1_1DenseBase.html#a4d4873e91be950c079f067fa97fd5c40">MatrixXf::Index</a> index;</div> <div class="line"> <span class="comment">// find nearest neighbour</span></div> <div class="line"> (m.<a class="code" href="classEigen_1_1DenseBase.html#a49a617f24129ca31a27fe8a67ec20370">colwise</a>() - v).colwise().<a class="code" href="classEigen_1_1MatrixBase.html#a229bd5cc6237359a1d85401743476ede">squaredNorm</a>().minCoeff(&index);</div> <div class="line"></div> <div class="line"> cout << <span class="stringliteral">"Nearest neighbour is column "</span> << index << <span class="stringliteral">":"</span> << endl;</div> <div class="line"> cout << m.<a class="code" href="classEigen_1_1DenseBase.html#a58c77695de3b33405f01f2fdf3dc389d">col</a>(index) << endl;</div> <div class="line">}</div> </div><!-- fragment --> </td><td><pre class="fragment">Nearest neighbour is column 0: 1 3 </pre> </td></tr> </table> <p>The line that does the job is </p> <div class="fragment"><div class="line">(m.<a class="code" href="classEigen_1_1DenseBase.html#a49a617f24129ca31a27fe8a67ec20370">colwise</a>() - v).colwise().squaredNorm().minCoeff(&index);</div> </div><!-- fragment --><p>We will go step by step to understand what is happening:</p> <ul> <li><code>m.colwise() - v</code> is a broadcasting operation, subtracting <code>v</code> from each column in <code>m</code>. The result of this operation is a new matrix whose size is the same as matrix <code>m</code>: <p class="formulaDsp"> <img class="formulaDsp" alt="\[ \mbox{m.colwise() - v} = \begin{bmatrix} -1 & 21 & 4 & 7 \\ 0 & 8 & 4 & -1 \end{bmatrix} \]" src="form_196.png"/> </p> </li> <li><code>(m.colwise() - v).colwise().squaredNorm()</code> is a partial reduction, computing the squared norm column-wise. The result of this operation is a row-vector where each coefficient is the squared Euclidean distance between each column in <code>m</code> and <code>v</code>: <p class="formulaDsp"> <img class="formulaDsp" alt="\[ \mbox{(m.colwise() - v).colwise().squaredNorm()} = \begin{bmatrix} 1 & 505 & 32 & 50 \end{bmatrix} \]" src="form_197.png"/> </p> </li> <li>Finally, <code>minCoeff(&index)</code> is used to obtain the index of the column in <code>m</code> that is closest to <code>v</code> in terms of Euclidean distance. </li> </ul> </div><!-- contents --> </div><!-- doc-content --> <!-- start footer part --> <div id="nav-path" class="navpath"><!-- id is needed for treeview function! --> <ul> <li class="footer">Generated on Mon Oct 28 2013 11:04:27 for Eigen by <a href="http://www.doxygen.org/index.html"> <img class="footer" src="doxygen.png" alt="doxygen"/></a> 1.8.5 </li> </ul> </div> <!-- Piwik --> <!-- <script type="text/javascript"> var pkBaseURL = (("https:" == document.location.protocol) ? "https://stats.sylphide-consulting.com/piwik/" : "http://stats.sylphide-consulting.com/piwik/"); document.write(unescape("%3Cscript src='" + pkBaseURL + "piwik.js' type='text/javascript'%3E%3C/script%3E")); </script><script type="text/javascript"> try { var piwikTracker = Piwik.getTracker(pkBaseURL + "piwik.php", 20); piwikTracker.trackPageView(); piwikTracker.enableLinkTracking(); } catch( err ) {} </script><noscript><p><img src="http://stats.sylphide-consulting.com/piwik/piwik.php?idsite=20" style="border:0" alt="" /></p></noscript> --> <!-- End Piwik Tracking Code --> </body> </html>