<!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: VectorwiseOp< ExpressionType, Direction > Class Template Reference</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('classEigen_1_1VectorwiseOp.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="summary"> <a href="classEigen_1_1VectorwiseOp-members.html">List of all members</a> | <a href="#pub-methods">Public Member Functions</a> </div> <div class="headertitle"> <div class="title">VectorwiseOp< ExpressionType, Direction > Class Template Reference<div class="ingroups"><a class="el" href="group__Core__Module.html">Core module</a></div></div> </div> </div><!--header--> <div class="contents"> <a name="details" id="details"></a><h2 class="groupheader">Detailed Description</h2> <div class="textblock"><h3>template<typename ExpressionType, int Direction><br/> class Eigen::VectorwiseOp< ExpressionType, Direction ></h3> <p>Pseudo expression providing partial reduction operations. </p> <dl class="params"><dt>Parameters</dt><dd> <table class="params"> <tr><td class="paramname">ExpressionType</td><td>the type of the object on which to do partial reductions </td></tr> <tr><td class="paramname">Direction</td><td>indicates the direction of the redux (<a class="el" href="group__enums.html#gga8ef30fa9c08e08c8706653571f9f5b81a2434cd8c1a594a4cdaa250f86639c600">Vertical</a> or <a class="el" href="group__enums.html#gga8ef30fa9c08e08c8706653571f9f5b81a7d5f78c516bedc0a066182a6fd606b8b">Horizontal</a>)</td></tr> </table> </dd> </dl> <p>This class represents a pseudo expression with partial reduction features. It is the return type of <a class="el" href="classEigen_1_1DenseBase.html#abe7ae69362c464b6721adbb47c655874">DenseBase::colwise()</a> and <a class="el" href="classEigen_1_1DenseBase.html#a8a7fd1e8004d4bd93a7ea36957aa8e99">DenseBase::rowwise()</a> and most of the time this is the only way it is used.</p> <p>Example: </p> <div class="fragment"><div class="line"><a class="code" href="group__matrixtypedefs.html#ga45a14b423c687c3e2e8325f148e27bf3">Matrix3d</a> m = <a class="code" href="classEigen_1_1DenseBase.html#a8e759dafdd9ecc446d397b7f5435f60a">Matrix3d::Random</a>();</div> <div class="line">cout << <span class="stringliteral">"Here is the matrix m:"</span> << endl << m << endl;</div> <div class="line">cout << <span class="stringliteral">"Here is the sum of each column:"</span> << endl << m.colwise().sum() << endl;</div> <div class="line">cout << <span class="stringliteral">"Here is the maximum absolute value of each column:"</span></div> <div class="line"> << endl << m.cwiseAbs().colwise().maxCoeff() << endl;</div> </div><!-- fragment --><p> Output: </p> <pre class="fragment">Here is the matrix m: 0.68 0.597 -0.33 -0.211 0.823 0.536 0.566 -0.605 -0.444 Here is the sum of each column: 1.04 0.815 -0.238 Here is the maximum absolute value of each column: 0.68 0.823 0.536 </pre><dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1DenseBase.html#abe7ae69362c464b6721adbb47c655874">DenseBase::colwise()</a>, <a class="el" href="classEigen_1_1DenseBase.html#a8a7fd1e8004d4bd93a7ea36957aa8e99">DenseBase::rowwise()</a>, class <a class="el" href="classEigen_1_1PartialReduxExpr.html" title="Generic expression of a partially reduxed matrix. ">PartialReduxExpr</a> </dd></dl> </div><table class="memberdecls"> <tr class="heading"><td colspan="2"><h2 class="groupheader"><a name="pub-methods"></a> Public Member Functions</h2></td></tr> <tr class="memitem:a3535582e42f802aac9b2f02316dbdd38"><td class="memItemLeft" align="right" valign="top">const ReturnType<br class="typebreak"/> < internal::member_all >::Type </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1VectorwiseOp.html#a3535582e42f802aac9b2f02316dbdd38">all</a> () const </td></tr> <tr class="separator:a3535582e42f802aac9b2f02316dbdd38"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a68ab4a7216ecfeeef6073e3879bd70c3"><td class="memItemLeft" align="right" valign="top">const ReturnType<br class="typebreak"/> < internal::member_any >::Type </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1VectorwiseOp.html#a68ab4a7216ecfeeef6073e3879bd70c3">any</a> () const </td></tr> <tr class="separator:a68ab4a7216ecfeeef6073e3879bd70c3"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a8804de1a54a51b252c394fb07f19d94a"><td class="memItemLeft" align="right" valign="top">const ReturnType<br class="typebreak"/> < internal::member_blueNorm, <br class="typebreak"/> RealScalar >::Type </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1VectorwiseOp.html#a8804de1a54a51b252c394fb07f19d94a">blueNorm</a> () const </td></tr> <tr class="separator:a8804de1a54a51b252c394fb07f19d94a"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a554a480bd3ca0f3a0113ecef4153b5b7"><td class="memItemLeft" align="right" valign="top">const <a class="el" href="classEigen_1_1PartialReduxExpr.html">PartialReduxExpr</a><br class="typebreak"/> < ExpressionType, <br class="typebreak"/> internal::member_count< Index ><br class="typebreak"/> , Direction > </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1VectorwiseOp.html#a554a480bd3ca0f3a0113ecef4153b5b7">count</a> () const </td></tr> <tr class="separator:a554a480bd3ca0f3a0113ecef4153b5b7"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a93a7a9d23fc28219f5dd82f3f90222ca"><td class="memTemplParams" colspan="2">template<typename OtherDerived > </td></tr> <tr class="memitem:a93a7a9d23fc28219f5dd82f3f90222ca"><td class="memTemplItemLeft" align="right" valign="top">const <a class="el" href="classEigen_1_1VectorwiseOp.html">VectorwiseOp</a><br class="typebreak"/> < ExpressionType, Direction ><br class="typebreak"/> ::CrossReturnType </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="classEigen_1_1VectorwiseOp.html#a93a7a9d23fc28219f5dd82f3f90222ca">cross</a> (const <a class="el" href="classEigen_1_1MatrixBase.html">MatrixBase</a>< OtherDerived > &other) const </td></tr> <tr class="separator:a93a7a9d23fc28219f5dd82f3f90222ca"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a1162b67a2744d546a6d22ecf56f47131"><td class="memItemLeft" align="right" valign="top">const <a class="el" href="classEigen_1_1CwiseBinaryOp.html">HNormalizedReturnType</a> </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1VectorwiseOp.html#a1162b67a2744d546a6d22ecf56f47131">hnormalized</a> () const </td></tr> <tr class="separator:a1162b67a2744d546a6d22ecf56f47131"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a635f8f9087821d4b3dcebd2997b5f9bd"><td class="memItemLeft" align="right" valign="top"><a class="el" href="classEigen_1_1Homogeneous.html">Homogeneous</a>< ExpressionType, <br class="typebreak"/> Direction > </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1VectorwiseOp.html#a635f8f9087821d4b3dcebd2997b5f9bd">homogeneous</a> () const </td></tr> <tr class="separator:a635f8f9087821d4b3dcebd2997b5f9bd"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:ac478d60f4ed6bbfb0833b704408999d7"><td class="memItemLeft" align="right" valign="top">const ReturnType<br class="typebreak"/> < internal::member_hypotNorm, <br class="typebreak"/> RealScalar >::Type </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1VectorwiseOp.html#ac478d60f4ed6bbfb0833b704408999d7">hypotNorm</a> () const </td></tr> <tr class="separator:ac478d60f4ed6bbfb0833b704408999d7"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:ae45bc1f13fd426ca7bea74ae4b8f99d6"><td class="memItemLeft" align="right" valign="top">const ReturnType<br class="typebreak"/> < internal::member_maxCoeff ><br class="typebreak"/> ::Type </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1VectorwiseOp.html#ae45bc1f13fd426ca7bea74ae4b8f99d6">maxCoeff</a> () const </td></tr> <tr class="separator:ae45bc1f13fd426ca7bea74ae4b8f99d6"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a3edf1e726d7e3f30d6c9655219dffbe5"><td class="memItemLeft" align="right" valign="top">const ReturnType<br class="typebreak"/> < internal::member_mean ><br class="typebreak"/> ::Type </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1VectorwiseOp.html#a3edf1e726d7e3f30d6c9655219dffbe5">mean</a> () const </td></tr> <tr class="separator:a3edf1e726d7e3f30d6c9655219dffbe5"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:af53dccc8ce1bcbfc641fc5d155e4bfe8"><td class="memItemLeft" align="right" valign="top">const ReturnType<br class="typebreak"/> < internal::member_minCoeff ><br class="typebreak"/> ::Type </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1VectorwiseOp.html#af53dccc8ce1bcbfc641fc5d155e4bfe8">minCoeff</a> () const </td></tr> <tr class="separator:af53dccc8ce1bcbfc641fc5d155e4bfe8"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a45bccb7e3e8031ebc99fdc4898a6806c"><td class="memItemLeft" align="right" valign="top">const ReturnType<br class="typebreak"/> < internal::member_norm, <br class="typebreak"/> RealScalar >::Type </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1VectorwiseOp.html#a45bccb7e3e8031ebc99fdc4898a6806c">norm</a> () const </td></tr> <tr class="separator:a45bccb7e3e8031ebc99fdc4898a6806c"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:acd0de676568888d848beb97dcc53ae47"><td class="memItemLeft" align="right" valign="top">void </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1VectorwiseOp.html#acd0de676568888d848beb97dcc53ae47">normalize</a> ()</td></tr> <tr class="separator:acd0de676568888d848beb97dcc53ae47"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:ac0d3937fc9cba6f251995fbd07110216"><td class="memItemLeft" align="right" valign="top"><a class="el" href="classEigen_1_1CwiseBinaryOp.html">CwiseBinaryOp</a><br class="typebreak"/> < internal::scalar_quotient_op<br class="typebreak"/> < Scalar >, const <br class="typebreak"/> ExpressionTypeNestedCleaned, <br class="typebreak"/> const typename <br class="typebreak"/> OppositeExtendedType< typename <br class="typebreak"/> ReturnType<br class="typebreak"/> < internal::member_norm, <br class="typebreak"/> RealScalar >::Type >::Type > </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1VectorwiseOp.html#ac0d3937fc9cba6f251995fbd07110216">normalized</a> () const </td></tr> <tr class="separator:ac0d3937fc9cba6f251995fbd07110216"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:ac154305f4f56b0f78d67c8aa1c52e87e"><td class="memTemplParams" colspan="2">template<typename OtherDerived > </td></tr> <tr class="memitem:ac154305f4f56b0f78d67c8aa1c52e87e"><td class="memTemplItemLeft" align="right" valign="top"><a class="el" href="classEigen_1_1CwiseBinaryOp.html">CwiseBinaryOp</a><br class="typebreak"/> < internal::scalar_product_op<br class="typebreak"/> < Scalar >, const <br class="typebreak"/> ExpressionTypeNestedCleaned, <br class="typebreak"/> const typename ExtendedType<br class="typebreak"/> < OtherDerived >::Type > </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="classEigen_1_1VectorwiseOp.html#ac154305f4f56b0f78d67c8aa1c52e87e">operator*</a> (const <a class="el" href="classEigen_1_1DenseBase.html">DenseBase</a>< OtherDerived > &other) const </td></tr> <tr class="separator:ac154305f4f56b0f78d67c8aa1c52e87e"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:abb151c90a593c62a6738e034f9218a23"><td class="memTemplParams" colspan="2">template<typename OtherDerived > </td></tr> <tr class="memitem:abb151c90a593c62a6738e034f9218a23"><td class="memTemplItemLeft" align="right" valign="top">ExpressionType & </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="classEigen_1_1VectorwiseOp.html#abb151c90a593c62a6738e034f9218a23">operator*=</a> (const <a class="el" href="classEigen_1_1DenseBase.html">DenseBase</a>< OtherDerived > &other)</td></tr> <tr class="separator:abb151c90a593c62a6738e034f9218a23"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:ade112b4c47a096790a047b6994d8d982"><td class="memTemplParams" colspan="2">template<typename OtherDerived > </td></tr> <tr class="memitem:ade112b4c47a096790a047b6994d8d982"><td class="memTemplItemLeft" align="right" valign="top"><a class="el" href="classEigen_1_1CwiseBinaryOp.html">CwiseBinaryOp</a><br class="typebreak"/> < internal::scalar_sum_op<br class="typebreak"/> < Scalar >, const <br class="typebreak"/> ExpressionTypeNestedCleaned, <br class="typebreak"/> const typename ExtendedType<br class="typebreak"/> < OtherDerived >::Type > </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="classEigen_1_1VectorwiseOp.html#ade112b4c47a096790a047b6994d8d982">operator+</a> (const <a class="el" href="classEigen_1_1DenseBase.html">DenseBase</a>< OtherDerived > &other) const </td></tr> <tr class="separator:ade112b4c47a096790a047b6994d8d982"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:af3ddaf5c89daf6fbe41ffb65f91fafdd"><td class="memTemplParams" colspan="2">template<typename OtherDerived > </td></tr> <tr class="memitem:af3ddaf5c89daf6fbe41ffb65f91fafdd"><td class="memTemplItemLeft" align="right" valign="top">ExpressionType & </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="classEigen_1_1VectorwiseOp.html#af3ddaf5c89daf6fbe41ffb65f91fafdd">operator+=</a> (const <a class="el" href="classEigen_1_1DenseBase.html">DenseBase</a>< OtherDerived > &other)</td></tr> <tr class="separator:af3ddaf5c89daf6fbe41ffb65f91fafdd"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a146098648dcf02be5c827938e9cfcc2d"><td class="memTemplParams" colspan="2">template<typename OtherDerived > </td></tr> <tr class="memitem:a146098648dcf02be5c827938e9cfcc2d"><td class="memTemplItemLeft" align="right" valign="top"><a class="el" href="classEigen_1_1CwiseBinaryOp.html">CwiseBinaryOp</a><br class="typebreak"/> < internal::scalar_difference_op<br class="typebreak"/> < Scalar >, const <br class="typebreak"/> ExpressionTypeNestedCleaned, <br class="typebreak"/> const typename ExtendedType<br class="typebreak"/> < OtherDerived >::Type > </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="classEigen_1_1VectorwiseOp.html#a146098648dcf02be5c827938e9cfcc2d">operator-</a> (const <a class="el" href="classEigen_1_1DenseBase.html">DenseBase</a>< OtherDerived > &other) const </td></tr> <tr class="separator:a146098648dcf02be5c827938e9cfcc2d"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a6d1e24810f6b86de45a8c2323aeafd19"><td class="memTemplParams" colspan="2">template<typename OtherDerived > </td></tr> <tr class="memitem:a6d1e24810f6b86de45a8c2323aeafd19"><td class="memTemplItemLeft" align="right" valign="top">ExpressionType & </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="classEigen_1_1VectorwiseOp.html#a6d1e24810f6b86de45a8c2323aeafd19">operator-=</a> (const <a class="el" href="classEigen_1_1DenseBase.html">DenseBase</a>< OtherDerived > &other)</td></tr> <tr class="separator:a6d1e24810f6b86de45a8c2323aeafd19"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a9234e61ac3a965ca03ce3a10ec5d2821"><td class="memTemplParams" colspan="2">template<typename OtherDerived > </td></tr> <tr class="memitem:a9234e61ac3a965ca03ce3a10ec5d2821"><td class="memTemplItemLeft" align="right" valign="top"><a class="el" href="classEigen_1_1CwiseBinaryOp.html">CwiseBinaryOp</a><br class="typebreak"/> < internal::scalar_quotient_op<br class="typebreak"/> < Scalar >, const <br class="typebreak"/> ExpressionTypeNestedCleaned, <br class="typebreak"/> const typename ExtendedType<br class="typebreak"/> < OtherDerived >::Type > </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="classEigen_1_1VectorwiseOp.html#a9234e61ac3a965ca03ce3a10ec5d2821">operator/</a> (const <a class="el" href="classEigen_1_1DenseBase.html">DenseBase</a>< OtherDerived > &other) const </td></tr> <tr class="separator:a9234e61ac3a965ca03ce3a10ec5d2821"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:abf8e26e47a2eb1a7d71e441a7c6cf455"><td class="memTemplParams" colspan="2">template<typename OtherDerived > </td></tr> <tr class="memitem:abf8e26e47a2eb1a7d71e441a7c6cf455"><td class="memTemplItemLeft" align="right" valign="top">ExpressionType & </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="classEigen_1_1VectorwiseOp.html#abf8e26e47a2eb1a7d71e441a7c6cf455">operator/=</a> (const <a class="el" href="classEigen_1_1DenseBase.html">DenseBase</a>< OtherDerived > &other)</td></tr> <tr class="separator:abf8e26e47a2eb1a7d71e441a7c6cf455"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a27b68a5f10afab40c9a6b63f0350a557"><td class="memTemplParams" colspan="2">template<typename OtherDerived > </td></tr> <tr class="memitem:a27b68a5f10afab40c9a6b63f0350a557"><td class="memTemplItemLeft" align="right" valign="top">ExpressionType & </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="classEigen_1_1VectorwiseOp.html#a27b68a5f10afab40c9a6b63f0350a557">operator=</a> (const <a class="el" href="classEigen_1_1DenseBase.html">DenseBase</a>< OtherDerived > &other)</td></tr> <tr class="separator:a27b68a5f10afab40c9a6b63f0350a557"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a64b778a6e321b207f62530fcc3ed690f"><td class="memItemLeft" align="right" valign="top">const ReturnType<br class="typebreak"/> < internal::member_prod ><br class="typebreak"/> ::Type </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1VectorwiseOp.html#a64b778a6e321b207f62530fcc3ed690f">prod</a> () const </td></tr> <tr class="separator:a64b778a6e321b207f62530fcc3ed690f"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a12ed8da7c7d83cdc68e9cbbde358a164"><td class="memTemplParams" colspan="2">template<typename BinaryOp > </td></tr> <tr class="memitem:a12ed8da7c7d83cdc68e9cbbde358a164"><td class="memTemplItemLeft" align="right" valign="top">const ReduxReturnType<br class="typebreak"/> < BinaryOp >::Type </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="classEigen_1_1VectorwiseOp.html#a12ed8da7c7d83cdc68e9cbbde358a164">redux</a> (const BinaryOp &func=BinaryOp()) const </td></tr> <tr class="separator:a12ed8da7c7d83cdc68e9cbbde358a164"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a9a078997c488c2d83fbd45108af26e1d"><td class="memItemLeft" align="right" valign="top">const <a class="el" href="classEigen_1_1Replicate.html">ReplicateReturnType</a> </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1VectorwiseOp.html#a9a078997c488c2d83fbd45108af26e1d">replicate</a> (Index factor) const </td></tr> <tr class="separator:a9a078997c488c2d83fbd45108af26e1d"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a5a64491d1ee839cf717328f17a1ffc45"><td class="memTemplParams" colspan="2">template<int Factor> </td></tr> <tr class="memitem:a5a64491d1ee839cf717328f17a1ffc45"><td class="memTemplItemLeft" align="right" valign="top">const <a class="el" href="classEigen_1_1Replicate.html">Replicate</a><br class="typebreak"/> < ExpressionType,(IsVertical?Factor:1),(IsHorizontal?Factor:1)> </td><td class="memTemplItemRight" valign="bottom"><a class="el" href="classEigen_1_1VectorwiseOp.html#a5a64491d1ee839cf717328f17a1ffc45">replicate</a> (Index factor=Factor) const </td></tr> <tr class="separator:a5a64491d1ee839cf717328f17a1ffc45"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:a04e38b3a68fe5ae6908594f3b1ec618c"><td class="memItemLeft" align="right" valign="top">const <a class="el" href="classEigen_1_1Reverse.html">Reverse</a>< ExpressionType, <br class="typebreak"/> Direction > </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1VectorwiseOp.html#a04e38b3a68fe5ae6908594f3b1ec618c">reverse</a> () const </td></tr> <tr class="separator:a04e38b3a68fe5ae6908594f3b1ec618c"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:ad56e72700d9ee45893c41533d8c46945"><td class="memItemLeft" align="right" valign="top">const ReturnType<br class="typebreak"/> < internal::member_squaredNorm, <br class="typebreak"/> RealScalar >::Type </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1VectorwiseOp.html#ad56e72700d9ee45893c41533d8c46945">squaredNorm</a> () const </td></tr> <tr class="separator:ad56e72700d9ee45893c41533d8c46945"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:ac32ca2926dfc6a21859c0aa9e8eb1bcb"><td class="memItemLeft" align="right" valign="top">const ReturnType<br class="typebreak"/> < internal::member_stableNorm, <br class="typebreak"/> RealScalar >::Type </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1VectorwiseOp.html#ac32ca2926dfc6a21859c0aa9e8eb1bcb">stableNorm</a> () const </td></tr> <tr class="separator:ac32ca2926dfc6a21859c0aa9e8eb1bcb"><td class="memSeparator" colspan="2"> </td></tr> <tr class="memitem:af27cb8c9043ed6408e64f94645439362"><td class="memItemLeft" align="right" valign="top">const ReturnType<br class="typebreak"/> < internal::member_sum >::Type </td><td class="memItemRight" valign="bottom"><a class="el" href="classEigen_1_1VectorwiseOp.html#af27cb8c9043ed6408e64f94645439362">sum</a> () const </td></tr> <tr class="separator:af27cb8c9043ed6408e64f94645439362"><td class="memSeparator" colspan="2"> </td></tr> </table> <h2 class="groupheader">Member Function Documentation</h2> <a class="anchor" id="a3535582e42f802aac9b2f02316dbdd38"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">const ReturnType<internal::member_all>::Type all </td> <td>(</td> <td class="paramtype">void </td> <td class="paramname"></td><td>)</td> <td> const</td> </tr> </table> </td> <td class="mlabels-right"> <span class="mlabels"><span class="mlabel">inline</span></span> </td> </tr> </table> </div><div class="memdoc"> <dl class="section return"><dt>Returns</dt><dd>a row (or column) vector expression representing whether <b>all</b> coefficients of each respective column (or row) are <code>true</code>.</dd></dl> <dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1DenseBase.html#aea914316af61df197f21629e14e7870a">DenseBase::all()</a> </dd></dl> </div> </div> <a class="anchor" id="a68ab4a7216ecfeeef6073e3879bd70c3"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">const ReturnType<internal::member_any>::Type any </td> <td>(</td> <td class="paramtype">void </td> <td class="paramname"></td><td>)</td> <td> const</td> </tr> </table> </td> <td class="mlabels-right"> <span class="mlabels"><span class="mlabel">inline</span></span> </td> </tr> </table> </div><div class="memdoc"> <dl class="section return"><dt>Returns</dt><dd>a row (or column) vector expression representing whether <b>at</b> <b>least</b> one coefficient of each respective column (or row) is <code>true</code>.</dd></dl> <dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1DenseBase.html#a42571e028736ca9103bac8b50f269824">DenseBase::any()</a> </dd></dl> </div> </div> <a class="anchor" id="a8804de1a54a51b252c394fb07f19d94a"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">const ReturnType<internal::member_blueNorm,RealScalar>::Type blueNorm </td> <td>(</td> <td class="paramname"></td><td>)</td> <td> const</td> </tr> </table> </td> <td class="mlabels-right"> <span class="mlabels"><span class="mlabel">inline</span></span> </td> </tr> </table> </div><div class="memdoc"> <dl class="section return"><dt>Returns</dt><dd>a row (or column) vector expression of the norm of each column (or row) of the referenced expression, using blue's algorithm.</dd></dl> <dl class="section see"><dt>See Also</dt><dd>DenseBase::blueNorm() </dd></dl> </div> </div> <a class="anchor" id="a554a480bd3ca0f3a0113ecef4153b5b7"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">const <a class="el" href="classEigen_1_1PartialReduxExpr.html">PartialReduxExpr</a><ExpressionType, internal::member_count<Index>, Direction> count </td> <td>(</td> <td class="paramname"></td><td>)</td> <td> const</td> </tr> </table> </td> <td class="mlabels-right"> <span class="mlabels"><span class="mlabel">inline</span></span> </td> </tr> </table> </div><div class="memdoc"> <dl class="section return"><dt>Returns</dt><dd>a row (or column) vector expression representing the number of <code>true</code> coefficients of each respective column (or row).</dd></dl> <p>Example: </p> <div class="fragment"><div class="line"><a class="code" href="group__matrixtypedefs.html#ga45a14b423c687c3e2e8325f148e27bf3">Matrix3d</a> m = <a class="code" href="classEigen_1_1DenseBase.html#a8e759dafdd9ecc446d397b7f5435f60a">Matrix3d::Random</a>();</div> <div class="line">cout << <span class="stringliteral">"Here is the matrix m:"</span> << endl << m << endl;</div> <div class="line">cout << <span class="stringliteral">"Here is the count of elements larger or equal than 0.5 of each row:"</span> << endl << (m.array() >= 0.5).rowwise().count() << endl;</div> </div><!-- fragment --><p> Output: </p> <pre class="fragment">Here is the matrix m: 0.68 0.597 -0.33 -0.211 0.823 0.536 0.566 -0.605 -0.444 Here is the count of elements larger or equal than 0.5 of each row: 2 2 1 </pre><dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1DenseBase.html#aa671b5ea336ba21a7644d3fa6577ee00">DenseBase::count()</a> </dd></dl> </div> </div> <a class="anchor" id="a93a7a9d23fc28219f5dd82f3f90222ca"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">const <a class="el" href="classEigen_1_1VectorwiseOp.html">VectorwiseOp</a><ExpressionType,Direction>::CrossReturnType cross </td> <td>(</td> <td class="paramtype">const <a class="el" href="classEigen_1_1MatrixBase.html">MatrixBase</a>< OtherDerived > & </td> <td class="paramname"><em>other</em></td><td>)</td> <td> const</td> </tr> </table> </div><div class="memdoc"> <dl class="section return"><dt>Returns</dt><dd>a matrix expression of the cross product of each column or row of the referenced expression with the <em>other</em> vector.</dd></dl> <p>The referenced matrix must have one dimension equal to 3. The result matrix has the same dimensions than the referenced one.</p> <p>This is defined in the Geometry module.</p> <div class="fragment"><div class="line"><span class="preprocessor">#include <Eigen/Geometry></span> </div> </div><!-- fragment --><dl class="section see"><dt>See Also</dt><dd>MatrixBase::cross() </dd></dl> <p>References <a class="el" href="classEigen_1_1DenseBase.html#a58c77695de3b33405f01f2fdf3dc389d">DenseBase< Derived >::col()</a>, <a class="el" href="classEigen_1_1DenseBase.html#aa8716d44f51321072ee5c88665c28813">DenseBase< Derived >::row()</a>, and <a class="el" href="group__enums.html#gga8ef30fa9c08e08c8706653571f9f5b81a2434cd8c1a594a4cdaa250f86639c600">Eigen::Vertical</a>.</p> </div> </div> <a class="anchor" id="a1162b67a2744d546a6d22ecf56f47131"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">const <a class="el" href="classEigen_1_1VectorwiseOp.html">VectorwiseOp</a>< ExpressionType, Direction >::<a class="el" href="classEigen_1_1CwiseBinaryOp.html">HNormalizedReturnType</a> hnormalized </td> <td>(</td> <td class="paramname"></td><td>)</td> <td> const</td> </tr> </table> </td> <td class="mlabels-right"> <span class="mlabels"><span class="mlabel">inline</span></span> </td> </tr> </table> </div><div class="memdoc"> <p>This is defined in the Geometry module.</p> <div class="fragment"><div class="line"><span class="preprocessor">#include <Eigen/Geometry></span> </div> </div><!-- fragment --><dl class="section return"><dt>Returns</dt><dd>an expression of the homogeneous normalized vector of <code>*this</code> </dd></dl> <p>Example: </p> <div class="fragment"></div><!-- fragment --><p> Output: </p> <pre class="fragment"></pre><dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1MatrixBase.html#adcf645f676d798fc19417d3ce43484d0">MatrixBase::hnormalized()</a> </dd></dl> <p>References <a class="el" href="group__enums.html#gga8ef30fa9c08e08c8706653571f9f5b81a7d5f78c516bedc0a066182a6fd606b8b">Eigen::Horizontal</a>, and <a class="el" href="group__enums.html#gga8ef30fa9c08e08c8706653571f9f5b81a2434cd8c1a594a4cdaa250f86639c600">Eigen::Vertical</a>.</p> </div> </div> <a class="anchor" id="a635f8f9087821d4b3dcebd2997b5f9bd"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname"><a class="el" href="classEigen_1_1Homogeneous.html">Homogeneous</a>< ExpressionType, Direction > homogeneous </td> <td>(</td> <td class="paramname"></td><td>)</td> <td> const</td> </tr> </table> </td> <td class="mlabels-right"> <span class="mlabels"><span class="mlabel">inline</span></span> </td> </tr> </table> </div><div class="memdoc"> <p>This is defined in the Geometry module.</p> <div class="fragment"><div class="line"><span class="preprocessor">#include <Eigen/Geometry></span> </div> </div><!-- fragment --><dl class="section return"><dt>Returns</dt><dd>a matrix expression of homogeneous column (or row) vectors</dd></dl> <p>Example: </p> <div class="fragment"></div><!-- fragment --><p> Output: </p> <pre class="fragment"></pre><dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1MatrixBase.html#af3513a0abfc3875e6d1912b8bb7b28ac">MatrixBase::homogeneous()</a> </dd></dl> </div> </div> <a class="anchor" id="ac478d60f4ed6bbfb0833b704408999d7"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">const ReturnType<internal::member_hypotNorm,RealScalar>::Type hypotNorm </td> <td>(</td> <td class="paramname"></td><td>)</td> <td> const</td> </tr> </table> </td> <td class="mlabels-right"> <span class="mlabels"><span class="mlabel">inline</span></span> </td> </tr> </table> </div><div class="memdoc"> <dl class="section return"><dt>Returns</dt><dd>a row (or column) vector expression of the norm of each column (or row) of the referenced expression, avoiding underflow and overflow using a concatenation of hypot() calls.</dd></dl> <dl class="section see"><dt>See Also</dt><dd>DenseBase::hypotNorm() </dd></dl> </div> </div> <a class="anchor" id="ae45bc1f13fd426ca7bea74ae4b8f99d6"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">const ReturnType<internal::member_maxCoeff>::Type maxCoeff </td> <td>(</td> <td class="paramname"></td><td>)</td> <td> const</td> </tr> </table> </td> <td class="mlabels-right"> <span class="mlabels"><span class="mlabel">inline</span></span> </td> </tr> </table> </div><div class="memdoc"> <dl class="section return"><dt>Returns</dt><dd>a row (or column) vector expression of the largest coefficient of each column (or row) of the referenced expression.</dd></dl> <dl class="section warning"><dt>Warning</dt><dd>the result is undefined if <code>*this</code> contains NaN.</dd></dl> <p>Example: </p> <div class="fragment"><div class="line"><a class="code" href="group__matrixtypedefs.html#ga45a14b423c687c3e2e8325f148e27bf3">Matrix3d</a> m = <a class="code" href="classEigen_1_1DenseBase.html#a8e759dafdd9ecc446d397b7f5435f60a">Matrix3d::Random</a>();</div> <div class="line">cout << <span class="stringliteral">"Here is the matrix m:"</span> << endl << m << endl;</div> <div class="line">cout << <span class="stringliteral">"Here is the maximum of each column:"</span> << endl << m.colwise().maxCoeff() << endl;</div> </div><!-- fragment --><p> Output: </p> <pre class="fragment">Here is the matrix m: 0.68 0.597 -0.33 -0.211 0.823 0.536 0.566 -0.605 -0.444 Here is the maximum of each column: 0.68 0.823 0.536 </pre><dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1DenseBase.html#a878f0dae18b28d8158c5f1c232edced2">DenseBase::maxCoeff()</a> </dd></dl> </div> </div> <a class="anchor" id="a3edf1e726d7e3f30d6c9655219dffbe5"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">const ReturnType<internal::member_mean>::Type mean </td> <td>(</td> <td class="paramname"></td><td>)</td> <td> const</td> </tr> </table> </td> <td class="mlabels-right"> <span class="mlabels"><span class="mlabel">inline</span></span> </td> </tr> </table> </div><div class="memdoc"> <dl class="section return"><dt>Returns</dt><dd>a row (or column) vector expression of the mean of each column (or row) of the referenced expression.</dd></dl> <dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1DenseBase.html#a0af2b3991862a079e3efaef3e4d17d96">DenseBase::mean()</a> </dd></dl> </div> </div> <a class="anchor" id="af53dccc8ce1bcbfc641fc5d155e4bfe8"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">const ReturnType<internal::member_minCoeff>::Type minCoeff </td> <td>(</td> <td class="paramname"></td><td>)</td> <td> const</td> </tr> </table> </td> <td class="mlabels-right"> <span class="mlabels"><span class="mlabel">inline</span></span> </td> </tr> </table> </div><div class="memdoc"> <dl class="section return"><dt>Returns</dt><dd>a row (or column) vector expression of the smallest coefficient of each column (or row) of the referenced expression.</dd></dl> <dl class="section warning"><dt>Warning</dt><dd>the result is undefined if <code>*this</code> contains NaN.</dd></dl> <p>Example: </p> <div class="fragment"><div class="line"><a class="code" href="group__matrixtypedefs.html#ga45a14b423c687c3e2e8325f148e27bf3">Matrix3d</a> m = <a class="code" href="classEigen_1_1DenseBase.html#a8e759dafdd9ecc446d397b7f5435f60a">Matrix3d::Random</a>();</div> <div class="line">cout << <span class="stringliteral">"Here is the matrix m:"</span> << endl << m << endl;</div> <div class="line">cout << <span class="stringliteral">"Here is the minimum of each column:"</span> << endl << m.colwise().minCoeff() << endl;</div> </div><!-- fragment --><p> Output: </p> <pre class="fragment">Here is the matrix m: 0.68 0.597 -0.33 -0.211 0.823 0.536 0.566 -0.605 -0.444 Here is the minimum of each column: -0.211 -0.605 -0.444 </pre><dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1DenseBase.html#add6cb2d85282829eb9adc9565ce784d6">DenseBase::minCoeff()</a> </dd></dl> </div> </div> <a class="anchor" id="a45bccb7e3e8031ebc99fdc4898a6806c"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">const ReturnType<internal::member_norm,RealScalar>::Type norm </td> <td>(</td> <td class="paramname"></td><td>)</td> <td> const</td> </tr> </table> </td> <td class="mlabels-right"> <span class="mlabels"><span class="mlabel">inline</span></span> </td> </tr> </table> </div><div class="memdoc"> <dl class="section return"><dt>Returns</dt><dd>a row (or column) vector expression of the norm of each column (or row) of the referenced expression.</dd></dl> <p>Example: </p> <div class="fragment"><div class="line"><a class="code" href="group__matrixtypedefs.html#ga45a14b423c687c3e2e8325f148e27bf3">Matrix3d</a> m = <a class="code" href="classEigen_1_1DenseBase.html#a8e759dafdd9ecc446d397b7f5435f60a">Matrix3d::Random</a>();</div> <div class="line">cout << <span class="stringliteral">"Here is the matrix m:"</span> << endl << m << endl;</div> <div class="line">cout << <span class="stringliteral">"Here is the norm of each column:"</span> << endl << m.colwise().norm() << endl;</div> </div><!-- fragment --><p> Output: </p> <pre class="fragment">Here is the matrix m: 0.68 0.597 -0.33 -0.211 0.823 0.536 0.566 -0.605 -0.444 Here is the norm of each column: 0.91 1.18 0.771 </pre><dl class="section see"><dt>See Also</dt><dd>DenseBase::norm() </dd></dl> </div> </div> <a class="anchor" id="acd0de676568888d848beb97dcc53ae47"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">void normalize </td> <td>(</td> <td class="paramname"></td><td>)</td> <td></td> </tr> </table> </td> <td class="mlabels-right"> <span class="mlabels"><span class="mlabel">inline</span></span> </td> </tr> </table> </div><div class="memdoc"> <p>Normalize in-place each row or columns of the referenced matrix. </p> <dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1MatrixBase.html#acd0de676568888d848beb97dcc53ae47">MatrixBase::normalize()</a>, <a class="el" href="classEigen_1_1VectorwiseOp.html#ac0d3937fc9cba6f251995fbd07110216">normalized()</a> </dd></dl> </div> </div> <a class="anchor" id="ac0d3937fc9cba6f251995fbd07110216"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname"><a class="el" href="classEigen_1_1CwiseBinaryOp.html">CwiseBinaryOp</a><internal::scalar_quotient_op<Scalar>, const ExpressionTypeNestedCleaned, const typename OppositeExtendedType<typename ReturnType<internal::member_norm,RealScalar>::Type>::Type> normalized </td> <td>(</td> <td class="paramname"></td><td>)</td> <td> const</td> </tr> </table> </td> <td class="mlabels-right"> <span class="mlabels"><span class="mlabel">inline</span></span> </td> </tr> </table> </div><div class="memdoc"> <dl class="section return"><dt>Returns</dt><dd>an expression where each column of row of the referenced matrix are normalized. The referenced matrix is <b>not</b> modified. </dd></dl> <dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1MatrixBase.html#a8ed1fb2e792b1079639a74e3581fbc74">MatrixBase::normalized()</a>, <a class="el" href="classEigen_1_1VectorwiseOp.html#acd0de676568888d848beb97dcc53ae47">normalize()</a> </dd></dl> </div> </div> <a class="anchor" id="ac154305f4f56b0f78d67c8aa1c52e87e"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname"><a class="el" href="classEigen_1_1CwiseBinaryOp.html">CwiseBinaryOp</a><internal::scalar_product_op<Scalar>, const ExpressionTypeNestedCleaned, const typename ExtendedType<OtherDerived>::Type> operator* </td> <td>(</td> <td class="paramtype">const <a class="el" href="classEigen_1_1DenseBase.html">DenseBase</a>< OtherDerived > & </td> <td class="paramname"><em>other</em></td><td>)</td> <td> const</td> </tr> </table> </td> <td class="mlabels-right"> <span class="mlabels"><span class="mlabel">inline</span></span> </td> </tr> </table> </div><div class="memdoc"> <p>Returns the expression where each subvector is the product of the vector <em>other</em> by the corresponding subvector of <code>*this</code> </p> </div> </div> <a class="anchor" id="abb151c90a593c62a6738e034f9218a23"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">ExpressionType& operator*= </td> <td>(</td> <td class="paramtype">const <a class="el" href="classEigen_1_1DenseBase.html">DenseBase</a>< OtherDerived > & </td> <td class="paramname"><em>other</em></td><td>)</td> <td></td> </tr> </table> </td> <td class="mlabels-right"> <span class="mlabels"><span class="mlabel">inline</span></span> </td> </tr> </table> </div><div class="memdoc"> <p>Multiples each subvector of <code>*this</code> by the vector <em>other</em> </p> </div> </div> <a class="anchor" id="ade112b4c47a096790a047b6994d8d982"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname"><a class="el" href="classEigen_1_1CwiseBinaryOp.html">CwiseBinaryOp</a><internal::scalar_sum_op<Scalar>, const ExpressionTypeNestedCleaned, const typename ExtendedType<OtherDerived>::Type> operator+ </td> <td>(</td> <td class="paramtype">const <a class="el" href="classEigen_1_1DenseBase.html">DenseBase</a>< OtherDerived > & </td> <td class="paramname"><em>other</em></td><td>)</td> <td> const</td> </tr> </table> </td> <td class="mlabels-right"> <span class="mlabels"><span class="mlabel">inline</span></span> </td> </tr> </table> </div><div class="memdoc"> <p>Returns the expression of the sum of the vector <em>other</em> to each subvector of <code>*this</code> </p> </div> </div> <a class="anchor" id="af3ddaf5c89daf6fbe41ffb65f91fafdd"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">ExpressionType& operator+= </td> <td>(</td> <td class="paramtype">const <a class="el" href="classEigen_1_1DenseBase.html">DenseBase</a>< OtherDerived > & </td> <td class="paramname"><em>other</em></td><td>)</td> <td></td> </tr> </table> </td> <td class="mlabels-right"> <span class="mlabels"><span class="mlabel">inline</span></span> </td> </tr> </table> </div><div class="memdoc"> <p>Adds the vector <em>other</em> to each subvector of <code>*this</code> </p> </div> </div> <a class="anchor" id="a146098648dcf02be5c827938e9cfcc2d"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname"><a class="el" href="classEigen_1_1CwiseBinaryOp.html">CwiseBinaryOp</a><internal::scalar_difference_op<Scalar>, const ExpressionTypeNestedCleaned, const typename ExtendedType<OtherDerived>::Type> operator- </td> <td>(</td> <td class="paramtype">const <a class="el" href="classEigen_1_1DenseBase.html">DenseBase</a>< OtherDerived > & </td> <td class="paramname"><em>other</em></td><td>)</td> <td> const</td> </tr> </table> </td> <td class="mlabels-right"> <span class="mlabels"><span class="mlabel">inline</span></span> </td> </tr> </table> </div><div class="memdoc"> <p>Returns the expression of the difference between each subvector of <code>*this</code> and the vector <em>other</em> </p> </div> </div> <a class="anchor" id="a6d1e24810f6b86de45a8c2323aeafd19"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">ExpressionType& operator-= </td> <td>(</td> <td class="paramtype">const <a class="el" href="classEigen_1_1DenseBase.html">DenseBase</a>< OtherDerived > & </td> <td class="paramname"><em>other</em></td><td>)</td> <td></td> </tr> </table> </td> <td class="mlabels-right"> <span class="mlabels"><span class="mlabel">inline</span></span> </td> </tr> </table> </div><div class="memdoc"> <p>Substracts the vector <em>other</em> to each subvector of <code>*this</code> </p> </div> </div> <a class="anchor" id="a9234e61ac3a965ca03ce3a10ec5d2821"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname"><a class="el" href="classEigen_1_1CwiseBinaryOp.html">CwiseBinaryOp</a><internal::scalar_quotient_op<Scalar>, const ExpressionTypeNestedCleaned, const typename ExtendedType<OtherDerived>::Type> operator/ </td> <td>(</td> <td class="paramtype">const <a class="el" href="classEigen_1_1DenseBase.html">DenseBase</a>< OtherDerived > & </td> <td class="paramname"><em>other</em></td><td>)</td> <td> const</td> </tr> </table> </td> <td class="mlabels-right"> <span class="mlabels"><span class="mlabel">inline</span></span> </td> </tr> </table> </div><div class="memdoc"> <p>Returns the expression where each subvector is the quotient of the corresponding subvector of <code>*this</code> by the vector <em>other</em> </p> </div> </div> <a class="anchor" id="abf8e26e47a2eb1a7d71e441a7c6cf455"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">ExpressionType& operator/= </td> <td>(</td> <td class="paramtype">const <a class="el" href="classEigen_1_1DenseBase.html">DenseBase</a>< OtherDerived > & </td> <td class="paramname"><em>other</em></td><td>)</td> <td></td> </tr> </table> </td> <td class="mlabels-right"> <span class="mlabels"><span class="mlabel">inline</span></span> </td> </tr> </table> </div><div class="memdoc"> <p>Divides each subvector of <code>*this</code> by the vector <em>other</em> </p> </div> </div> <a class="anchor" id="a27b68a5f10afab40c9a6b63f0350a557"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">ExpressionType& operator= </td> <td>(</td> <td class="paramtype">const <a class="el" href="classEigen_1_1DenseBase.html">DenseBase</a>< OtherDerived > & </td> <td class="paramname"><em>other</em></td><td>)</td> <td></td> </tr> </table> </td> <td class="mlabels-right"> <span class="mlabels"><span class="mlabel">inline</span></span> </td> </tr> </table> </div><div class="memdoc"> <p>Copies the vector <em>other</em> to each subvector of <code>*this</code> </p> </div> </div> <a class="anchor" id="a64b778a6e321b207f62530fcc3ed690f"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">const ReturnType<internal::member_prod>::Type prod </td> <td>(</td> <td class="paramname"></td><td>)</td> <td> const</td> </tr> </table> </td> <td class="mlabels-right"> <span class="mlabels"><span class="mlabel">inline</span></span> </td> </tr> </table> </div><div class="memdoc"> <dl class="section return"><dt>Returns</dt><dd>a row (or column) vector expression of the product of each column (or row) of the referenced expression.</dd></dl> <p>Example: </p> <div class="fragment"><div class="line"><a class="code" href="group__matrixtypedefs.html#ga45a14b423c687c3e2e8325f148e27bf3">Matrix3d</a> m = <a class="code" href="classEigen_1_1DenseBase.html#a8e759dafdd9ecc446d397b7f5435f60a">Matrix3d::Random</a>();</div> <div class="line">cout << <span class="stringliteral">"Here is the matrix m:"</span> << endl << m << endl;</div> <div class="line">cout << <span class="stringliteral">"Here is the product of each row:"</span> << endl << m.rowwise().prod() << endl;</div> </div><!-- fragment --><p> Output: </p> <pre class="fragment">Here is the matrix m: 0.68 0.597 -0.33 -0.211 0.823 0.536 0.566 -0.605 -0.444 Here is the product of each row: -0.134 -0.0933 0.152 </pre><dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1DenseBase.html#a6bdcbfa7e3b07d3246ad80de7170b0f5">DenseBase::prod()</a> </dd></dl> </div> </div> <a class="anchor" id="a12ed8da7c7d83cdc68e9cbbde358a164"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">const ReduxReturnType<BinaryOp>::Type redux </td> <td>(</td> <td class="paramtype">const BinaryOp & </td> <td class="paramname"><em>func</em> = <code>BinaryOp()</code></td><td>)</td> <td> const</td> </tr> </table> </td> <td class="mlabels-right"> <span class="mlabels"><span class="mlabel">inline</span></span> </td> </tr> </table> </div><div class="memdoc"> <dl class="section return"><dt>Returns</dt><dd>a row or column vector expression of <code>*this</code> reduxed by <em>func</em> </dd></dl> <p>The template parameter <em>BinaryOp</em> is the type of the functor of the custom redux operator. Note that func must be an associative operator.</p> <dl class="section see"><dt>See Also</dt><dd>class <a class="el" href="classEigen_1_1VectorwiseOp.html" title="Pseudo expression providing partial reduction operations. ">VectorwiseOp</a>, <a class="el" href="classEigen_1_1DenseBase.html#abe7ae69362c464b6721adbb47c655874">DenseBase::colwise()</a>, <a class="el" href="classEigen_1_1DenseBase.html#a8a7fd1e8004d4bd93a7ea36957aa8e99">DenseBase::rowwise()</a> </dd></dl> </div> </div> <a class="anchor" id="a9a078997c488c2d83fbd45108af26e1d"></a> <div class="memitem"> <div class="memproto"> <table class="memname"> <tr> <td class="memname">const <a class="el" href="classEigen_1_1VectorwiseOp.html">VectorwiseOp</a>< ExpressionType, Direction >::<a class="el" href="classEigen_1_1Replicate.html">ReplicateReturnType</a> replicate </td> <td>(</td> <td class="paramtype">Index </td> <td class="paramname"><em>factor</em></td><td>)</td> <td> const</td> </tr> </table> </div><div class="memdoc"> <dl class="section return"><dt>Returns</dt><dd>an expression of the replication of each column (or row) of <code>*this</code> </dd></dl> <p>Example: </p> <div class="fragment"><div class="line"><a class="code" href="group__matrixtypedefs.html#gaa4931fe5bb599038466be823fdfadd04">Vector3i</a> v = <a class="code" href="classEigen_1_1DenseBase.html#a8e759dafdd9ecc446d397b7f5435f60a">Vector3i::Random</a>();</div> <div class="line">cout << <span class="stringliteral">"Here is the vector v:"</span> << endl << v << endl;</div> <div class="line">cout << <span class="stringliteral">"v.rowwise().replicate(5) = ..."</span> << endl;</div> <div class="line">cout << v.rowwise().replicate(5) << endl;</div> </div><!-- fragment --><p> Output: </p> <pre class="fragment">Here is the vector v: 7 -2 6 v.rowwise().replicate(5) = ... 7 7 7 7 7 -2 -2 -2 -2 -2 6 6 6 6 6 </pre><dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1VectorwiseOp.html#a9a078997c488c2d83fbd45108af26e1d">VectorwiseOp::replicate()</a>, <a class="el" href="classEigen_1_1DenseBase.html#afca0e8ff7921ee0e3ab4422818ecb214">DenseBase::replicate()</a>, class <a class="el" href="classEigen_1_1Replicate.html" title="Expression of the multiple replication of a matrix or vector. ">Replicate</a> </dd></dl> <p>References <a class="el" href="group__enums.html#gga8ef30fa9c08e08c8706653571f9f5b81a7d5f78c516bedc0a066182a6fd606b8b">Eigen::Horizontal</a>, and <a class="el" href="group__enums.html#gga8ef30fa9c08e08c8706653571f9f5b81a2434cd8c1a594a4cdaa250f86639c600">Eigen::Vertical</a>.</p> </div> </div> <a class="anchor" id="a5a64491d1ee839cf717328f17a1ffc45"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">const <a class="el" href="classEigen_1_1Replicate.html">Replicate</a><ExpressionType,(IsVertical?Factor:1),(IsHorizontal?Factor:1)> replicate </td> <td>(</td> <td class="paramtype">Index </td> <td class="paramname"><em>factor</em> = <code>Factor</code></td><td>)</td> <td> const</td> </tr> </table> </td> <td class="mlabels-right"> <span class="mlabels"><span class="mlabel">inline</span></span> </td> </tr> </table> </div><div class="memdoc"> <dl class="section return"><dt>Returns</dt><dd>an expression of the replication of each column (or row) of <code>*this</code> </dd></dl> <p>Example: </p> <div class="fragment"><div class="line"><a class="code" href="group__matrixtypedefs.html#ga5b9d0d11e78fd355d23221154b7620e4">MatrixXi</a> m = <a class="code" href="classEigen_1_1DenseBase.html#a8e759dafdd9ecc446d397b7f5435f60a">MatrixXi::Random</a>(2,3);</div> <div class="line">cout << <span class="stringliteral">"Here is the matrix m:"</span> << endl << m << endl;</div> <div class="line">cout << <span class="stringliteral">"m.colwise().replicate<3>() = ..."</span> << endl;</div> <div class="line">cout << m.colwise().replicate<3>() << endl;</div> </div><!-- fragment --><p> Output: </p> <pre class="fragment">Here is the matrix m: 7 6 9 -2 6 -6 m.colwise().replicate<3>() = ... 7 6 9 -2 6 -6 7 6 9 -2 6 -6 7 6 9 -2 6 -6 </pre><dl class="section see"><dt>See Also</dt><dd>VectorwiseOp::replicate(Index), <a class="el" href="classEigen_1_1DenseBase.html#afca0e8ff7921ee0e3ab4422818ecb214">DenseBase::replicate()</a>, class <a class="el" href="classEigen_1_1Replicate.html" title="Expression of the multiple replication of a matrix or vector. ">Replicate</a> </dd></dl> <p>References <a class="el" href="group__enums.html#gga8ef30fa9c08e08c8706653571f9f5b81a7d5f78c516bedc0a066182a6fd606b8b">Eigen::Horizontal</a>, and <a class="el" href="group__enums.html#gga8ef30fa9c08e08c8706653571f9f5b81a2434cd8c1a594a4cdaa250f86639c600">Eigen::Vertical</a>.</p> </div> </div> <a class="anchor" id="a04e38b3a68fe5ae6908594f3b1ec618c"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">const <a class="el" href="classEigen_1_1Reverse.html">Reverse</a><ExpressionType, Direction> reverse </td> <td>(</td> <td class="paramname"></td><td>)</td> <td> const</td> </tr> </table> </td> <td class="mlabels-right"> <span class="mlabels"><span class="mlabel">inline</span></span> </td> </tr> </table> </div><div class="memdoc"> <dl class="section return"><dt>Returns</dt><dd>a matrix expression where each column (or row) are reversed.</dd></dl> <p>Example: </p> <div class="fragment"><div class="line"><a class="code" href="group__matrixtypedefs.html#ga5b9d0d11e78fd355d23221154b7620e4">MatrixXi</a> m = <a class="code" href="classEigen_1_1DenseBase.html#a8e759dafdd9ecc446d397b7f5435f60a">MatrixXi::Random</a>(3,4);</div> <div class="line">cout << <span class="stringliteral">"Here is the matrix m:"</span> << endl << m << endl;</div> <div class="line">cout << <span class="stringliteral">"Here is the rowwise reverse of m:"</span> << endl << m.rowwise().reverse() << endl;</div> <div class="line">cout << <span class="stringliteral">"Here is the colwise reverse of m:"</span> << endl << m.colwise().reverse() << endl;</div> <div class="line"></div> <div class="line">cout << <span class="stringliteral">"Here is the coefficient (1,0) in the rowise reverse of m:"</span> << endl</div> <div class="line"><< m.rowwise().reverse()(1,0) << endl;</div> <div class="line">cout << <span class="stringliteral">"Let us overwrite this coefficient with the value 4."</span> << endl;</div> <div class="line"><span class="comment">//m.colwise().reverse()(1,0) = 4;</span></div> <div class="line">cout << <span class="stringliteral">"Now the matrix m is:"</span> << endl << m << endl;</div> </div><!-- fragment --><p> Output: </p> <pre class="fragment">Here is the matrix m: 7 6 -3 1 -2 9 6 0 6 -6 -5 3 Here is the rowwise reverse of m: 1 -3 6 7 0 6 9 -2 3 -5 -6 6 Here is the colwise reverse of m: 6 -6 -5 3 -2 9 6 0 7 6 -3 1 Here is the coefficient (1,0) in the rowise reverse of m: 0 Let us overwrite this coefficient with the value 4. Now the matrix m is: 7 6 -3 1 -2 9 6 0 6 -6 -5 3 </pre><dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1DenseBase.html#a6e354bb81f0c7b6888c6a0ce4b4649e2">DenseBase::reverse()</a> </dd></dl> </div> </div> <a class="anchor" id="ad56e72700d9ee45893c41533d8c46945"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">const ReturnType<internal::member_squaredNorm,RealScalar>::Type squaredNorm </td> <td>(</td> <td class="paramname"></td><td>)</td> <td> const</td> </tr> </table> </td> <td class="mlabels-right"> <span class="mlabels"><span class="mlabel">inline</span></span> </td> </tr> </table> </div><div class="memdoc"> <dl class="section return"><dt>Returns</dt><dd>a row (or column) vector expression of the squared norm of each column (or row) of the referenced expression.</dd></dl> <p>Example: </p> <div class="fragment"><div class="line"><a class="code" href="group__matrixtypedefs.html#ga45a14b423c687c3e2e8325f148e27bf3">Matrix3d</a> m = <a class="code" href="classEigen_1_1DenseBase.html#a8e759dafdd9ecc446d397b7f5435f60a">Matrix3d::Random</a>();</div> <div class="line">cout << <span class="stringliteral">"Here is the matrix m:"</span> << endl << m << endl;</div> <div class="line">cout << <span class="stringliteral">"Here is the square norm of each row:"</span> << endl << m.rowwise().squaredNorm() << endl;</div> </div><!-- fragment --><p> Output: </p> <pre class="fragment">Here is the matrix m: 0.68 0.597 -0.33 -0.211 0.823 0.536 0.566 -0.605 -0.444 Here is the square norm of each row: 0.928 1.01 0.884 </pre><dl class="section see"><dt>See Also</dt><dd>DenseBase::squaredNorm() </dd></dl> </div> </div> <a class="anchor" id="ac32ca2926dfc6a21859c0aa9e8eb1bcb"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">const ReturnType<internal::member_stableNorm,RealScalar>::Type stableNorm </td> <td>(</td> <td class="paramname"></td><td>)</td> <td> const</td> </tr> </table> </td> <td class="mlabels-right"> <span class="mlabels"><span class="mlabel">inline</span></span> </td> </tr> </table> </div><div class="memdoc"> <dl class="section return"><dt>Returns</dt><dd>a row (or column) vector expression of the norm of each column (or row) of the referenced expression, avoiding underflow and overflow.</dd></dl> <dl class="section see"><dt>See Also</dt><dd>DenseBase::stableNorm() </dd></dl> </div> </div> <a class="anchor" id="af27cb8c9043ed6408e64f94645439362"></a> <div class="memitem"> <div class="memproto"> <table class="mlabels"> <tr> <td class="mlabels-left"> <table class="memname"> <tr> <td class="memname">const ReturnType<internal::member_sum>::Type sum </td> <td>(</td> <td class="paramname"></td><td>)</td> <td> const</td> </tr> </table> </td> <td class="mlabels-right"> <span class="mlabels"><span class="mlabel">inline</span></span> </td> </tr> </table> </div><div class="memdoc"> <dl class="section return"><dt>Returns</dt><dd>a row (or column) vector expression of the sum of each column (or row) of the referenced expression.</dd></dl> <p>Example: </p> <div class="fragment"><div class="line"><a class="code" href="group__matrixtypedefs.html#ga45a14b423c687c3e2e8325f148e27bf3">Matrix3d</a> m = <a class="code" href="classEigen_1_1DenseBase.html#a8e759dafdd9ecc446d397b7f5435f60a">Matrix3d::Random</a>();</div> <div class="line">cout << <span class="stringliteral">"Here is the matrix m:"</span> << endl << m << endl;</div> <div class="line">cout << <span class="stringliteral">"Here is the sum of each row:"</span> << endl << m.rowwise().sum() << endl;</div> </div><!-- fragment --><p> Output: </p> <pre class="fragment">Here is the matrix m: 0.68 0.597 -0.33 -0.211 0.823 0.536 0.566 -0.605 -0.444 Here is the sum of each row: 0.948 1.15 -0.483 </pre><dl class="section see"><dt>See Also</dt><dd><a class="el" href="classEigen_1_1DenseBase.html#a3a3b3fb530d3364ecef0bf9c9daf0983">DenseBase::sum()</a> </dd></dl> </div> </div> <hr/>The documentation for this class was generated from the following files:<ul> <li><a class="el" href="ForwardDeclarations_8h_source.html">ForwardDeclarations.h</a></li> <li><a class="el" href="VectorwiseOp_8h_source.html">VectorwiseOp.h</a></li> <li><a class="el" href="Replicate_8h_source.html">Replicate.h</a></li> <li><a class="el" href="Homogeneous_8h_source.html">Homogeneous.h</a></li> <li><a class="el" href="OrthoMethods_8h_source.html">OrthoMethods.h</a></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="navelem"><a class="el" href="namespaceEigen.html">Eigen</a></li><li class="navelem"><a class="el" href="classEigen_1_1VectorwiseOp.html">VectorwiseOp</a></li> <li class="footer">Generated on Mon Oct 28 2013 11:04:31 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>