-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathadaptRefinementThb_example.html
543 lines (541 loc) · 57.6 KB
/
adaptRefinementThb_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
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
<!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: adaptRefinementThb_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('adaptRefinementThb_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">adaptRefinementThb_example.cpp</div></div>
</div><!--header-->
<div class="contents">
<div class="textblock"><p>This tutorial, shows how to use some functions for adaptive refinement. It is based on <a class="el" href="poisson_example.html">poisson_example.cpp</a>, see there for details.</p>
<p>More details on the hierarchical splines can be found in <a class="el" href="thbSplineBasis_example.html">thbSplineBasis_example.cpp</a> as well as in the <a class="el" href="group__HSplines.html">Hierarchical splines module</a> description.</p>
<h1><a class="anchor" id="sSetup"></a>
Setting up the problem</h1>
<p>As in <a class="el" href="poisson_example.html">poisson_example.cpp</a>, we include the gismo namespace...</p>
<div class="fragment"><div class="line"><span class="preprocessor"># include <<a class="code" href="gismo_8h.html">gismo.h</a>></span></div>
<div class="line"><span class="preprocessor"># include <<a class="code" href="gsAdaptiveRefUtils_8h.html">gsAssembler/gsAdaptiveRefUtils.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="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="agsAdaptiveRefUtils_8h_html"><div class="ttname"><a href="gsAdaptiveRefUtils_8h.html">gsAdaptiveRefUtils.h</a></div><div class="ttdoc">Provides class for adaptive refinement.</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 --><p>...and a command-line argument which allows the user to create a Paraview visualization file.</p>
<div class="fragment"><div class="line"> <span class="keywordtype">bool</span> plot = <span class="keyword">false</span>;</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Number of refinement loops to be done</span></div>
<div class="line"> <a class="code hl_define" href="gsConfig_8h.html#ad722078677c063a09661059674fb996c">index_t</a> numRefinementLoops = 4;</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 on solving a Poisson problem."</span>);</div>
<div class="line"> cmd.addSwitch(<span class="stringliteral">"plot"</span>, <span class="stringliteral">"Create a ParaView visualization file with the solution"</span>, plot);</div>
<div class="line"> cmd.addInt(<span class="stringliteral">"r"</span>, <span class="stringliteral">"numLoop"</span>, <span class="stringliteral">"number of refinement loops"</span>, numRefinementLoops);</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="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="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><!-- fragment --><h2><a class="anchor" id="sProblemDef"></a>
Problem definition</h2>
<p>We consider the Laplace problem on an L-shaped domain \( (-1,1)^2 \setminus [0,1]^2 \) which is represented by three squares.</p>
<div class="image">
<img src="tutAdaRefTHB_geo.png" alt=""/>
<div class="caption">
L-shaped domain, 3 patches.</div></div>
<p>We solve the Laplace problem </p><p class="formulaDsp">
\[ -\Delta u = 0~\mathrm{in}~\Omega,\quad u = u_{ex} ~\mathrm{on}~\partial \Omega, \]
</p>
<p> where the exact solution \(u_{ex}\) is given in polar coordinates by </p><p class="formulaDsp">
\[ u_{ex}(r,\varphi) = r^{\frac{2}{3}} \sin\Big( ( 2\varphi - \pi)/3 \Big). \]
</p>
<p>We define exact solution (in Cartesian coordinates) and the (homogenous) right-hand-side of the PDE: </p><div class="fragment"><div class="line"> <span class="comment">// Define exact solution (will be used for specifying Dirichlet boundary conditions</span></div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsFunctionExpr.html">gsFunctionExpr<></a> g(<span class="stringliteral">"if( y>0, ( (x^2+y^2)^(1.0/3.0) )*sin( (2*atan2(y,x) - pi)/3.0 ), ( (x^2+y^2)^(1.0/3.0) )*sin( (2*atan2(y,x)+3*pi)/3.0 ) )"</span>, 2);</div>
<div class="line"> <span class="comment">// Define source function</span></div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsFunctionExpr.html">gsFunctionExpr<></a> f(<span class="stringliteral">"0"</span>,2);</div>
<div class="ttc" id="aclassgismo_1_1gsFunctionExpr_html"><div class="ttname"><a href="classgismo_1_1gsFunctionExpr.html">gismo::gsFunctionExpr</a></div><div class="ttdoc">Class defining a multivariate (real or vector) function given by a string mathematical expression.</div><div class="ttdef"><b>Definition</b> gsFunctionExpr.h:52</div></div>
</div><!-- fragment --><h2><a class="anchor" id="sGetGeoBasis"></a>
Getting geometry and basis</h2>
<p>The <a class="el" href="classgismo_1_1gsMultiPatch.html" title="Container class for a set of geometry patches and their topology, that is, the interface connections ...">gsMultiPatch</a> object describing the geometry is read in from an .xml-file:</p>
<div class="fragment"><div class="line"> <span class="comment">// Read xml and create gsMultiPatch</span></div>
<div class="line"> std::string fileSrc( <span class="stringliteral">"planar/lshape2d_3patches_thb.xml"</span> );</div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsMultiPatch.html">gsMultiPatch<real_t></a> patches;</div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsReadFile.html">gsReadFile<real_t></a>( fileSrc, patches);</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="aclassgismo_1_1gsReadFile_html"><div class="ttname"><a href="classgismo_1_1gsReadFile.html">gismo::gsReadFile</a></div><div class="ttdoc">Reads an object from a data file, if such the requested object exists in the file.</div><div class="ttdef"><b>Definition</b> gsReadFile.h:43</div></div>
</div><!-- fragment --><p>The information on patch-interfaces and/or boundaries is not necessarily contained in the .xml-file. To make sure that interfaces and boundaries are correctly set up, call <a class="el" href="classgismo_1_1gsMultiPatch.html#a64586846e869ce2e515ff4255b528832" title="Attempt to compute interfaces and boundaries automatically.">gsMultiPatch::computeTopology()</a>.</p>
<div class="fragment"><div class="line"> <span class="comment">// Get all interfaces and boundaries:</span></div>
<div class="line"> patches.<a class="code hl_function" href="classgismo_1_1gsMultiPatch.html#a64586846e869ce2e515ff4255b528832">computeTopology</a>();</div>
<div class="ttc" id="aclassgismo_1_1gsMultiPatch_html_a64586846e869ce2e515ff4255b528832"><div class="ttname"><a href="classgismo_1_1gsMultiPatch.html#a64586846e869ce2e515ff4255b528832">gismo::gsMultiPatch::computeTopology</a></div><div class="ttdeci">bool computeTopology(T tol=1e-4, bool cornersOnly=false, bool tjunctions=false)</div><div class="ttdoc">Attempt to compute interfaces and boundaries automatically.</div><div class="ttdef"><b>Definition</b> gsMultiPatch.hpp:377</div></div>
</div><!-- fragment --><p>To prescribe Dirichlet boundary conditions on all boundaries, we use the gsMultiPatch::const_biterator to iterate through all boundaries:</p>
<div class="fragment"><div class="line"> <a class="code hl_class" href="classgismo_1_1gsBoundaryConditions.html">gsBoundaryConditions<></a> bcInfo;</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// For simplicity, set Dirichlet boundary conditions</span></div>
<div class="line"> <span class="comment">// given by exact solution g on all boundaries:</span></div>
<div class="line"> <span class="keywordflow">for</span> ( gsMultiPatch<>::const_biterator</div>
<div class="line"> bit = patches.<a class="code hl_function" href="classgismo_1_1gsBoxTopology.html#a53808683176e02e2ad93426825006759">bBegin</a>(); bit != patches.<a class="code hl_function" href="classgismo_1_1gsBoxTopology.html#a6151815505965a3f1e4e181b4db227e9">bEnd</a>(); ++bit)</div>
<div class="line"> {</div>
<div class="line"> bcInfo.<a class="code hl_function" href="classgismo_1_1gsBoundaryConditions.html#a1dd0c4f4cdae8082d37db999fc3932d2">addCondition</a>( *bit, <a class="code hl_enumvalue" href="structgismo_1_1condition__type.html#a7aead736a07eaf25623ad7bfa1f0ee2da1b4f535d1f53363637f933002efc2b4b">condition_type::dirichlet</a>, &g );</div>
<div class="line"> }</div>
<div class="ttc" id="aclassgismo_1_1gsBoundaryConditions_html"><div class="ttname"><a href="classgismo_1_1gsBoundaryConditions.html">gismo::gsBoundaryConditions</a></div><div class="ttdoc">Class containing a set of boundary conditions.</div><div class="ttdef"><b>Definition</b> gsBoundaryConditions.h:342</div></div>
<div class="ttc" id="aclassgismo_1_1gsBoundaryConditions_html_a1dd0c4f4cdae8082d37db999fc3932d2"><div class="ttname"><a href="classgismo_1_1gsBoundaryConditions.html#a1dd0c4f4cdae8082d37db999fc3932d2">gismo::gsBoundaryConditions::addCondition</a></div><div class="ttdeci">void addCondition(int p, boxSide s, condition_type::type t, gsFunctionSet< T > *f, short_t unknown=0, bool parametric=false, int comp=-1)</div><div class="ttdoc">Adds another boundary condition.</div><div class="ttdef"><b>Definition</b> gsBoundaryConditions.h:650</div></div>
<div class="ttc" id="aclassgismo_1_1gsBoxTopology_html_a53808683176e02e2ad93426825006759"><div class="ttname"><a href="classgismo_1_1gsBoxTopology.html#a53808683176e02e2ad93426825006759">gismo::gsBoxTopology::bBegin</a></div><div class="ttdeci">const_biterator bBegin() const</div><div class="ttdef"><b>Definition</b> gsBoxTopology.h:139</div></div>
<div class="ttc" id="aclassgismo_1_1gsBoxTopology_html_a6151815505965a3f1e4e181b4db227e9"><div class="ttname"><a href="classgismo_1_1gsBoxTopology.html#a6151815505965a3f1e4e181b4db227e9">gismo::gsBoxTopology::bEnd</a></div><div class="ttdeci">const_biterator bEnd() const</div><div class="ttdef"><b>Definition</b> gsBoxTopology.h:144</div></div>
<div class="ttc" id="astructgismo_1_1condition__type_html_a7aead736a07eaf25623ad7bfa1f0ee2da1b4f535d1f53363637f933002efc2b4b"><div class="ttname"><a href="structgismo_1_1condition__type.html#a7aead736a07eaf25623ad7bfa1f0ee2da1b4f535d1f53363637f933002efc2b4b">gismo::condition_type::dirichlet</a></div><div class="ttdeci">@ dirichlet</div><div class="ttdoc">Dirichlet type.</div><div class="ttdef"><b>Definition</b> gsBoundaryConditions.h:31</div></div>
</div><!-- fragment --><p>Obtaining the <a class="el" href="classgismo_1_1gsMultiBasis.html" title="Holds a set of patch-wise bases and their topology information.">gsMultiBasis</a> from the <a class="el" href="classgismo_1_1gsMultiPatch.html" title="Container class for a set of geometry patches and their topology, that is, the interface connections ...">gsMultiPatch</a> is straightforward, just call the corresponding constructor:</p>
<div class="fragment"><div class="line"> <span class="comment">// Copy basis from the geometry</span></div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsMultiBasis.html">gsMultiBasis<></a> bases( patches );</div>
<div class="ttc" id="aclassgismo_1_1gsMultiBasis_html"><div class="ttname"><a href="classgismo_1_1gsMultiBasis.html">gismo::gsMultiBasis</a></div><div class="ttdoc">Holds a set of patch-wise bases and their topology information.</div><div class="ttdef"><b>Definition</b> gsMultiBasis.h:37</div></div>
</div><!-- fragment --><p>Note that the resulting <a class="el" href="classgismo_1_1gsMultiBasis.html" title="Holds a set of patch-wise bases and their topology information.">gsMultiBasis</a> will be composed of the bases defining the geometry. In this example, the geometry is given as a <a class="el" href="classgismo_1_1gsTHBSpline.html" title="A truncated hierarchical B-Spline function, in d dimensions.">gsTHBSpline</a>. Accordingly, the created <a class="el" href="classgismo_1_1gsMultiBasis.html" title="Holds a set of patch-wise bases and their topology information.">gsMultiBasis</a> is composed of <a class="el" href="classgismo_1_1gsTHBSplineBasis.html" title="Truncated hierarchical B-spline basis.">gsTHBSplineBasis</a>.</p>
<h3><a class="anchor" id="sGeoFromTens"></a>
Get a gsTHBSplineBasis from a gsTensorBSpline</h3>
<p>If the geometry is not given as a <a class="el" href="classgismo_1_1gsTHBSpline.html" title="A truncated hierarchical B-Spline function, in d dimensions.">gsTHBSpline</a>, but as a <a class="el" href="classgismo_1_1gsTensorBSpline.html" title="A tensor product of d B-spline functions, with arbitrary target dimension.">gsTensorBSpline</a>...</p>
<div class="fragment"><div class="line"> std::string fileSrcTens( <span class="stringliteral">"planar/lshape2d_3patches_tens.xml"</span> );</div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsMultiPatch.html">gsMultiPatch<real_t></a> patchesTens;</div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsReadFile.html">gsReadFile<real_t></a>( fileSrcTens, patchesTens);</div>
<div class="line"> patchesTens.<a class="code hl_function" href="classgismo_1_1gsMultiPatch.html#a64586846e869ce2e515ff4255b528832">computeTopology</a>();</div>
</div><!-- fragment --><p>...one can create a <a class="el" href="classgismo_1_1gsTHBSplineBasis.html" title="Truncated hierarchical B-spline basis.">gsTHBSplineBasis</a> as follows:</p>
<div class="fragment"><div class="line"> <span class="comment">// Copy tensor basis</span></div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsMultiBasis.html">gsMultiBasis<real_t></a> basesTens( patchesTens );</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Create a "basisContainer"</span></div>
<div class="line"> std::vector< gsBasis<real_t>* > basisContainer;</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// fill the "basisContainer" with patch-wise...</span></div>
<div class="line"> <span class="keywordflow">for</span> ( <span class="keywordtype">size_t</span> i = 0; i < basesTens.nBases(); i++)</div>
<div class="line"> basisContainer.push_back(<span class="keyword">new</span> <a class="code hl_class" href="classgismo_1_1gsTHBSplineBasis.html">gsTHBSplineBasis<2,real_t></a>( basesTens.basis(i) ));</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// finally, create the gsMultiBasis containing gsTHBSpline ...</span></div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsMultiBasis.html">gsMultiBasis<real_t></a> basesFromTens( basisContainer, patches );</div>
<div class="ttc" id="aclassgismo_1_1gsTHBSplineBasis_html"><div class="ttname"><a href="classgismo_1_1gsTHBSplineBasis.html">gismo::gsTHBSplineBasis</a></div><div class="ttdoc">Truncated hierarchical B-spline basis.</div><div class="ttdef"><b>Definition</b> gsTHBSplineBasis.h:36</div></div>
</div><!-- fragment --><h1><a class="anchor" id="sMainLoop"></a>
Loop</h1>
<p>Instead of just solving the problem once, we wrap a loop around the solving process. For each computed solution, we compute the distribution of the error. Based on this error indication, we mark a certain number of cells for refinement.</p>
<h2><a class="anchor" id="sAdaRefSettings"></a>
Setting up adaptive refinement parameters</h2>
<p>For adaptive refinement, we need to specify a criterion for marking cells, and a parameter. Currently implemented are PUCA, GARU, errorFraction (a.k.a. Doerfel-marking), see <a class="el" href="group__Assembler.html#ga01b4eb9039518258e67fa541adbc5fd1" title="Marks elements/cells for refinement.">gsMarkElementsForRef()</a>.<br />
Setting the parameter to 0 will result in global refinement, setting it to 1 will result in (almost) no refinement.</p>
<div class="fragment"><div class="line"> </div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Specify cell-marking strategy...</span></div>
<div class="line"> MarkingStrategy adaptRefCrit = PUCA;</div>
<div class="line"> <span class="comment">//MarkingStrategy adaptRefCrit = GARU;</span></div>
<div class="line"> <span class="comment">//MarkingStrategy adaptRefCrit = errorFraction;</span></div>
<div class="line"> </div>
<div class="line"> <span class="comment">// ... and parameter.</span></div>
<div class="line"> <span class="keyword">const</span> real_t adaptRefParam = 0.9;</div>
</div><!-- fragment --><h2><a class="anchor" id="solving"></a>
Solving the problem</h2>
<p>If desired or necessary (e.g., because the geometry is given by a very coarse representation), apply initial uniform (i.e., global) refinement steps:</p>
<div class="fragment"><div class="line"> <span class="comment">// Number of initial uniform refinement steps:</span></div>
<div class="line"> <span class="keywordtype">int</span> numInitUniformRefine = 2;</div>
<div class="line"> </div>
<div class="line"> <span class="keywordflow">for</span> (<span class="keywordtype">int</span> i = 0; i < numInitUniformRefine; ++i)</div>
<div class="line"> bases.uniformRefine();</div>
</div><!-- fragment --><p>Wrap a loop around the solution process</p>
<div class="fragment"><div class="line"> <span class="keywordflow">for</span>( <span class="keywordtype">int</span> refLoop = 0; refLoop <= numRefinementLoops; refLoop++)</div>
<div class="line"> {</div>
</div><!-- fragment --><p>The part re solving the equation is the same as in <a class="el" href="poisson_example.html">poisson_example.cpp</a>.</p>
<div class="fragment"><div class="line"> <span class="comment">// Construct assembler</span></div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsPoissonAssembler.html">gsPoissonAssembler<real_t></a> PoissonAssembler(patches,bases,bcInfo,f);</div>
<div class="line"> PoissonAssembler.options().setInt(<span class="stringliteral">"DirichletValues"</span>, dirichlet::l2Projection);</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Generate system matrix and load vector</span></div>
<div class="line"> PoissonAssembler.assemble();</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Initialize the conjugate gradient solver</span></div>
<div class="line"> gsSparseSolver<>::CGDiagonal solver( PoissonAssembler.matrix() );</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Solve the linear system</span></div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsMatrix.html">gsMatrix<></a> solVector = solver.solve( PoissonAssembler.rhs() );</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Construct the isogeometric solution</span></div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsMultiPatch.html">gsMultiPatch<></a> sol;</div>
<div class="line"> PoissonAssembler.constructSolution(solVector, sol);</div>
<div class="line"> <span class="comment">// Associate the solution to the patches (isogeometric field)</span></div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsField.html">gsField<></a> solField(patches, sol);</div>
<div class="ttc" id="aclassgismo_1_1gsField_html"><div class="ttname"><a href="classgismo_1_1gsField.html">gismo::gsField</a></div><div class="ttdoc">A scalar of vector field defined on a m_parametric geometry.</div><div class="ttdef"><b>Definition</b> gsField.h:55</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_1gsPoissonAssembler_html"><div class="ttname"><a href="classgismo_1_1gsPoissonAssembler.html">gismo::gsPoissonAssembler</a></div><div class="ttdoc">Implementation of an (multiple right-hand side) Poisson assembler.</div><div class="ttdef"><b>Definition</b> gsPoissonAssembler.h:37</div></div>
</div><!-- fragment --><h2><a class="anchor" id="sAdaRefExecute"></a>
Adaptive refinement</h2>
<p>After having computed the solution, we compute the local errors. Either by some a posteriori error estimation technique, or, as in this academic example, by computation of the exact error.</p>
<div class="fragment"><div class="line"> <span class="comment">// Compute the error in the H1-seminorm ( = energy norm in this example )</span></div>
<div class="line"> <span class="comment">// using the known exact solution.</span></div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsExprEvaluator.html">gsExprEvaluator<></a> ev;</div>
<div class="line"> ev.<a class="code hl_function" href="classgismo_1_1gsExprEvaluator.html#a0ec66d534b5822bf6f68607c69af9c33">setIntegrationElements</a>(PoissonAssembler.multiBasis());</div>
<div class="line"> gsExprEvaluator<>::geometryMap Gm = ev.<a class="code hl_function" href="classgismo_1_1gsExprEvaluator.html#acd8fe15286df0845bdce1e8689908304">getMap</a>(patches);</div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1expr_1_1gsFeVariable.html">gsExprEvaluator<>::variable</a> is = ev.<a class="code hl_function" href="classgismo_1_1gsExprEvaluator.html#ac146357502c87ba4fe61a7fddf63c86d">getVariable</a>(sol);</div>
<div class="line"> <span class="keyword">auto</span> ms = ev.<a class="code hl_function" href="classgismo_1_1gsExprEvaluator.html#ac146357502c87ba4fe61a7fddf63c86d">getVariable</a>(g, Gm);</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Get the element-wise norms.</span></div>
<div class="line"> ev.<a class="code hl_function" href="classgismo_1_1gsExprEvaluator.html#a465bfc4ca998bdb2717d09465a3bfb4e">integralElWise</a>( ( igrad(is,Gm) - igrad(ms)).sqNorm()*meas(Gm) );</div>
<div class="line"> <span class="keyword">const</span> std::vector<real_t> & eltErrs = ev.<a class="code hl_function" href="classgismo_1_1gsExprEvaluator.html#a33f9c828d5c44f05ecffe0f3eaea31da">elementwise</a>();</div>
<div class="ttc" id="aclassgismo_1_1expr_1_1gsFeVariable_html"><div class="ttname"><a href="classgismo_1_1expr_1_1gsFeVariable.html">gismo::expr::gsFeVariable</a></div><div class="ttdef"><b>Definition</b> gsExpressions.h:928</div></div>
<div class="ttc" id="aclassgismo_1_1gsExprEvaluator_html"><div class="ttname"><a href="classgismo_1_1gsExprEvaluator.html">gismo::gsExprEvaluator</a></div><div class="ttdoc">Generic evaluator of isogeometric expressions.</div><div class="ttdef"><b>Definition</b> gsExprEvaluator.h:39</div></div>
<div class="ttc" id="aclassgismo_1_1gsExprEvaluator_html_a0ec66d534b5822bf6f68607c69af9c33"><div class="ttname"><a href="classgismo_1_1gsExprEvaluator.html#a0ec66d534b5822bf6f68607c69af9c33">gismo::gsExprEvaluator::setIntegrationElements</a></div><div class="ttdeci">void setIntegrationElements(const gsMultiBasis< T > &mesh)</div><div class="ttdoc">Sets the domain of integration.</div><div class="ttdef"><b>Definition</b> gsExprEvaluator.h:110</div></div>
<div class="ttc" id="aclassgismo_1_1gsExprEvaluator_html_a33f9c828d5c44f05ecffe0f3eaea31da"><div class="ttname"><a href="classgismo_1_1gsExprEvaluator.html#a33f9c828d5c44f05ecffe0f3eaea31da">gismo::gsExprEvaluator::elementwise</a></div><div class="ttdeci">const std::vector< T > & elementwise() const</div><div class="ttdoc">Returns an std::vector containing the last computed values per element.</div><div class="ttdef"><b>Definition</b> gsExprEvaluator.h:99</div></div>
<div class="ttc" id="aclassgismo_1_1gsExprEvaluator_html_a465bfc4ca998bdb2717d09465a3bfb4e"><div class="ttname"><a href="classgismo_1_1gsExprEvaluator.html#a465bfc4ca998bdb2717d09465a3bfb4e">gismo::gsExprEvaluator::integralElWise</a></div><div class="ttdeci">T integralElWise(const expr::_expr< E > &expr)</div><div class="ttdoc">Calculates the integral of the expression expr on each element.</div><div class="ttdef"><b>Definition</b> gsExprEvaluator.h:159</div></div>
<div class="ttc" id="aclassgismo_1_1gsExprEvaluator_html_ac146357502c87ba4fe61a7fddf63c86d"><div class="ttname"><a href="classgismo_1_1gsExprEvaluator.html#ac146357502c87ba4fe61a7fddf63c86d">gismo::gsExprEvaluator::getVariable</a></div><div class="ttdeci">variable getVariable(const gsFunctionSet< T > &func, index_t dim=1)</div><div class="ttdoc">Registers func as a variable and returns a handle to it.</div><div class="ttdef"><b>Definition</b> gsExprEvaluator.h:124</div></div>
<div class="ttc" id="aclassgismo_1_1gsExprEvaluator_html_acd8fe15286df0845bdce1e8689908304"><div class="ttname"><a href="classgismo_1_1gsExprEvaluator.html#acd8fe15286df0845bdce1e8689908304">gismo::gsExprEvaluator::getMap</a></div><div class="ttdeci">geometryMap getMap(const gsMultiPatch< T > &mp)</div><div class="ttdoc">Registers mp as an isogeometric geometry map and return a handle to it.</div><div class="ttdef"><b>Definition</b> gsExprEvaluator.h:116</div></div>
</div><!-- fragment --><p>Having computed the element-wise errors, we select those that need to be marked by calling <a class="el" href="group__Assembler.html#ga01b4eb9039518258e67fa541adbc5fd1" title="Marks elements/cells for refinement.">gsMarkElementsForRef()</a>, where the previously defined parameters <em>adaptRefCrit</em> and <em>adaptRefParam</em> are used. The actual refinement is done by <a class="el" href="group__Assembler.html#ga33738199285ae0ea097c96b9c572d786" title="Refine a gsMultiBasis, based on a vector of element-markings.">gsRefineMarkedElements()</a>.</p>
<div class="fragment"><div class="line"> <span class="comment">// Mark elements for refinement, based on the computed local errors and</span></div>
<div class="line"> <span class="comment">// the refinement-criterion and -parameter.</span></div>
<div class="line"> std::vector<bool> elMarked;</div>
<div class="line"> <a class="code hl_function" href="group__Assembler.html#ga01b4eb9039518258e67fa541adbc5fd1">gsMarkElementsForRef</a>( eltErrs, adaptRefCrit, adaptRefParam, elMarked);</div>
<div class="line"> <span class="keywordflow">for</span> (<span class="keywordtype">size_t</span> k=0; k!=elMarked.size(); k++) <a class="code hl_define" href="gsDebug_8h.html#af6177e1e282640c3cb9175565f84a432">gsInfo</a><<<span class="stringliteral">" "</span><<elMarked[k];</div>
<div class="line"> <a class="code hl_define" href="gsDebug_8h.html#af6177e1e282640c3cb9175565f84a432">gsInfo</a><<<span class="stringliteral">"\n"</span>;</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Refine the marked elements with a 1-ring of cells around marked elements</span></div>
<div class="line"> <a class="code hl_function" href="group__Assembler.html#ga33738199285ae0ea097c96b9c572d786">gsRefineMarkedElements</a>( bases, elMarked, 1 );</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// std::vector<bool> elCMarked;</span></div>
<div class="line"> <span class="comment">// for (size_t k=0; k!=eltErrs.size(); k++) eltErrs[k] = -eltErrs[k];</span></div>
<div class="line"> <span class="comment">//gsMarkElementsForRef( eltErrs, adaptRefCrit, adaptRefParam, elCMarked);</span></div>
<div class="line"> <span class="comment">//gsUnrefineMarkedElements( bases, elCMarked, 1 );</span></div>
<div class="line"> </div>
<div class="ttc" id="agroup__Assembler_html_ga01b4eb9039518258e67fa541adbc5fd1"><div class="ttname"><a href="group__Assembler.html#ga01b4eb9039518258e67fa541adbc5fd1">gismo::gsMarkElementsForRef</a></div><div class="ttdeci">void gsMarkElementsForRef(const std::vector< T > &elError, int refCriterion, T refParameter, std::vector< bool > &elMarked)</div><div class="ttdoc">Marks elements/cells for refinement.</div><div class="ttdef"><b>Definition</b> gsAdaptiveRefUtils.h:196</div></div>
<div class="ttc" id="agroup__Assembler_html_ga33738199285ae0ea097c96b9c572d786"><div class="ttname"><a href="group__Assembler.html#ga33738199285ae0ea097c96b9c572d786">gismo::gsRefineMarkedElements</a></div><div class="ttdeci">void gsRefineMarkedElements(gsMultiBasis< T > &basis, const std::vector< bool > &elMarked, index_t refExtension=0)</div><div class="ttdoc">Refine a gsMultiBasis, based on a vector of element-markings.</div><div class="ttdef"><b>Definition</b> gsAdaptiveRefUtils.h:294</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><!-- fragment --><p><b>Note/Warning:</b><br />
As of now, we cannot access a single element of a <a class="el" href="classgismo_1_1gsHTensorBasis.html" title="Class representing a (scalar) hierarchical tensor basis of functions .">gsHTensorBasis</a> directly (e.g., by some index). <a class="el" href="group__Assembler.html#ga01b4eb9039518258e67fa541adbc5fd1" title="Marks elements/cells for refinement.">gsMarkElementsForRef()</a> and <a class="el" href="group__Assembler.html#ga33738199285ae0ea097c96b9c572d786" title="Refine a gsMultiBasis, based on a vector of element-markings.">gsRefineMarkedElements()</a> rely on the assumption that the <a class="el" href="classgismo_1_1gsDomainIterator.html" title="Class which enables iteration over all elements of a parameter domain.">gsDomainIterator</a> used in these two functions will iterate over the elements in the same order. <br />
This is also the reason why <em>elMarked</em> is a vector of booleans indicating "refine!" or "don't refine!" for each element (instead of, e.g., a list of inidces of elements).</p>
<p>Recall from <a class="el" href="thbSplineBasis_example.html">thbSplineBasis_example.cpp</a> that the refined area must contain the support of at least one basis function. Due to this, we also refine the 1-ring of cells around marked cells.</p>
<p>Adaptive refinement on a <a class="el" href="classgismo_1_1gsMultiPatch.html" title="Container class for a set of geometry patches and their topology, that is, the interface connections ...">gsMultiPatch</a> geometry can result in interfaces which are no longer fully matching. Thus, it is important to call <a class="el" href="classgismo_1_1gsMultiBasis.html#ad30dc690b202454b3dbc31efe58fe244" title="Checks if the interfaces bivec are fully matching, and if not, repairs them, i.e.,...">gsMultiBasis::repairInterfaces()</a>.</p>
<div class="fragment"><div class="line"> <span class="comment">// Call repair interfaces to make sure that the new meshes</span></div>
<div class="line"> <span class="comment">// match along patch interfaces.</span></div>
<div class="line"> bases.repairInterfaces( patches.<a class="code hl_function" href="classgismo_1_1gsBoxTopology.html#a0a094ffe7da8592546292486cea117da">interfaces</a>() );</div>
<div class="ttc" id="aclassgismo_1_1gsBoxTopology_html_a0a094ffe7da8592546292486cea117da"><div class="ttname"><a href="classgismo_1_1gsBoxTopology.html#a0a094ffe7da8592546292486cea117da">gismo::gsBoxTopology::interfaces</a></div><div class="ttdeci">const ifContainer & interfaces() const</div><div class="ttdoc">Return the vector of interfaces.</div><div class="ttdef"><b>Definition</b> gsBoxTopology.h:252</div></div>
</div><!-- fragment --><h1><a class="anchor" id="sPlot"></a>
Plotting final solution</h1>
<p>In the last iteration, we export the solution to paraview files</p>
<div class="fragment"><div class="line"> <span class="comment">// Export the final solution</span></div>
<div class="line"> <span class="keywordflow">if</span>( plot && refLoop == numRefinementLoops )</div>
<div class="line"> {</div>
<div class="line"> <span class="comment">// Write the computed solution to paraview files</span></div>
<div class="line"> gsWriteParaview<>(solField, <span class="stringliteral">"adaptRef"</span>, 1000, <span class="keyword">true</span>);</div>
<div class="line"> }</div>
</div><!-- fragment --><p>Plot the solution in Paraview</p>
<div class="fragment"><div class="line"> <span class="keywordflow">if</span>( plot )</div>
<div class="line"> {</div>
<div class="line"> <span class="comment">// Run paraview</span></div>
<div class="line"> <a class="code hl_function" href="classgismo_1_1gsFileManager.html#a0f725eecb78d3bc1f482326cb9e85eee">gsFileManager::open</a>(<span class="stringliteral">"adaptRef.pvd"</span>);</div>
<div class="line"> }</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><!-- fragment --><div class="image">
<img src="tutAdaRefTHB_LSres4.png" alt=""/>
<div class="caption">
Mesh (left, coloured by patch-index) and solution after 4 adaptive refinement steps.</div></div>
<h1><a class="anchor" id="adaptRefinementThb_exampleAnnotated"></a>
Annotated source file</h1>
<p>Here is the full file <code>examples/adaptRefinementThb_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 <<a class="code" href="gismo_8h.html">gismo.h</a>></span></div>
<div class="line"><span class="preprocessor"># include <<a class="code" href="gsAdaptiveRefUtils_8h.html">gsAssembler/gsAdaptiveRefUtils.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>;</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Number of refinement loops to be done</span></div>
<div class="line"> <a class="code hl_define" href="gsConfig_8h.html#ad722078677c063a09661059674fb996c">index_t</a> numRefinementLoops = 4;</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 on solving a Poisson problem."</span>);</div>
<div class="line"> cmd.addSwitch(<span class="stringliteral">"plot"</span>, <span class="stringliteral">"Create a ParaView visualization file with the solution"</span>, plot);</div>
<div class="line"> cmd.addInt(<span class="stringliteral">"r"</span>, <span class="stringliteral">"numLoop"</span>, <span class="stringliteral">"number of refinement loops"</span>, numRefinementLoops);</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">// --------------- specify exact solution and right-hand-side ---------------</span></div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Define exact solution (will be used for specifying Dirichlet boundary conditions</span></div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsFunctionExpr.html">gsFunctionExpr<></a> g(<span class="stringliteral">"if( y>0, ( (x^2+y^2)^(1.0/3.0) )*sin( (2*atan2(y,x) - pi)/3.0 ), ( (x^2+y^2)^(1.0/3.0) )*sin( (2*atan2(y,x)+3*pi)/3.0 ) )"</span>, 2);</div>
<div class="line"> <span class="comment">// Define source function</span></div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsFunctionExpr.html">gsFunctionExpr<></a> f(<span class="stringliteral">"0"</span>,2);</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Print out source function and solution</span></div>
<div class="line"> <a class="code hl_define" href="gsDebug_8h.html#af6177e1e282640c3cb9175565f84a432">gsInfo</a><<<span class="stringliteral">"Source function "</span> << f << <span class="stringliteral">"\n"</span>;</div>
<div class="line"> <a class="code hl_define" href="gsDebug_8h.html#af6177e1e282640c3cb9175565f84a432">gsInfo</a><<<span class="stringliteral">"Exact solution "</span> << g << <span class="stringliteral">"\n\n"</span>;</div>
<div class="line"> </div>
<div class="line"> </div>
<div class="line"> <span class="comment">// --------------- read geometry from file ---------------</span></div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Read xml and create gsMultiPatch</span></div>
<div class="line"> std::string fileSrc( <span class="stringliteral">"planar/lshape2d_3patches_thb.xml"</span> );</div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsMultiPatch.html">gsMultiPatch<real_t></a> patches;</div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsReadFile.html">gsReadFile<real_t></a>( fileSrc, patches);</div>
<div class="line"> <a class="code hl_define" href="gsDebug_8h.html#af6177e1e282640c3cb9175565f84a432">gsInfo</a> << <span class="stringliteral">"The domain is a "</span><< patches <<<span class="stringliteral">"\n"</span>;</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Get all interfaces and boundaries:</span></div>
<div class="line"> patches.<a class="code hl_function" href="classgismo_1_1gsMultiPatch.html#a64586846e869ce2e515ff4255b528832">computeTopology</a>();</div>
<div class="line"> </div>
<div class="line"> </div>
<div class="line"> std::string fileSrcTens( <span class="stringliteral">"planar/lshape2d_3patches_tens.xml"</span> );</div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsMultiPatch.html">gsMultiPatch<real_t></a> patchesTens;</div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsReadFile.html">gsReadFile<real_t></a>( fileSrcTens, patchesTens);</div>
<div class="line"> patchesTens.<a class="code hl_function" href="classgismo_1_1gsMultiPatch.html#a64586846e869ce2e515ff4255b528832">computeTopology</a>();</div>
<div class="line"> </div>
<div class="line"> </div>
<div class="line"> <span class="comment">// --------------- add bonudary conditions ---------------</span></div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsBoundaryConditions.html">gsBoundaryConditions<></a> bcInfo;</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// For simplicity, set Dirichlet boundary conditions</span></div>
<div class="line"> <span class="comment">// given by exact solution g on all boundaries:</span></div>
<div class="line"> <span class="keywordflow">for</span> ( gsMultiPatch<>::const_biterator</div>
<div class="line"> bit = patches.<a class="code hl_function" href="classgismo_1_1gsBoxTopology.html#a53808683176e02e2ad93426825006759">bBegin</a>(); bit != patches.<a class="code hl_function" href="classgismo_1_1gsBoxTopology.html#a6151815505965a3f1e4e181b4db227e9">bEnd</a>(); ++bit)</div>
<div class="line"> {</div>
<div class="line"> bcInfo.<a class="code hl_function" href="classgismo_1_1gsBoundaryConditions.html#a1dd0c4f4cdae8082d37db999fc3932d2">addCondition</a>( *bit, <a class="code hl_enumvalue" href="structgismo_1_1condition__type.html#a7aead736a07eaf25623ad7bfa1f0ee2da1b4f535d1f53363637f933002efc2b4b">condition_type::dirichlet</a>, &g );</div>
<div class="line"> }</div>
<div class="line"> </div>
<div class="line"> </div>
<div class="line"> <span class="comment">// --------------- set up basis ---------------</span></div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Copy basis from the geometry</span></div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsMultiBasis.html">gsMultiBasis<></a> bases( patches );</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Copy tensor basis</span></div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsMultiBasis.html">gsMultiBasis<real_t></a> basesTens( patchesTens );</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Create a "basisContainer"</span></div>
<div class="line"> std::vector< gsBasis<real_t>* > basisContainer;</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// fill the "basisContainer" with patch-wise...</span></div>
<div class="line"> <span class="keywordflow">for</span> ( <span class="keywordtype">size_t</span> i = 0; i < basesTens.nBases(); i++)</div>
<div class="line"> basisContainer.push_back(<span class="keyword">new</span> <a class="code hl_class" href="classgismo_1_1gsTHBSplineBasis.html">gsTHBSplineBasis<2,real_t></a>( basesTens.basis(i) ));</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// finally, create the gsMultiBasis containing gsTHBSpline ...</span></div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsMultiBasis.html">gsMultiBasis<real_t></a> basesFromTens( basisContainer, patches );</div>
<div class="line"> </div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Number of initial uniform refinement steps:</span></div>
<div class="line"> <span class="keywordtype">int</span> numInitUniformRefine = 2;</div>
<div class="line"> </div>
<div class="line"> <span class="keywordflow">for</span> (<span class="keywordtype">int</span> i = 0; i < numInitUniformRefine; ++i)</div>
<div class="line"> bases.uniformRefine();</div>
<div class="line"> </div>
<div class="line"> </div>
<div class="line"> <span class="comment">// --------------- set up adaptive refinement loop ---------------</span></div>
<div class="line"> </div>
<div class="line"> </div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Specify cell-marking strategy...</span></div>
<div class="line"> MarkingStrategy adaptRefCrit = PUCA;</div>
<div class="line"> <span class="comment">//MarkingStrategy adaptRefCrit = GARU;</span></div>
<div class="line"> <span class="comment">//MarkingStrategy adaptRefCrit = errorFraction;</span></div>
<div class="line"> </div>
<div class="line"> <span class="comment">// ... and parameter.</span></div>
<div class="line"> <span class="keyword">const</span> real_t adaptRefParam = 0.9;</div>
<div class="line"> </div>
<div class="line"> </div>
<div class="line"> <span class="comment">// --------------- adaptive refinement loop ---------------</span></div>
<div class="line"> </div>
<div class="line"> <span class="keywordflow">for</span>( <span class="keywordtype">int</span> refLoop = 0; refLoop <= numRefinementLoops; refLoop++)</div>
<div class="line"> {</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// --------------- solving ---------------</span></div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Construct assembler</span></div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsPoissonAssembler.html">gsPoissonAssembler<real_t></a> PoissonAssembler(patches,bases,bcInfo,f);</div>
<div class="line"> PoissonAssembler.options().setInt(<span class="stringliteral">"DirichletValues"</span>, dirichlet::l2Projection);</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Generate system matrix and load vector</span></div>
<div class="line"> PoissonAssembler.assemble();</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Initialize the conjugate gradient solver</span></div>
<div class="line"> gsSparseSolver<>::CGDiagonal solver( PoissonAssembler.matrix() );</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Solve the linear system</span></div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsMatrix.html">gsMatrix<></a> solVector = solver.solve( PoissonAssembler.rhs() );</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Construct the isogeometric solution</span></div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsMultiPatch.html">gsMultiPatch<></a> sol;</div>
<div class="line"> PoissonAssembler.constructSolution(solVector, sol);</div>
<div class="line"> <span class="comment">// Associate the solution to the patches (isogeometric field)</span></div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsField.html">gsField<></a> solField(patches, sol);</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// --------------- error estimation/computation ---------------</span></div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Compute the error in the H1-seminorm ( = energy norm in this example )</span></div>
<div class="line"> <span class="comment">// using the known exact solution.</span></div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1gsExprEvaluator.html">gsExprEvaluator<></a> ev;</div>
<div class="line"> ev.<a class="code hl_function" href="classgismo_1_1gsExprEvaluator.html#a0ec66d534b5822bf6f68607c69af9c33">setIntegrationElements</a>(PoissonAssembler.multiBasis());</div>
<div class="line"> gsExprEvaluator<>::geometryMap Gm = ev.<a class="code hl_function" href="classgismo_1_1gsExprEvaluator.html#acd8fe15286df0845bdce1e8689908304">getMap</a>(patches);</div>
<div class="line"> <a class="code hl_class" href="classgismo_1_1expr_1_1gsFeVariable.html">gsExprEvaluator<>::variable</a> is = ev.<a class="code hl_function" href="classgismo_1_1gsExprEvaluator.html#ac146357502c87ba4fe61a7fddf63c86d">getVariable</a>(sol);</div>
<div class="line"> <span class="keyword">auto</span> ms = ev.<a class="code hl_function" href="classgismo_1_1gsExprEvaluator.html#ac146357502c87ba4fe61a7fddf63c86d">getVariable</a>(g, Gm);</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Get the element-wise norms.</span></div>
<div class="line"> ev.<a class="code hl_function" href="classgismo_1_1gsExprEvaluator.html#a465bfc4ca998bdb2717d09465a3bfb4e">integralElWise</a>( ( igrad(is,Gm) - igrad(ms)).sqNorm()*meas(Gm) );</div>
<div class="line"> <span class="keyword">const</span> std::vector<real_t> & eltErrs = ev.<a class="code hl_function" href="classgismo_1_1gsExprEvaluator.html#a33f9c828d5c44f05ecffe0f3eaea31da">elementwise</a>();</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// --------------- adaptive refinement ---------------</span></div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Mark elements for refinement, based on the computed local errors and</span></div>
<div class="line"> <span class="comment">// the refinement-criterion and -parameter.</span></div>
<div class="line"> std::vector<bool> elMarked;</div>
<div class="line"> <a class="code hl_function" href="group__Assembler.html#ga01b4eb9039518258e67fa541adbc5fd1">gsMarkElementsForRef</a>( eltErrs, adaptRefCrit, adaptRefParam, elMarked);</div>
<div class="line"> <span class="keywordflow">for</span> (<span class="keywordtype">size_t</span> k=0; k!=elMarked.size(); k++) <a class="code hl_define" href="gsDebug_8h.html#af6177e1e282640c3cb9175565f84a432">gsInfo</a><<<span class="stringliteral">" "</span><<elMarked[k];</div>
<div class="line"> <a class="code hl_define" href="gsDebug_8h.html#af6177e1e282640c3cb9175565f84a432">gsInfo</a><<<span class="stringliteral">"\n"</span>;</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Refine the marked elements with a 1-ring of cells around marked elements</span></div>
<div class="line"> <a class="code hl_function" href="group__Assembler.html#ga33738199285ae0ea097c96b9c572d786">gsRefineMarkedElements</a>( bases, elMarked, 1 );</div>
<div class="line"> </div>
<div class="line"> <span class="comment">// std::vector<bool> elCMarked;</span></div>
<div class="line"> <span class="comment">// for (size_t k=0; k!=eltErrs.size(); k++) eltErrs[k] = -eltErrs[k];</span></div>
<div class="line"> <span class="comment">//gsMarkElementsForRef( eltErrs, adaptRefCrit, adaptRefParam, elCMarked);</span></div>
<div class="line"> <span class="comment">//gsUnrefineMarkedElements( bases, elCMarked, 1 );</span></div>
<div class="line"> </div>
<div class="line"> </div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Call repair interfaces to make sure that the new meshes</span></div>
<div class="line"> <span class="comment">// match along patch interfaces.</span></div>
<div class="line"> bases.repairInterfaces( patches.<a class="code hl_function" href="classgismo_1_1gsBoxTopology.html#a0a094ffe7da8592546292486cea117da">interfaces</a>() );</div>
<div class="line"> </div>
<div class="line"> </div>
<div class="line"> <span class="comment">// Export the final solution</span></div>
<div class="line"> <span class="keywordflow">if</span>( plot && refLoop == numRefinementLoops )</div>
<div class="line"> {</div>
<div class="line"> <span class="comment">// Write the computed solution to paraview files</span></div>
<div class="line"> gsWriteParaview<>(solField, <span class="stringliteral">"adaptRef"</span>, 1000, <span class="keyword">true</span>);</div>
<div class="line"> }</div>
<div class="line"> </div>
<div class="line"> }</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">// Run paraview</span></div>
<div class="line"> <a class="code hl_function" href="classgismo_1_1gsFileManager.html#a0f725eecb78d3bc1f482326cb9e85eee">gsFileManager::open</a>(<span class="stringliteral">"adaptRef.pvd"</span>);</div>
<div class="line"> }</div>
<div class="line"> <span class="keywordflow">else</span></div>
<div class="line"> {</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 Plotting image data.\n"</span>;</div>
<div class="line"> }</div>
<div class="line"> </div>
<div class="line"> <span class="keywordflow">return</span> EXIT_SUCCESS;</div>
<div class="line"> </div>
<div class="line">}<span class="comment">// end main</span></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>