-
Notifications
You must be signed in to change notification settings - Fork 3
/
Copy path9_mixed_problems.html
499 lines (481 loc) · 52.7 KB
/
9_mixed_problems.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
<!DOCTYPE html>
<html lang="en" data-content_root="./">
<head>
<meta charset="utf-8" />
<meta name="viewport" content="width=device-width, initial-scale=1.0" /><meta name="viewport" content="width=device-width, initial-scale=1" />
<title>9. Mixed problems — Finite element course 2024.0 documentation</title>
<link rel="stylesheet" type="text/css" href="_static/pygments.css?v=03e43079" />
<link rel="stylesheet" type="text/css" href="_static/fenics.css?v=7793b76c" />
<link rel="stylesheet" type="text/css" href="_static/proof.css" />
<link rel="stylesheet" type="text/css" href="_static/sphinx-design.min.css?v=95c83b7e" />
<script src="_static/documentation_options.js?v=7ff0cb77"></script>
<script src="_static/doctools.js?v=9a2dae69"></script>
<script src="_static/sphinx_highlight.js?v=dc90522c"></script>
<script src="_static/proof.js"></script>
<script src="_static/design-tabs.js?v=36754332"></script>
<script async="async" src="https://cdn.mathjax.org/mathjax/latest/MathJax.js?config=TeX-AMS-MML_HTMLorMML"></script>
<link rel="index" title="Index" href="genindex.html" />
<link rel="search" title="Search" href="search.html" />
<link rel="next" title="10. References" href="zbibliography.html" />
<link rel="prev" title="8. Nonlinear problems" href="8_nonlinear_problems.html" />
<!--[if lte IE 6]>
<link rel="stylesheet" href="_static/ie6.css" type="text/css" media="screen" charset="utf-8" />
<![endif]-->
<link rel="stylesheet" href="_static/featured.css">
<link rel="shortcut icon" href="_static/icon.ico" />
</head><body>
<div class="wrapper">
<a href="index.html"><img src="_static/banner.svg" width="900px" alt="FInAT Project Banner" /></a>
<div id="access">
<div class="menu">
<ul>
<li class="page_item"><a href="https://github.com/finite-element/finite-element-course" title="GitHub">GitHub</a></li>
</ul>
</div><!-- .menu -->
</div><!-- #access -->
</div><!-- #wrapper -->
<div class="document">
<div class="documentwrapper">
<div class="bodywrapper">
<div class="body" role="main">
<section id="mixed-problems">
<span id="mixed"></span><h1><span class="section-number">9. </span>Mixed problems<a class="headerlink" href="#mixed-problems" title="Link to this heading">¶</a></h1>
<div class="admonition note">
<p class="admonition-title">Note</p>
<p>This section is the mastery exercise for this module. This exercise
is explicitly intended to test whether you can bring together what
has been learned in the rest of the module in order to go beyond
what has been covered in lectures and labs.</p>
<p>This exercise is not a part of the third year version of this module.</p>
</div>
<p>As an example of a mixed problem, let’s take the Stokes problem presented in
<a class="reference internal" href="L6_stokes.html#stokes"><span class="std std-numref">Section 6</span></a> of the analysis part of the course. The weak form of the
Stokes problem presented in <a class="reference internal" href="L6_stokes.html#weak-stokes"><span class="std std-numref">Definition 6.1</span></a> is find
<span class="math notranslate nohighlight">\((u,p)\in V\times Q\)</span> such that:</p>
<div class="math notranslate nohighlight" id="equation-stokes-ch9">
<span class="eqno">(9.1)<a class="headerlink" href="#equation-stokes-ch9" title="Link to this equation">¶</a></span>\[ \begin{align}\begin{aligned}a(u,v) + b(v, p) & = \int_\Omega f\cdot v\,\mathrm{d}\, x,\\b(u,q) & = 0, \quad \forall (v,q) \in V\times Q,\end{aligned}\end{align} \]</div>
<p>where</p>
<div class="math notranslate nohighlight" id="equation-9-mixed-problems-0">
<span class="eqno">(9.2)<a class="headerlink" href="#equation-9-mixed-problems-0" title="Link to this equation">¶</a></span>\[ \begin{align}\begin{aligned}a(u,v) & = \mu\int_\Omega \epsilon(u):\epsilon(v)\,\mathrm{d}\, x,\\b(v,q) & = \int_\Omega q \nabla\cdot v\,\mathrm{d}\, x,\\V & =(\mathring{H}^1(\Omega))^d,\\Q & =\mathring{L}^2(\Omega),\end{aligned}\end{align} \]</div>
<p>and where <span class="math notranslate nohighlight">\((\mathring{H}^1(\Omega))^d\)</span> is the subspace of <span class="math notranslate nohighlight">\(H^1(\Omega)^d\)</span> for
which all components vanish on the boundary, and <span class="math notranslate nohighlight">\(\mathring{L}^2(\Omega)\)</span> is
the subspace of <span class="math notranslate nohighlight">\(L^2(\Omega)\)</span> for which the integral of the function over the
domain vanishes. This last constraint was introduced to remove the null space
of constant pressure functions from the system. This constraint introduces a
little complexity into the implementation. So instead, we will redefine
<span class="math notranslate nohighlight">\(\mathring{L}^2(\Omega)\)</span> to be the subspace of <span class="math notranslate nohighlight">\(L^2(\Omega)\)</span> constrained to
take the value 0 at some arbitrary but specified point. For example, one might
choose to require the pressure at the origin to vanish. This is also an
effective way to remove the nullspace, but it is simpler to implement. We will
implement the two-dimensional case (<span class="math notranslate nohighlight">\(d=2\)</span>) and, for simplicity, we will assume <span class="math notranslate nohighlight">\(\mu=1\)</span>.</p>
<p>The colon (<span class="math notranslate nohighlight">\(:\)</span>) indicates an inner product so:</p>
<div class="math notranslate nohighlight" id="equation-9-mixed-problems-1">
<span class="eqno">(9.3)<a class="headerlink" href="#equation-9-mixed-problems-1" title="Link to this equation">¶</a></span>\[\epsilon(u):\epsilon(v) = \sum_{\alpha\beta} \epsilon(u)_{\alpha\beta}\epsilon(v)_{\alpha\beta}\]</div>
<p>In choosing a finite element subspace of <span class="math notranslate nohighlight">\(V \times Q\)</span> we will similarly choose
a simpler to implement, yet still stable, space than was chosen in
<a class="reference internal" href="L6_stokes.html#stokes"><span class="std std-numref">Analysis Section 6</span></a>. The space we will use is the lowest order
Taylor-Hood space <span id="id1">[<a class="reference internal" href="zbibliography.html#id6" title="C. Taylor and P. Hood. A numerical solution of the navier-stokes equations using the finite element technique. Computers & Fluids, 1(1):73-100, 1973. URL: https://www.sciencedirect.com/science/article/pii/0045793073900273, doi:https://doi.org/10.1016/0045-7930(73)90027-3.">TH73</a>]</span>. This has:</p>
<div class="math notranslate nohighlight" id="equation-9-mixed-problems-2">
<span class="eqno">(9.4)<a class="headerlink" href="#equation-9-mixed-problems-2" title="Link to this equation">¶</a></span>\[ \begin{align}\begin{aligned}V^h & = P2(\Omega)^2\\Q^h & = P1(\Omega)\end{aligned}\end{align} \]</div>
<p>i.e. quadratic velocity and linear pressure on each cell. We note that
<span class="math notranslate nohighlight">\(Q^h\in H^1(\Omega)\)</span> but since <span class="math notranslate nohighlight">\(H^1(\Omega) \subset L^2(\Omega)\)</span>, this is not
itself a problem. We will impose further constraints on the spaces in the
form of Dirichlet boundary conditions to ensure that the solutions found are in
fact in <span class="math notranslate nohighlight">\(V \times Q\)</span>.</p>
<p>In addition to the finite element functionality we have already implemented,
there are two further challenges we need to address. First, the implementation
of the vector-valued space <span class="math notranslate nohighlight">\(P2(\Omega)^2`m and second, the implementation of
functions and matrices defined over the mixed space `V^h \times Q^h\)</span>.</p>
<section id="vector-valued-finite-elements">
<h2><span class="section-number">9.1. </span>Vector-valued finite elements<a class="headerlink" href="#vector-valued-finite-elements" title="Link to this heading">¶</a></h2>
<p>Consider the local representation of <span class="math notranslate nohighlight">\(P2(\Omega)^2\)</span> on one element. The scalar
<span class="math notranslate nohighlight">\(P2\)</span> element on one triangle has 6 nodes, one at each vertex and one at each
edge. If we write <span class="math notranslate nohighlight">\(\{\Phi_i\}_{i=0}^{5}\)</span> for the scalar basis, then a basis
<span class="math notranslate nohighlight">\(\{\mathbf{\Phi}_i\}_{i=0}^{11}\)</span> for the vector-valued space is given by:</p>
<div class="math notranslate nohighlight" id="equation-vector-basis">
<span class="eqno">(9.5)<a class="headerlink" href="#equation-vector-basis" title="Link to this equation">¶</a></span>\[\mathbf{\Phi}_i(X) = \Phi_{i\,//\,2}(X)\,\mathbf{e}_{i\,\%\,2}\]</div>
<p>where <span class="math notranslate nohighlight">\(//\)</span> is the integer division operator, <span class="math notranslate nohighlight">\(\%\)</span> is the modulus operator, and
<span class="math notranslate nohighlight">\({\mathbf{e}_0, \mathbf{e}_1}\)</span> is the standard basis for <span class="math notranslate nohighlight">\(\mathbb{R}^2\)</span>. That is to say, we
interleave <span class="math notranslate nohighlight">\(x\)</span> and <span class="math notranslate nohighlight">\(y\)</span> component basis functions.</p>
<figure class="align-default" id="id2">
<img alt="_images/p2vec.svg" src="_images/p2vec.svg" />
<figcaption>
<p><span class="caption-number">Fig. 9.1 </span><span class="caption-text">The local numbering of vector degrees of freedom.</span><a class="headerlink" href="#id2" title="Link to this image">¶</a></p>
</figcaption>
</figure>
<p>We can practically implement vector function spaces by implementing a new class
<code class="xref py py-class docutils literal notranslate"><span class="pre">fe_utils.finite_elements.VectorFiniteElement</span></code>. The constructor
(<a class="reference external" href="https://docs.python.org/3/reference/datamodel.html#object.__init__" title="(in Python v3.13)"><code class="xref py py-meth docutils literal notranslate"><span class="pre">__init__()</span></code></a>) of this new class should take a
<a class="reference internal" href="fe_utils.html#fe_utils.finite_elements.FiniteElement" title="fe_utils.finite_elements.FiniteElement"><code class="xref py py-class docutils literal notranslate"><span class="pre">FiniteElement</span></code></a> and construct the
corresponding vector element. For current purposes we can assume that the
vector element will always have a vector dimension equal to the element
geometric and topological dimension (i.e. 2 if the element is defined on a
triangle). We’ll refer to this dimension as <span class="math notranslate nohighlight">\(d\)</span>.</p>
<section id="implementing-vectorfiniteelement">
<h3><span class="section-number">9.1.1. </span>Implementing <code class="xref py py-class docutils literal notranslate"><span class="pre">VectorFiniteElement</span></code><a class="headerlink" href="#implementing-vectorfiniteelement" title="Link to this heading">¶</a></h3>
<p><code class="xref py py-class docutils literal notranslate"><span class="pre">VectorFiniteElement</span></code> needs to implement as far as possible the same
interface as <a class="reference internal" href="fe_utils.html#fe_utils.finite_elements.FiniteElement" title="fe_utils.finite_elements.FiniteElement"><code class="xref py py-class docutils literal notranslate"><span class="pre">FiniteElement</span></code></a>. Let’s think
about how to do that.</p>
<dl>
<dt><code class="xref py py-data docutils literal notranslate"><span class="pre">cell</span></code>, <code class="xref py py-data docutils literal notranslate"><span class="pre">degree</span></code></dt><dd><p>Same as for the input <a class="reference internal" href="fe_utils.html#fe_utils.finite_elements.FiniteElement" title="fe_utils.finite_elements.FiniteElement"><code class="xref py py-class docutils literal notranslate"><span class="pre">FiniteElement</span></code></a>.</p>
</dd>
<dt><code class="xref py py-data docutils literal notranslate"><span class="pre">entity_nodes</span></code></dt><dd><p>There will be twice as many nodes, and the node ordering is such that each
node is replaced by a <span class="math notranslate nohighlight">\(d\)</span>-tuple. For example the scalar triangle P1
entity-node list is:</p>
<div class="highlight-python3 notranslate"><div class="highlight"><pre><span></span><span class="p">{</span>
<span class="mi">0</span> <span class="p">:</span> <span class="p">{</span><span class="mi">0</span> <span class="p">:</span> <span class="p">[</span><span class="mi">0</span><span class="p">],</span> <span class="mi">1</span> <span class="p">:</span> <span class="p">[</span><span class="mi">1</span><span class="p">],</span> <span class="mi">2</span> <span class="p">:</span> <span class="p">[</span><span class="mi">2</span><span class="p">]},</span>
<span class="mi">1</span> <span class="p">:</span> <span class="p">{</span><span class="mi">0</span> <span class="p">:</span> <span class="p">[],</span> <span class="mi">1</span> <span class="p">:</span> <span class="p">[],</span> <span class="mi">2</span> <span class="p">:</span> <span class="p">[]},</span>
<span class="mi">2</span> <span class="p">:</span> <span class="p">{</span><span class="mi">0</span> <span class="p">:</span> <span class="p">[]}</span>
<span class="p">}</span>
</pre></div>
</div>
<p>The vector version is achieved by looping over the scalar version and
returning a mapping with the pair <span class="math notranslate nohighlight">\(2n, 2(n+1)\)</span> in place of node <span class="math notranslate nohighlight">\(n\)</span>:</p>
<div class="highlight-python3 notranslate"><div class="highlight"><pre><span></span><span class="p">{</span>
<span class="mi">0</span> <span class="p">:</span> <span class="p">{</span><span class="mi">0</span> <span class="p">:</span> <span class="p">[</span><span class="mi">0</span><span class="p">,</span> <span class="mi">1</span><span class="p">],</span> <span class="mi">1</span> <span class="p">:</span> <span class="p">[</span><span class="mi">2</span><span class="p">,</span> <span class="mi">3</span><span class="p">],</span> <span class="mi">2</span> <span class="p">:</span> <span class="p">[</span><span class="mi">4</span><span class="p">,</span> <span class="mi">5</span><span class="p">]},</span>
<span class="mi">1</span> <span class="p">:</span> <span class="p">{</span><span class="mi">0</span> <span class="p">:</span> <span class="p">[],</span> <span class="mi">1</span> <span class="p">:</span> <span class="p">[],</span> <span class="mi">2</span> <span class="p">:</span> <span class="p">[]},</span>
<span class="mi">2</span> <span class="p">:</span> <span class="p">{</span><span class="mi">0</span> <span class="p">:</span> <span class="p">[]}</span>
<span class="p">}</span>
</pre></div>
</div>
</dd>
<dt><code class="xref py py-data docutils literal notranslate"><span class="pre">nodes_per_entity</span></code>:</dt><dd><p>Each entry will be <span class="math notranslate nohighlight">\(d\)</span> times that on the input
<a class="reference internal" href="fe_utils.html#fe_utils.finite_elements.FiniteElement" title="fe_utils.finite_elements.FiniteElement"><code class="xref py py-class docutils literal notranslate"><span class="pre">FiniteElement</span></code></a>.</p>
</dd>
</dl>
</section>
<section id="tabulation">
<h3><span class="section-number">9.1.2. </span>Tabulation<a class="headerlink" href="#tabulation" title="Link to this heading">¶</a></h3>
<p>In order to tabulate the element, we can use <a class="reference internal" href="#equation-vector-basis">(9.5)</a>. We first
call the tabulate method from the input
<a class="reference internal" href="fe_utils.html#fe_utils.finite_elements.FiniteElement" title="fe_utils.finite_elements.FiniteElement"><code class="xref py py-class docutils literal notranslate"><span class="pre">FiniteElement</span></code></a>, and we use this and
<a class="reference internal" href="#equation-vector-basis">(9.5)</a> to produce the array to return. Notice that the array
will both have a basis functions dimension which is <span class="math notranslate nohighlight">\(d\)</span> times longer than the
input element, and will also have an extra dimension to account for the
multiplication by <span class="math notranslate nohighlight">\(\mathbf{e}_{i\,\%\,d}\)</span>. This means that the tabulation array
with <code class="xref py py-data docutils literal notranslate"><span class="pre">grad=False</span></code> will now be rank 3, and that with <code class="xref py py-data docutils literal notranslate"><span class="pre">grad=True</span></code>
will be rank 4. Make sure you keep track of which rank is which!
The <code class="xref py py-class docutils literal notranslate"><span class="pre">VectorFiniteElement</span></code> will need to keep a reference to the
input <a class="reference internal" href="fe_utils.html#fe_utils.finite_elements.FiniteElement" title="fe_utils.finite_elements.FiniteElement"><code class="xref py py-class docutils literal notranslate"><span class="pre">FiniteElement</span></code></a> in order to facilitate
tabulation.</p>
</section>
<section id="nodes">
<h3><span class="section-number">9.1.3. </span>Nodes<a class="headerlink" href="#nodes" title="Link to this heading">¶</a></h3>
<p>Even though we didn’t need the nodes of the <code class="xref py py-class docutils literal notranslate"><span class="pre">VectorFiniteElement</span></code> to
construct its basis, we will need them to implement interpolation. In
<a class="reference internal" href="2_finite_elements.html#nodalbasis"><span class="std std-numref">Definition 2.44</span></a> we learned that
the nodes of a finite element are related to the corresponding nodal basis by:</p>
<div class="math notranslate nohighlight" id="equation-9-mixed-problems-3">
<span class="eqno">(9.6)<a class="headerlink" href="#equation-9-mixed-problems-3" title="Link to this equation">¶</a></span>\[\mathbf{\Phi}^*_i(\mathbf{\Phi}_j) = \delta_{ij}\]</div>
<p>From <a class="reference internal" href="#equation-vector-basis">(9.5)</a> and assuming, as we have throughout the course,
that the scalar finite element has point evaluation nodes given by:</p>
<div class="math notranslate nohighlight" id="equation-9-mixed-problems-4">
<span class="eqno">(9.7)<a class="headerlink" href="#equation-9-mixed-problems-4" title="Link to this equation">¶</a></span>\[\Phi_i(v) = v(X_i),\]</div>
<p>it follows that:</p>
<div class="math notranslate nohighlight" id="equation-vectornodes">
<span class="eqno">(9.8)<a class="headerlink" href="#equation-vectornodes" title="Link to this equation">¶</a></span>\[ \begin{align}\begin{aligned}\mathbf{\Phi}^*_i(v) & = \Phi^*_{i\,//\,d}(\mathbf{e}_{i\,\%\,d}\cdot v)\\& = \mathbf{e}_{i\,\%\,d}\cdot v(X_{i\,//\,d})\end{aligned}\end{align} \]</div>
<div class="admonition hint">
<p class="admonition-title">Hint</p>
<p>To see that this is the correct nodal basis, choose
<span class="math notranslate nohighlight">\(v=\mathbf{\Phi}_j\)</span> in <a class="reference internal" href="#equation-vectornodes">(9.8)</a> and substitute
<a class="reference internal" href="#equation-vector-basis">(9.5)</a> for <span class="math notranslate nohighlight">\(\mathbf{\Phi}_j\)</span>.</p>
</div>
<p>This means that the correct <code class="xref py py-data docutils literal notranslate"><span class="pre">VectorFiniteElement.nodes</span></code> attribute is
the list of nodal points from the input
<a class="reference internal" href="fe_utils.html#fe_utils.finite_elements.FiniteElement" title="fe_utils.finite_elements.FiniteElement"><code class="xref py py-class docutils literal notranslate"><span class="pre">FiniteElement</span></code></a> but with each point repeated
<span class="math notranslate nohighlight">\(d\)</span> times. It will also be necessary to add another attribute, perhaps
<code class="xref py py-data docutils literal notranslate"><span class="pre">node_weights</span></code> which is a rank 2 array whose <span class="math notranslate nohighlight">\(i\)</span>-th row is the correct
canonical basis vector to contract with the function value at the <span class="math notranslate nohighlight">\(i\)</span>-th node (<span class="math notranslate nohighlight">\(\mathbf{e}_{i\,\%\,d}\)</span>).</p>
</section>
</section>
<section id="vector-valued-function-spaces">
<h2><span class="section-number">9.2. </span>Vector-valued function spaces<a class="headerlink" href="#vector-valued-function-spaces" title="Link to this heading">¶</a></h2>
<p>Assuming we correctly implement <code class="xref py py-class docutils literal notranslate"><span class="pre">VectorFiniteElement</span></code>,
<a class="reference internal" href="fe_utils.html#fe_utils.function_spaces.FunctionSpace" title="fe_utils.function_spaces.FunctionSpace"><code class="xref py py-class docutils literal notranslate"><span class="pre">FunctionSpace</span></code></a> should work out of the box.
In particular, the global numbering algorithm only depends on having a correct
local numbering so this should work unaltered.</p>
</section>
<section id="functions-in-vector-valued-spaces">
<h2><span class="section-number">9.3. </span>Functions in vector-valued spaces<a class="headerlink" href="#functions-in-vector-valued-spaces" title="Link to this heading">¶</a></h2>
<p>The general form of a function in a vector-valued function space is:</p>
<div class="math notranslate nohighlight" id="equation-9-mixed-problems-5">
<span class="eqno">(9.9)<a class="headerlink" href="#equation-9-mixed-problems-5" title="Link to this equation">¶</a></span>\[f = f_i \mathbf{\Phi}_i(X).\]</div>
<p>That is to say, the basis functions are vector valued and their coefficients
are still scalar. This means that if the <code class="xref py py-class docutils literal notranslate"><span class="pre">VectorFiniteElement</span></code> had a
correct entity-node list then the core functionality of the existing
<a class="reference internal" href="fe_utils.html#fe_utils.function_spaces.Function" title="fe_utils.function_spaces.Function"><code class="xref py py-class docutils literal notranslate"><span class="pre">Function</span></code></a> will automatically be correct. In
particular, the array of values will have the correct extent. However,
interpolation and plotting of vector valued fields will require some
adjustment.</p>
<section id="interpolating-into-vector-valued-spaces">
<h3><span class="section-number">9.3.1. </span>Interpolating into vector-valued spaces<a class="headerlink" href="#interpolating-into-vector-valued-spaces" title="Link to this heading">¶</a></h3>
<p>Since the form of the nodes of a <code class="xref py py-class docutils literal notranslate"><span class="pre">VectorFiniteElement</span></code> is different from
that of a scalar element, there will be some changes required in the
<a class="reference internal" href="fe_utils.html#fe_utils.function_spaces.Function.interpolate" title="fe_utils.function_spaces.Function.interpolate"><code class="xref py py-meth docutils literal notranslate"><span class="pre">interpolate()</span></code></a> method. Specifically:</p>
<div class="highlight-python3 notranslate"><div class="highlight"><pre><span></span><span class="bp">self</span><span class="o">.</span><span class="n">values</span><span class="p">[</span><span class="n">fs</span><span class="o">.</span><span class="n">cell_nodes</span><span class="p">[</span><span class="n">c</span><span class="p">,</span> <span class="p">:]]</span> <span class="o">=</span> <span class="p">[</span><span class="n">fn</span><span class="p">(</span><span class="n">x</span><span class="p">)</span> <span class="k">for</span> <span class="n">x</span> <span class="ow">in</span> <span class="n">node_coords</span><span class="p">]</span>
</pre></div>
</div>
<p>This line will need to take into account the dot product with the
canonical basis from <a class="reference internal" href="#equation-vectornodes">(9.8)</a>, which you have implemented as
<code class="xref py py-data docutils literal notranslate"><span class="pre">VectorFiniteElement.node_weights</span></code>. This change will need to be made
conditional on the class of finite element passed in, so that the code doesn’t
break in the scalar element case.</p>
</section>
<section id="plotting-functions-in-vector-valued-spaces">
<h3><span class="section-number">9.3.2. </span>Plotting functions in vector-valued spaces<a class="headerlink" href="#plotting-functions-in-vector-valued-spaces" title="Link to this heading">¶</a></h3>
<p>The coloured surface plots that we’ve used thus far for two-dimensional scalar
functions don’t extend easily to vector quantities. Instead, a frequently used
visualisation technique is the quiver plot. This draws a set of
arrows representing the function value at a set of points. For our purposes,
the nodes of the function space in question are a good choice of evaluation
points. <a class="reference internal" href="#qplot"><span class="std std-numref">Listing 9.1</span></a> provides the code you will need. Notice that at line 3
we interpolated the function <span class="math notranslate nohighlight">\(f(x)=x\)</span> into the function space in order to
obtain a list of the global coordinates of the node locations. At lines 6 and 7
we use what we know about the node ordering to recover vector values from the
list of basis function coefficients.</p>
<div class="literal-block-wrapper docutils container" id="id3">
<span id="qplot"></span><div class="code-block-caption"><span class="caption-number">Listing 9.1 </span><span class="caption-text">Code implementing quiver plots to visualise functions in vector
function spaces. This code should be added to
<a class="reference internal" href="fe_utils.html#fe_utils.function_spaces.Function.plot" title="fe_utils.function_spaces.Function.plot"><code class="xref py py-meth docutils literal notranslate"><span class="pre">plot()</span></code></a> immediately after the
definition of <code class="xref py py-data docutils literal notranslate"><span class="pre">fs</span></code>.</span><a class="headerlink" href="#id3" title="Link to this code">¶</a></div>
<div class="highlight-python3 notranslate"><div class="highlight"><pre><span></span><span class="linenos"> 1</span><span class="k">if</span> <span class="nb">isinstance</span><span class="p">(</span><span class="n">fs</span><span class="o">.</span><span class="n">element</span><span class="p">,</span> <span class="n">VectorFiniteElement</span><span class="p">):</span>
<span class="linenos"> 2</span> <span class="n">coords</span> <span class="o">=</span> <span class="n">Function</span><span class="p">(</span><span class="n">fs</span><span class="p">)</span>
<span class="linenos"> 3</span> <span class="n">coords</span><span class="o">.</span><span class="n">interpolate</span><span class="p">(</span><span class="k">lambda</span> <span class="n">x</span><span class="p">:</span> <span class="n">x</span><span class="p">)</span>
<span class="linenos"> 4</span> <span class="n">fig</span> <span class="o">=</span> <span class="n">plt</span><span class="o">.</span><span class="n">figure</span><span class="p">()</span>
<span class="linenos"> 5</span> <span class="n">ax</span> <span class="o">=</span> <span class="n">fig</span><span class="o">.</span><span class="n">add_subplot</span><span class="p">()</span>
<span class="linenos"> 6</span> <span class="n">x</span> <span class="o">=</span> <span class="n">coords</span><span class="o">.</span><span class="n">values</span><span class="o">.</span><span class="n">reshape</span><span class="p">(</span><span class="o">-</span><span class="mi">1</span><span class="p">,</span> <span class="mi">2</span><span class="p">)</span>
<span class="linenos"> 7</span> <span class="n">v</span> <span class="o">=</span> <span class="bp">self</span><span class="o">.</span><span class="n">values</span><span class="o">.</span><span class="n">reshape</span><span class="p">(</span><span class="o">-</span><span class="mi">1</span><span class="p">,</span> <span class="mi">2</span><span class="p">)</span>
<span class="linenos"> 8</span> <span class="n">plt</span><span class="o">.</span><span class="n">quiver</span><span class="p">(</span><span class="n">x</span><span class="p">[:,</span> <span class="mi">0</span><span class="p">],</span> <span class="n">x</span><span class="p">[:,</span> <span class="mi">1</span><span class="p">],</span> <span class="n">v</span><span class="p">[:,</span> <span class="mi">0</span><span class="p">],</span> <span class="n">v</span><span class="p">[:,</span> <span class="mi">1</span><span class="p">])</span>
<span class="linenos"> 9</span> <span class="n">plt</span><span class="o">.</span><span class="n">show</span><span class="p">()</span>
<span class="linenos">10</span> <span class="k">return</span>
</pre></div>
</div>
</div>
<p>Once this code has been inserted, then running the code in
<a class="reference internal" href="#quiverplotcode"><span class="std std-numref">Listing 9.2</span></a> will result in a plot rather like
<a class="reference internal" href="#quiverplot"><span class="std std-numref">Fig. 9.2</span></a>.</p>
<div class="literal-block-wrapper docutils container" id="id4">
<span id="quiverplotcode"></span><div class="code-block-caption"><span class="caption-number">Listing 9.2 </span><span class="caption-text">Creation of a vector function space, interpolation of a given
function into it, and subsequent plot creation.</span><a class="headerlink" href="#id4" title="Link to this code">¶</a></div>
<div class="highlight-python3 notranslate"><div class="highlight"><pre><span></span><span class="linenos"> 1</span><span class="kn">from</span><span class="w"> </span><span class="nn">fe_utils</span><span class="w"> </span><span class="kn">import</span> <span class="o">*</span>
<span class="linenos"> 2</span><span class="kn">from</span><span class="w"> </span><span class="nn">math</span><span class="w"> </span><span class="kn">import</span> <span class="n">cos</span><span class="p">,</span> <span class="n">sin</span><span class="p">,</span> <span class="n">pi</span>
<span class="linenos"> 3</span>
<span class="linenos"> 4</span><span class="n">se</span> <span class="o">=</span> <span class="n">LagrangeElement</span><span class="p">(</span><span class="n">ReferenceTriangle</span><span class="p">,</span> <span class="mi">2</span><span class="p">)</span>
<span class="linenos"> 5</span><span class="n">ve</span> <span class="o">=</span> <span class="n">VectorFiniteElement</span><span class="p">(</span><span class="n">se</span><span class="p">)</span>
<span class="linenos"> 6</span><span class="n">m</span> <span class="o">=</span> <span class="n">UnitSquareMesh</span><span class="p">(</span><span class="mi">10</span><span class="p">,</span><span class="mi">10</span><span class="p">)</span>
<span class="linenos"> 7</span><span class="n">fs</span> <span class="o">=</span> <span class="n">FunctionSpace</span><span class="p">(</span><span class="n">m</span><span class="p">,</span> <span class="n">ve</span><span class="p">)</span>
<span class="linenos"> 8</span><span class="n">f</span> <span class="o">=</span> <span class="n">Function</span><span class="p">(</span><span class="n">fs</span><span class="p">)</span>
<span class="linenos"> 9</span><span class="n">f</span><span class="o">.</span><span class="n">interpolate</span><span class="p">(</span><span class="k">lambda</span> <span class="n">x</span><span class="p">:</span> <span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">pi</span><span class="o">*</span><span class="p">(</span><span class="mi">1</span> <span class="o">-</span> <span class="n">cos</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">pi</span><span class="o">*</span><span class="n">x</span><span class="p">[</span><span class="mi">0</span><span class="p">]))</span><span class="o">*</span><span class="n">sin</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">pi</span><span class="o">*</span><span class="n">x</span><span class="p">[</span><span class="mi">1</span><span class="p">]),</span>
<span class="linenos">10</span> <span class="o">-</span><span class="mi">2</span><span class="o">*</span><span class="n">pi</span><span class="o">*</span><span class="p">(</span><span class="mi">1</span> <span class="o">-</span> <span class="n">cos</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">pi</span><span class="o">*</span><span class="n">x</span><span class="p">[</span><span class="mi">1</span><span class="p">]))</span><span class="o">*</span><span class="n">sin</span><span class="p">(</span><span class="mi">2</span><span class="o">*</span><span class="n">pi</span><span class="o">*</span><span class="n">x</span><span class="p">[</span><span class="mi">0</span><span class="p">])))</span>
<span class="linenos">11</span><span class="n">f</span><span class="o">.</span><span class="n">plot</span><span class="p">()</span>
</pre></div>
</div>
</div>
<figure class="align-default" id="id5">
<span id="quiverplot"></span><a class="reference internal image-reference" href="_images/quiver.png"><img alt="_images/quiver.png" src="_images/quiver.png" style="width: 70%;" />
</a>
<figcaption>
<p><span class="caption-number">Fig. 9.2 </span><span class="caption-text">The quiver plot resulting from <a class="reference internal" href="#quiverplotcode"><span class="std std-numref">Listing 9.2</span></a>.</span><a class="headerlink" href="#id5" title="Link to this image">¶</a></p>
</figcaption>
</figure>
</section>
<section id="solving-vector-valued-systems">
<h3><span class="section-number">9.3.3. </span>Solving vector-valued systems<a class="headerlink" href="#solving-vector-valued-systems" title="Link to this heading">¶</a></h3>
<p>Solving a finite element problem in a vector-valued space is essentially
similar to the scalar problems you have already solved. It does, nonetheless,
provide a useful check on the correctness of your code before adding in the
additional complications of mixed systems. As a very simple example, consider
computing vector-valued field which is the gradient of a known function. For
some suitable finite element space <span class="math notranslate nohighlight">\(V\subset H^1(\Omega)^2\)</span> and
<span class="math notranslate nohighlight">\(f:\Omega\rightarrow \mathbb{R}\)</span>, find <span class="math notranslate nohighlight">\(u\in V\)</span> such that:</p>
<div class="math notranslate nohighlight" id="equation-9-mixed-problems-6">
<span class="eqno">(9.10)<a class="headerlink" href="#equation-9-mixed-problems-6" title="Link to this equation">¶</a></span>\[\int_\Omega u\cdot v\,\mathrm{d}x = \int_\Omega \nabla f\cdot v\,\mathrm{d}x\quad \forall v\in V.\]</div>
<p>If <span class="math notranslate nohighlight">\(f\)</span> is chosen such that <span class="math notranslate nohighlight">\(\nabla f\in V\)</span> then this projection is exact up to
roundoff, and the following calculation should result in a good approximation
to zero:</p>
<div class="math notranslate nohighlight" id="equation-9-mixed-problems-7">
<span class="eqno">(9.11)<a class="headerlink" href="#equation-9-mixed-problems-7" title="Link to this equation">¶</a></span>\[e = \int_\Omega (u -\nabla f)\cdot(u -\nabla f)\,\mathrm{d}x\]</div>
<div class="admonition note">
<p class="admonition-title">Note</p>
<p>The computations in this subsection are not required to complete the
mastery exercise. They are, nonetheless, strongly recommended as a
mechanism for checking your implementation thus far.</p>
</div>
</section>
</section>
<section id="mixed-function-spaces">
<h2><span class="section-number">9.4. </span>Mixed function spaces<a class="headerlink" href="#mixed-function-spaces" title="Link to this heading">¶</a></h2>
<p>The Stokes equations are defined over the mixed function space <span class="math notranslate nohighlight">\(W = V \times Q\)</span>.
Here “mixed” simply means that there are two solution variables, and therefore
two solution spaces. Functions in <span class="math notranslate nohighlight">\(W\)</span> are pairs <span class="math notranslate nohighlight">\((u, p)\)</span> where <span class="math notranslate nohighlight">\(u\in V\)</span> and
<span class="math notranslate nohighlight">\(p\in Q\)</span>. If <span class="math notranslate nohighlight">\(\{\phi_i\}_{i=0}^{m-1}\)</span> is a basis for <span class="math notranslate nohighlight">\(V\)</span> and
<span class="math notranslate nohighlight">\(\{\psi_j\}_{j=0}^{n-1}\)</span> then a basis for <span class="math notranslate nohighlight">\(W\)</span> is given by:</p>
<div class="math notranslate nohighlight" id="equation-9-mixed-problems-8">
<span class="eqno">(9.12)<a class="headerlink" href="#equation-9-mixed-problems-8" title="Link to this equation">¶</a></span>\[\{\omega_i\}_{i=0}^{m+n-1}=\{(\phi_i, 0)\}_{i=0}^{m-1} \cup \{(0,
\psi_{j-n})\}_{j=m}^{m+n-1}.\]</div>
<p>This in turn enables us to write <span class="math notranslate nohighlight">\(w\in W\)</span> in the form <span class="math notranslate nohighlight">\(w=w_i\omega_i\)</span> as we
would expect for a function in a finite element space. The Cartesian product
structure of the mixed space <span class="math notranslate nohighlight">\(W\)</span> means that the first <span class="math notranslate nohighlight">\(m\)</span> coefficients are
simply the coefficients of the <span class="math notranslate nohighlight">\(V\)</span> basis functions, and the latter <span class="math notranslate nohighlight">\(n\)</span>
coefficients correspond to the <span class="math notranslate nohighlight">\(Q\)</span> basis functions. This means that our full
mixed finite element system is simply a linear system of block matrices and
block vectors. If we disregard boundary conditions, including the pressure
constraint, this system has the following form:</p>
<div class="math notranslate nohighlight" id="equation-9-mixed-problems-9">
<span class="eqno">(9.13)<a class="headerlink" href="#equation-9-mixed-problems-9" title="Link to this equation">¶</a></span>\[\begin{split}\begin{bmatrix}
A & B^\mathrm{T} \\
B & 0
\end{bmatrix}
\begin{bmatrix}
U \\
P
\end{bmatrix}
=
\begin{bmatrix}
F \\
0
\end{bmatrix}\end{split}\]</div>
<p>where:</p>
<div class="math notranslate nohighlight" id="equation-blocks">
<span class="eqno">(9.14)<a class="headerlink" href="#equation-blocks" title="Link to this equation">¶</a></span>\[ \begin{align}\begin{aligned}A_{ij} = a(\phi_j, \phi_i),\\B_{ij} = b(\phi_j, \psi_i),\\F_i = \int_\Omega f\cdot v\, d\, x,\\U_i = u_i = w_i,\\P_i = p_i = w_{i+m}.\end{aligned}\end{align} \]</div>
<p>This means that the assembly of the mixed problem comes down to the assembly of
several finite operators of the form that we have already encountered. These
then need to be assembled into the full block matrix and right hand side
vector, before the system is solved and the resulting solution vector pulled
appart and interpreted as the coefficients of <span class="math notranslate nohighlight">\(u\)</span> and <span class="math notranslate nohighlight">\(p\)</span>. Observe in
<a class="reference internal" href="#equation-blocks">(9.14)</a> that the order of the indices <span class="math notranslate nohighlight">\(i\)</span> and <span class="math notranslate nohighlight">\(j\)</span> is reversed on the right
hand side of the equations. This reflects the differing conventions for matrix
indices and bilinear form arguments, and is a source of unending confusion in
this field.</p>
<section id="assembling-block-systems">
<h3><span class="section-number">9.4.1. </span>Assembling block systems<a class="headerlink" href="#assembling-block-systems" title="Link to this heading">¶</a></h3>
<p>The procedure for assembling the individual blocks of the block matrix and the
block vectors is the one you are familiar with, but we will need to do
something new to assemble the block structures. What is required differs
slightly between the matrix and the vectors.</p>
<p>In the case of the vectors, then it is sufficient to know that a slice into a
<a class="reference external" href="https://numpy.org/doc/stable/reference/generated/numpy.ndarray.html#numpy.ndarray" title="(in NumPy v2.2)"><code class="xref py py-class docutils literal notranslate"><span class="pre">numpy.ndarray</span></code></a> returns a view on the same memory as the full vector.
This is most easily understood through an example:</p>
<div class="highlight-ipython3 notranslate"><div class="highlight"><pre><span></span><span class="n">In</span> <span class="p">[</span><span class="mi">1</span><span class="p">]:</span> <span class="kn">import</span><span class="w"> </span><span class="nn">numpy</span><span class="w"> </span><span class="k">as</span><span class="w"> </span><span class="nn">np</span>
<span class="n">In</span> <span class="p">[</span><span class="mi">2</span><span class="p">]:</span> <span class="n">a</span> <span class="o">=</span> <span class="n">np</span><span class="o">.</span><span class="n">zeros</span><span class="p">(</span><span class="mi">10</span><span class="p">)</span>
<span class="n">In</span> <span class="p">[</span><span class="mi">3</span><span class="p">]:</span> <span class="n">b</span> <span class="o">=</span> <span class="n">a</span><span class="p">[:</span><span class="mi">5</span><span class="p">]</span>
<span class="n">In</span> <span class="p">[</span><span class="mi">4</span><span class="p">]:</span> <span class="n">b</span><span class="p">[</span><span class="mi">2</span><span class="p">]</span> <span class="o">=</span> <span class="mi">1</span>
<span class="n">In</span> <span class="p">[</span><span class="mi">5</span><span class="p">]:</span> <span class="n">a</span>
<span class="n">Out</span><span class="p">[</span><span class="mi">5</span><span class="p">]:</span> <span class="n">array</span><span class="p">([</span><span class="mf">0.</span><span class="p">,</span> <span class="mf">0.</span><span class="p">,</span> <span class="mf">1.</span><span class="p">,</span> <span class="mf">0.</span><span class="p">,</span> <span class="mf">0.</span><span class="p">,</span> <span class="mf">0.</span><span class="p">,</span> <span class="mf">0.</span><span class="p">,</span> <span class="mf">0.</span><span class="p">,</span> <span class="mf">0.</span><span class="p">,</span> <span class="mf">0.</span><span class="p">])</span>
</pre></div>
</div>
<p>This means that one can first create a full vector of length <span class="math notranslate nohighlight">\(n+m\)</span> and then
slice it to create subvectors that can be used for assembly.</p>
<p>Conversely, <a class="reference external" href="http://scipy.github.io/devdocs/reference/sparse.html#module-scipy.sparse" title="(in SciPy v1.16.0.dev)"><code class="xref py py-mod docutils literal notranslate"><span class="pre">scipy.sparse</span></code></a> provides the <a class="reference external" href="http://scipy.github.io/devdocs/reference/generated/scipy.sparse.bmat.html#scipy.sparse.bmat" title="(in SciPy v1.16.0.dev)"><code class="xref py py-func docutils literal notranslate"><span class="pre">bmat()</span></code></a>
function which will stitch together a larger sparse matrix from blocks. In
order to have the full indexing options you are likely to want for imposing the
boundary conditions, you will probably want to specify that the resulting
matrix is in <code class="xref py py-data docutils literal notranslate"><span class="pre">"lil"</span></code> format.</p>
</section>
<section id="boundary-conditions">
<h3><span class="section-number">9.4.2. </span>Boundary conditions<a class="headerlink" href="#boundary-conditions" title="Link to this heading">¶</a></h3>
<p>The imposition of the constraint in <span class="math notranslate nohighlight">\((\mathring{H}^1(\Omega))^2\)</span> that solutions
vanish on the boundary is a Dirichlet condition of the type that you have
encountered before. Observe that the condition changes the test space, which
affects whole rows of the block system, so you will want to impose the boundary
condition <em>after</em> assembling the block matrix. You will also need to ensure
that the constraint is applied to both the <span class="math notranslate nohighlight">\(x\)</span> and <span class="math notranslate nohighlight">\(y\)</span> components of the space.</p>
<p>The imposition of the constraint in <span class="math notranslate nohighlight">\(\mathring{L}^2(\Omega)\)</span> that the solution
is zero at some prescribed point can be achieved by selecting an arbitrary
basis function and applying a zero Dirichlet condition for that degree of
freedom. In this regard we can observe that there is nothing about the
implementation of Dirichlet conditions that constrains them to lie on the
boundary. Rather, they should be understood as specifying a subspace on which
the solution is prescribed rather than solved for. In this particular case,
that subspace is one-dimensional.</p>
</section>
<section id="solving-the-matrix-system">
<h3><span class="section-number">9.4.3. </span>Solving the matrix system<a class="headerlink" href="#solving-the-matrix-system" title="Link to this heading">¶</a></h3>
<p>The block matrix system that you eventually produce will be larger than many of
those we have previously encountered, and will have non-zero entries further
from the diagonal. This can cause the matrix solver to become expensive in both
time and memory. Fortunately, <a class="reference external" href="http://scipy.github.io/devdocs/reference/sparse.linalg.html#module-scipy.sparse.linalg" title="(in SciPy v1.16.0.dev)"><code class="xref py py-mod docutils literal notranslate"><span class="pre">scipy.sparse.linalg</span></code></a> now incorporates an
interface to <a class="reference external" href="https://portal.nersc.gov/project/sparse/superlu/">SuperLU</a>,
which is a high-performance direct sparse solver. The recommended solution
strategy is therefore:</p>
<ol class="arabic simple">
<li><p>Convert your block matrix to <a class="reference external" href="http://scipy.github.io/devdocs/reference/generated/scipy.sparse.csc_matrix.html#scipy.sparse.csc_matrix" title="(in SciPy v1.16.0.dev)"><code class="xref py py-class docutils literal notranslate"><span class="pre">scipy.sparse.csc_matrix</span></code></a>, which is the
format that SuperLU requires.</p></li>
<li><p>Factorise the matrix using <a class="reference external" href="http://scipy.github.io/devdocs/reference/generated/scipy.sparse.linalg.splu.html#scipy.sparse.linalg.splu" title="(in SciPy v1.16.0.dev)"><code class="xref py py-func docutils literal notranslate"><span class="pre">scipy.sparse.linalg.splu()</span></code></a>.</p></li>
<li><p>Use the resulting <a class="reference external" href="http://scipy.github.io/devdocs/reference/generated/scipy.sparse.linalg.SuperLU.html#scipy.sparse.linalg.SuperLU" title="(in SciPy v1.16.0.dev)"><code class="xref py py-class docutils literal notranslate"><span class="pre">SuperLU</span></code></a> object to finally solve
the system.</p></li>
</ol>
</section>
<section id="computing-the-error">
<h3><span class="section-number">9.4.4. </span>Computing the error<a class="headerlink" href="#computing-the-error" title="Link to this heading">¶</a></h3>
<p>We will wish to compute the convergence of our solution in the <span class="math notranslate nohighlight">\(L^2\)</span> norm. For
<span class="math notranslate nohighlight">\(w\in W\)</span>, this is given by:</p>
<div class="math notranslate nohighlight" id="equation-9-mixed-problems-10">
<span class="eqno">(9.15)<a class="headerlink" href="#equation-9-mixed-problems-10" title="Link to this equation">¶</a></span>\[\|w\|_{L^2} = \sqrt{\int_\Omega w\cdot w\,\mathrm{dx}}\]</div>
<p>When we expand this in terms of the basis of <span class="math notranslate nohighlight">\(W\)</span>, it will be useful to note
that basis functions from the different component spaces are orthogonal. That
is to say:</p>
<div class="math notranslate nohighlight" id="equation-9-mixed-problems-11">
<span class="eqno">(9.16)<a class="headerlink" href="#equation-9-mixed-problems-11" title="Link to this equation">¶</a></span>\[(\phi, 0) \cdot (0, \psi) = 0 \quad \forall \phi\in V,\, \psi \in Q.\]</div>
<p>The direct result of this is that if <span class="math notranslate nohighlight">\(w = (u, p)\)</span> then:</p>
<div class="math notranslate nohighlight" id="equation-9-mixed-problems-12">
<span class="eqno">(9.17)<a class="headerlink" href="#equation-9-mixed-problems-12" title="Link to this equation">¶</a></span>\[\|w\|_{L^2}^2 = \|u\|_{L^2}^2 + \|p\|_{L^2}^2.\]</div>
</section>
</section>
<section id="manufacturing-a-solution-to-the-stokes-equations">
<h2><span class="section-number">9.5. </span>Manufacturing a solution to the Stokes equations<a class="headerlink" href="#manufacturing-a-solution-to-the-stokes-equations" title="Link to this heading">¶</a></h2>
<p>As previously, we will wish to check our code using the method of manufactured
solutions. The Stokes equations represent a form of incompressible fluid
mechanics, so it is usually preferable to select a target solution for which
<span class="math notranslate nohighlight">\(\nabla\cdot u = 0\)</span>. The straightforward way to do this is to choose a scalar
field <span class="math notranslate nohighlight">\(\gamma: \Omega\rightarrow \mathbb{R}\)</span> to use as a streamfunction. We can
then define <span class="math notranslate nohighlight">\(u = \nabla^{\perp}\gamma\)</span> and rely on the vector calculus identity
<span class="math notranslate nohighlight">\(\nabla\cdot\nabla^{\perp} \gamma = 0\)</span> to guarantee that the velocity field is
divergence-free. We also need to ensure that <span class="math notranslate nohighlight">\(u\)</span> satisfies the boundary
conditions, which amounts to choosing <span class="math notranslate nohighlight">\(\gamma\)</span> such that its gradient vanishes
on the domain boundary. The following function is a suitable choice on a unit
square domain:</p>
<div class="math notranslate nohighlight" id="equation-stream">
<span class="eqno">(9.18)<a class="headerlink" href="#equation-stream" title="Link to this equation">¶</a></span>\[\gamma(x,y) = \big(1-\cos(2\pi x)\big)\big(1-\cos(2\pi y)\big)\]</div>
</section>
<section id="implementing-the-stokes-problem">
<h2><span class="section-number">9.6. </span>Implementing the Stokes problem<a class="headerlink" href="#implementing-the-stokes-problem" title="Link to this heading">¶</a></h2>
<div class="proof proof-type-exercise" id="id6">
<div class="proof-title">
<span class="proof-type">Exercise 9.1</span>
</div><div class="proof-content">
<p>The goal of this exercise is to implement a solver for the Stokes
equations, on a unit square. Implement
<a class="reference internal" href="fe_utils.solvers.html#fe_utils.solvers.mastery.solve_mastery" title="fe_utils.solvers.mastery.solve_mastery"><code class="xref py py-func docutils literal notranslate"><span class="pre">solve_mastery()</span></code></a> so that it solves
<a class="reference internal" href="#equation-stokes-ch9">(9.1)</a> using the forcing function derived from <a class="reference internal" href="#equation-stream">(9.18)</a>.</p>
<p>Your full solution should:</p>
<ol class="arabic simple">
<li><p>Implement <code class="xref py py-class docutils literal notranslate"><span class="pre">VectorFiniteElement</span></code>.</p></li>
<li><p>Make the consequential changes to
<a class="reference internal" href="fe_utils.html#fe_utils.function_spaces.Function" title="fe_utils.function_spaces.Function"><code class="xref py py-class docutils literal notranslate"><span class="pre">Function</span></code></a> to enable values
to be interpolated into vector-valued functions, and to create quiver plots.</p></li>
<li><p>Assemble and solve the required mixed system.</p></li>
<li><p>Compute the <span class="math notranslate nohighlight">\(L^2\)</span> error of the mixed solution from the analytic solution.</p></li>
</ol>
<p>A convergence test for your code is provided in
<code class="docutils literal notranslate"><span class="pre">test/test_13_mastery_convergence.py</span></code>. In order to be compatible with
this code, your implementation of
<a class="reference internal" href="fe_utils.solvers.html#fe_utils.solvers.mastery.solve_mastery" title="fe_utils.solvers.mastery.solve_mastery"><code class="xref py py-func docutils literal notranslate"><span class="pre">solve_mastery()</span></code></a> should return its results
as a tuple of the form <code class="xref py py-data docutils literal notranslate"><span class="pre">(u,</span> <span class="pre">p),</span> <span class="pre">error</span></code>. This is a slight change from
the comment in the code which takes into account that the problem is mixed.
The obvious consequential change will be needed at the end of
<a class="reference internal" href="fe_utils.solvers.html#module-fe_utils.solvers.mastery" title="fe_utils.solvers.mastery"><code class="xref py py-mod docutils literal notranslate"><span class="pre">fe_utils.solvers.mastery</span></code></a>.</p>
</div></div></section>
</section>
<div class="clearer"></div>
</div>
</div>
</div>
<div class="clearer"></div>
</div>
<div class="footer" role="contentinfo">
© Copyright 2014-2024, David A. Ham and Colin J. Cotter.
Created using <a href="https://www.sphinx-doc.org/">Sphinx</a> 7.4.7.
</div>
</body>
</html>