small_gicp/doc_cpp/projection_8hpp_source.html

171 lines
24 KiB
HTML

<!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>small_gicp: /home/runner/work/small_gicp/small_gicp/include/small_gicp/ann/projection.hpp Source File</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="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>
<link href="doxygen.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">small_gicp
</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>
<!-- 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 id="nav-path" class="navpath">
<ul>
<li class="navelem"><a class="el" href="dir_d44c64559bbebec7f509842c48db8b23.html">include</a></li><li class="navelem"><a class="el" href="dir_8f93768a68ffc2367e7100c40b4e616a.html">small_gicp</a></li><li class="navelem"><a class="el" href="dir_c59a592f06e82849d331f7a80f9235b4.html">ann</a></li> </ul>
</div>
</div><!-- top -->
<div class="header">
<div class="headertitle">
<div class="title">projection.hpp</div> </div>
</div><!--header-->
<div class="contents">
<a href="projection_8hpp.html">Go to the documentation of this file.</a><div class="fragment"><div class="line"><a name="l00001"></a><span class="lineno"> 1</span>&#160;<span class="comment">// SPDX-FileCopyrightText: Copyright 2024 Kenji Koide</span></div>
<div class="line"><a name="l00002"></a><span class="lineno"> 2</span>&#160;<span class="comment">// SPDX-License-Identifier: MIT</span></div>
<div class="line"><a name="l00003"></a><span class="lineno"> 3</span>&#160;<span class="preprocessor">#pragma once</span></div>
<div class="line"><a name="l00004"></a><span class="lineno"> 4</span>&#160; </div>
<div class="line"><a name="l00005"></a><span class="lineno"> 5</span>&#160;<span class="preprocessor">#include &lt;array&gt;</span></div>
<div class="line"><a name="l00006"></a><span class="lineno"> 6</span>&#160;<span class="preprocessor">#include &lt;Eigen/Core&gt;</span></div>
<div class="line"><a name="l00007"></a><span class="lineno"> 7</span>&#160;<span class="preprocessor">#include &lt;Eigen/Eigen&gt;</span></div>
<div class="line"><a name="l00008"></a><span class="lineno"> 8</span>&#160; </div>
<div class="line"><a name="l00009"></a><span class="lineno"> 9</span>&#160;<span class="keyword">namespace </span><a class="code" href="namespacesmall__gicp.html">small_gicp</a> {</div>
<div class="line"><a name="l00010"></a><span class="lineno"> 10</span>&#160; </div>
<div class="line"><a name="l00012"></a><span class="lineno"><a class="line" href="structsmall__gicp_1_1ProjectionSetting.html"> 12</a></span>&#160;<span class="keyword">struct </span><a class="code" href="structsmall__gicp_1_1ProjectionSetting.html">ProjectionSetting</a> {</div>
<div class="line"><a name="l00013"></a><span class="lineno"><a class="line" href="structsmall__gicp_1_1ProjectionSetting.html#ae8935cacacd4c725389e0879ae27dd85"> 13</a></span>&#160; <span class="keywordtype">int</span> <a class="code" href="structsmall__gicp_1_1ProjectionSetting.html#ae8935cacacd4c725389e0879ae27dd85">max_scan_count</a> = 128; </div>
<div class="line"><a name="l00014"></a><span class="lineno"> 14</span>&#160;};</div>
<div class="line"><a name="l00015"></a><span class="lineno"> 15</span>&#160; </div>
<div class="line"><a name="l00018"></a><span class="lineno"><a class="line" href="structsmall__gicp_1_1AxisAlignedProjection.html"> 18</a></span>&#160;<span class="keyword">struct </span><a class="code" href="structsmall__gicp_1_1AxisAlignedProjection.html">AxisAlignedProjection</a> {</div>
<div class="line"><a name="l00019"></a><span class="lineno"> 19</span>&#160;<span class="keyword">public</span>:</div>
<div class="line"><a name="l00023"></a><span class="lineno"><a class="line" href="structsmall__gicp_1_1AxisAlignedProjection.html#ae930ad0caeb1dc98980f6f0b58c4cef0"> 23</a></span>&#160; <span class="keywordtype">double</span> <a class="code" href="structsmall__gicp_1_1AxisAlignedProjection.html#ae930ad0caeb1dc98980f6f0b58c4cef0">operator()</a>(<span class="keyword">const</span> Eigen::Vector4d&amp; pt)<span class="keyword"> const </span>{ <span class="keywordflow">return</span> pt[<a class="code" href="structsmall__gicp_1_1AxisAlignedProjection.html#a662ba94e134976ecdc877f263974afcd">axis</a>]; }</div>
<div class="line"><a name="l00024"></a><span class="lineno"> 24</span>&#160; </div>
<div class="line"><a name="l00031"></a><span class="lineno"> 31</span>&#160; <span class="keyword">template</span> &lt;<span class="keyword">typename</span> Po<span class="keywordtype">int</span>Cloud, <span class="keyword">typename</span> IndexConstIterator&gt;</div>
<div class="line"><a name="l00032"></a><span class="lineno"><a class="line" href="structsmall__gicp_1_1AxisAlignedProjection.html#aae18b5b687d7413dd82dfd79e23d8d40"> 32</a></span>&#160; <span class="keyword">static</span> <a class="code" href="structsmall__gicp_1_1AxisAlignedProjection.html">AxisAlignedProjection</a> <a class="code" href="structsmall__gicp_1_1AxisAlignedProjection.html#aae18b5b687d7413dd82dfd79e23d8d40">find_axis</a>(<span class="keyword">const</span> <a class="code" href="structsmall__gicp_1_1PointCloud.html">PointCloud</a>&amp; points, IndexConstIterator first, IndexConstIterator last, <span class="keyword">const</span> <a class="code" href="structsmall__gicp_1_1ProjectionSetting.html">ProjectionSetting</a>&amp; setting) {</div>
<div class="line"><a name="l00033"></a><span class="lineno"> 33</span>&#160; <span class="keyword">const</span> <span class="keywordtype">size_t</span> N = std::distance(first, last);</div>
<div class="line"><a name="l00034"></a><span class="lineno"> 34</span>&#160; Eigen::Vector4d sum_pt = Eigen::Vector4d::Zero();</div>
<div class="line"><a name="l00035"></a><span class="lineno"> 35</span>&#160; Eigen::Vector4d sum_sq = Eigen::Vector4d::Zero();</div>
<div class="line"><a name="l00036"></a><span class="lineno"> 36</span>&#160; </div>
<div class="line"><a name="l00037"></a><span class="lineno"> 37</span>&#160; <span class="keyword">const</span> <span class="keywordtype">size_t</span> step = N &lt; setting.<a class="code" href="structsmall__gicp_1_1ProjectionSetting.html#ae8935cacacd4c725389e0879ae27dd85">max_scan_count</a> ? 1 : N / setting.<a class="code" href="structsmall__gicp_1_1ProjectionSetting.html#ae8935cacacd4c725389e0879ae27dd85">max_scan_count</a>;</div>
<div class="line"><a name="l00038"></a><span class="lineno"> 38</span>&#160; <span class="keyword">const</span> <span class="keywordtype">size_t</span> num_steps = N / step;</div>
<div class="line"><a name="l00039"></a><span class="lineno"> 39</span>&#160; <span class="keywordflow">for</span> (<span class="keywordtype">int</span> i = 0; i &lt; num_steps; i++) {</div>
<div class="line"><a name="l00040"></a><span class="lineno"> 40</span>&#160; <span class="keyword">const</span> <span class="keyword">auto</span> itr = first + step * i;</div>
<div class="line"><a name="l00041"></a><span class="lineno"> 41</span>&#160; <span class="keyword">const</span> Eigen::Vector4d pt = <a class="code" href="namespacesmall__gicp_1_1traits.html#a66b660d3a0c2e81cdf06abc8efebdad1">traits::point</a>(points, *itr);</div>
<div class="line"><a name="l00042"></a><span class="lineno"> 42</span>&#160; sum_pt += pt;</div>
<div class="line"><a name="l00043"></a><span class="lineno"> 43</span>&#160; sum_sq += pt.cwiseProduct(pt);</div>
<div class="line"><a name="l00044"></a><span class="lineno"> 44</span>&#160; }</div>
<div class="line"><a name="l00045"></a><span class="lineno"> 45</span>&#160; </div>
<div class="line"><a name="l00046"></a><span class="lineno"> 46</span>&#160; <span class="keyword">const</span> Eigen::Vector4d mean = sum_pt / sum_pt.w();</div>
<div class="line"><a name="l00047"></a><span class="lineno"> 47</span>&#160; <span class="keyword">const</span> Eigen::Vector4d var = (sum_sq - mean.cwiseProduct(sum_pt));</div>
<div class="line"><a name="l00048"></a><span class="lineno"> 48</span>&#160; </div>
<div class="line"><a name="l00049"></a><span class="lineno"> 49</span>&#160; <span class="keywordflow">return</span> <a class="code" href="structsmall__gicp_1_1AxisAlignedProjection.html">AxisAlignedProjection</a>{var[0] &gt; var[1] ? (var[0] &gt; var[2] ? 0 : 2) : (var[1] &gt; var[2] ? 1 : 2)};</div>
<div class="line"><a name="l00050"></a><span class="lineno"> 50</span>&#160; }</div>
<div class="line"><a name="l00051"></a><span class="lineno"> 51</span>&#160; </div>
<div class="line"><a name="l00052"></a><span class="lineno"> 52</span>&#160;<span class="keyword">public</span>:</div>
<div class="line"><a name="l00053"></a><span class="lineno"><a class="line" href="structsmall__gicp_1_1AxisAlignedProjection.html#a662ba94e134976ecdc877f263974afcd"> 53</a></span>&#160; <span class="keywordtype">int</span> <a class="code" href="structsmall__gicp_1_1AxisAlignedProjection.html#a662ba94e134976ecdc877f263974afcd">axis</a>; </div>
<div class="line"><a name="l00054"></a><span class="lineno"> 54</span>&#160;};</div>
<div class="line"><a name="l00055"></a><span class="lineno"> 55</span>&#160; </div>
<div class="line"><a name="l00058"></a><span class="lineno"><a class="line" href="structsmall__gicp_1_1NormalProjection.html"> 58</a></span>&#160;<span class="keyword">struct </span><a class="code" href="structsmall__gicp_1_1NormalProjection.html">NormalProjection</a> {</div>
<div class="line"><a name="l00059"></a><span class="lineno"> 59</span>&#160;<span class="keyword">public</span>:</div>
<div class="line"><a name="l00063"></a><span class="lineno"><a class="line" href="structsmall__gicp_1_1NormalProjection.html#ab249d80db1af4c8c093bebc8baefa801"> 63</a></span>&#160; <span class="keywordtype">double</span> <a class="code" href="structsmall__gicp_1_1NormalProjection.html#ab249d80db1af4c8c093bebc8baefa801">operator()</a>(<span class="keyword">const</span> Eigen::Vector4d&amp; pt)<span class="keyword"> const </span>{ <span class="keywordflow">return</span> Eigen::Map&lt;const Eigen::Vector3d&gt;(<a class="code" href="structsmall__gicp_1_1NormalProjection.html#acda76271fab04ef9a84e8676e2937e8b">normal</a>.data()).dot(pt.head&lt;3&gt;()); }</div>
<div class="line"><a name="l00064"></a><span class="lineno"> 64</span>&#160; </div>
<div class="line"><a name="l00071"></a><span class="lineno"> 71</span>&#160; <span class="keyword">template</span> &lt;<span class="keyword">typename</span> Po<span class="keywordtype">int</span>Cloud, <span class="keyword">typename</span> IndexConstIterator&gt;</div>
<div class="line"><a name="l00072"></a><span class="lineno"><a class="line" href="structsmall__gicp_1_1NormalProjection.html#a656395ece294173196c38999f5409c0d"> 72</a></span>&#160; <span class="keyword">static</span> <a class="code" href="structsmall__gicp_1_1NormalProjection.html">NormalProjection</a> <a class="code" href="structsmall__gicp_1_1NormalProjection.html#a656395ece294173196c38999f5409c0d">find_axis</a>(<span class="keyword">const</span> <a class="code" href="structsmall__gicp_1_1PointCloud.html">PointCloud</a>&amp; points, IndexConstIterator first, IndexConstIterator last, <span class="keyword">const</span> <a class="code" href="structsmall__gicp_1_1ProjectionSetting.html">ProjectionSetting</a>&amp; setting) {</div>
<div class="line"><a name="l00073"></a><span class="lineno"> 73</span>&#160; <span class="keyword">const</span> <span class="keywordtype">size_t</span> N = std::distance(first, last);</div>
<div class="line"><a name="l00074"></a><span class="lineno"> 74</span>&#160; Eigen::Vector4d sum_pt = Eigen::Vector4d::Zero();</div>
<div class="line"><a name="l00075"></a><span class="lineno"> 75</span>&#160; Eigen::Matrix4d sum_sq = Eigen::Matrix4d::Zero();</div>
<div class="line"><a name="l00076"></a><span class="lineno"> 76</span>&#160; </div>
<div class="line"><a name="l00077"></a><span class="lineno"> 77</span>&#160; <span class="keyword">const</span> <span class="keywordtype">size_t</span> step = N &lt; setting.<a class="code" href="structsmall__gicp_1_1ProjectionSetting.html#ae8935cacacd4c725389e0879ae27dd85">max_scan_count</a> ? 1 : N / setting.<a class="code" href="structsmall__gicp_1_1ProjectionSetting.html#ae8935cacacd4c725389e0879ae27dd85">max_scan_count</a>;</div>
<div class="line"><a name="l00078"></a><span class="lineno"> 78</span>&#160; <span class="keyword">const</span> <span class="keywordtype">size_t</span> num_steps = N / step;</div>
<div class="line"><a name="l00079"></a><span class="lineno"> 79</span>&#160; <span class="keywordflow">for</span> (<span class="keywordtype">int</span> i = 0; i &lt; num_steps; i++) {</div>
<div class="line"><a name="l00080"></a><span class="lineno"> 80</span>&#160; <span class="keyword">const</span> <span class="keyword">auto</span> itr = first + step * i;</div>
<div class="line"><a name="l00081"></a><span class="lineno"> 81</span>&#160; <span class="keyword">const</span> Eigen::Vector4d pt = <a class="code" href="namespacesmall__gicp_1_1traits.html#a66b660d3a0c2e81cdf06abc8efebdad1">traits::point</a>(points, *itr);</div>
<div class="line"><a name="l00082"></a><span class="lineno"> 82</span>&#160; sum_pt += pt;</div>
<div class="line"><a name="l00083"></a><span class="lineno"> 83</span>&#160; sum_sq += pt * pt.transpose();</div>
<div class="line"><a name="l00084"></a><span class="lineno"> 84</span>&#160; }</div>
<div class="line"><a name="l00085"></a><span class="lineno"> 85</span>&#160; </div>
<div class="line"><a name="l00086"></a><span class="lineno"> 86</span>&#160; <span class="keyword">const</span> Eigen::Vector4d mean = sum_pt / sum_pt.w();</div>
<div class="line"><a name="l00087"></a><span class="lineno"> 87</span>&#160; <span class="keyword">const</span> Eigen::Matrix4d <a class="code" href="namespacesmall__gicp_1_1traits.html#ac9146681bcef16a4c73d3fd05ad295c7">cov</a> = (sum_sq - mean * sum_pt.transpose()) / sum_pt.w();</div>
<div class="line"><a name="l00088"></a><span class="lineno"> 88</span>&#160; </div>
<div class="line"><a name="l00089"></a><span class="lineno"> 89</span>&#160; Eigen::SelfAdjointEigenSolver&lt;Eigen::Matrix3d&gt; eig;</div>
<div class="line"><a name="l00090"></a><span class="lineno"> 90</span>&#160; eig.computeDirect(<a class="code" href="namespacesmall__gicp_1_1traits.html#ac9146681bcef16a4c73d3fd05ad295c7">cov</a>.block&lt;3, 3&gt;(0, 0));</div>
<div class="line"><a name="l00091"></a><span class="lineno"> 91</span>&#160; </div>
<div class="line"><a name="l00092"></a><span class="lineno"> 92</span>&#160; <a class="code" href="structsmall__gicp_1_1NormalProjection.html">NormalProjection</a> p;</div>
<div class="line"><a name="l00093"></a><span class="lineno"> 93</span>&#160; Eigen::Map&lt;Eigen::Vector3d&gt;(p.<a class="code" href="structsmall__gicp_1_1NormalProjection.html#acda76271fab04ef9a84e8676e2937e8b">normal</a>.data()) = eig.eigenvectors().col(2);</div>
<div class="line"><a name="l00094"></a><span class="lineno"> 94</span>&#160; <span class="keywordflow">return</span> p;</div>
<div class="line"><a name="l00095"></a><span class="lineno"> 95</span>&#160; }</div>
<div class="line"><a name="l00096"></a><span class="lineno"> 96</span>&#160; </div>
<div class="line"><a name="l00097"></a><span class="lineno"> 97</span>&#160;<span class="keyword">public</span>:</div>
<div class="line"><a name="l00098"></a><span class="lineno"><a class="line" href="structsmall__gicp_1_1NormalProjection.html#acda76271fab04ef9a84e8676e2937e8b"> 98</a></span>&#160; std::array&lt;double, 3&gt; <a class="code" href="structsmall__gicp_1_1NormalProjection.html#acda76271fab04ef9a84e8676e2937e8b">normal</a>; </div>
<div class="line"><a name="l00099"></a><span class="lineno"> 99</span>&#160;};</div>
<div class="line"><a name="l00100"></a><span class="lineno"> 100</span>&#160; </div>
<div class="line"><a name="l00101"></a><span class="lineno"> 101</span>&#160;} <span class="comment">// namespace small_gicp</span></div>
<div class="ttc" id="anamespacesmall__gicp_1_1traits_html_a66b660d3a0c2e81cdf06abc8efebdad1"><div class="ttname"><a href="namespacesmall__gicp_1_1traits.html#a66b660d3a0c2e81cdf06abc8efebdad1">small_gicp::traits::point</a></div><div class="ttdeci">auto point(const T &amp;points, size_t i)</div><div class="ttdoc">Get i-th point. 4D vector is used to take advantage of SIMD intrinsics. The last element must be fill...</div><div class="ttdef"><b>Definition:</b> traits.hpp:40</div></div>
<div class="ttc" id="anamespacesmall__gicp_1_1traits_html_ac9146681bcef16a4c73d3fd05ad295c7"><div class="ttname"><a href="namespacesmall__gicp_1_1traits.html#ac9146681bcef16a4c73d3fd05ad295c7">small_gicp::traits::cov</a></div><div class="ttdeci">auto cov(const T &amp;points, size_t i)</div><div class="ttdoc">Get i-th covariance. Only the top-left 3x3 matrix is filled, and the bottom row and the right col mus...</div><div class="ttdef"><b>Definition:</b> traits.hpp:52</div></div>
<div class="ttc" id="anamespacesmall__gicp_html"><div class="ttname"><a href="namespacesmall__gicp.html">small_gicp</a></div><div class="ttdef"><b>Definition:</b> flat_container.hpp:12</div></div>
<div class="ttc" id="astructsmall__gicp_1_1AxisAlignedProjection_html"><div class="ttname"><a href="structsmall__gicp_1_1AxisAlignedProjection.html">small_gicp::AxisAlignedProjection</a></div><div class="ttdoc">Conventional axis-aligned projection (i.e., selecting any of XYZ axes with the largest variance).</div><div class="ttdef"><b>Definition:</b> projection.hpp:18</div></div>
<div class="ttc" id="astructsmall__gicp_1_1AxisAlignedProjection_html_a662ba94e134976ecdc877f263974afcd"><div class="ttname"><a href="structsmall__gicp_1_1AxisAlignedProjection.html#a662ba94e134976ecdc877f263974afcd">small_gicp::AxisAlignedProjection::axis</a></div><div class="ttdeci">int axis</div><div class="ttdoc">Axis index (0: X, 1: Y, 2: Z)</div><div class="ttdef"><b>Definition:</b> projection.hpp:53</div></div>
<div class="ttc" id="astructsmall__gicp_1_1AxisAlignedProjection_html_aae18b5b687d7413dd82dfd79e23d8d40"><div class="ttname"><a href="structsmall__gicp_1_1AxisAlignedProjection.html#aae18b5b687d7413dd82dfd79e23d8d40">small_gicp::AxisAlignedProjection::find_axis</a></div><div class="ttdeci">static AxisAlignedProjection find_axis(const PointCloud &amp;points, IndexConstIterator first, IndexConstIterator last, const ProjectionSetting &amp;setting)</div><div class="ttdoc">Find the axis with the largest variance.</div><div class="ttdef"><b>Definition:</b> projection.hpp:32</div></div>
<div class="ttc" id="astructsmall__gicp_1_1AxisAlignedProjection_html_ae930ad0caeb1dc98980f6f0b58c4cef0"><div class="ttname"><a href="structsmall__gicp_1_1AxisAlignedProjection.html#ae930ad0caeb1dc98980f6f0b58c4cef0">small_gicp::AxisAlignedProjection::operator()</a></div><div class="ttdeci">double operator()(const Eigen::Vector4d &amp;pt) const</div><div class="ttdoc">Project the point to the selected axis.</div><div class="ttdef"><b>Definition:</b> projection.hpp:23</div></div>
<div class="ttc" id="astructsmall__gicp_1_1NormalProjection_html"><div class="ttname"><a href="structsmall__gicp_1_1NormalProjection.html">small_gicp::NormalProjection</a></div><div class="ttdoc">Normal projection (i.e., selecting the 3D direction with the largest variance of the points).</div><div class="ttdef"><b>Definition:</b> projection.hpp:58</div></div>
<div class="ttc" id="astructsmall__gicp_1_1NormalProjection_html_a656395ece294173196c38999f5409c0d"><div class="ttname"><a href="structsmall__gicp_1_1NormalProjection.html#a656395ece294173196c38999f5409c0d">small_gicp::NormalProjection::find_axis</a></div><div class="ttdeci">static NormalProjection find_axis(const PointCloud &amp;points, IndexConstIterator first, IndexConstIterator last, const ProjectionSetting &amp;setting)</div><div class="ttdoc">Find the direction with the largest variance.</div><div class="ttdef"><b>Definition:</b> projection.hpp:72</div></div>
<div class="ttc" id="astructsmall__gicp_1_1NormalProjection_html_ab249d80db1af4c8c093bebc8baefa801"><div class="ttname"><a href="structsmall__gicp_1_1NormalProjection.html#ab249d80db1af4c8c093bebc8baefa801">small_gicp::NormalProjection::operator()</a></div><div class="ttdeci">double operator()(const Eigen::Vector4d &amp;pt) const</div><div class="ttdoc">Project the point to the normal direction.</div><div class="ttdef"><b>Definition:</b> projection.hpp:63</div></div>
<div class="ttc" id="astructsmall__gicp_1_1NormalProjection_html_acda76271fab04ef9a84e8676e2937e8b"><div class="ttname"><a href="structsmall__gicp_1_1NormalProjection.html#acda76271fab04ef9a84e8676e2937e8b">small_gicp::NormalProjection::normal</a></div><div class="ttdeci">std::array&lt; double, 3 &gt; normal</div><div class="ttdoc">Projection direction.</div><div class="ttdef"><b>Definition:</b> projection.hpp:98</div></div>
<div class="ttc" id="astructsmall__gicp_1_1PointCloud_html"><div class="ttname"><a href="structsmall__gicp_1_1PointCloud.html">small_gicp::PointCloud</a></div><div class="ttdoc">Point cloud.</div><div class="ttdef"><b>Definition:</b> point_cloud.hpp:15</div></div>
<div class="ttc" id="astructsmall__gicp_1_1ProjectionSetting_html"><div class="ttname"><a href="structsmall__gicp_1_1ProjectionSetting.html">small_gicp::ProjectionSetting</a></div><div class="ttdoc">Parameters to control the projection axis search.</div><div class="ttdef"><b>Definition:</b> projection.hpp:12</div></div>
<div class="ttc" id="astructsmall__gicp_1_1ProjectionSetting_html_ae8935cacacd4c725389e0879ae27dd85"><div class="ttname"><a href="structsmall__gicp_1_1ProjectionSetting.html#ae8935cacacd4c725389e0879ae27dd85">small_gicp::ProjectionSetting::max_scan_count</a></div><div class="ttdeci">int max_scan_count</div><div class="ttdoc">Maximum number of points to use for the axis search.</div><div class="ttdef"><b>Definition:</b> projection.hpp:13</div></div>
</div><!-- fragment --></div><!-- contents -->
<!-- start footer part -->
<hr class="footer"/><address class="footer"><small>
Generated by&#160;<a href="https://www.doxygen.org/index.html"><img class="footer" src="doxygen.svg" width="104" height="31" alt="doxygen"/></a> 1.9.1
</small></address>
</body>
</html>