<!DOCTYPE html PUBLIC "-//W3C//DTD XHTML 1.0 Transitional//EN" "https://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.9.1"/>
<meta name="viewport" content="width=device-width, initial-scale=1"/>
<title>DGtal: DGtal::functors::IIPrincipalCurvatures3DFunctor&lt; TSpace, TMatrix &gt; Struct 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="navtreedata.js"></script>
<script type="text/javascript" src="navtree.js"></script>
<link href="search/search.css" rel="stylesheet" type="text/css"/>
<script type="text/javascript" src="search/searchdata.js"></script>
<script type="text/javascript" src="search/search.js"></script>
<script type="text/x-mathjax-config">
  MathJax.Hub.Config({
    extensions: ["tex2jax.js", "TeX/AMSmath.js", "TeX/AMSsymbols.js"],
    jax: ["input/TeX","output/HTML-CSS"],
});
</script>
<script type="text/javascript" async="async" src="https://cdnjs.cloudflare.com/ajax/libs/mathjax/2.7.3/MathJax.js?config=TeX-MML-AM_CHTML/MathJax.js"></script>
<link href="doxygen.css" rel="stylesheet" type="text/css" />
<link href="doxygen-awesome.css" rel="stylesheet" type="text/css"/>
</head>
<body>
<div id="top"><!-- do not remove this div, it is closed by doxygen! -->
<div id="titlearea">
<table cellspacing="0" cellpadding="0">
 <tbody>
 <tr style="height: 56px;">
  <td id="projectalign" style="padding-left: 0.5em;">
   <div id="projectname">DGtal
   &#160;<span id="projectnumber">1.4.2</span>
   </div>
  </td>
 </tr>
 </tbody>
</table>
</div>
<!-- end header part -->
<!-- Generated by Doxygen 1.9.1 -->
<script type="text/javascript">
/* @license magnet:?xt=urn:btih:cf05388f2679ee054f2beb29a391d25f4e673ac3&amp;dn=gpl-2.0.txt GPL-v2 */
var searchBox = new SearchBox("searchBox", "search",false,'Search','.html');
/* @license-end */
</script>
<script type="text/javascript" src="menudata.js"></script>
<script type="text/javascript" src="menu.js"></script>
<script type="text/javascript">
/* @license magnet:?xt=urn:btih:cf05388f2679ee054f2beb29a391d25f4e673ac3&amp;dn=gpl-2.0.txt GPL-v2 */
$(function() {
  initMenu('',true,false,'search.php','Search');
  $(document).ready(function() { init_search(); });
});
/* @license-end */</script>
<div id="main-nav"></div>
</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">
/* @license magnet:?xt=urn:btih:cf05388f2679ee054f2beb29a391d25f4e673ac3&amp;dn=gpl-2.0.txt GPL-v2 */
$(document).ready(function(){initNavTree('structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html',''); initResizable(); });
/* @license-end */
</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)">
</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="#pub-types">Public Types</a> &#124;
<a href="#pub-methods">Public Member Functions</a> &#124;
<a href="#pri-attribs">Private Attributes</a>  </div>
  <div class="headertitle">
<div class="title">DGtal::functors::IIPrincipalCurvatures3DFunctor&lt; TSpace, TMatrix &gt; Struct Template Reference</div>  </div>
</div><!--header-->
<div class="contents">

<p>Aim: A functor Matrix -&gt; std::pair&lt;Real,Real&gt; that returns the first and the second principal curvature value by diagonalizing the given covariance matrix. This functor is valid starting from 3D space. Note that by first we mean the value with first greatest curvature in absolute value.  
 <a href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#details">More...</a></p>

<p><code>#include &lt;<a class="el" href="IIGeometricFunctors_8h_source.html">DGtal/geometry/surfaces/estimation/IIGeometricFunctors.h</a>&gt;</code></p>
<table class="memberdecls">
<tr class="heading"><td colspan="2"><h2 class="groupheader"><a name="pub-types"></a>
Public Types</h2></td></tr>
<tr class="memitem:a8179253dfd6e269395da349910ba9e9e"><td class="memItemLeft" align="right" valign="top">typedef <a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html">IIPrincipalCurvatures3DFunctor</a>&lt; TSpace &gt;&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a8179253dfd6e269395da349910ba9e9e">Self</a></td></tr>
<tr class="separator:a8179253dfd6e269395da349910ba9e9e"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a929f8df064f1e3b66db6dd9adaa97b38"><td class="memItemLeft" align="right" valign="top">typedef TSpace&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a929f8df064f1e3b66db6dd9adaa97b38">Space</a></td></tr>
<tr class="separator:a929f8df064f1e3b66db6dd9adaa97b38"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a53becf0c31439d565f071bc6f88ff8a4"><td class="memItemLeft" align="right" valign="top">typedef <a class="el" href="classDGtal_1_1SpaceND.html#ac103823b62d88adef786537c9c2a77fb">Space::RealVector</a>&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a53becf0c31439d565f071bc6f88ff8a4">RealVector</a></td></tr>
<tr class="separator:a53becf0c31439d565f071bc6f88ff8a4"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a64da36486bc805c95fd34ffd2d0988dd"><td class="memItemLeft" align="right" valign="top">typedef <a class="el" href="classDGtal_1_1PointVector.html#ab5f7548d4edc0d8edb4b64cc332db8f2">RealVector::Component</a>&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a64da36486bc805c95fd34ffd2d0988dd">Component</a></td></tr>
<tr class="separator:a64da36486bc805c95fd34ffd2d0988dd"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:af678b7274eb6d3609df360bca9bd1c0d"><td class="memItemLeft" align="right" valign="top">typedef TMatrix&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#af678b7274eb6d3609df360bca9bd1c0d">Matrix</a></td></tr>
<tr class="separator:af678b7274eb6d3609df360bca9bd1c0d"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:ad3e58211c10a4891e2fd2d0ff473c96a"><td class="memItemLeft" align="right" valign="top">typedef <a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#af678b7274eb6d3609df360bca9bd1c0d">Matrix</a>&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#ad3e58211c10a4891e2fd2d0ff473c96a">Argument</a></td></tr>
<tr class="separator:ad3e58211c10a4891e2fd2d0ff473c96a"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a03e5d38b68f9e82cea9305ab1202073e"><td class="memItemLeft" align="right" valign="top">typedef std::pair&lt; <a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a64da36486bc805c95fd34ffd2d0988dd">Component</a>, <a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a64da36486bc805c95fd34ffd2d0988dd">Component</a> &gt;&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a03e5d38b68f9e82cea9305ab1202073e">Quantity</a></td></tr>
<tr class="separator:a03e5d38b68f9e82cea9305ab1202073e"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:ae62a9a04b43e44c90857e587e33141e8"><td class="memItemLeft" align="right" valign="top">typedef <a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a03e5d38b68f9e82cea9305ab1202073e">Quantity</a>&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#ae62a9a04b43e44c90857e587e33141e8">Value</a></td></tr>
<tr class="separator:ae62a9a04b43e44c90857e587e33141e8"><td class="memSeparator" colspan="2">&#160;</td></tr>
</table><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:a016cf3b748a1a8036d7e0a0eac70bdf0"><td class="memItemLeft" align="right" valign="top">&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a016cf3b748a1a8036d7e0a0eac70bdf0">BOOST_CONCEPT_ASSERT</a> ((<a class="el" href="structDGtal_1_1concepts_1_1CMatrix.html">concepts::CMatrix</a>&lt; <a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#af678b7274eb6d3609df360bca9bd1c0d">Matrix</a> &gt;))</td></tr>
<tr class="separator:a016cf3b748a1a8036d7e0a0eac70bdf0"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:af3ace363ac47ae6435cd14355937c880"><td class="memItemLeft" align="right" valign="top">&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#af3ace363ac47ae6435cd14355937c880">BOOST_CONCEPT_ASSERT</a> ((<a class="el" href="structDGtal_1_1concepts_1_1CSpace.html">concepts::CSpace</a>&lt; TSpace &gt;))</td></tr>
<tr class="separator:af3ace363ac47ae6435cd14355937c880"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a3f904ea85339434c375bba97ca8ab580"><td class="memItemLeft" align="right" valign="top">&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a3f904ea85339434c375bba97ca8ab580">BOOST_STATIC_ASSERT</a> ((Space::dimension==3))</td></tr>
<tr class="separator:a3f904ea85339434c375bba97ca8ab580"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:ac4fee375eecdfaa026008f7000bb56af"><td class="memItemLeft" align="right" valign="top"><a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#ae62a9a04b43e44c90857e587e33141e8">Value</a>&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#ac4fee375eecdfaa026008f7000bb56af">operator()</a> (const <a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#ad3e58211c10a4891e2fd2d0ff473c96a">Argument</a> &amp;arg) const</td></tr>
<tr class="separator:ac4fee375eecdfaa026008f7000bb56af"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a10e807e1271fa86ef64d5563a90a2c60"><td class="memItemLeft" align="right" valign="top">void&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a10e807e1271fa86ef64d5563a90a2c60">init</a> (<a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a64da36486bc805c95fd34ffd2d0988dd">Component</a> h, <a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a64da36486bc805c95fd34ffd2d0988dd">Component</a> r)</td></tr>
<tr class="separator:a10e807e1271fa86ef64d5563a90a2c60"><td class="memSeparator" colspan="2">&#160;</td></tr>
</table><table class="memberdecls">
<tr class="heading"><td colspan="2"><h2 class="groupheader"><a name="pri-attribs"></a>
Private Attributes</h2></td></tr>
<tr class="memitem:a220caf9d2e46f781c24ab05f8312f462"><td class="memItemLeft" align="right" valign="top">double&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a220caf9d2e46f781c24ab05f8312f462">dh5</a></td></tr>
<tr class="separator:a220caf9d2e46f781c24ab05f8312f462"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:ada43c69f9fa35d7d30ef858c83abbfce"><td class="memItemLeft" align="right" valign="top">double&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#ada43c69f9fa35d7d30ef858c83abbfce">d6_PIr6</a></td></tr>
<tr class="separator:ada43c69f9fa35d7d30ef858c83abbfce"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a0ee4ebb4e1ec5c07648a009c242a7171"><td class="memItemLeft" align="right" valign="top">double&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a0ee4ebb4e1ec5c07648a009c242a7171">d8_5r</a></td></tr>
<tr class="separator:a0ee4ebb4e1ec5c07648a009c242a7171"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:a2308a8adc544387faa793be87290ee39"><td class="memItemLeft" align="right" valign="top"><a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#af678b7274eb6d3609df360bca9bd1c0d">Matrix</a>&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a2308a8adc544387faa793be87290ee39">eigenVectors</a></td></tr>
<tr class="memdesc:a2308a8adc544387faa793be87290ee39"><td class="mdescLeft">&#160;</td><td class="mdescRight">A data member only used for temporary calculations.  <a href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a2308a8adc544387faa793be87290ee39">More...</a><br /></td></tr>
<tr class="separator:a2308a8adc544387faa793be87290ee39"><td class="memSeparator" colspan="2">&#160;</td></tr>
<tr class="memitem:ae0fe492a01b65ccd59b90bb9d5aad1f5"><td class="memItemLeft" align="right" valign="top"><a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a53becf0c31439d565f071bc6f88ff8a4">RealVector</a>&#160;</td><td class="memItemRight" valign="bottom"><a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#ae0fe492a01b65ccd59b90bb9d5aad1f5">eigenValues</a></td></tr>
<tr class="memdesc:ae0fe492a01b65ccd59b90bb9d5aad1f5"><td class="mdescLeft">&#160;</td><td class="mdescRight">A data member only used for temporary calculations.  <a href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#ae0fe492a01b65ccd59b90bb9d5aad1f5">More...</a><br /></td></tr>
<tr class="separator:ae0fe492a01b65ccd59b90bb9d5aad1f5"><td class="memSeparator" colspan="2">&#160;</td></tr>
</table>
<a name="details" id="details"></a><h2 class="groupheader">Detailed Description</h2>
<div class="textblock"><h3>template&lt;typename TSpace, typename TMatrix = SimpleMatrix&lt; typename TSpace::RealVector::Component, TSpace::dimension, TSpace::dimension&gt;&gt;<br />
struct DGtal::functors::IIPrincipalCurvatures3DFunctor&lt; TSpace, TMatrix &gt;</h3>

<p>Aim: A functor Matrix -&gt; std::pair&lt;Real,Real&gt; that returns the first and the second principal curvature value by diagonalizing the given covariance matrix. This functor is valid starting from 3D space. Note that by first we mean the value with first greatest curvature in absolute value. </p>
<p>Description of template class '<a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html" title="Aim: A functor Matrix -&gt; std::pair&lt;Real,Real&gt; that returns the first and the second principal curvatu...">IIPrincipalCurvatures3DFunctor</a>' </p>
<dl class="tparams"><dt>Template Parameters</dt><dd>
  <table class="tparams">
    <tr><td class="paramname">TSpace</td><td>a model of CSpace, for instance <a class="el" href="classDGtal_1_1SpaceND.html">SpaceND</a>. </td></tr>
    <tr><td class="paramname">TMatrix</td><td>a model of CMatrix, for instance <a class="el" href="classDGtal_1_1SimpleMatrix.html" title="Aim: implements basic MxN Matrix services (M,N&gt;=1).">SimpleMatrix</a>.</td></tr>
  </table>
  </dd>
</dl>
<dl class="section see"><dt>See also</dt><dd><a class="el" href="classDGtal_1_1IntegralInvariantCovarianceEstimator.html" title="Aim: This class implement an Integral Invariant estimator which computes for each surfel the covarian...">IntegralInvariantCovarianceEstimator</a> </dd></dl>

<p class="definition">Definition at line <a class="el" href="IIGeometricFunctors_8h_source.html#l00912">912</a> of file <a class="el" href="IIGeometricFunctors_8h_source.html">IIGeometricFunctors.h</a>.</p>
</div><h2 class="groupheader">Member Typedef Documentation</h2>
<a id="ad3e58211c10a4891e2fd2d0ff473c96a"></a>
<h2 class="memtitle"><span class="permalink"><a href="#ad3e58211c10a4891e2fd2d0ff473c96a">&#9670;&nbsp;</a></span>Argument</h2>

<div class="memitem">
<div class="memproto">
<div class="memtemplate">
template&lt;typename TSpace , typename TMatrix  = SimpleMatrix&lt; typename TSpace::RealVector::Component, TSpace::dimension, TSpace::dimension&gt;&gt; </div>
      <table class="memname">
        <tr>
          <td class="memname">typedef <a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#af678b7274eb6d3609df360bca9bd1c0d">Matrix</a> <a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html">DGtal::functors::IIPrincipalCurvatures3DFunctor</a>&lt; TSpace, TMatrix &gt;::<a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#ad3e58211c10a4891e2fd2d0ff473c96a">Argument</a></td>
        </tr>
      </table>
</div><div class="memdoc">

<p class="definition">Definition at line <a class="el" href="IIGeometricFunctors_8h_source.html#l00921">921</a> of file <a class="el" href="IIGeometricFunctors_8h_source.html">IIGeometricFunctors.h</a>.</p>

</div>
</div>
<a id="a64da36486bc805c95fd34ffd2d0988dd"></a>
<h2 class="memtitle"><span class="permalink"><a href="#a64da36486bc805c95fd34ffd2d0988dd">&#9670;&nbsp;</a></span>Component</h2>

<div class="memitem">
<div class="memproto">
<div class="memtemplate">
template&lt;typename TSpace , typename TMatrix  = SimpleMatrix&lt; typename TSpace::RealVector::Component, TSpace::dimension, TSpace::dimension&gt;&gt; </div>
      <table class="memname">
        <tr>
          <td class="memname">typedef <a class="el" href="classDGtal_1_1PointVector.html#ab5f7548d4edc0d8edb4b64cc332db8f2">RealVector::Component</a> <a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html">DGtal::functors::IIPrincipalCurvatures3DFunctor</a>&lt; TSpace, TMatrix &gt;::<a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a64da36486bc805c95fd34ffd2d0988dd">Component</a></td>
        </tr>
      </table>
</div><div class="memdoc">

<p class="definition">Definition at line <a class="el" href="IIGeometricFunctors_8h_source.html#l00919">919</a> of file <a class="el" href="IIGeometricFunctors_8h_source.html">IIGeometricFunctors.h</a>.</p>

</div>
</div>
<a id="af678b7274eb6d3609df360bca9bd1c0d"></a>
<h2 class="memtitle"><span class="permalink"><a href="#af678b7274eb6d3609df360bca9bd1c0d">&#9670;&nbsp;</a></span>Matrix</h2>

<div class="memitem">
<div class="memproto">
<div class="memtemplate">
template&lt;typename TSpace , typename TMatrix  = SimpleMatrix&lt; typename TSpace::RealVector::Component, TSpace::dimension, TSpace::dimension&gt;&gt; </div>
      <table class="memname">
        <tr>
          <td class="memname">typedef TMatrix <a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html">DGtal::functors::IIPrincipalCurvatures3DFunctor</a>&lt; TSpace, TMatrix &gt;::<a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#af678b7274eb6d3609df360bca9bd1c0d">Matrix</a></td>
        </tr>
      </table>
</div><div class="memdoc">

<p class="definition">Definition at line <a class="el" href="IIGeometricFunctors_8h_source.html#l00920">920</a> of file <a class="el" href="IIGeometricFunctors_8h_source.html">IIGeometricFunctors.h</a>.</p>

</div>
</div>
<a id="a03e5d38b68f9e82cea9305ab1202073e"></a>
<h2 class="memtitle"><span class="permalink"><a href="#a03e5d38b68f9e82cea9305ab1202073e">&#9670;&nbsp;</a></span>Quantity</h2>

<div class="memitem">
<div class="memproto">
<div class="memtemplate">
template&lt;typename TSpace , typename TMatrix  = SimpleMatrix&lt; typename TSpace::RealVector::Component, TSpace::dimension, TSpace::dimension&gt;&gt; </div>
      <table class="memname">
        <tr>
          <td class="memname">typedef std::pair&lt;<a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a64da36486bc805c95fd34ffd2d0988dd">Component</a>, <a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a64da36486bc805c95fd34ffd2d0988dd">Component</a>&gt; <a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html">DGtal::functors::IIPrincipalCurvatures3DFunctor</a>&lt; TSpace, TMatrix &gt;::<a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a03e5d38b68f9e82cea9305ab1202073e">Quantity</a></td>
        </tr>
      </table>
</div><div class="memdoc">

<p class="definition">Definition at line <a class="el" href="IIGeometricFunctors_8h_source.html#l00922">922</a> of file <a class="el" href="IIGeometricFunctors_8h_source.html">IIGeometricFunctors.h</a>.</p>

</div>
</div>
<a id="a53becf0c31439d565f071bc6f88ff8a4"></a>
<h2 class="memtitle"><span class="permalink"><a href="#a53becf0c31439d565f071bc6f88ff8a4">&#9670;&nbsp;</a></span>RealVector</h2>

<div class="memitem">
<div class="memproto">
<div class="memtemplate">
template&lt;typename TSpace , typename TMatrix  = SimpleMatrix&lt; typename TSpace::RealVector::Component, TSpace::dimension, TSpace::dimension&gt;&gt; </div>
      <table class="memname">
        <tr>
          <td class="memname">typedef <a class="el" href="classDGtal_1_1SpaceND.html#ac103823b62d88adef786537c9c2a77fb">Space::RealVector</a> <a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html">DGtal::functors::IIPrincipalCurvatures3DFunctor</a>&lt; TSpace, TMatrix &gt;::<a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a53becf0c31439d565f071bc6f88ff8a4">RealVector</a></td>
        </tr>
      </table>
</div><div class="memdoc">

<p class="definition">Definition at line <a class="el" href="IIGeometricFunctors_8h_source.html#l00918">918</a> of file <a class="el" href="IIGeometricFunctors_8h_source.html">IIGeometricFunctors.h</a>.</p>

</div>
</div>
<a id="a8179253dfd6e269395da349910ba9e9e"></a>
<h2 class="memtitle"><span class="permalink"><a href="#a8179253dfd6e269395da349910ba9e9e">&#9670;&nbsp;</a></span>Self</h2>

<div class="memitem">
<div class="memproto">
<div class="memtemplate">
template&lt;typename TSpace , typename TMatrix  = SimpleMatrix&lt; typename TSpace::RealVector::Component, TSpace::dimension, TSpace::dimension&gt;&gt; </div>
      <table class="memname">
        <tr>
          <td class="memname">typedef <a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html">IIPrincipalCurvatures3DFunctor</a>&lt;TSpace&gt; <a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html">DGtal::functors::IIPrincipalCurvatures3DFunctor</a>&lt; TSpace, TMatrix &gt;::<a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a8179253dfd6e269395da349910ba9e9e">Self</a></td>
        </tr>
      </table>
</div><div class="memdoc">

<p class="definition">Definition at line <a class="el" href="IIGeometricFunctors_8h_source.html#l00916">916</a> of file <a class="el" href="IIGeometricFunctors_8h_source.html">IIGeometricFunctors.h</a>.</p>

</div>
</div>
<a id="a929f8df064f1e3b66db6dd9adaa97b38"></a>
<h2 class="memtitle"><span class="permalink"><a href="#a929f8df064f1e3b66db6dd9adaa97b38">&#9670;&nbsp;</a></span>Space</h2>

<div class="memitem">
<div class="memproto">
<div class="memtemplate">
template&lt;typename TSpace , typename TMatrix  = SimpleMatrix&lt; typename TSpace::RealVector::Component, TSpace::dimension, TSpace::dimension&gt;&gt; </div>
      <table class="memname">
        <tr>
          <td class="memname">typedef TSpace <a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html">DGtal::functors::IIPrincipalCurvatures3DFunctor</a>&lt; TSpace, TMatrix &gt;::<a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a929f8df064f1e3b66db6dd9adaa97b38">Space</a></td>
        </tr>
      </table>
</div><div class="memdoc">

<p class="definition">Definition at line <a class="el" href="IIGeometricFunctors_8h_source.html#l00917">917</a> of file <a class="el" href="IIGeometricFunctors_8h_source.html">IIGeometricFunctors.h</a>.</p>

</div>
</div>
<a id="ae62a9a04b43e44c90857e587e33141e8"></a>
<h2 class="memtitle"><span class="permalink"><a href="#ae62a9a04b43e44c90857e587e33141e8">&#9670;&nbsp;</a></span>Value</h2>

<div class="memitem">
<div class="memproto">
<div class="memtemplate">
template&lt;typename TSpace , typename TMatrix  = SimpleMatrix&lt; typename TSpace::RealVector::Component, TSpace::dimension, TSpace::dimension&gt;&gt; </div>
      <table class="memname">
        <tr>
          <td class="memname">typedef <a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a03e5d38b68f9e82cea9305ab1202073e">Quantity</a> <a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html">DGtal::functors::IIPrincipalCurvatures3DFunctor</a>&lt; TSpace, TMatrix &gt;::<a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#ae62a9a04b43e44c90857e587e33141e8">Value</a></td>
        </tr>
      </table>
</div><div class="memdoc">

<p class="definition">Definition at line <a class="el" href="IIGeometricFunctors_8h_source.html#l00923">923</a> of file <a class="el" href="IIGeometricFunctors_8h_source.html">IIGeometricFunctors.h</a>.</p>

</div>
</div>
<h2 class="groupheader">Member Function Documentation</h2>
<a id="a016cf3b748a1a8036d7e0a0eac70bdf0"></a>
<h2 class="memtitle"><span class="permalink"><a href="#a016cf3b748a1a8036d7e0a0eac70bdf0">&#9670;&nbsp;</a></span>BOOST_CONCEPT_ASSERT() <span class="overload">[1/2]</span></h2>

<div class="memitem">
<div class="memproto">
<div class="memtemplate">
template&lt;typename TSpace , typename TMatrix  = SimpleMatrix&lt; typename TSpace::RealVector::Component, TSpace::dimension, TSpace::dimension&gt;&gt; </div>
      <table class="memname">
        <tr>
          <td class="memname"><a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html">DGtal::functors::IIPrincipalCurvatures3DFunctor</a>&lt; TSpace, TMatrix &gt;::BOOST_CONCEPT_ASSERT </td>
          <td>(</td>
          <td class="paramtype">(<a class="el" href="structDGtal_1_1concepts_1_1CMatrix.html">concepts::CMatrix</a>&lt; <a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#af678b7274eb6d3609df360bca9bd1c0d">Matrix</a> &gt;)&#160;</td>
          <td class="paramname"></td><td>)</td>
          <td></td>
        </tr>
      </table>
</div><div class="memdoc">

</div>
</div>
<a id="af3ace363ac47ae6435cd14355937c880"></a>
<h2 class="memtitle"><span class="permalink"><a href="#af3ace363ac47ae6435cd14355937c880">&#9670;&nbsp;</a></span>BOOST_CONCEPT_ASSERT() <span class="overload">[2/2]</span></h2>

<div class="memitem">
<div class="memproto">
<div class="memtemplate">
template&lt;typename TSpace , typename TMatrix  = SimpleMatrix&lt; typename TSpace::RealVector::Component, TSpace::dimension, TSpace::dimension&gt;&gt; </div>
      <table class="memname">
        <tr>
          <td class="memname"><a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html">DGtal::functors::IIPrincipalCurvatures3DFunctor</a>&lt; TSpace, TMatrix &gt;::BOOST_CONCEPT_ASSERT </td>
          <td>(</td>
          <td class="paramtype">(<a class="el" href="structDGtal_1_1concepts_1_1CSpace.html">concepts::CSpace</a>&lt; TSpace &gt;)&#160;</td>
          <td class="paramname"></td><td>)</td>
          <td></td>
        </tr>
      </table>
</div><div class="memdoc">

</div>
</div>
<a id="a3f904ea85339434c375bba97ca8ab580"></a>
<h2 class="memtitle"><span class="permalink"><a href="#a3f904ea85339434c375bba97ca8ab580">&#9670;&nbsp;</a></span>BOOST_STATIC_ASSERT()</h2>

<div class="memitem">
<div class="memproto">
<div class="memtemplate">
template&lt;typename TSpace , typename TMatrix  = SimpleMatrix&lt; typename TSpace::RealVector::Component, TSpace::dimension, TSpace::dimension&gt;&gt; </div>
      <table class="memname">
        <tr>
          <td class="memname"><a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html">DGtal::functors::IIPrincipalCurvatures3DFunctor</a>&lt; TSpace, TMatrix &gt;::BOOST_STATIC_ASSERT </td>
          <td>(</td>
          <td class="paramtype">(Space::dimension==3)&#160;</td>
          <td class="paramname"></td><td>)</td>
          <td></td>
        </tr>
      </table>
</div><div class="memdoc">

</div>
</div>
<a id="a10e807e1271fa86ef64d5563a90a2c60"></a>
<h2 class="memtitle"><span class="permalink"><a href="#a10e807e1271fa86ef64d5563a90a2c60">&#9670;&nbsp;</a></span>init()</h2>

<div class="memitem">
<div class="memproto">
<div class="memtemplate">
template&lt;typename TSpace , typename TMatrix  = SimpleMatrix&lt; typename TSpace::RealVector::Component, TSpace::dimension, TSpace::dimension&gt;&gt; </div>
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">void <a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html">DGtal::functors::IIPrincipalCurvatures3DFunctor</a>&lt; TSpace, TMatrix &gt;::init </td>
          <td>(</td>
          <td class="paramtype"><a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a64da36486bc805c95fd34ffd2d0988dd">Component</a>&#160;</td>
          <td class="paramname"><em>h</em>, </td>
        </tr>
        <tr>
          <td class="paramkey"></td>
          <td></td>
          <td class="paramtype"><a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a64da36486bc805c95fd34ffd2d0988dd">Component</a>&#160;</td>
          <td class="paramname"><em>r</em>&#160;</td>
        </tr>
        <tr>
          <td></td>
          <td>)</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>Initializes the functor with the gridstep and the ball Euclidean radius.</p>
<dl class="params"><dt>Parameters</dt><dd>
  <table class="params">
    <tr><td class="paramname">h</td><td>the gridstep </td></tr>
    <tr><td class="paramname">r</td><td>the ball radius </td></tr>
  </table>
  </dd>
</dl>

<p class="definition">Definition at line <a class="el" href="IIGeometricFunctors_8h_source.html#l00960">960</a> of file <a class="el" href="IIGeometricFunctors_8h_source.html">IIGeometricFunctors.h</a>.</p>
<div class="fragment"><div class="line"><a name="l00961"></a><span class="lineno">  961</span>&#160;      {</div>
<div class="line"><a name="l00962"></a><span class="lineno">  962</span>&#160;        <span class="keywordtype">double</span> r3 = r * r * r;</div>
<div class="line"><a name="l00963"></a><span class="lineno">  963</span>&#160;        <span class="keywordtype">double</span> r6 = r3 * r3;</div>
<div class="line"><a name="l00964"></a><span class="lineno">  964</span>&#160;        <a class="code" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#ada43c69f9fa35d7d30ef858c83abbfce">d6_PIr6</a> = 6.0 / ( M_PI * r6 );</div>
<div class="line"><a name="l00965"></a><span class="lineno">  965</span>&#160;        <a class="code" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a0ee4ebb4e1ec5c07648a009c242a7171">d8_5r</a> = 8.0 / ( 5.0 * r );</div>
<div class="line"><a name="l00966"></a><span class="lineno">  966</span>&#160;        <span class="keywordtype">double</span> h2 = h * h;</div>
<div class="line"><a name="l00967"></a><span class="lineno">  967</span>&#160;        <a class="code" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a220caf9d2e46f781c24ab05f8312f462">dh5</a> = h2 * h2 * h;</div>
<div class="line"><a name="l00968"></a><span class="lineno">  968</span>&#160;      }</div>
<div class="ttc" id="astructDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor_html_a0ee4ebb4e1ec5c07648a009c242a7171"><div class="ttname"><a href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a0ee4ebb4e1ec5c07648a009c242a7171">DGtal::functors::IIPrincipalCurvatures3DFunctor::d8_5r</a></div><div class="ttdeci">double d8_5r</div><div class="ttdef"><b>Definition:</b> <a href="IIGeometricFunctors_8h_source.html#l00973">IIGeometricFunctors.h:973</a></div></div>
<div class="ttc" id="astructDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor_html_a220caf9d2e46f781c24ab05f8312f462"><div class="ttname"><a href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a220caf9d2e46f781c24ab05f8312f462">DGtal::functors::IIPrincipalCurvatures3DFunctor::dh5</a></div><div class="ttdeci">double dh5</div><div class="ttdef"><b>Definition:</b> <a href="IIGeometricFunctors_8h_source.html#l00971">IIGeometricFunctors.h:971</a></div></div>
<div class="ttc" id="astructDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor_html_ada43c69f9fa35d7d30ef858c83abbfce"><div class="ttname"><a href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#ada43c69f9fa35d7d30ef858c83abbfce">DGtal::functors::IIPrincipalCurvatures3DFunctor::d6_PIr6</a></div><div class="ttdeci">double d6_PIr6</div><div class="ttdef"><b>Definition:</b> <a href="IIGeometricFunctors_8h_source.html#l00972">IIGeometricFunctors.h:972</a></div></div>
</div><!-- fragment -->
<p class="reference">References <a class="el" href="IIGeometricFunctors_8h_source.html#l00972">DGtal::functors::IIPrincipalCurvatures3DFunctor&lt; TSpace, TMatrix &gt;::d6_PIr6</a>, <a class="el" href="IIGeometricFunctors_8h_source.html#l00973">DGtal::functors::IIPrincipalCurvatures3DFunctor&lt; TSpace, TMatrix &gt;::d8_5r</a>, and <a class="el" href="IIGeometricFunctors_8h_source.html#l00971">DGtal::functors::IIPrincipalCurvatures3DFunctor&lt; TSpace, TMatrix &gt;::dh5</a>.</p>

</div>
</div>
<a id="ac4fee375eecdfaa026008f7000bb56af"></a>
<h2 class="memtitle"><span class="permalink"><a href="#ac4fee375eecdfaa026008f7000bb56af">&#9670;&nbsp;</a></span>operator()()</h2>

<div class="memitem">
<div class="memproto">
<div class="memtemplate">
template&lt;typename TSpace , typename TMatrix  = SimpleMatrix&lt; typename TSpace::RealVector::Component, TSpace::dimension, TSpace::dimension&gt;&gt; </div>
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname"><a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#ae62a9a04b43e44c90857e587e33141e8">Value</a> <a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html">DGtal::functors::IIPrincipalCurvatures3DFunctor</a>&lt; TSpace, TMatrix &gt;::operator() </td>
          <td>(</td>
          <td class="paramtype">const <a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#ad3e58211c10a4891e2fd2d0ff473c96a">Argument</a> &amp;&#160;</td>
          <td class="paramname"><em>arg</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>Apply operator. </p><dl class="params"><dt>Parameters</dt><dd>
  <table class="params">
    <tr><td class="paramname">arg</td><td>any symmetric positive matrix (covariance matrix</td></tr>
  </table>
  </dd>
</dl>
<dl class="section return"><dt>Returns</dt><dd>the first and the second principal curvature value in a std::pair for the II covariance matrix, which are the first and the second highest eigenvalue. </dd></dl>

<p class="definition">Definition at line <a class="el" href="IIGeometricFunctors_8h_source.html#l00937">937</a> of file <a class="el" href="IIGeometricFunctors_8h_source.html">IIGeometricFunctors.h</a>.</p>
<div class="fragment"><div class="line"><a name="l00938"></a><span class="lineno">  938</span>&#160;      {</div>
<div class="line"><a name="l00939"></a><span class="lineno">  939</span>&#160;        <a class="code" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#ad3e58211c10a4891e2fd2d0ff473c96a">Argument</a> cp_arg = arg;</div>
<div class="line"><a name="l00940"></a><span class="lineno">  940</span>&#160;        cp_arg *= <a class="code" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a220caf9d2e46f781c24ab05f8312f462">dh5</a>;</div>
<div class="line"><a name="l00941"></a><span class="lineno">  941</span>&#160;        <a class="code" href="classDGtal_1_1EigenDecomposition.html#a062ae590027b118c107dac78f0e9773c">EigenDecomposition&lt;Space::dimension, Component, Matrix&gt;</a></div>
<div class="line"><a name="l00942"></a><span class="lineno">  942</span>&#160;<a class="code" href="classDGtal_1_1EigenDecomposition.html#a062ae590027b118c107dac78f0e9773c">          ::getEigenDecomposition</a>( cp_arg, <a class="code" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a2308a8adc544387faa793be87290ee39">eigenVectors</a>, <a class="code" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#ae0fe492a01b65ccd59b90bb9d5aad1f5">eigenValues</a> );</div>
<div class="line"><a name="l00943"></a><span class="lineno">  943</span>&#160; </div>
<div class="line"><a name="l00944"></a><span class="lineno">  944</span>&#160;        ASSERT ( !std::isnan(<a class="code" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#ae0fe492a01b65ccd59b90bb9d5aad1f5">eigenValues</a>[0]) ); <span class="comment">// NaN</span></div>
<div class="line"><a name="l00945"></a><span class="lineno">  945</span>&#160;        ASSERT ( (std::abs(<a class="code" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#ae0fe492a01b65ccd59b90bb9d5aad1f5">eigenValues</a>[0]) &lt;= std::abs(<a class="code" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#ae0fe492a01b65ccd59b90bb9d5aad1f5">eigenValues</a>[1]))</div>
<div class="line"><a name="l00946"></a><span class="lineno">  946</span>&#160;              &amp;&amp; (std::abs(<a class="code" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#ae0fe492a01b65ccd59b90bb9d5aad1f5">eigenValues</a>[1]) &lt;= std::abs(<a class="code" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#ae0fe492a01b65ccd59b90bb9d5aad1f5">eigenValues</a>[2])) );</div>
<div class="line"><a name="l00947"></a><span class="lineno">  947</span>&#160; </div>
<div class="line"><a name="l00948"></a><span class="lineno">  948</span>&#160;        <span class="keywordflow">return</span> <a class="code" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#ae62a9a04b43e44c90857e587e33141e8">Value</a>(</div>
<div class="line"><a name="l00949"></a><span class="lineno">  949</span>&#160;                <a class="code" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#ada43c69f9fa35d7d30ef858c83abbfce">d6_PIr6</a> * ( <a class="code" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#ae0fe492a01b65ccd59b90bb9d5aad1f5">eigenValues</a>[2] - ( 3.0 * <a class="code" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#ae0fe492a01b65ccd59b90bb9d5aad1f5">eigenValues</a>[1] )) + <a class="code" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a0ee4ebb4e1ec5c07648a009c242a7171">d8_5r</a>,</div>
<div class="line"><a name="l00950"></a><span class="lineno">  950</span>&#160;                <a class="code" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#ada43c69f9fa35d7d30ef858c83abbfce">d6_PIr6</a> * ( <a class="code" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#ae0fe492a01b65ccd59b90bb9d5aad1f5">eigenValues</a>[1] - ( 3.0 * <a class="code" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#ae0fe492a01b65ccd59b90bb9d5aad1f5">eigenValues</a>[2] )) + <a class="code" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a0ee4ebb4e1ec5c07648a009c242a7171">d8_5r</a></div>
<div class="line"><a name="l00951"></a><span class="lineno">  951</span>&#160;               );</div>
<div class="line"><a name="l00952"></a><span class="lineno">  952</span>&#160;      }</div>
<div class="ttc" id="aclassDGtal_1_1EigenDecomposition_html_a062ae590027b118c107dac78f0e9773c"><div class="ttname"><a href="classDGtal_1_1EigenDecomposition.html#a062ae590027b118c107dac78f0e9773c">DGtal::EigenDecomposition::getEigenDecomposition</a></div><div class="ttdeci">static void getEigenDecomposition(const Matrix &amp;matrix, Matrix &amp;eigenVectors, Vector &amp;eigenValues)</div><div class="ttdoc">Compute both eigen vectors and eigen values from an input matrix.</div></div>
<div class="ttc" id="astructDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor_html_a2308a8adc544387faa793be87290ee39"><div class="ttname"><a href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a2308a8adc544387faa793be87290ee39">DGtal::functors::IIPrincipalCurvatures3DFunctor::eigenVectors</a></div><div class="ttdeci">Matrix eigenVectors</div><div class="ttdoc">A data member only used for temporary calculations.</div><div class="ttdef"><b>Definition:</b> <a href="IIGeometricFunctors_8h_source.html#l00976">IIGeometricFunctors.h:976</a></div></div>
<div class="ttc" id="astructDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor_html_ad3e58211c10a4891e2fd2d0ff473c96a"><div class="ttname"><a href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#ad3e58211c10a4891e2fd2d0ff473c96a">DGtal::functors::IIPrincipalCurvatures3DFunctor::Argument</a></div><div class="ttdeci">Matrix Argument</div><div class="ttdef"><b>Definition:</b> <a href="IIGeometricFunctors_8h_source.html#l00921">IIGeometricFunctors.h:921</a></div></div>
<div class="ttc" id="astructDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor_html_ae0fe492a01b65ccd59b90bb9d5aad1f5"><div class="ttname"><a href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#ae0fe492a01b65ccd59b90bb9d5aad1f5">DGtal::functors::IIPrincipalCurvatures3DFunctor::eigenValues</a></div><div class="ttdeci">RealVector eigenValues</div><div class="ttdoc">A data member only used for temporary calculations.</div><div class="ttdef"><b>Definition:</b> <a href="IIGeometricFunctors_8h_source.html#l00978">IIGeometricFunctors.h:978</a></div></div>
<div class="ttc" id="astructDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor_html_ae62a9a04b43e44c90857e587e33141e8"><div class="ttname"><a href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#ae62a9a04b43e44c90857e587e33141e8">DGtal::functors::IIPrincipalCurvatures3DFunctor::Value</a></div><div class="ttdeci">Quantity Value</div><div class="ttdef"><b>Definition:</b> <a href="IIGeometricFunctors_8h_source.html#l00923">IIGeometricFunctors.h:923</a></div></div>
</div><!-- fragment -->
<p class="reference">References <a class="el" href="IIGeometricFunctors_8h_source.html#l00972">DGtal::functors::IIPrincipalCurvatures3DFunctor&lt; TSpace, TMatrix &gt;::d6_PIr6</a>, <a class="el" href="IIGeometricFunctors_8h_source.html#l00973">DGtal::functors::IIPrincipalCurvatures3DFunctor&lt; TSpace, TMatrix &gt;::d8_5r</a>, <a class="el" href="IIGeometricFunctors_8h_source.html#l00971">DGtal::functors::IIPrincipalCurvatures3DFunctor&lt; TSpace, TMatrix &gt;::dh5</a>, <a class="el" href="IIGeometricFunctors_8h_source.html#l00978">DGtal::functors::IIPrincipalCurvatures3DFunctor&lt; TSpace, TMatrix &gt;::eigenValues</a>, <a class="el" href="IIGeometricFunctors_8h_source.html#l00976">DGtal::functors::IIPrincipalCurvatures3DFunctor&lt; TSpace, TMatrix &gt;::eigenVectors</a>, and <a class="el" href="classDGtal_1_1EigenDecomposition.html#a062ae590027b118c107dac78f0e9773c">DGtal::EigenDecomposition&lt; TN, TComponent, TMatrix &gt;::getEigenDecomposition()</a>.</p>

</div>
</div>
<h2 class="groupheader">Field Documentation</h2>
<a id="ada43c69f9fa35d7d30ef858c83abbfce"></a>
<h2 class="memtitle"><span class="permalink"><a href="#ada43c69f9fa35d7d30ef858c83abbfce">&#9670;&nbsp;</a></span>d6_PIr6</h2>

<div class="memitem">
<div class="memproto">
<div class="memtemplate">
template&lt;typename TSpace , typename TMatrix  = SimpleMatrix&lt; typename TSpace::RealVector::Component, TSpace::dimension, TSpace::dimension&gt;&gt; </div>
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">double <a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html">DGtal::functors::IIPrincipalCurvatures3DFunctor</a>&lt; TSpace, TMatrix &gt;::d6_PIr6</td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">private</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">

<p class="definition">Definition at line <a class="el" href="IIGeometricFunctors_8h_source.html#l00972">972</a> of file <a class="el" href="IIGeometricFunctors_8h_source.html">IIGeometricFunctors.h</a>.</p>

<p class="reference">Referenced by <a class="el" href="IIGeometricFunctors_8h_source.html#l00960">DGtal::functors::IIPrincipalCurvatures3DFunctor&lt; TSpace, TMatrix &gt;::init()</a>, and <a class="el" href="IIGeometricFunctors_8h_source.html#l00937">DGtal::functors::IIPrincipalCurvatures3DFunctor&lt; TSpace, TMatrix &gt;::operator()()</a>.</p>

</div>
</div>
<a id="a0ee4ebb4e1ec5c07648a009c242a7171"></a>
<h2 class="memtitle"><span class="permalink"><a href="#a0ee4ebb4e1ec5c07648a009c242a7171">&#9670;&nbsp;</a></span>d8_5r</h2>

<div class="memitem">
<div class="memproto">
<div class="memtemplate">
template&lt;typename TSpace , typename TMatrix  = SimpleMatrix&lt; typename TSpace::RealVector::Component, TSpace::dimension, TSpace::dimension&gt;&gt; </div>
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">double <a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html">DGtal::functors::IIPrincipalCurvatures3DFunctor</a>&lt; TSpace, TMatrix &gt;::d8_5r</td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">private</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">

<p class="definition">Definition at line <a class="el" href="IIGeometricFunctors_8h_source.html#l00973">973</a> of file <a class="el" href="IIGeometricFunctors_8h_source.html">IIGeometricFunctors.h</a>.</p>

<p class="reference">Referenced by <a class="el" href="IIGeometricFunctors_8h_source.html#l00960">DGtal::functors::IIPrincipalCurvatures3DFunctor&lt; TSpace, TMatrix &gt;::init()</a>, and <a class="el" href="IIGeometricFunctors_8h_source.html#l00937">DGtal::functors::IIPrincipalCurvatures3DFunctor&lt; TSpace, TMatrix &gt;::operator()()</a>.</p>

</div>
</div>
<a id="a220caf9d2e46f781c24ab05f8312f462"></a>
<h2 class="memtitle"><span class="permalink"><a href="#a220caf9d2e46f781c24ab05f8312f462">&#9670;&nbsp;</a></span>dh5</h2>

<div class="memitem">
<div class="memproto">
<div class="memtemplate">
template&lt;typename TSpace , typename TMatrix  = SimpleMatrix&lt; typename TSpace::RealVector::Component, TSpace::dimension, TSpace::dimension&gt;&gt; </div>
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname">double <a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html">DGtal::functors::IIPrincipalCurvatures3DFunctor</a>&lt; TSpace, TMatrix &gt;::dh5</td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">private</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">

<p class="definition">Definition at line <a class="el" href="IIGeometricFunctors_8h_source.html#l00971">971</a> of file <a class="el" href="IIGeometricFunctors_8h_source.html">IIGeometricFunctors.h</a>.</p>

<p class="reference">Referenced by <a class="el" href="IIGeometricFunctors_8h_source.html#l00960">DGtal::functors::IIPrincipalCurvatures3DFunctor&lt; TSpace, TMatrix &gt;::init()</a>, and <a class="el" href="IIGeometricFunctors_8h_source.html#l00937">DGtal::functors::IIPrincipalCurvatures3DFunctor&lt; TSpace, TMatrix &gt;::operator()()</a>.</p>

</div>
</div>
<a id="ae0fe492a01b65ccd59b90bb9d5aad1f5"></a>
<h2 class="memtitle"><span class="permalink"><a href="#ae0fe492a01b65ccd59b90bb9d5aad1f5">&#9670;&nbsp;</a></span>eigenValues</h2>

<div class="memitem">
<div class="memproto">
<div class="memtemplate">
template&lt;typename TSpace , typename TMatrix  = SimpleMatrix&lt; typename TSpace::RealVector::Component, TSpace::dimension, TSpace::dimension&gt;&gt; </div>
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname"><a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#a53becf0c31439d565f071bc6f88ff8a4">RealVector</a> <a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html">DGtal::functors::IIPrincipalCurvatures3DFunctor</a>&lt; TSpace, TMatrix &gt;::eigenValues</td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">mutable</span><span class="mlabel">private</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">

<p>A data member only used for temporary calculations. </p>

<p class="definition">Definition at line <a class="el" href="IIGeometricFunctors_8h_source.html#l00978">978</a> of file <a class="el" href="IIGeometricFunctors_8h_source.html">IIGeometricFunctors.h</a>.</p>

<p class="reference">Referenced by <a class="el" href="IIGeometricFunctors_8h_source.html#l00937">DGtal::functors::IIPrincipalCurvatures3DFunctor&lt; TSpace, TMatrix &gt;::operator()()</a>.</p>

</div>
</div>
<a id="a2308a8adc544387faa793be87290ee39"></a>
<h2 class="memtitle"><span class="permalink"><a href="#a2308a8adc544387faa793be87290ee39">&#9670;&nbsp;</a></span>eigenVectors</h2>

<div class="memitem">
<div class="memproto">
<div class="memtemplate">
template&lt;typename TSpace , typename TMatrix  = SimpleMatrix&lt; typename TSpace::RealVector::Component, TSpace::dimension, TSpace::dimension&gt;&gt; </div>
<table class="mlabels">
  <tr>
  <td class="mlabels-left">
      <table class="memname">
        <tr>
          <td class="memname"><a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html#af678b7274eb6d3609df360bca9bd1c0d">Matrix</a> <a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html">DGtal::functors::IIPrincipalCurvatures3DFunctor</a>&lt; TSpace, TMatrix &gt;::eigenVectors</td>
        </tr>
      </table>
  </td>
  <td class="mlabels-right">
<span class="mlabels"><span class="mlabel">mutable</span><span class="mlabel">private</span></span>  </td>
  </tr>
</table>
</div><div class="memdoc">

<p>A data member only used for temporary calculations. </p>

<p class="definition">Definition at line <a class="el" href="IIGeometricFunctors_8h_source.html#l00976">976</a> of file <a class="el" href="IIGeometricFunctors_8h_source.html">IIGeometricFunctors.h</a>.</p>

<p class="reference">Referenced by <a class="el" href="IIGeometricFunctors_8h_source.html#l00937">DGtal::functors::IIPrincipalCurvatures3DFunctor&lt; TSpace, TMatrix &gt;::operator()()</a>.</p>

</div>
</div>
<hr/>The documentation for this struct was generated from the following file:<ul>
<li><a class="el" href="IIGeometricFunctors_8h_source.html">IIGeometricFunctors.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="namespaceDGtal.html">DGtal</a></li><li class="navelem"><a class="el" href="namespaceDGtal_1_1functors.html">functors</a></li><li class="navelem"><a class="el" href="structDGtal_1_1functors_1_1IIPrincipalCurvatures3DFunctor.html">IIPrincipalCurvatures3DFunctor</a></li>
    <li class="footer">Generated on Mon Dec 23 2024 13:19:09 for DGtal by <a href="https://www.doxygen.org/index.html"><img class="footer" src="doxygen.svg" width="104" height="31" alt="doxygen"/></a> 1.9.1 </li>
  </ul>
</div>
</body>
</html>