-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathbSplineCurve_example.html
267 lines (265 loc) · 22.4 KB
/
bSplineCurve_example.html
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
<!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.9.8"/>
<title>G+Smo: bSplineCurve_example.cpp</title>
<link href="gismodoxy_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/javascript">
/* @license magnet:?xt=urn:btih:d3d9a9a6595521f9666a5e94cc830dab83b65699&dn=expat.txt MIT */
$(document).ready(function() { init_search(); });
/* @license-end */
</script>
<script type="text/x-mathjax-config">
MathJax.Hub.Config({
extensions: ["tex2jax.js"],
jax: ["input/TeX","output/HTML-CSS"],
});
</script>
<script type="text/javascript" async="async" src="https://people.ricam.oeaw.ac.at/gismo/mj/MathJax.js"></script>
<link href="doxygen.css" rel="stylesheet" type="text/css" />
<link href="gismodoxy.css" rel="stylesheet" type="text/css">
<!-- -->
<!-- <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="gslogo-doxy.jpg"/></td>
<td style="padding-left: 0.5em;">
<div id="projectname"><a href="index.html"><font style="font-variant:small-caps;">G+S</font><font style="font-variant:small-caps;"
color="#000000">mo</font></a>
 <span id="projectnumber">25.01.0</span>
</div>
<div id="projectbrief">Geometry + Simulation Modules</div>
</td>
<td> <div id="MSearchBox" class="MSearchBoxInactive">
<span class="left">
<span id="MSearchSelect" onmouseover="return searchBox.OnSearchSelectShow()" onmouseout="return searchBox.OnSearchSelectHide()"> </span>
<input type="text" id="MSearchField" value="" placeholder="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.svg" alt=""/></a>
</span>
</div>
</td>
</tr>
</tbody>
</table>
</div>
<!-- end header part -->
<!-- Generated by Doxygen 1.9.8 -->
<script type="text/javascript">
/* @license magnet:?xt=urn:btih:d3d9a9a6595521f9666a5e94cc830dab83b65699&dn=expat.txt MIT */
var searchBox = new SearchBox("searchBox", "search/",'.html');
/* @license-end */
</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">
/* @license magnet:?xt=urn:btih:d3d9a9a6595521f9666a5e94cc830dab83b65699&dn=expat.txt MIT */
$(document).ready(function(){initNavTree('bSplineCurve_example.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">
<div id="MSearchResults">
<div class="SRPage">
<div id="SRIndex">
<div id="SRResults"></div>
<div class="SRStatus" id="Loading">Loading...</div>
<div class="SRStatus" id="Searching">Searching...</div>
<div class="SRStatus" id="NoMatches">No Matches</div>
</div>
</div>
</div>
</div>
<div><div class="header">
<div class="headertitle"><div class="title">bSplineCurve_example.cpp</div></div>
</div><!--header-->
<div class="contents">
<div class="textblock"><h1><a class="anchor" id="bSplineCurve_exampleAnnotated"></a>
Annotated source file</h1>
<p>Here is the full file <code>examples/bSplineCurve_example.cpp</code>. Clicking on a function or class name will lead you to its reference documentation.</p>
<div class="fragment"><div class="line"> </div>
<div class="line"><span class="preprocessor">#include <iostream></span></div>
<div class="line"> </div>
<div class="line"><span class="preprocessor">#include <<a class="code" href="gismo_8h.html">gismo.h</a>></span></div>
<div class="line"> </div>
<div class="line"><span class="keyword">using namespace </span><a class="code hl_namespace" href="namespacegismo.html">gismo</a>;</div>
<div class="line"> </div>
<div class="line"><span class="keywordtype">int</span> main(<span class="keywordtype">int</span> argc, <span class="keywordtype">char</span> *argv[])</div>
<div class="line">{</div>
<div class="line"> <span class="keywordtype">bool</span> plot = <span class="keyword">false</span>; <span class="comment">// If set to true, paraview file is generated and launched on exit</span></div>
<div class="line"> <span class="keywordtype">bool</span> trim = <span class="keyword">false</span>; <span class="comment">// If set to true, trim/merge operations are displayed</span></div>
<div class="line"> <span class="keywordtype">bool</span> intersect = <span class="keyword">false</span>; <span class="comment">// If set to true, intersection example is displayed</span></div>
<div class="line"> </div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsCmdLine.html">gsCmdLine</a> cmd(<span class="stringliteral">"Tutorial 01 shows the use of BSpline curves."</span>);</div>
<div class="line"> cmd.addSwitch(<span class="stringliteral">"plot"</span>, <span class="stringliteral">"Plot result in ParaView format"</span>, plot);</div>
<div class="line"> cmd.addSwitch(<span class="stringliteral">"trim"</span>, <span class="stringliteral">"Basic trim/merge operations"</span>, trim);</div>
<div class="line"> cmd.addSwitch(<span class="stringliteral">"intersect"</span>, <span class="stringliteral">"Intersection operations"</span>, intersect);</div>
<div class="line"> </div>
<div class="line"> <span class="keywordflow">try</span> { cmd.getValues(argc,argv); } <span class="keywordflow">catch</span> (<span class="keywordtype">int</span> rv) { <span class="keywordflow">return</span> rv; }</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Make a BSpline curve</span></div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsKnotVector.html">gsKnotVector<></a> kv(0, 1, 1, 3);<span class="comment">//start,end,interior knots, start/end multiplicites of knots</span></div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsMatrix.html">gsMatrix<></a> coefs(4, 3);</div>
<div class="line"> coefs << 0, 0, 0,</div>
<div class="line"> 1, 2, 3,</div>
<div class="line"> 2, 1, 4,</div>
<div class="line"> 4, 4, 4;</div>
<div class="line"> </div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsBSpline.html">gsBSpline<></a> curve(kv, coefs);</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Print the Bspline curve</span></div>
<div class="line"> <a class="code hl_define" href="gsDebug_8h.html#af6177e1e282640c3cb9175565f84a432">gsInfo</a> << <span class="stringliteral">"I am a "</span> << curve << <span class="stringliteral">"\n"</span>;</div>
<div class="line"> </div>
<div class="line"> <span class="keywordflow">if</span> (plot)</div>
<div class="line"> { </div>
<div class="line"> <span class="comment">// Output a paraview file</span></div>
<div class="line"> coefs.transposeInPlace();</div>
<div class="line"> gsWriteParaview( curve, <span class="stringliteral">"bsplinecurve0"</span>, 100);</div>
<div class="line"> gsWriteParaview( curve, <span class="stringliteral">"bsplinecurve"</span>, 100, <span class="keyword">true</span>, <span class="keyword">true</span>);</div>
<div class="line"> <a class="code hl_function" href="namespacegismo.html#ae547d2af24b488450296ced53b2418c7">gsWriteParaviewPoints</a>( coefs, <span class="stringliteral">"coefficients"</span>);</div>
<div class="line"> <a class="code hl_function" href="classgismo_1_1gsFileManager.html#a0f725eecb78d3bc1f482326cb9e85eee">gsFileManager::open</a>(<span class="stringliteral">"bsplinecurve0.pvd"</span>);</div>
<div class="line"> }</div>
<div class="line"> <span class="keywordflow">else</span></div>
<div class="line"> <a class="code hl_define" href="gsDebug_8h.html#af6177e1e282640c3cb9175565f84a432">gsInfo</a> << <span class="stringliteral">"Done. No output created, re-run with --plot to get a ParaView "</span></div>
<div class="line"> <span class="stringliteral">"file containing the solution.\n"</span>;</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Basic trim/merge operations on BSpline curves - @Ye</span></div>
<div class="line"> <span class="keywordflow">if</span> (trim)</div>
<div class="line"> {</div>
<div class="line"> <a class="code hl_define" href="gsDebug_8h.html#af6177e1e282640c3cb9175565f84a432">gsInfo</a> << <span class="stringliteral">"Original BSpline curve: "</span> << curve << <span class="stringliteral">"\n"</span>;</div>
<div class="line"> gsWriteParaview( curve, <span class="stringliteral">"originalCurve"</span>, 100); <span class="comment">// Output the original curve</span></div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Segment this BSpline curve between parameters 0.3 and 0.8</span></div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsBSpline.html">gsBSpline<></a> segment = curve.<a class="code hl_function" href="classgismo_1_1gsBSpline.html#ac0edfc609c3f30eee881fa7927a409c4">segmentFromTo</a>(0.3, 0.8);</div>
<div class="line"> <a class="code hl_define" href="gsDebug_8h.html#af6177e1e282640c3cb9175565f84a432">gsInfo</a> << <span class="stringliteral">"Curve segment from u0 = 0.3 to u1 = 0.8: "</span> << segment << <span class="stringliteral">"\n"</span>;</div>
<div class="line"> gsWriteParaview(segment, <span class="stringliteral">"segment"</span>, 100); <span class="comment">// Output the curve segment</span></div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Split the curve at parameter 0.4 into two parts</span></div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsBSpline.html">gsBSpline<></a> segmentLeft, segmentRight;</div>
<div class="line"> curve.<a class="code hl_function" href="classgismo_1_1gsBSpline.html#a7ea4fd12ea1541fc2105cbbb05daa2b7">splitAt</a>(0.4, segmentLeft, segmentRight);</div>
<div class="line"> <a class="code hl_define" href="gsDebug_8h.html#af6177e1e282640c3cb9175565f84a432">gsInfo</a> << <span class="stringliteral">"Curve segment from u0 = 0.0 to u1 = 0.4: "</span> << segmentLeft << <span class="stringliteral">"\n"</span>;</div>
<div class="line"> <a class="code hl_define" href="gsDebug_8h.html#af6177e1e282640c3cb9175565f84a432">gsInfo</a> << <span class="stringliteral">"Curve segment from u0 = 0.4 to u1 = 1.0: "</span> << segmentRight << <span class="stringliteral">"\n"</span>;</div>
<div class="line"> gsWriteParaview( segmentLeft, <span class="stringliteral">"segmentLeft"</span>, 100);</div>
<div class="line"> gsWriteParaview( segmentRight, <span class="stringliteral">"segmentRight"</span>, 100);</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Merge the left and right segments back to the original curve</span></div>
<div class="line"> <span class="comment">// Note: Due to the segmentation, an inner knot value of 0.4 is introduced, while</span></div>
<div class="line"> <span class="comment">// the geometry remains exactly the same as the original one</span></div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsBSpline.html">gsBSpline<></a> mergedCurve = segmentLeft;</div>
<div class="line"> mergedCurve.<a class="code hl_function" href="classgismo_1_1gsBSpline.html#a10f68d8b1f6660bfb809d908418b7905">merge</a>(&segmentRight);</div>
<div class="line"> <a class="code hl_define" href="gsDebug_8h.html#af6177e1e282640c3cb9175565f84a432">gsInfo</a> << <span class="stringliteral">"The merged curve: "</span> << mergedCurve << <span class="stringliteral">"\n"</span>;</div>
<div class="line"> gsWriteParaview( mergedCurve, <span class="stringliteral">"mergedCurve"</span>, 100);</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// convert it into bezier segments</span></div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsMultiPatch.html">gsMultiPatch<></a> bezSegments = mergedCurve.toBezier();</div>
<div class="line"> gsWriteParaview(bezSegments, <span class="stringliteral">"bezierContainer"</span>, 100);</div>
<div class="line"> }</div>
<div class="line"> <span class="keywordflow">else</span></div>
<div class="line"> <a class="code hl_define" href="gsDebug_8h.html#af6177e1e282640c3cb9175565f84a432">gsInfo</a> << <span class="stringliteral">"Done. Re-run with --trim to learn basic trim/merge operations\n"</span>;</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Basic intersection operations between two BSpline curves - @Ye</span></div>
<div class="line"> <span class="keywordflow">if</span> (intersect)</div>
<div class="line"> {</div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsMatrix.html">gsMatrix<real_t></a> ctrPts1(4, 2);</div>
<div class="line"> ctrPts1 << 0,0, 1,1, 2,1, 3,1;</div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsBSpline.html">gsBSpline<real_t></a> bsp1(0, 1, 0, 3, ctrPts1);</div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsMatrix.html">gsMatrix<real_t></a> ctrPts2(4, 2);</div>
<div class="line"> ctrPts2 << 0,0, 1,2, 2,2, 3,0;</div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsBSpline.html">gsBSpline<real_t></a> bsp2(0, 1, 0, 3, ctrPts2);</div>
<div class="line"> </div>
<div class="line"> <span class="keyword">auto</span> intersectPts = bsp1.intersect(bsp2, 1e-5);</div>
<div class="line"> </div>
<div class="line"> <a class="code hl_define" href="gsDebug_8h.html#af6177e1e282640c3cb9175565f84a432">gsInfo</a> << intersectPts.size() << <span class="stringliteral">" intersections are found!"</span> << <span class="stringliteral">"\n"</span>;</div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsMatrix.html">gsMatrix<></a> iPts(bsp1.geoDim(), intersectPts.size());</div>
<div class="line"> <span class="keywordflow">for</span> (<span class="keywordtype">size_t</span> j = 0; j < intersectPts.size(); ++j)</div>
<div class="line"> {</div>
<div class="line"> iPts.col(j) = intersectPts[j].getPoint();</div>
<div class="line"> }</div>
<div class="line"> <span class="keywordflow">if</span> (!intersectPts.empty())</div>
<div class="line"> {</div>
<div class="line"> <a class="code hl_function" href="namespacegismo.html#ae547d2af24b488450296ced53b2418c7">gsWriteParaviewPoints</a>(iPts, <span class="stringliteral">"intersect"</span>);</div>
<div class="line"> }</div>
<div class="line"> </div>
<div class="line"> gsWriteParaview(bsp1, <span class="stringliteral">"bsp1"</span>, 2000);</div>
<div class="line"> gsWriteParaview(bsp2, <span class="stringliteral">"bsp2"</span>, 2000);</div>
<div class="line"> }</div>
<div class="line"> <span class="keywordflow">else</span></div>
<div class="line"> <a class="code hl_define" href="gsDebug_8h.html#af6177e1e282640c3cb9175565f84a432">gsInfo</a> << <span class="stringliteral">"Done. Re-run with --intersect to learn intersection operations\n"</span>;</div>
<div class="line"> </div>
<div class="line"> <span class="keywordflow">return</span> 0;</div>
<div class="line">}</div>
<div class="ttc" id="aclassgismo_1_1gsBSpline_html"><div class="ttname"><a href="classgismo_1_1gsBSpline.html">gismo::gsBSpline</a></div><div class="ttdoc">A B-spline function of one argument, with arbitrary target dimension.</div><div class="ttdef"><b>Definition</b> gsBSpline.h:51</div></div>
<div class="ttc" id="aclassgismo_1_1gsBSpline_html_a10f68d8b1f6660bfb809d908418b7905"><div class="ttname"><a href="classgismo_1_1gsBSpline.html#a10f68d8b1f6660bfb809d908418b7905">gismo::gsBSpline::merge</a></div><div class="ttdeci">void merge(gsGeometry< T > *otherG)</div><div class="ttdoc">Merge other B-spline into this one.</div><div class="ttdef"><b>Definition</b> gsBSpline.hpp:29</div></div>
<div class="ttc" id="aclassgismo_1_1gsBSpline_html_a7ea4fd12ea1541fc2105cbbb05daa2b7"><div class="ttname"><a href="classgismo_1_1gsBSpline.html#a7ea4fd12ea1541fc2105cbbb05daa2b7">gismo::gsBSpline::splitAt</a></div><div class="ttdeci">void splitAt(T u0, gsBSpline< T > &left, gsBSpline< T > &right, T tolerance=1e-15) const</div><div class="ttdef"><b>Definition</b> gsBSpline.hpp:255</div></div>
<div class="ttc" id="aclassgismo_1_1gsBSpline_html_ac0edfc609c3f30eee881fa7927a409c4"><div class="ttname"><a href="classgismo_1_1gsBSpline.html#ac0edfc609c3f30eee881fa7927a409c4">gismo::gsBSpline::segmentFromTo</a></div><div class="ttdeci">gsBSpline< T > segmentFromTo(T u0, T u1, T tolerance=1e-15) const</div><div class="ttdef"><b>Definition</b> gsBSpline.hpp:87</div></div>
<div class="ttc" id="aclassgismo_1_1gsCmdLine_html"><div class="ttname"><a href="classgismo_1_1gsCmdLine.html">gismo::gsCmdLine</a></div><div class="ttdoc">Class for command-line argument parsing.</div><div class="ttdef"><b>Definition</b> gsCmdLine.h:57</div></div>
<div class="ttc" id="aclassgismo_1_1gsFileManager_html_a0f725eecb78d3bc1f482326cb9e85eee"><div class="ttname"><a href="classgismo_1_1gsFileManager.html#a0f725eecb78d3bc1f482326cb9e85eee">gismo::gsFileManager::open</a></div><div class="ttdeci">static void open(const std::string &fn)</div><div class="ttdoc">Opens the file fn using the preferred application of the OS.</div><div class="ttdef"><b>Definition</b> gsFileManager.cpp:688</div></div>
<div class="ttc" id="aclassgismo_1_1gsKnotVector_html"><div class="ttname"><a href="classgismo_1_1gsKnotVector.html">gismo::gsKnotVector</a></div><div class="ttdoc">Class for representing a knot vector.</div><div class="ttdef"><b>Definition</b> gsKnotVector.h:80</div></div>
<div class="ttc" id="aclassgismo_1_1gsMatrix_html"><div class="ttname"><a href="classgismo_1_1gsMatrix.html">gismo::gsMatrix</a></div><div class="ttdoc">A matrix with arbitrary coefficient type and fixed or dynamic size.</div><div class="ttdef"><b>Definition</b> gsMatrix.h:41</div></div>
<div class="ttc" id="aclassgismo_1_1gsMultiPatch_html"><div class="ttname"><a href="classgismo_1_1gsMultiPatch.html">gismo::gsMultiPatch</a></div><div class="ttdoc">Container class for a set of geometry patches and their topology, that is, the interface connections ...</div><div class="ttdef"><b>Definition</b> gsMultiPatch.h:100</div></div>
<div class="ttc" id="agismo_8h_html"><div class="ttname"><a href="gismo_8h.html">gismo.h</a></div><div class="ttdoc">Main header to be included by clients using the G+Smo library.</div></div>
<div class="ttc" id="agsDebug_8h_html_af6177e1e282640c3cb9175565f84a432"><div class="ttname"><a href="gsDebug_8h.html#af6177e1e282640c3cb9175565f84a432">gsInfo</a></div><div class="ttdeci">#define gsInfo</div><div class="ttdef"><b>Definition</b> gsDebug.h:43</div></div>
<div class="ttc" id="anamespacegismo_html"><div class="ttname"><a href="namespacegismo.html">gismo</a></div><div class="ttdoc">The G+Smo namespace, containing all definitions for the library.</div></div>
<div class="ttc" id="anamespacegismo_html_ae547d2af24b488450296ced53b2418c7"><div class="ttname"><a href="namespacegismo.html#ae547d2af24b488450296ced53b2418c7">gismo::gsWriteParaviewPoints</a></div><div class="ttdeci">void gsWriteParaviewPoints(gsMatrix< T > const &X, gsMatrix< T > const &Y, std::string const &fn)</div><div class="ttdoc">Export 2D Point set to Paraview file.</div><div class="ttdef"><b>Definition</b> gsWriteParaview.hpp:1383</div></div>
</div><!-- fragment --> </div></div><!-- contents -->
</div><!-- PageDoc -->
</div><!-- doc-content -->
<!-- start footer part -->
<div id="nav-path" class="navpath"><!-- id is needed for treeview function! -->
<ul>
<li class="footer">Generated on Sun Jan 12 2025 14:02:07 for G+Smo by Doxygen v1.9.8 </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>