-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathNumGy_LV05.qmd
424 lines (325 loc) · 11.3 KB
/
NumGy_LV05.qmd
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
---
title: "Numerical Simulation Methods in Geophysics, Part 5: Timestepping"
subtitle: "1. MGPY+MGIN"
author: "[email protected]"
title-slide-attributes:
data-background-image: pics/tubaf-logo.png
data-background-size: 20%
data-background-position: 10% 95%
format:
tubaf-revealjs:
html-math-method: mathjax
chalkboard: true
include-in-header:
- text: |
<script>
window.MathJax = {
loader: {
load: ['[tex]/physics']
},
tex: {
packages: {'[+]': ['physics']}
}
};
</script>
slide-number: c/t
transition: slide
transition-speed: fast
menu:
side: left
# jupyter: python3
---
# Recap last lessons & exercises
## The general case
$\Delta x \ne 1$ & $a \ne 1$ $\Rightarrow$ $a \pdv{u}{x} \approx a_i\frac{u_{i+1}-u_i}{x_{i+1}-x_i}$
```{python}
#| eval: true
#| echo: false
import numpy as np
import matplotlib.pyplot as plt
import pygimli as pg
import pygimli.meshtools as mt
nx = 6
grid = mt.createGrid(nx, 2)
ax, cb = pg.show(grid)
ax.set_yticks([])
for i in range(nx):
ax.plot(i, 0, "ko", ms=10)
for i in range(nx-1):
ax.text(i+0.5, 0.5, rf"$a_{i}$", ha="center", va="center", fontsize=18)
ax.set_xticklabels([f"$u_{i}$" for i in range(6)])
for i in range(5):
ax.text(i+0.5, -0.1, "$u'_{"+f"{i}"+".5}$", va="top");
ax.set_xlabel("")
ax.set_ylabel("");
```
$\dv{x}(a \pdv{u}{x}) \approx (a_i\frac{u_{i+1}-u_i}{x_{i+1}-x_i} - a_{i-1}\frac{u_{i}-u_{i-1}}{x_{i}-x_{i-1}}) / (x_{i+1}-x_{i-1}) \cdot 2$
$$ A_{i, i-1} = a_{i-1} / (x_{i}-x_{i-1}) / (x_{i+1}-x_{i-1}) \cdot 2 $$
## The coupling coefficients
$$ A_{i,i-1} = C_i^{left} = a_{i-1} / (x_{i}-x_{i-1}) / (x_{i+1}-x_{i-1}) \cdot 2 $$
$$ A_{i,i+1} = C_i^{right} = a_i / (x_{i+1}-x_{i}) / (x_{i+1}-x_{i-1}) \cdot 2 $$
$$ \begin{bmatrix}
+1 & 0 & 0 & \ldots & & \\
C_1^L & -(C_1^L+C_1^R) & C_1^R & 0 & \ldots & \\
\vdots & \vdots & \ddots & \vdots & \\
\ldots & \ldots & 0 & C_N^L & -(C_N^L+C_N^R) & C_N^R \\
\ldots & \ldots & & 0 & -1 & +1
\end{bmatrix} \cdot\vb{u} =
\begin{bmatrix} u_B\\ f_1 \\ \vdots \\ f_N \\ g_B \Delta x_N \end{bmatrix} $$
## Symmetry
$$ A_{i+1,i} = C_{i+1}^{left} = a_i / (x_{i+1}-x_{i}) / (x_{i+2}-x_{i}) \cdot 2 $$
$$ A_{i,i+1} = C_i^{right} = a_i / (x_{i+1}-x_{i}) / (x_{i+1}-x_{i-1}) \cdot 2 $$
only symmetric if $\Delta x$ is constant around $x_i$, better take $a_i/(\Delta x_i)^2$
$$ A_{i,i-1} = C_i^{left} = a_{i-1} / (x_{i}-x_{i-1})^2 $$
$$ A_{i,i+1} = C_i^{right} = a_i / (x_{i+1}-x_{i})^2 $$
$\Rightarrow$ inaccuracies expected for non-equidistant discretization
## A closer look at the Dirichlet boundary
$$ \begin{bmatrix}
+1 & 0 & 0 & \ldots & & \\
C_1^L & -(C_1^L+C_1^R) & C_1^R & 0 & \ldots &
\end{bmatrix} \begin{bmatrix} u_B\\ f_1 \end{bmatrix}$$
1. $C$ can be differently scaled from 1 $\Rightarrow$ multiply with $C_i^L$
1. Matrix is non-symmetric
<!-- $\Rightarrow$ add $C_1^L$ -->
$$ \begin{bmatrix}
C_1^L & 0 & 0 & \ldots & & \\
C_1^L & -(C_1^L+C_1^R) & C_1^R & 0 & \ldots &
\end{bmatrix} \begin{bmatrix} u_B C_i^L\\ f_1 \end{bmatrix} $$
## A closer look at the Neumann boundary
$$ \begin{bmatrix}
\ldots & \ldots & 0 & C_N^L & -(C_N^L+C_N^R) & C_N^R \\
\ldots & \ldots & & 0 & -1 & +1
\end{bmatrix} \cdot\vb{u} =
\begin{bmatrix} f_N \\ g_B \Delta x_N \end{bmatrix} $$
1. $C$ can be differently scaled from 1 $\Rightarrow$ multiply with $C_N^R$
1. Matrix is non-symmetric $\Rightarrow$ multiply with -1
$$ \begin{bmatrix}
\ldots & \ldots & 0 & C_N^L & -(C_N^L+C_N^R) & C_N^R \\
\ldots & \ldots & & 0 & C_N^R & -C_N^R
\end{bmatrix} \cdot\vb{u} =
\begin{bmatrix} f_N \\ -g_B \Delta x_N C_N^R \end{bmatrix} $$
## Accuracy
How can we prove the accuracy of our solution?
* compare with analytical solutions
* single (Point) $f$ lead to piece-wise linear $u$ (correct)
* what about continuous source terms?<br> (e.g. radioactive elements in the Earth's crust)
## Analytical solution
$a$=1, $f(x)$=1 $\Rightarrow$ double integration $\Rightarrow$ quadratic function
$$u(x) = C_0 + C_1 x -1/2 x^2 $$
Left BC | Right BC | $C_0$ | $C_1$
----------|-----------|-----------------------|-------------------------
Dirichlet | Dirichlet | $u_L$ | $L/2 + (u_R-u_L)/ L$
Dirichlet | Neumann | $u_L$ | $g_R + L$
Neumann | Dirichlet | $u_R - g_L L + L^2/2$ | $g_L$
with $L=(x_N-x_0)$
## Tasks stationary heat equation
* finish your implementation so that it can work with any $x$, $a$ and $f$, and at least D-D or D-N boundary conditions
* make sure you achieve the same (phenomenological) results like in the "collection notebook" (variation of, a, x, f)
* simulate the two different cases (D-D, D-N) with some choices of uL and uR or gR
* make sure the N-D case is analog to D-N
* compute analytical solution and plot it with numerical solution
* change the discretization and improve the solution
## Next spatial dimension
:::: {.columns}
::: {.column}
::: {.r-stack}
![Simple 2D conductivity grid](notebooks/grid32sigma.svg)
![Simple 2D conductivity grid with FD stencil](notebooks/grid32sigmaStencil.svg){.fragment}
:::
:::
::: {.column}
:::: {.fragment}
$$ C_{i,j}^{right} = a_{i,j-1/2} / (x_{i+1}-x_{i})^2 $$
$a_{i,j-1/2}=(a_{i,j-1}+a_{i,j})/2$ ?
harmonic, geometric? weighting?
$a_{i,j-1/2}=\frac{a_{i,j-1}\Delta y_{j-1}+a_{i,j}\Delta y_{j}}{y_j+1-y_{j-1}}$?
::::
:::
::::
# Parabolic PDEs
## Instationary heat flow in 3D
$$ \pdv{T}{t} - \div (a\grad{T}) = \div q_s $$
* $a=\frac{k}{\rho c_p}$ [m²/s] thermal diffusivity - measure of heat transfer
* $k$ [W/m/K] thermal conductivity - measure of temperature transfer
* $c_p$ [J/kg/K] - heat capacity - measure of heat storage per mass
* $\rho$ (kg/m³) density
Water $k$=0.6 W/m/K, $\rho$=1000 kg/m³, $c$=4180 J/kg/K $\Rightarrow$ $a$=1.43e-7 m²/s
## Periodic boundary conditions
:::: {.columns}
::: {.column}
Upper boundary: daily/yearly variation
$$ T_0(t) = T_m + \Delta T \sin \omega t $$
$T_m$ mean temperature (e.g. 12°C), $\Delta T$ variation, e.g. 8°C
$\omega=2\pi/t_P$ daily ($T_d$=3600*24s) or yearly ($T_y=365 T_d$) cycle
:::
::: {.column}
```{python}
#| eval: true
#| echo: false
from math import sqrt, pi
import numpy as np
import matplotlib.pyplot as plt
```
```{python}
#| eval: true
#| echo: true
day = 3600 * 24
T0, dT = 12, 8
t = np.arange(100) / 50 * day
T = T0+dT*np.sin(t/day*2*np.pi)
plt.plot(t/day*24, T)
plt.xlabel("t (h)")
plt.ylabel("T (°C)")
plt.grid()
```
:::
::::
## Separation of variables
$$ \pdv{T}{t} = a\pdv[2]{T}{z} $$
$T(t, z)/\Delta T - T_0 = \theta(t) Z(z)$
$$ Z \pdv{\theta}{t} = a \theta \pdv[2]{Z}{z} $$
$$ \frac1\theta \pdv{\theta}{t} = C = a \frac{1}{Z}\pdv[2]{Z}{z} $$
## Solution
regarding the BC $e^{\imath\omega t}$ leads to $C=\imath\omega$ and thus $\theta=\theta_0 e^{\imath\omega t}$
$$ \pdv[2]{Z}{z} - \imath \frac{\omega}{a} Z = \pdv[2]{Z}{z} + n^2 Z = 0 $$
Helmholtz equation with solution $Z=Z_0 e^{\imath n z}$ ($n^2=\imath\omega/a$)
$$ Z=Z_0 e^{\imath n z} = Z_0 e^{\sqrt{\imath\omega/a}z} =
Z_0 e^{\sqrt{\omega/2a}(1+\imath)z} $$
$$ T(t, z)/\Delta T + T_0 = Z(z)\theta(t)=
Z_0\theta_0 e^{-\sqrt{\omega/2a}z} e^{i(\omega t-\sqrt{\omega/2a}z)} $$
## Interpretation
replacing the term $\sqrt{2a/\omega}=\sqrt{a t_P/\pi}=d$ leads to
$$ T(z, t) = T_0 + \Delta T e^{-z/d} \sin(\omega t - z/d) $$
* exponential damping of the temperature variation with decay depth $d$
* phase lag $z/d$ increases with depth, $z_\pi=\sqrt{2 a/\omega}\pi=\sqrt{a t_P \pi}$
1. Daily cycle: decay depth $d$=6cm, minimum depth=20cm
2. Yearly cycle: decay depth $d$=1.2m, minimum depth=4m
## Depth profiles
```{python}
#| output-location: column
#| code-line-numbers: false
#| echo: true
a = 1.5e-7
year = day*365
d = sqrt(a*year/pi)
t = np.arange(0, 1, 0.1) * year
z = np.arange(0, 6, 0.1)
fig, ax = plt.subplots()
for ti in t:
Tz = np.exp(-z/d)*np.sin(ti*2*pi/year-
z/d) * dT + T0
ax.plot(Tz, z, label="t={:.1f}".format(
ti/year))
ax.invert_yaxis()
ax.legend()
ax.grid()
```
## Temporal behaviour
```{python}
#| output-location: column
#| code-line-numbers: false
#| echo: true
t = np.arange(0, 1.01, 0.01) * year
z = np.arange(0, 3, 0.4)
fig, ax = plt.subplots()
for zi in z:
Tt = np.exp(-zi/d)*np.sin(t*2*pi/year-
zi/d) * dT + T0
ax.plot(t/year*12, Tt, label=f"z={zi:.1f}")
ax.set_xlim(0, 12)
ax.set_xlabel("t (months)")
ax.set_ylabel("T (°C)")
ax.legend()
ax.grid()
```
<!-- ## Experimental data from North Sea beach
![](pics/temperature2022-2023.svg) -->
## Time stepping - explicit method
$$ \pdv{T}{t} - a \pdv[2]{T}{z} = 0 $$
Finite-difference approximation
$$ \pdv{T}{t}^n \approx \frac{T^{n+1}-T^n}{\Delta t} = a \pdv[2]{T^n}{z} $$
## Explicit
:::: {.columns}
::: {.column}
Start $T^0$ with initial condition (e.g. 0)
Update field by
$$ T^{n+1} = T^n + a \pdv[2]{T^n}{z} \cdot \Delta t$$
E.g. by using the matrix $A$
:::
::: {.column}
```{python}
#| eval: true
#| echo: false
import numpy as np
import matplotlib.pyplot as plt
import pygimli as pg
import pygimli.meshtools as mt
grid = mt.createGrid(5, 4)
ax, cb = pg.show(grid)
ax.set_yticks([])
ax.plot([1, 2, 3], [1, 1, 1], "bo-", ms=12, lw=4)
ax.plot([2, 2], [1, 2], "g-", lw=4)
ax.plot(2, 2, "ro", ms=12)
```
:::
::::
## Implicit methods
:::: {.columns}
::: {.column}
$$ \pdv{T}{t}^{n+1} \approx \frac{T^{n+1}-T^n}{\Delta t} = a \pdv[2]{T}{z} ^{n+1} $$
$$ \frac{1}{\Delta t} T^{n+1} - a \pdv[2]{T}{z} ^{n+1} = \frac{1}{\Delta t} T^n $$
$$ (\vb M - \vb A) \vb u^{n+1} = \vb M \vb u^n $$
$\vb M$ - mass matrix
:::
::: {.column}
```{python}
#| eval: true
#| echo: false
import numpy as np
import matplotlib.pyplot as plt
import pygimli as pg
import pygimli.meshtools as mt
grid = mt.createGrid(5, 4)
ax, cb = pg.show(grid)
ax.set_yticks([])
ax.plot([1, 2, 3], [2, 2, 2], "bo-", ms=12, lw=4)
ax.plot([2, 2], [1, 2], "go-", ms=12, lw=4)
ax.plot(2, 2, "ro", ms=12)
```
:::
::::
## Mixed - Crank-Nicholson method
$$ \pdv{T}{t}^{n+1/2} \approx \frac{T^{n+1}-T^n}{\Delta t} = \frac12 a \pdv[2]{T}{z} ^n + \frac12 a \pdv[2]{T}{z} ^{n+1} $$
:::: {.columns}
::: {.column}
$$ \frac{2}{\Delta t} T^{n+1} - a \pdv[2]{T^{n+1}}{z} =
\frac{2}{\Delta t} T^n + a \pdv[2]{T^n}{z} $$
$$ (2\vb M - \vb A) \vb u^{n+1} = (2\vb M + \vb A) \vb u^n $$
<!-- $$ \pdv{T}{t}^n \approx \frac{T^{n+1}-T^n}{\Delta t} = \frac12 a \pdv[2]{T}{z} ^n + \frac12 a \pdv[2]{T}{z} ^{n+1} $$ -->
:::
::: {.column}
```{python}
#| eval: true
#| echo: false
import numpy as np
import matplotlib.pyplot as plt
import pygimli as pg
import pygimli.meshtools as mt
grid = mt.createGrid(5, 4)
ax, cb = pg.show(grid)
ax.set_yticks([])
ax.plot([1, 2, 3], [1, 1, 1], "bo-", ms=12, lw=4)
ax.plot([1, 2, 3], [2, 2, 2], "bo-", ms=12, lw=4)
ax.plot([2, 2], [1, 2], "g-", ms=12, lw=4)
ax.plot(2, 2, "ro", ms=12)
```
:::
::::
## Tasks instationary heat equation
* setup a discretization and compute the stiffness matrix $A$ for some $a$
* choose an initial condition (e.g. homogeneous)
* choose a time step $\Delta t$ and perform the explicit method using the surface temperature
* change the spatial/temporal discretization and observe the solution
* setup mass matrix and implement the implicit method for diff. $\Delta t$
* implement the Crank-Nicholson method and compare all three
* compare the solutions with the analytical solution