-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathbSplineBasis_example.html
256 lines (254 loc) · 19.7 KB
/
bSplineBasis_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
<!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: bSplineBasis_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('bSplineBasis_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">bSplineBasis_example.cpp</div></div>
</div><!--header-->
<div class="contents">
<div class="textblock"><h1><a class="anchor" id="bSplineBasis_exampleAnnotated"></a>
Annotated source file</h1>
<p>Here is the full file <code>examples/bSplineBasis_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="comment">// Look also in tuturialBasis for more functionality of gsBSplineBasis.</span></div>
<div class="line"> </div>
<div class="line"><span class="preprocessor">#include <iostream></span></div>
<div class="line"><span class="preprocessor">#include <string></span></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="comment">// forward declaration of some utility functions</span></div>
<div class="line"><span class="keywordtype">void</span> print(<span class="keyword">const</span> <a class="code hl_class" href="classgismo_1_1gsBSplineBasis.html">gsBSplineBasis<></a>& bsb, <span class="keyword">const</span> std::string& name);</div>
<div class="line"><span class="keywordtype">void</span> printToParaview(<span class="keyword">const</span> <a class="code hl_class" href="classgismo_1_1gsBSplineBasis.html">gsBSplineBasis<></a>& bsb, <span class="keyword">const</span> std::string& name);</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="comment">// ======================================================================</span></div>
<div class="line"> <span class="comment">// different construction of a knot vector</span></div>
<div class="line"> <span class="comment">// ======================================================================</span></div>
<div class="line"> </div>
<div class="line"> </div>
<div class="line"> real_t a = 0; <span class="comment">// starting knot</span></div>
<div class="line"> real_t b = 1; <span class="comment">// ending knot</span></div>
<div class="line"> <a class="code hl_define" href="gsConfig_8h.html#ad722078677c063a09661059674fb996c">index_t</a> interior = 4; <span class="comment">// number of interior knots</span></div>
<div class="line"> <a class="code hl_define" href="gsConfig_8h.html#ad722078677c063a09661059674fb996c">index_t</a> multEnd = 3; <span class="comment">// multiplicity at the two end knots</span></div>
<div class="line"> <span class="keywordtype">bool</span> paraview = <span class="keyword">false</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">"This is a tutorial on the gsBSplineBasis class."</span>);</div>
<div class="line"> cmd.addReal(<span class="stringliteral">""</span>,<span class="stringliteral">"starting"</span>,<span class="stringliteral">"Starting knot"</span>,a);</div>
<div class="line"> cmd.addReal(<span class="stringliteral">""</span>,<span class="stringliteral">"ending"</span>,<span class="stringliteral">"Ending knot"</span>,b);</div>
<div class="line"> cmd.addInt(<span class="stringliteral">"n"</span>,<span class="stringliteral">"interior"</span>,<span class="stringliteral">"Number of interior knots"</span>,interior);</div>
<div class="line"> cmd.addInt(<span class="stringliteral">"m"</span>,<span class="stringliteral">"mult"</span>,<span class="stringliteral">"Multiplicity at the two end knots"</span>,multEnd);</div>
<div class="line"> cmd.addSwitch(<span class="stringliteral">"plot"</span>,<span class="stringliteral">"Plot with paraview"</span>,paraview);</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"> <a class="code hl_define" href="gsDebug_8h.html#af6177e1e282640c3cb9175565f84a432">gsInfo</a> << <span class="stringliteral">"------------- Constructions -------------\n"</span>;</div>
<div class="line"> </div>
<div class="line"> <span class="keywordtype">int</span> degree = multEnd - 1;</div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsKnotVector.html">gsKnotVector<></a> kv(a, b, interior, multEnd);</div>
<div class="line"> </div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsBSplineBasis.html">gsBSplineBasis<></a> bsb0(kv);</div>
<div class="line"> print(bsb0, <span class="stringliteral">"bsb0"</span>);</div>
<div class="line"> </div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsBSplineBasis.html">gsBSplineBasis<></a> bsb1(a, b, interior, degree);</div>
<div class="line"> print(bsb1, <span class="stringliteral">"bsb1"</span>);</div>
<div class="line"> </div>
<div class="line"> </div>
<div class="line"> <span class="comment">// ======================================================================</span></div>
<div class="line"> <span class="comment">// some properties</span></div>
<div class="line"> <span class="comment">// ======================================================================</span></div>
<div class="line"> </div>
<div class="line"> </div>
<div class="line"> <a class="code hl_define" href="gsDebug_8h.html#af6177e1e282640c3cb9175565f84a432">gsInfo</a> << <span class="stringliteral">"------------- Some properties -------------\n\n"</span>;</div>
<div class="line"> </div>
<div class="line"> <a class="code hl_define" href="gsDebug_8h.html#af6177e1e282640c3cb9175565f84a432">gsInfo</a> << <span class="stringliteral">"bsb0.size(): "</span> << bsb0.size() << <span class="stringliteral">"\n\n"</span></div>
<div class="line"> << <span class="stringliteral">"bsb0.numElements(): "</span> << bsb0.numElements() << <span class="stringliteral">"\n\n"</span></div>
<div class="line"> << <span class="stringliteral">"bsb0.degree(): "</span> << bsb0.degree() << <span class="stringliteral">"\n\n"</span>;</div>
<div class="line"> </div>
<div class="line"> </div>
<div class="line"> <span class="comment">// ======================================================================</span></div>
<div class="line"> <span class="comment">// some operations</span></div>
<div class="line"> <span class="comment">// ======================================================================</span></div>
<div class="line"> </div>
<div class="line"> <a class="code hl_define" href="gsDebug_8h.html#af6177e1e282640c3cb9175565f84a432">gsInfo</a> << <span class="stringliteral">"------------- Some operations -------------\n\n"</span>;</div>
<div class="line"> </div>
<div class="line"> <span class="keyword">const</span> <a class="code hl_class" href="classgismo_1_1gsKnotVector.html">gsKnotVector<></a>& knots = bsb0.knots();</div>
<div class="line"> <a class="code hl_define" href="gsDebug_8h.html#af6177e1e282640c3cb9175565f84a432">gsInfo</a> << <span class="stringliteral">"Knots: \n"</span>;</div>
<div class="line"> knots.<a class="code hl_function" href="classgismo_1_1gsKnotVector.html#a24c5a17355cbf7272ea902189b747839">print</a>(<a class="code hl_define" href="gsDebug_8h.html#af6177e1e282640c3cb9175565f84a432">gsInfo</a>);</div>
<div class="line"> <a class="code hl_define" href="gsDebug_8h.html#af6177e1e282640c3cb9175565f84a432">gsInfo</a> << <span class="stringliteral">"\n\n"</span>;</div>
<div class="line"> </div>
<div class="line"> <span class="keywordflow">if</span> (paraview)</div>
<div class="line"> printToParaview(bsb0, <span class="stringliteral">"basis"</span>);</div>
<div class="line"> </div>
<div class="line"> <a class="code hl_define" href="gsDebug_8h.html#af6177e1e282640c3cb9175565f84a432">gsInfo</a> << <span class="stringliteral">"bsb0.uniformRefine()\n"</span>;</div>
<div class="line"> bsb0.uniformRefine();</div>
<div class="line"> <span class="keywordflow">if</span> (paraview)</div>
<div class="line"> printToParaview(bsb0, <span class="stringliteral">"basisRefined"</span>);</div>
<div class="line"> </div>
<div class="line"> <a class="code hl_define" href="gsDebug_8h.html#af6177e1e282640c3cb9175565f84a432">gsInfo</a> << <span class="stringliteral">"bsb0.degreeElevate()\n"</span>;</div>
<div class="line"> bsb0.degreeElevate();</div>
<div class="line"> <span class="keywordflow">if</span> (paraview)</div>
<div class="line"> printToParaview(bsb0, <span class="stringliteral">"basisElevated"</span>);</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">"files containing the solution.\n"</span>;</div>
<div class="line"> </div>
<div class="line"> <span class="keywordflow">return</span> 0;</div>
<div class="line">}</div>
<div class="line"> </div>
<div class="line"><span class="keywordtype">void</span> printToParaview(<span class="keyword">const</span> <a class="code hl_class" href="classgismo_1_1gsBSplineBasis.html">gsBSplineBasis<></a>& bsb,</div>
<div class="line"> <span class="keyword">const</span> std::string& name)</div>
<div class="line">{</div>
<div class="line"> <a class="code hl_define" href="gsDebug_8h.html#af6177e1e282640c3cb9175565f84a432">gsInfo</a> << <span class="stringliteral">"Writing bsb0 to paraview in a file: "</span> << name << <span class="stringliteral">"\n\n"</span>;</div>
<div class="line"> gsWriteParaview(bsb, name);</div>
<div class="line">}</div>
<div class="line"> </div>
<div class="line"><span class="keywordtype">void</span> print(<span class="keyword">const</span> <a class="code hl_class" href="classgismo_1_1gsBSplineBasis.html">gsBSplineBasis<></a>& bsb,</div>
<div class="line"> <span class="keyword">const</span> std::string& name)</div>
<div class="line">{</div>
<div class="line"> <a class="code hl_define" href="gsDebug_8h.html#af6177e1e282640c3cb9175565f84a432">gsInfo</a> << name << <span class="stringliteral">": \n"</span>;</div>
<div class="line"> bsb.<a class="code hl_function" href="classgismo_1_1gsTensorBSplineBasis_3_011_00_01T_01_4.html#ada405ef534eb608b95aedda5ca0170b6">print</a>(<a class="code hl_define" href="gsDebug_8h.html#af6177e1e282640c3cb9175565f84a432">gsInfo</a>);</div>
<div class="line"> <a class="code hl_define" href="gsDebug_8h.html#af6177e1e282640c3cb9175565f84a432">gsInfo</a> << <span class="stringliteral">"\n\n"</span>;</div>
<div class="line">}</div>
<div class="line"> </div>
<div class="line"> </div>
<div class="ttc" id="aclassgismo_1_1gsBSplineBasis_html"><div class="ttname"><a href="classgismo_1_1gsBSplineBasis.html">gismo::gsBSplineBasis</a></div><div class="ttdoc">A univariate B-spline basis.</div><div class="ttdef"><b>Definition</b> gsBSplineBasis.h:700</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_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_1gsKnotVector_html_a24c5a17355cbf7272ea902189b747839"><div class="ttname"><a href="classgismo_1_1gsKnotVector.html#a24c5a17355cbf7272ea902189b747839">gismo::gsKnotVector::print</a></div><div class="ttdeci">std::ostream & print(std::ostream &os=gsInfo) const</div><div class="ttdef"><b>Definition</b> gsKnotVector.hpp:493</div></div>
<div class="ttc" id="aclassgismo_1_1gsTensorBSplineBasis_3_011_00_01T_01_4_html_ada405ef534eb608b95aedda5ca0170b6"><div class="ttname"><a href="classgismo_1_1gsTensorBSplineBasis_3_011_00_01T_01_4.html#ada405ef534eb608b95aedda5ca0170b6">gismo::gsTensorBSplineBasis< 1, T >::print</a></div><div class="ttdeci">std::ostream & print(std::ostream &os) const</div><div class="ttdoc">Prints the object as a string.</div><div class="ttdef"><b>Definition</b> gsBSplineBasis.hpp:1100</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="agsConfig_8h_html_ad722078677c063a09661059674fb996c"><div class="ttname"><a href="gsConfig_8h.html#ad722078677c063a09661059674fb996c">index_t</a></div><div class="ttdeci">#define index_t</div><div class="ttdef"><b>Definition</b> gsConfig.h:32</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><!-- 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>