Sophie

Sophie

distrib > Mageia > 4 > x86_64 > by-pkgid > 99cb5ede6a5329071fbeecc8218deb35 > files > 219

eigen3-doc-3.2-3.mga4.noarch.rpm

<!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: Writing Functions Taking %Eigen Types as Parameters</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>
   &#160;<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('TopicFunctionTakingEigenTypes.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">&#160;</span>All</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(1)"><span class="SelectionMark">&#160;</span>Classes</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(2)"><span class="SelectionMark">&#160;</span>Namespaces</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(3)"><span class="SelectionMark">&#160;</span>Functions</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(4)"><span class="SelectionMark">&#160;</span>Variables</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(5)"><span class="SelectionMark">&#160;</span>Typedefs</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(6)"><span class="SelectionMark">&#160;</span>Enumerations</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(7)"><span class="SelectionMark">&#160;</span>Enumerator</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(8)"><span class="SelectionMark">&#160;</span>Friends</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(9)"><span class="SelectionMark">&#160;</span>Groups</a><a class="SelectItem" href="javascript:void(0)" onclick="searchBox.OnSelectItem(10)"><span class="SelectionMark">&#160;</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">Writing Functions Taking Eigen Types as Parameters </div>  </div>
</div><!--header-->
<div class="contents">
<div class="textblock"><p>Eigen's use of expression templates results in potentially every expression being of a different type. If you pass such an expression to a function taking a parameter of type <a class="el" href="classEigen_1_1Matrix.html" title="The matrix class, also used for vectors and row-vectors. ">Matrix</a>, your expression will implicitly be evaluated into a temporary <a class="el" href="classEigen_1_1Matrix.html" title="The matrix class, also used for vectors and row-vectors. ">Matrix</a>, which will then be passed to the function. This means that you lose the benefit of expression templates. Concretely, this has two drawbacks: </p>
<ul>
<li>The evaluation into a temporary may be useless and inefficient; </li>
<li>This only allows the function to read from the expression, not to write to it.</li>
</ul>
<p>Fortunately, all this myriad of expression types have in common that they all inherit a few common, templated base classes. By letting your function take templated parameters of these base types, you can let them play nicely with Eigen's expression templates.</p>
<h1><a class="anchor" id="TopicFirstExamples"></a>
Some First Examples</h1>
<p>This section will provide simple examples for different types of objects Eigen is offering. Before starting with the actual examples, we need to recapitulate which base objects we can work with (see also <a class="el" href="TopicClassHierarchy.html">The class hierarchy</a>).</p>
<ul>
<li><a class="el" href="classEigen_1_1MatrixBase.html" title="Base class for all dense matrices, vectors, and expressions. ">MatrixBase</a>: The common base class for all dense matrix expressions (as opposed to array expressions, as opposed to sparse and special matrix classes). Use it in functions that are meant to work only on dense matrices. </li>
<li><a class="el" href="classEigen_1_1ArrayBase.html" title="Base class for all 1D and 2D array, and related expressions. ">ArrayBase</a>: The common base class for all dense array expressions (as opposed to matrix expressions, etc). Use it in functions that are meant to work only on arrays. </li>
<li><a class="el" href="classEigen_1_1DenseBase.html" title="Base class for all dense matrices, vectors, and arrays. ">DenseBase</a>: The common base class for all dense matrix expression, that is, the base class for both <code><a class="el" href="classEigen_1_1MatrixBase.html" title="Base class for all dense matrices, vectors, and expressions. ">MatrixBase</a></code> and <code><a class="el" href="classEigen_1_1ArrayBase.html" title="Base class for all 1D and 2D array, and related expressions. ">ArrayBase</a></code>. It can be used in functions that are meant to work on both matrices and arrays. </li>
<li><a class="el" href="structEigen_1_1EigenBase.html">EigenBase</a>: The base class unifying all types of objects that can be evaluated into dense matrices or arrays, for example special matrix classes such as diagonal matrices, permutation matrices, etc. It can be used in functions that are meant to work on any such general type.</li>
</ul>
<p><b> EigenBase Example </b><br/>
<br/>
 Prints the dimensions of the most generic object present in Eigen. It could be any matrix expressions, any dense or sparse matrix and any 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 &lt;iostream&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;Eigen/Core&gt;</span></div>
<div class="line"><span class="keyword">using namespace </span>Eigen;</div>
<div class="line"></div>
<div class="line"><span class="keyword">template</span> &lt;<span class="keyword">typename</span> Derived&gt;</div>
<div class="line"><span class="keywordtype">void</span> print_size(<span class="keyword">const</span> <a class="code" href="structEigen_1_1EigenBase.html">EigenBase&lt;Derived&gt;</a>&amp; b)</div>
<div class="line">{</div>
<div class="line">  std::cout &lt;&lt; <span class="stringliteral">&quot;size (rows, cols): &quot;</span> &lt;&lt; b.<a class="code" href="structEigen_1_1EigenBase.html#a76f5bc8a03ec105ab4be1f2b91b7b5d5">size</a>() &lt;&lt; <span class="stringliteral">&quot; (&quot;</span> &lt;&lt; b.<a class="code" href="structEigen_1_1EigenBase.html#a5552abd83dbd03c85cea6d61fd8875a5">rows</a>()</div>
<div class="line">            &lt;&lt; <span class="stringliteral">&quot;, &quot;</span> &lt;&lt; b.<a class="code" href="structEigen_1_1EigenBase.html#aaca1908a5ec508a25ff0a8bca803e5f3">cols</a>() &lt;&lt; <span class="stringliteral">&quot;)&quot;</span> &lt;&lt; std::endl;</div>
<div class="line">}</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">Vector3f</a> v;</div>
<div class="line">    print_size(v);</div>
<div class="line">    <span class="comment">// v.asDiagonal() returns a 3x3 diagonal matrix pseudo-expression</span></div>
<div class="line">    print_size(v.<a class="code" href="classEigen_1_1MatrixBase.html#adaf22d3a2069ec2c0df912cb87329e9c">asDiagonal</a>());</div>
<div class="line">}</div>
</div><!-- fragment -->  </td><td><pre class="fragment">size (rows, cols): 3 (3, 1)
size (rows, cols): 9 (3, 3)
</pre> </td></tr>
</table>
<p><b> DenseBase Example </b><br/>
<br/>
 Prints a sub-block of the dense expression. Accepts any dense matrix or array expression, but no sparse objects and no special matrix classes such as <a class="el" href="classEigen_1_1DiagonalMatrix.html" title="Represents a diagonal matrix with its storage. ">DiagonalMatrix</a>. </p>
<div class="fragment"><div class="line"><span class="keyword">template</span> &lt;<span class="keyword">typename</span> Derived&gt;</div>
<div class="line"><span class="keywordtype">void</span> print_block(<span class="keyword">const</span> DenseBase&lt;Derived&gt;&amp; b, <span class="keywordtype">int</span> x, <span class="keywordtype">int</span> y, <span class="keywordtype">int</span> r, <span class="keywordtype">int</span> c)</div>
<div class="line">{</div>
<div class="line">  std::cout &lt;&lt; <span class="stringliteral">&quot;block: &quot;</span> &lt;&lt; b.block(x,y,r,c) &lt;&lt; std::endl;</div>
<div class="line">}</div>
</div><!-- fragment --><p> <b> ArrayBase Example </b><br/>
<br/>
 Prints the maximum coefficient of the array or array-expression. </p>
<div class="fragment"><div class="line"><span class="keyword">template</span> &lt;<span class="keyword">typename</span> Derived&gt;</div>
<div class="line"><span class="keywordtype">void</span> print_max_coeff(<span class="keyword">const</span> ArrayBase&lt;Derived&gt; &amp;a)</div>
<div class="line">{</div>
<div class="line">  std::cout &lt;&lt; <span class="stringliteral">&quot;max: &quot;</span> &lt;&lt; a.maxCoeff() &lt;&lt; std::endl;</div>
<div class="line">}</div>
</div><!-- fragment --><p> <b> MatrixBase Example </b><br/>
<br/>
 Prints the inverse condition number of the given matrix or matrix-expression. </p>
<div class="fragment"><div class="line"><span class="keyword">template</span> &lt;<span class="keyword">typename</span> Derived&gt;</div>
<div class="line"><span class="keywordtype">void</span> print_inv_cond(<span class="keyword">const</span> MatrixBase&lt;Derived&gt;&amp; a)</div>
<div class="line">{</div>
<div class="line">  <span class="keyword">const</span> <span class="keyword">typename</span> JacobiSVD&lt;typename Derived::PlainObject&gt;::SingularValuesType&amp;</div>
<div class="line">    sing_vals = a.jacobiSvd().singularValues();</div>
<div class="line">  std::cout &lt;&lt; <span class="stringliteral">&quot;inv cond: &quot;</span> &lt;&lt; sing_vals(sing_vals.size()-1) / sing_vals(0) &lt;&lt; std::endl;</div>
<div class="line">}</div>
</div><!-- fragment --><p> <b> Multiple templated arguments example </b><br/>
<br/>
 Calculate the Euclidean distance between two points. </p>
<div class="fragment"><div class="line"><span class="keyword">template</span> &lt;<span class="keyword">typename</span> DerivedA,<span class="keyword">typename</span> DerivedB&gt;</div>
<div class="line"><span class="keyword">typename</span> DerivedA::Scalar squaredist(<span class="keyword">const</span> MatrixBase&lt;DerivedA&gt;&amp; p1,<span class="keyword">const</span> MatrixBase&lt;DerivedB&gt;&amp; p2)</div>
<div class="line">{</div>
<div class="line">  <span class="keywordflow">return</span> (p1-p2).squaredNorm();</div>
<div class="line">}</div>
</div><!-- fragment --><p> Notice that we used two template parameters, one per argument. This permits the function to handle inputs of different types, e.g., </p>
<div class="fragment"><div class="line">squaredist(v1,2*v2)</div>
</div><!-- fragment --><p> where the first argument <code>v1</code> is a vector and the second argument <code>2*v2</code> is an expression. <br/>
<br/>
</p>
<p>These examples are just intended to give the reader a first impression of how functions can be written which take a plain and constant <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> argument. They are also intended to give the reader an idea about the most common base classes being the optimal candidates for functions. In the next section we will look in more detail at an example and the different ways it can be implemented, while discussing each implementation's problems and advantages. For the discussion below, <a class="el" href="classEigen_1_1Matrix.html" title="The matrix class, also used for vectors and row-vectors. ">Matrix</a> and <a class="el" href="classEigen_1_1Array.html" title="General-purpose arrays with easy API for coefficient-wise operations. ">Array</a> as well as <a class="el" href="classEigen_1_1MatrixBase.html" title="Base class for all dense matrices, vectors, and expressions. ">MatrixBase</a> and <a class="el" href="classEigen_1_1ArrayBase.html" title="Base class for all 1D and 2D array, and related expressions. ">ArrayBase</a> can be exchanged and all arguments still hold.</p>
<h1><a class="anchor" id="TopicUsingRefClass"></a>
How to write generic, but non-templated function?</h1>
<p>In all the previous examples, the functions had to be template functions. This approach allows to write very generic code, but it is often desirable to write non templated function and still keep some level of genericity to avoid stupid copies of the arguments. The typical example is to write functions accepting both a MatrixXf or a block of a MatrixXf. This exactly the purpose of the <a class="el" href="classEigen_1_1Ref.html" title="A matrix or vector expression mapping an existing expressions. ">Ref</a> class. Here is a simple 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 &lt;iostream&gt;</span></div>
<div class="line"><span class="preprocessor">#include &lt;Eigen/SVD&gt;</span></div>
<div class="line"><span class="keyword">using namespace </span>Eigen;</div>
<div class="line"><span class="keyword">using namespace </span>std;</div>
<div class="line"></div>
<div class="line"><span class="keywordtype">float</span> inv_cond(<span class="keyword">const</span> <a class="code" href="classEigen_1_1Ref.html">Ref&lt;const MatrixXf&gt;</a>&amp; a)</div>
<div class="line">{</div>
<div class="line">  <span class="keyword">const</span> <a class="code" href="classEigen_1_1Matrix.html">VectorXf</a> sing_vals = a.jacobiSvd().singularValues();</div>
<div class="line">  <span class="keywordflow">return</span> sing_vals(sing_vals.size()-1) / sing_vals(0);</div>
<div class="line">}</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">Matrix4f</a> m = <a class="code" href="classEigen_1_1DenseBase.html#a8e759dafdd9ecc446d397b7f5435f60a">Matrix4f::Random</a>();</div>
<div class="line">  cout &lt;&lt; <span class="stringliteral">&quot;matrix m:&quot;</span> &lt;&lt; endl &lt;&lt; m &lt;&lt; endl &lt;&lt; endl;</div>
<div class="line">  cout &lt;&lt; <span class="stringliteral">&quot;inv_cond(m):          &quot;</span> &lt;&lt; inv_cond(m)                      &lt;&lt; endl;</div>
<div class="line">  cout &lt;&lt; <span class="stringliteral">&quot;inv_cond(m(1:3,1:3)): &quot;</span> &lt;&lt; inv_cond(m.<a class="code" href="classEigen_1_1DenseBase.html#a6f5fc5fe9d3fb70e62d4a9b1795704a8">topLeftCorner</a>(3,3))   &lt;&lt; endl;</div>
<div class="line">  cout &lt;&lt; <span class="stringliteral">&quot;inv_cond(m+I):        &quot;</span> &lt;&lt; inv_cond(m+<a class="code" href="classEigen_1_1MatrixBase.html#a0650b65c6ae6c3d19a138b72a6d68568">Matrix4f::Identity</a>()) &lt;&lt; endl;</div>
<div class="line">}</div>
</div><!-- fragment -->  </td><td><pre class="fragment">matrix m:
   0.68   0.823  -0.444   -0.27
 -0.211  -0.605   0.108  0.0268
  0.566   -0.33 -0.0452   0.904
  0.597   0.536   0.258   0.832

inv_cond(m):          0.0562343
inv_cond(m(1:3,1:3)): 0.0836819
inv_cond(m+I):        0.160203
</pre> </td></tr>
</table>
<p>In the first two calls to inv_cond, no copy occur because the memory layout of the arguments matches the memory layout accepted by Ref&lt;MatrixXf&gt;. However, in the last call, we have a generic expression that will be automatically evaluated into a temporary MatrixXf by the Ref&lt;&gt; object.</p>
<p>A <a class="el" href="classEigen_1_1Ref.html" title="A matrix or vector expression mapping an existing expressions. ">Ref</a> object can also be writable. Here is an example of a function computing the covariance matrix of two input matrices where each row is an observation: </p>
<div class="fragment"><div class="line"><span class="keywordtype">void</span> cov(<span class="keyword">const</span> Ref&lt;const MatrixXf&gt; x, <span class="keyword">const</span> Ref&lt;const MatrixXf&gt; y, Ref&lt;MatrixXf&gt; C)</div>
<div class="line">{</div>
<div class="line">  <span class="keyword">const</span> <span class="keywordtype">float</span> num_observations = <span class="keyword">static_cast&lt;</span><span class="keywordtype">float</span><span class="keyword">&gt;</span>(x.rows());</div>
<div class="line">  <span class="keyword">const</span> <a class="code" href="group__matrixtypedefs.html#gad950c8ebb50f405c8eabecb68e873ff4">RowVectorXf</a> x_mean = x.colwise().sum() / num_observations;</div>
<div class="line">  <span class="keyword">const</span> <a class="code" href="group__matrixtypedefs.html#gad950c8ebb50f405c8eabecb68e873ff4">RowVectorXf</a> y_mean = y.colwise().sum() / num_observations;</div>
<div class="line">  C = (x.rowwise() - x_mean).transpose() * (y.rowwise() - y_mean) / num_observations;</div>
<div class="line">}</div>
</div><!-- fragment --><p> and here are two examples calling cov without any copy: </p>
<div class="fragment"><div class="line"><a class="code" href="group__matrixtypedefs.html#gabab09c32e96cfa9829a88400627af162">MatrixXf</a> m1, m2, m3</div>
<div class="line">cov(m1, m2, m3);</div>
<div class="line">cov(m1.leftCols&lt;3&gt;(), m2.leftCols&lt;3&gt;(), m3.topLeftCorner&lt;3,3&gt;());</div>
</div><!-- fragment --><p> The Ref&lt;&gt; class has two other optional template arguments allowing to control the kind of memory layout that can be accepted without any copy. See the class <a class="el" href="classEigen_1_1Ref.html" title="A matrix or vector expression mapping an existing expressions. ">Ref</a> documentation for the details.</p>
<h1><a class="anchor" id="TopicPlainFunctionsWorking"></a>
In which cases do functions taking plain Matrix or Array arguments work?</h1>
<p>Without using template functions, and without the <a class="el" href="classEigen_1_1Ref.html" title="A matrix or vector expression mapping an existing expressions. ">Ref</a> class, a naive implementation of the previous cov function might look like this </p>
<div class="fragment"><div class="line"><a class="code" href="group__matrixtypedefs.html#gabab09c32e96cfa9829a88400627af162">MatrixXf</a> cov(<span class="keyword">const</span> <a class="code" href="group__matrixtypedefs.html#gabab09c32e96cfa9829a88400627af162">MatrixXf</a>&amp; x, <span class="keyword">const</span> <a class="code" href="group__matrixtypedefs.html#gabab09c32e96cfa9829a88400627af162">MatrixXf</a>&amp; y)</div>
<div class="line">{</div>
<div class="line">  <span class="keyword">const</span> <span class="keywordtype">float</span> num_observations = <span class="keyword">static_cast&lt;</span><span class="keywordtype">float</span><span class="keyword">&gt;</span>(x.rows());</div>
<div class="line">  <span class="keyword">const</span> <a class="code" href="group__matrixtypedefs.html#gad950c8ebb50f405c8eabecb68e873ff4">RowVectorXf</a> x_mean = x.colwise().sum() / num_observations;</div>
<div class="line">  <span class="keyword">const</span> <a class="code" href="group__matrixtypedefs.html#gad950c8ebb50f405c8eabecb68e873ff4">RowVectorXf</a> y_mean = y.colwise().sum() / num_observations;</div>
<div class="line">  <span class="keywordflow">return</span> (x.rowwise() - x_mean).transpose() * (y.rowwise() - y_mean) / num_observations;</div>
<div class="line">}</div>
</div><!-- fragment --><p> and contrary to what one might think at first, this implementation is fine unless you require a generic implementation that works with double matrices too and unless you do not care about temporary objects. Why is that the case? Where are temporaries involved? How can code as given below compile? </p>
<div class="fragment"><div class="line"><a class="code" href="group__matrixtypedefs.html#gabab09c32e96cfa9829a88400627af162">MatrixXf</a> x,y,z;</div>
<div class="line"><a class="code" href="group__matrixtypedefs.html#gabab09c32e96cfa9829a88400627af162">MatrixXf</a> C = cov(x,y+z);</div>
</div><!-- fragment --><p> In this special case, the example is fine and will be working because both parameters are declared as <em>const</em> references. The compiler creates a temporary and evaluates the expression x+z into this temporary. Once the function is processed, the temporary is released and the result is assigned to C.</p>
<p><b>Note:</b> Functions taking <em>const</em> references to <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>) can process expressions at the cost of temporaries.</p>
<h1><a class="anchor" id="TopicPlainFunctionsFailing"></a>
In which cases do functions taking a plain Matrix or Array argument fail?</h1>
<p>Here, we consider a slightly modified version of the function given above. This time, we do not want to return the result but pass an additional non-const paramter which allows us to store the result. A first naive implementation might look as follows. </p>
<div class="fragment"><div class="line"><span class="comment">// Note: This code is flawed!</span></div>
<div class="line"><span class="keywordtype">void</span> cov(<span class="keyword">const</span> <a class="code" href="group__matrixtypedefs.html#gabab09c32e96cfa9829a88400627af162">MatrixXf</a>&amp; x, <span class="keyword">const</span> <a class="code" href="group__matrixtypedefs.html#gabab09c32e96cfa9829a88400627af162">MatrixXf</a>&amp; y, <a class="code" href="group__matrixtypedefs.html#gabab09c32e96cfa9829a88400627af162">MatrixXf</a>&amp; C)</div>
<div class="line">{</div>
<div class="line">  <span class="keyword">const</span> <span class="keywordtype">float</span> num_observations = <span class="keyword">static_cast&lt;</span><span class="keywordtype">float</span><span class="keyword">&gt;</span>(x.rows());</div>
<div class="line">  <span class="keyword">const</span> <a class="code" href="group__matrixtypedefs.html#gad950c8ebb50f405c8eabecb68e873ff4">RowVectorXf</a> x_mean = x.colwise().sum() / num_observations;</div>
<div class="line">  <span class="keyword">const</span> <a class="code" href="group__matrixtypedefs.html#gad950c8ebb50f405c8eabecb68e873ff4">RowVectorXf</a> y_mean = y.colwise().sum() / num_observations;</div>
<div class="line">  C = (x.rowwise() - x_mean).transpose() * (y.rowwise() - y_mean) / num_observations;</div>
<div class="line">}</div>
</div><!-- fragment --><p> When trying to execute the following code </p>
<div class="fragment"><div class="line"><a class="code" href="group__matrixtypedefs.html#gabab09c32e96cfa9829a88400627af162">MatrixXf</a> C = <a class="code" href="classEigen_1_1DenseBase.html#a2e36721b4833498b713e94a7948c6a29">MatrixXf::Zero</a>(3,6);</div>
<div class="line">cov(x,y, C.block(0,0,3,3));</div>
</div><!-- fragment --><p> the compiler will fail, because it is not possible to convert the expression returned by <code><a class="el" href="classEigen_1_1DenseBase.html#a1dbaa2fc7b809720407130f48dfacf8f">MatrixXf::block()</a></code> into a non-const <code>MatrixXf&amp;</code>. This is the case because the compiler wants to protect you from writing your result to a temporary object. In this special case this protection is not intended &ndash; we want to write to a temporary object. So how can we overcome this problem?</p>
<p>The solution which is preferred at the moment is based on a little <em>hack</em>. One needs to pass a const reference to the matrix and internally the constness needs to be cast away. The correct implementation for C98 compliant compilers would be </p>
<div class="fragment"><div class="line"><span class="keyword">template</span> &lt;<span class="keyword">typename</span> Derived, <span class="keyword">typename</span> OtherDerived&gt;</div>
<div class="line"><span class="keywordtype">void</span> cov(<span class="keyword">const</span> MatrixBase&lt;Derived&gt;&amp; x, <span class="keyword">const</span> MatrixBase&lt;Derived&gt;&amp; y, MatrixBase&lt;OtherDerived&gt; <span class="keyword">const</span> &amp; C)</div>
<div class="line">{</div>
<div class="line">  <span class="keyword">typedef</span> <span class="keyword">typename</span> Derived::Scalar Scalar;</div>
<div class="line">  <span class="keyword">typedef</span> <span class="keyword">typename</span> internal::plain_row_type&lt;Derived&gt;::type RowVectorType;</div>
<div class="line"></div>
<div class="line">  <span class="keyword">const</span> Scalar num_observations = <span class="keyword">static_cast&lt;</span>Scalar<span class="keyword">&gt;</span>(x.rows());</div>
<div class="line"></div>
<div class="line">  <span class="keyword">const</span> RowVectorType x_mean = x.colwise().sum() / num_observations;</div>
<div class="line">  <span class="keyword">const</span> RowVectorType y_mean = y.colwise().sum() / num_observations;</div>
<div class="line"></div>
<div class="line">  <span class="keyword">const_cast&lt;</span> MatrixBase&lt;OtherDerived&gt;&amp; <span class="keyword">&gt;</span>(C) =</div>
<div class="line">    (x.rowwise() - x_mean).transpose() * (y.rowwise() - y_mean) / num_observations;</div>
<div class="line">}</div>
</div><!-- fragment --><p> The implementation above does now not only work with temporary expressions but it also allows to use the function with matrices of arbitrary floating point scalar types.</p>
<p><b>Note:</b> The const cast hack will only work with templated functions. It will not work with the MatrixXf implementation because it is not possible to cast a <a class="el" href="classEigen_1_1Block.html" title="Expression of a fixed-size or dynamic-size block. ">Block</a> expression to a <a class="el" href="classEigen_1_1Matrix.html" title="The matrix class, also used for vectors and row-vectors. ">Matrix</a> reference!</p>
<h1><a class="anchor" id="TopicResizingInGenericImplementations"></a>
How to resize matrices in generic implementations?</h1>
<p>One might think we are done now, right? This is not completely true because in order for our covariance function to be generically applicable, we want the follwing code to work </p>
<div class="fragment"><div class="line"><a class="code" href="group__matrixtypedefs.html#gabab09c32e96cfa9829a88400627af162">MatrixXf</a> x = <a class="code" href="classEigen_1_1DenseBase.html#a8e759dafdd9ecc446d397b7f5435f60a">MatrixXf::Random</a>(100,3);</div>
<div class="line"><a class="code" href="group__matrixtypedefs.html#gabab09c32e96cfa9829a88400627af162">MatrixXf</a> y = <a class="code" href="classEigen_1_1DenseBase.html#a8e759dafdd9ecc446d397b7f5435f60a">MatrixXf::Random</a>(100,3);</div>
<div class="line"><a class="code" href="group__matrixtypedefs.html#gabab09c32e96cfa9829a88400627af162">MatrixXf</a> C;</div>
<div class="line">cov(x, y, C);</div>
</div><!-- fragment --><p> This is not the case anymore, when we are using an implementation taking <a class="el" href="classEigen_1_1MatrixBase.html" title="Base class for all dense matrices, vectors, and expressions. ">MatrixBase</a> as a parameter. In general, Eigen supports automatic resizing but it is not possible to do so on expressions. Why should resizing of a matrix <a class="el" href="classEigen_1_1Block.html" title="Expression of a fixed-size or dynamic-size block. ">Block</a> be allowed? It is a reference to a sub-matrix and we definitely don't want to resize that. So how can we incorporate resizing if we cannot resize on <a class="el" href="classEigen_1_1MatrixBase.html" title="Base class for all dense matrices, vectors, and expressions. ">MatrixBase</a>? The solution is to resize the derived object as in this implementation. </p>
<div class="fragment"><div class="line"><span class="keyword">template</span> &lt;<span class="keyword">typename</span> Derived, <span class="keyword">typename</span> OtherDerived&gt;</div>
<div class="line"><span class="keywordtype">void</span> cov(<span class="keyword">const</span> MatrixBase&lt;Derived&gt;&amp; x, <span class="keyword">const</span> MatrixBase&lt;Derived&gt;&amp; y, MatrixBase&lt;OtherDerived&gt; <span class="keyword">const</span> &amp; C_)</div>
<div class="line">{</div>
<div class="line">  <span class="keyword">typedef</span> <span class="keyword">typename</span> Derived::Scalar Scalar;</div>
<div class="line">  <span class="keyword">typedef</span> <span class="keyword">typename</span> internal::plain_row_type&lt;Derived&gt;::type RowVectorType;</div>
<div class="line"></div>
<div class="line">  <span class="keyword">const</span> Scalar num_observations = <span class="keyword">static_cast&lt;</span>Scalar<span class="keyword">&gt;</span>(x.rows());</div>
<div class="line"></div>
<div class="line">  <span class="keyword">const</span> RowVectorType x_mean = x.colwise().sum() / num_observations;</div>
<div class="line">  <span class="keyword">const</span> RowVectorType y_mean = y.colwise().sum() / num_observations;</div>
<div class="line"></div>
<div class="line">  MatrixBase&lt;OtherDerived&gt;&amp; C = <span class="keyword">const_cast&lt;</span> MatrixBase&lt;OtherDerived&gt;&amp; <span class="keyword">&gt;</span>(C_);</div>
<div class="line">  </div>
<div class="line">  C.derived().resize(x.cols(),x.cols()); <span class="comment">// resize the derived object</span></div>
<div class="line">  C = (x.rowwise() - x_mean).transpose() * (y.rowwise() - y_mean) / num_observations;</div>
<div class="line">}</div>
</div><!-- fragment --><p> This implementation is now working for parameters being expressions and for parameters being matrices and having the wrong size. Resizing the expressions does not do any harm in this case unless they actually require resizing. That means, passing an expression with the wrong dimensions will result in a run-time error (in debug mode only) while passing expressions of the correct size will just work fine.</p>
<p><b>Note:</b> In the above discussion the terms <a class="el" href="classEigen_1_1Matrix.html" title="The matrix class, also used for vectors and row-vectors. ">Matrix</a> and <a class="el" href="classEigen_1_1Array.html" title="General-purpose arrays with easy API for coefficient-wise operations. ">Array</a> and <a class="el" href="classEigen_1_1MatrixBase.html" title="Base class for all dense matrices, vectors, and expressions. ">MatrixBase</a> and <a class="el" href="classEigen_1_1ArrayBase.html" title="Base class for all 1D and 2D array, and related expressions. ">ArrayBase</a> can be exchanged and all arguments still hold.</p>
<h1><a class="anchor" id="TopicSummary"></a>
Summary</h1>
<ul>
<li>To summarize, the implementation of functions taking non-writable (const referenced) objects is not a big issue and does not lead to problematic situations in terms of compiling and running your program. However, a naive implementation is likely to introduce unnecessary temporary objects in your code. In order to avoid evaluating parameters into temporaries, pass them as (const) references to <a class="el" href="classEigen_1_1MatrixBase.html" title="Base class for all dense matrices, vectors, and expressions. ">MatrixBase</a> or <a class="el" href="classEigen_1_1ArrayBase.html" title="Base class for all 1D and 2D array, and related expressions. ">ArrayBase</a> (so templatize your function).</li>
<li>Functions taking writable (non-const) parameters must take const references and cast away constness within the function body.</li>
<li>Functions that take as parameters <a class="el" href="classEigen_1_1MatrixBase.html" title="Base class for all dense matrices, vectors, and expressions. ">MatrixBase</a> (or <a class="el" href="classEigen_1_1ArrayBase.html" title="Base class for all 1D and 2D array, and related expressions. ">ArrayBase</a>) objects, and potentially need to resize them (in the case where they are resizable), must call resize() on the derived class, as returned by derived(). </li>
</ul>
</div></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:26 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>