-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathBidomain_1D.jl
2211 lines (1817 loc) · 73.3 KB
/
Bidomain_1D.jl
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
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
### A Pluto.jl notebook ###
# v0.17.2
using Markdown
using InteractiveUtils
# This Pluto notebook uses @bind for interactivity. When running this notebook outside of Pluto, the following 'mock version' of @bind gives bound variables a default value (instead of an error).
macro bind(def, element)
quote
local iv = try Base.loaded_modules[Base.PkgId(Base.UUID("6e696c72-6542-2067-7265-42206c756150"), "AbstractPlutoDingetjes")].Bonds.initial_value catch; b -> missing; end
local el = $(esc(element))
global $(esc(def)) = Core.applicable(Base.get, el) ? Base.get(el) : iv(el)
el
end
end
# ╔═╡ 2595b432-a4bc-4f64-856b-c44851122a14
begin
using PlutoUI, Plots, Roots , PyPlot, PlutoUI,HypertextLiteral,
ExtendableGrids,VoronoiFVM, PlutoVista,GridVisualize,LinearAlgebra, JLD2
# using DifferentialEquations
default_plotter!(PlutoVista)
pyplot()
TableOfContents() # TODO: use
end;
# ╔═╡ 735ba655-139a-4e74-ac87-2755daaf373f
begin
using ShortCodes
references=[
DOI("10.1137/070680503"),
DOI("10.5281/zenodo.3529808"),
DOI("10.21105/joss.01848")
]
end
# ╔═╡ e9d7e1b3-cd57-4e6d-ad95-94b5de379dab
md"""
# Bidomain Model
## Introduction
This is the report of Implementation of Bidomain Model as Formulated in [1], Usng **VoronoiFVM.jl** [2], as the Final project of the course "Scienfic Computing" at Techinical University Of Belin in winter semster of 21/22, under Supervision of Dr. Juergen Fuhrmann.
Presented and implemented by Alexander Quinlan and Erfan Baradarantohidi
## Bidomain Model
The bidomain is a system of partial differential equations used to model the propagation of electrical potential waves in myocardium. It is composed of coupled parabolic and eliptic PDEs, as well as at least one ordinary differential equation to model the ion activity through the cardic cell membrance [1].
The electrical properties of the myocardium are generally described by the bidomain equations, a set of coupled parabolic and elliptic partial differential equations (PDEs) that represents the tissue as two separate, distinct continua - one intracellular and the other extracellular. The intracellular and the extracellular media are connected via the cell membrane, and thus the two PDEs are coupled at each point in space through a set of complex, non-linear ordinary differential equations (ODEs), which describe the ionic transport across the cell membrane. Certain modelling environments use the monodomain representation of cardiac activity, which involves solving a single parabolic PDE, by assuming either that the extracellular potentials are negligible, or that the anisotropy ratios are equal in the intracellular and the extracellular domains[3].
"""
# ╔═╡ 0652ee6a-7705-40e1-ae7b-9dc6b04ce0e6
md"""
# Biodomain Problem Modeled as a system of PDEs
$u= u_i-u_e, \ u_e \; and \; \ v \;on \: \ Q\times\Omega= [0,\ T=30] \times [0, \ L=70]$
$\frac{\partial u}{\partial t} -
\nabla\cdot(\sigma_i(\nabla u+ \nabla u_e))
- \frac{1}{\epsilon} f(u,v) = 0$
$\nabla\cdot(\sigma_i\nabla u+(\sigma_i +\sigma_e)\nabla u_e)=0$
$\frac{\partial v}{\partial t} - \epsilon g(u,v)=0$
We take the initial and boundary conditions $u$ and $v$ constant at the equilibrium value of the system $f(u,v)=g(u,v)=0$ and $u_e$ constant at $0$ , except on the interval $[0,L/20]$ where we fix $u$ at a supercritical value, $u(x)=2$
So we write :
$u_e(0,x)=0 \:\; \forall x\in \Omega$
$u(0,x)=2 \;\;\forall x \in [0,L/20]$
$f(u(0,x_1),v(0,x))=g(u(0,x_1),v(0,x))=0 \;\;\forall x_1 \in [L/20,L] \;\; and\;\; \forall x \in [0,L]$
Where $f(u,v) = u - \frac{u^3}{3} - v$ and $g(u,v) = u + \beta -\gamma v$
and Neumann boundary conditions:
$\partial_x u(0)=\partial_x u(L)=0 \;\; and \;\; \partial_x u_e(0)=\partial_x u_e(L)=0$
Since pure Neumann boundary conditions are ill-posed, we avoid solving a singular system with:
$u_e(0)=0.$
"""
# ╔═╡ b14082d3-2140-48f5-80f7-f539600a5dcc
md"""
### What do these variables describe?
At each point in the computational domain, two electrical potential:
- ``u_i`` the **interacellular** potentionial
- ``u_e`` the **extracellular** potentionial
are recovered, representing the average of the electrical potential over the extracellular and the interacellular space, recpectively, in the neighborhoods of that point.
- ``u=u_i - u_e`` is the **transmembrance** potential
- ``\sigma_i`` and ``\sigma_e`` are second order tensor, representing **electrical conductivity** in each spatial direcation (in this implementation, it is assumed that conductivity of heart tissue is same in each direction, hence it is constant)
- ``v`` is a lumped **ionic variable** (the ODE as last equation is not properly part of the bidomain model but rather models the ionic activity across the cellular membrance which is responsible foe electrical activation of the tissue)
- ``\epsilon`` is a paarmeter linled to the **ratio between the repolarization rate and the tissue excitation rate**
- function ``f`` and ``g`` :
$f(u,v) = u - \frac{u^3}{3} - v$
$g(u,v) = u + \beta -\gamma v$
- ``\gamma `` control **ion transport** , ``\beta`` is linked to **cell excitabillity**
"""
# ╔═╡ 4e838ed5-9373-4298-82f5-ca82864545b1
md"""
### Bidomain as a reaction-diffusion system
In order to define bidomain problem as a system of $n$ coupled PDEs so that we can handle it to VoronoiFVM we define "**reaction**" $r(u)$, ""**storage"**" $s(u)$, ""**source"**" $f$ and ""**flux"**" $\vec{j}$ in vector form as follows:
$\partial_t s(u) + \nabla \cdot \vec{j}(u) + r(u) = f$
$u = \begin{bmatrix} u\\u_e\\v\end{bmatrix}$
$s(u) = \begin{bmatrix} u\\0\\v\end{bmatrix}$
$\vec{j}(u) = \begin{bmatrix} -\sigma_i(\nabla u+\nabla u_e)\\
\sigma_i\nabla u + (\sigma_i+\sigma_e)\nabla u_e
\\0\end{bmatrix}$
$r(u) = \begin{bmatrix} -f(u,v)/\epsilon\\0\\ -\epsilon g(u,v)\end{bmatrix}$
$f = 0$.
storage, recation and flux terms in our system are not dependent on $\vec{x}$ and $t$, but in general they can be.
"""
# ╔═╡ b1b889d3-e505-438a-b1f9-1c105614ab74
md"""
# Discretization over Finite Volumes $\omega_k$:
#### Constructing control volumes
We construct and formulate excatly as it was done in the Lecture
"""
# ╔═╡ 16399344-0b4f-4f8e-a364-3221448d4c17
md"""
Assume $\Omega\subset \mathbb{R}^d$ is a polygonal domain such that
$\partial\Omega=\bigcup_{m\in\mathcal{G}} \Gamma_m$, where $\Gamma_m$ are planar such that $\vec{n}|_{\Gamma_m}=\vec{n}_m$.
Subdivide $\Omega$ into into a finite number of
_control volumes_
$\bar\Omega= \bigcup_{k\in \mathcal N} \bar \omega_k$
such that
- $\omega_k$ are open convex domains such that $\omega_k\cap \omega_l = \emptyset$ if
$\omega_k\neq \omega_l$
- $\sigma_{kl}=\bar \omega_k \cap \bar \omega_l$ are either empty,
points or straight lines.
If $|\sigma_{kl}|>0$ we say that $\omega_k$, $\omega_l$
are neighbours.
- $\vec{n}_{kl}\perp \sigma_{kl}$: normal of $\partial\omega_k$ at $\sigma_{kl}$
- $\mathcal{N}_k = \{l \in \mathcal{N}: |\sigma_{kl}|>0\}$: set of neighbours of $\omega_k$
- $\gamma_{km}=\partial\omega_k \cap \Gamma_m$: boundary part of $\partial\omega_k$
- $\mathcal{G}_k= \{ m \in \mathcal{G}: |\gamma_{km}|>0\}$: set of non-empty boundary parts of $\partial\omega_k$.
$\Rightarrow$ $\partial\omega_k= \left(\cup_{l\in \mathcal{N}_k} \sigma_{kl}\right)\bigcup\left(\cup_{m\in \mathcal{G}_k} \gamma_{km}\right)$
"""
# ╔═╡ 1bddacdb-ae8c-4ca5-8165-246d9a41bc5f
md"""
To each control volume $\omega_k$ assign a _collocation
point_: $\vec{x}_k \in \bar \omega_k$ such that\\
- _Admissibility condition_:if $l\in \mathcal N_k$ then the
line $\vec{x}_k\vec{x}_l$ is orthogonal to $\sigma_{kl}$
- For a given function $u:\Omega \to \mathbb{R}$ this will allow to associate its value $u_k=u(\vec{x}_k)$
as the value of an unknown at $\vec{x}_k$.
- For two neigboring control volumes $\omega_k, \omega_l$ , this will allow to approximate
$\vec\nabla u \cdot \vec{n}_{kl} \approx \frac{u_l - u_k}{h_{kl}}$
- _Placement of boundary unknowns at the boundary_: if
$\omega_k$ is situated at the boundary, i.e. for
$|\partial \omega_k \cap \partial \Omega| >0$,
then $\vec{x}_k \in \partial\Omega$
- This will allow to apply boundary conditions in a direct manner
"""
# ╔═╡ bf447d84-7ad8-4e9f-bfaa-30789c0532f1
md"""
### First Equation
For $k\in\mathcal{N}$ Integrate over each control volume $\omega_k$:
```math
\newcommand{\eps}{\varepsilon}
\newcommand{\pth}[1]{\left(#1\right)}
\begin{equation}
\int_{\omega_k}\frac{\partial u}{\partial t}=\int_{\omega_k}\frac{1}{\eps}f(u,v)d\omega +\\
\int_{\omega_k}\nabla \cdot (\sigma_i \nabla u)d\omega + \int_{\omega_k}\nabla \cdot (\sigma_i \nabla u_e)d\omega
\end{equation}
```
By use of Gauss's thereom so the integral of divergence of flux over volume become integral of flux multiply by normal vector over bundary:
```math
\begin{align*}
\int_{\omega_k}\frac{\partial u}{\partial t} =\int_{\omega_k}\frac{1}{\eps}f(u,v)d\omega +
\int_{\partial\omega_k}\sigma_i \nabla u \cdot \vec{n}ds +
\int_{\partial\omega_k}\sigma_i \nabla u_e \cdot \vec{n}ds
\end{align*}
```
Since the $\partial\Omega_k$ is either edges in domain boundary $\mathcal{G}_k$ or neighbours $\mathcal{N}_k$
```math
\begin{align*}
\int_{\omega_k}\frac{\partial u}{\partial t} &=\int_{\omega_k}\frac{1}{\eps}f(u,v)d\omega + \sum_{l \in N_k}\int_{\sigma_{kl}} \sigma_i \nabla u \cdot \vec{n}_{kl}ds +
\sum_{m \in \mathcal{G}_k}\int_{\gamma_{kl}} \sigma_i \nabla u \cdot \vec{n}_{m}ds \\
&+ \sum_{l \in N_k}\int_{\sigma_{kl}} \sigma_i \nabla u_e \cdot \vec{n}_{kl}ds +
\sum_{m \in \mathcal{G}_k}\int_{\gamma_{kl}} \sigma_i \nabla u_e \cdot \vec{n}_{m}ds
\end{align*}
```
Finite difference approximation of normal derivative:
```math
\begin{align*}
h_{kl} = |x_k - x_l|\\
\sigma_i \nabla u \cdot \vec{n} \approx \sigma_i \frac{u_k - u_l}{h_{kl}}\\
\sigma_i \nabla u_e \cdot \vec{n} \approx \sigma_i \frac{u_{e_k} - u_{e_l}}{h_{kl}}
\end{align*}
```
We also approximate integrals by the length of the edge times approximated flux:
```math
\begin{align*}
|\omega_k|\frac{\partial u_k}{\partial t} &= \int_{w_k}\frac{1}{\eps}(u - \frac{u^3}{3} - v)d\omega + \sum_{l \in N_k} |\sigma_{kl}| \sigma_i \frac{u_k - u_l}{h_{kl}} +
\sum_{m \in \mathcal{G}_k} |\gamma_{km}| \sigma_i \frac{u_k - u_l}{h_{kl}} \\
&+ \sum_{l \in N_k} |\sigma_{kl}| \sigma_i \frac{u_{e_k} - u_{e_l}}{h_{kl}} +
\sum_{m \in \mathcal{G}_k} |\gamma_{km}| \sigma_i \frac{u_{e_k} - u_{e_l}}{h_{kl}}
\end{align*}
```
Approximate reaction integral by multiplying $|\omega_k|$ and combining $u$ and $u_e$ sum:
```math
\begin{align*}
|\omega_k|\frac{\partial u_k}{\partial t} &= \frac{|\omega_k|}{\eps}\pth{u_k - \frac{u_k^3}{3} - v_k}+ \sum_{l \in N_k} |\sigma_{kl}| \sigma_i \frac{u_k - u_l + u_{e_k} - u_{e_l}}{h_{kl}} \\
&+ \sum_{m \in \mathcal{G}_k} |\gamma_{km}| \sigma_i \frac{u_k - u_l + u_{e_k} - u_{e_l}}{h_{kl}}
\end{align*}
```
Then discretizing in time in forward Euler method yields:
```math
\begin{align}
\frac{|\omega_k|}{\Delta t}(u_k^{n+1} - u_k^n) = \sum_{l \in N_k} |\sigma_{kl}| \sigma_i \frac{u_k^{n} - u_l^{n} + u_{e_k}^n - u_{e_l}^{n}}{h_{kl}}
&+ \frac{|\omega_k|}{\eps}\pth{u_k^{n} - \frac{(u_k^{n})^3}{3} - v_k^{n}}
\\+
\sum_{m \in \mathcal{G}_k} |\gamma_{km}| \sigma_i \frac{u_k^{n} - u_l^{n} + u_{e_k}^{n} - u_{e_l}^{n}}{h_{kl}}
\end{align}
```
As $u_k^n = u(\vec{x_k},n\Delta t)$
"""
# ╔═╡ 616e63fa-e7d6-42f1-b332-408a0734b531
md"""
### Second Equation
After following same procedure as the discretization of Equation 1:
Since there is no stoarge term in this equation
$u_k=u_k^n = u(\vec{x_k},n\Delta t)$
```math
\begin{align}
\sum_{l \in N_k} |\sigma_{kl}| \sigma_i \frac{u_k - u_l}{h_{kl}}
+ \sum_{l \in N_k} |\sigma_{kl}| \pth{\sigma_i+\sigma_e} \frac{u_{e_k} - u_{e_l}}{h_{kl}} \\
+\sum_{m \in \mathcal{G}_k} |\gamma_{km}| \sigma_i \frac{u_k - u_l}{h_{kl}}
+\sum_{m \in \mathcal{G}_k} |\gamma_{km}| (\sigma_i+\sigma_e) \frac{u_{e_k} - u_{e_l}}{h_{kl}}=0
\end{align}
```
"""
# ╔═╡ 3b4d813c-eca4-4f36-8d89-29b2f9ebbb2b
md"""
### Last Equation
```math
\begin{align*}
\int_{\omega_k}\frac{\partial v}{\partial t} d\omega = \int_{\omega_k} \eps(u + \beta - \gamma v)d\omega
\end{align*}
```
Approximation of Integrals:
```math
\begin{align*}
|\omega_k|\frac{\partial v_k}{\partial t} = \eps |\omega_k| (u_k + \beta - \gamma v_k)
\end{align*}
```
Then forward discretizing in time :
```math
\begin{align*}
\frac{|\omega_k|}{\Delta t}(v_k^{n+1} - v_k^{n}) - \eps |\omega_k| (u_k^{n} + \beta - \gamma v_k^{n})
\end{align*}
```
"""
# ╔═╡ b5aad722-99cb-424d-a2af-a37cefac1a59
md"""
### Observations
We successfully recreated the 1d problem without much difficulty, but had a lot of trouble recreating the spiral pattern observed in [1]. It turns out that the solution is very sensitive to spatial step-size-- less than roughly 100 cells per L=70 in each dimension would result in an incorrect solution that did not display conductivity in that direction. To avoid this, we ran our spiral solution overnight with N=120x120 and ``\Delta t = 10^{-2}``.
"""
# ╔═╡ daf3c9ee-01ff-494b-a0b0-a6a6560b75cc
md"""
# Simulations
"""
# ╔═╡ 62b11f3a-7d3d-4790-bacf-df10f297888d
md"""
### Parameters
"""
# ╔═╡ 768f4e8d-a8f4-4f20-8b68-e180521edcff
begin
# shared parameters
ϵ = 0.1
β = 1
γ = 0.5
species = ["u", "u_e", "v"]
σᵢ_anisotropic = 25*[0.263 0; 0 0.0263]
σₑ_anisotropic = 25*[0.263 0; 0 0.1087]
σᵢ = 1
σₑ = 1
# L = 70
L = 70
tf = 70
Δt = 1e-1;
# Storage and replay
overwrite_storage = false
compression = 100
end;
# ╔═╡ fb6aa962-a7d6-479b-bf11-52b8be2483f8
# 1D Parameters
begin
# 1d
N = 500
end
# ╔═╡ df3a31b9-2a70-44d5-8a7e-5d427dc327cc
# 2D Parameters
begin
Δt₂ = 1e-1
end
# ╔═╡ 6dcd9399-8bd2-48a6-8e06-5af201d4f730
md"""
### Initial Condition
"""
# ╔═╡ 54b8a5a0-2ca5-4973-b7b4-35d420e675c2
begin
ut(u0)=u0^3 + 3*u0 + 6
u₀=find_zero(ut,0)
v₀= 2*u₀ + 2
end;
# ╔═╡ a0249b52-5f13-4d37-9709-a069b065278d
function u_init(x)
u = u₀; v = v₀
if 0 ≤ x ≤ L/20
u = 2
end
uₑ = 0
return [u, uₑ, v]
end
# ╔═╡ 9b06c02d-ec46-4590-acac-821733809c6f
md"""
### Problem Implementation in VoronoiFVM
#### Physics Defintion
"""
# ╔═╡ bb37d88a-13e4-4a8d-a852-175d34a06afd
function storage(f,u,node)
f[1]= u[1]
f[2]= 0
f[3]= u[3]
end
# ╔═╡ adbc6875-1d0f-4e49-9abd-f12eced388ef
begin
f(u,v)= u - u^3/3 - v
g(u,v)= u + β - γ*v
end;
# ╔═╡ 7a367571-9e02-4449-808a-c319997c3df9
function reaction(y,u,node)
y[1]= -f(u[1],u[3])/ϵ
y[2]= 0
y[3]= -ϵ*g(u[1],u[3])
end
# ╔═╡ b484a050-2e8a-4a0f-b975-62692be5cea9
function flux(f,u,edge)
f[1] = -σᵢ*((u[1,2] - u[1,1]) + (u[2,2] - u[2,1]))
f[2] = σᵢ*(u[1,2] - u[1,1]) + (σᵢ+σₑ) * (u[2,2] - u[2,1])
f[3] = 0
end
# ╔═╡ 7528709e-39f3-462a-8845-5ee2c78dec7e
"""
From the paper: "To avoid solving a singular system we add the condition u_e(0) = 0"
"""
function breaction(f,u,node)
if node.coord[:,node.index] == [0, 0, 0]
f[2] = 0
end
end
# ╔═╡ 593f1fa3-ce33-41b1-aaa6-492d13a6584a
md"""
### Boundary Condition
"""
# ╔═╡ 64a70e08-dea2-40bc-8db9-cc3ced84ac2f
function bcondition_1d(f)
boundary_neumann!(f,species=1,region=1,value=0)
boundary_neumann!(f,species=1,region=2,value=0)
boundary_neumann!(f,species=2,region=1,value=0)
boundary_neumann!(f,species=2,region=2,value=0)
boundary_dirichlet!(f,species=2,region=1,value=0)
end
# ╔═╡ 089b0d84-3c9f-499f-8e41-9e6e71b74036
function grid_1d(N)
X = LinRange(0,L,N)
return grid=simplexgrid(X)
end
# ╔═╡ b38b57f6-ad80-4a50-b8bd-952d4d3e4866
function bidomain_1D(;N=100, Δt=1e-1, tf=30, Plotter=Plots)
grid = grid_1d(N)
#system Def
physics = VoronoiFVM.Physics(
storage=storage, flux=flux, reaction=reaction, breaction=breaction)
system=VoronoiFVM.System(
grid, physics, unknown_storage=:sparse)
enable_species!(system, species=[1,2,3])
bcondition_1d(system)
system
# Solving ODE
init = unknowns(system)
inival = map(u_init, grid)
init .= [tuple[k] for tuple in inival, k in 1:3]'
U = unknowns(system)
SolArray = copy(init)
tgrid = 0:Δt:tf
for t ∈ tgrid[2:end]
VoronoiFVM.solve!(U, init, system; tstep=Δt)
init .= U
SolArray = cat(SolArray, copy(U), dims=3)
end
vis = GridVisualizer(resolution=(400,300), dim=1, Plotter=Plots)
return grid, tgrid, SolArray, vis, system
end
# ╔═╡ 6e2fab73-30a8-4209-ba55-529f9c77cbc3
begin
xgrid, tgrid, sol, vis, system= bidomain_1D(N=N, Δt=Δt);
system
end
# ╔═╡ 28dcd70a-952d-492a-84f8-ee04eb83e360
function plot_at_t(t,vis,xgrid,sol)
species = ["u", "u_e", "v"]
tₛ = Int16(round(t/Δt))+1
scalarplot!(vis, xgrid, sol[1,:,tₛ], linestyle=:solid)
scalarplot!(vis, xgrid, sol[2,:,tₛ], linestyle=:dash)
plotted_sol = scalarplot!(
vis, xgrid, sol[3,:,tₛ], linestyle=:dot, legend=:best, show=true)
for (i,spec) ∈ zip(2 .* (1:length(species)), species)
plotted_sol[1][i][:label] = spec
end
Plots.plot!(labels=species)
Plots.plot!(ylims=(-2.5,2.5))
plotted_sol
end
# ╔═╡ e7cc8d45-4816-4623-9d4e-0d84d78c8fd2
begin
anim_tf = 30
step_size =3*Δt
anim = @animate for t in 0:step_size:anim_tf
plot_at_t(t,vis,xgrid,sol)
end
gif(anim)
end
# ╔═╡ 03fa9a1b-2e6a-4189-8f1b-8a702c9dcbf0
function contour_plot(spec)
PyPlot.clf()
Xgrid = xgrid[Coordinates][:]
Tgrid = tgrid
Zgrid = sol[spec,:,:]
PyPlot.suptitle("Space-time contour plot of "*species[spec])
PyPlot.contourf(Xgrid,Tgrid,Zgrid',cmap="viridis",levels=100)
axes=PyPlot.gca()
axes.set_aspect(2)
PyPlot.colorbar()
PyPlot.xlabel(L"x")
PyPlot.ylabel(L"t")
figure=PyPlot.gcf()
figure.set_size_inches(5,5)
figure
end
# ╔═╡ ad2b0f17-789c-47b8-8775-45331823036e
function contour_2d_at_t(spec, t, Δt, xgrid, sol)
tₛ = Int16(round(t/Δt))+1
p = scalarplot(
xgrid,sol[spec,:,tₛ], Plotter=PyPlot, colormap=:viridis,
title=species[spec]*
" at t="*string(t))
PyPlot.xlabel(L"x")
PyPlot.ylabel(L"y")
figure=PyPlot.gcf()
figure.set_size_inches(5,5)
figure
end
# ╔═╡ df344c5f-06e6-41cf-adf8-63a184445388
@bind species_1d PlutoUI.Select([1, 2, 3])
# ╔═╡ 4828996d-6084-4279-b4d9-aa52254682c2
contour_plot(species_1d)
# ╔═╡ c7b6354b-398f-43cb-b3f1-ab410ecbea02
md"""
### 2D Problem
"""
# ╔═╡ 65194fa2-a30c-4423-b5c6-eab8507af235
""""
Create the 2d Grid
"""
function grid_2d(N)
X = LinRange(0,L,N[1])
Y = LinRange(0,L,N[2])
grid2d = simplexgrid(X,Y)
return grid2d
end
# ╔═╡ 6a2270d4-cb7c-41a1-a300-d87079666145
"""
Initialize the 2d problem
"""
function u_2d(x,y)
u = u₀
v = v₀
uₑ = 0
# Rectangle in middle
mid = L/2
# if mid ≤ x ≤ 5 + mid && mid ≤ y ≤ 5 + mid
# u = 2
# end
# Classic
if 0 ≤ x ≤ 3.5 && 0 ≤ y ≤ 70
u = 2
end
if 31 ≤ x ≤ 39 && 0 ≤ y ≤ 35
v = 2
end
return [u, uₑ, v]
end
# ╔═╡ 54c521f5-2165-442d-baee-64b9f413d99f
function u_2d_as_1d(x,y)
return u_init(x)
end
# ╔═╡ 51d646e9-a066-4068-9598-b9673bdda8f8
"""
2D Boundaries
"""
function bcondition_2d(f)
for species ∈ [1 2] # Species 2 as well??
for region ∈ [1 2 3 4]
boundary_neumann!(f,species=species,region=region,value=0)
end
end
end
# ╔═╡ cbcbcfac-9e76-4a44-b2a9-40c4e8d089d6
function bidomain_2d(;N=100, Δt=1e-1, tf=70, Plotter=Plots, setup="1d")
xgrid = grid_2d(N)
#system Def
physics = VoronoiFVM.Physics(
storage=storage, flux=flux, reaction=reaction, breaction=breaction)
system=VoronoiFVM.System(
xgrid, physics, unknown_storage=:sparse)
enable_species!(system, species=[1,2,3])
bcondition_2d(system)
init = unknowns(system)
U = unknowns(system)
if setup == "1d"
inival = map(u_2d_as_1d, xgrid)
else
inival = map(u_2d, xgrid)
end
init .= [tuple[k] for tuple in inival, k in 1:3]'
SolArray = copy(init)
vis = GridVisualizer(resolution=(400,300), dim=2, Plotter=Plots)
# tgrid = 0:Δt:T
tgrid = 0:Δt:tf
for t ∈ tgrid[2:end]
if t % 1 == 0
println(string("At t=", t))
# contour_2d_at_t(species_select,t,Δt₂,xgrid,SolArray)
end
VoronoiFVM.solve!(U, init, system; tstep=Δt)
init .= U
SolArray = cat(SolArray, copy(U), dims=3)
end
return xgrid, tgrid, SolArray, vis, system
end
# ╔═╡ af5a7305-7154-4053-8f8c-375c07f141e2
md"""
### 2D Problem with 1D Setup
"""
# ╔═╡ 4d6b1dc5-4a62-4b86-81de-c9a5e4c76d60
begin
N₂ = (100,25)
xgrid₂, tgrid₂, sol₂, vis₂, sys₂ = bidomain_2d(N=N₂,tf=tf, Δt=Δt₂, Plotter=PyPlot);
end
# ╔═╡ 4cfec74c-46de-4281-a8c1-7b5a76ea851f
@bind species_select PlutoUI.Select([1, 2, 3])
# ╔═╡ 3f00ecee-2003-4c14-aaa5-d525f3b11607
@bind time₂ PlutoUI.Slider(0:tf,show_value=true)
# ╔═╡ 9f834292-40d4-4e70-a31b-8bb37884176d
contour_2d_at_t(species_select,time₂,Δt₂,xgrid₂,sol₂)
# ╔═╡ 07009ae2-b8fa-4f6c-ad76-7a9901d3c304
md"""
## 2D Problem with 2D Setup
#### Using precompiled data
"""
# ╔═╡ 49c346dc-1a3f-4c3e-ab28-68f00c2b7ca1
@bind species_replay PlutoUI.Select([1, 2, 3])
# ╔═╡ fe52e576-d533-4a79-9707-ee867ea9df1c
@bind time_replay PlutoUI.Slider(0:tf,show_value=true)
# ╔═╡ 7695a4ff-15d1-4978-8c84-f22bc8cb9ccc
contour_2d_at_t(
species_replay,
time_replay,
Δt₂ * compression / 10, # This was run at 1/10 the timestep so dividing by another factor of 10
load_object("xgrid.jld2"),
load_object("sol.jld2")
)
# ╔═╡ ceebe18c-2e85-4e5a-bcad-437bc60453b7
begin
if overwrite_storage
xgrid3, tgrid3, sol3, vis3, sys3 = bidomain_2d(N=N₂,T=tf, Δt=Δt₂, Plotter=PyPlot, setup="2d");
end
end
# ╔═╡ f22b494e-f0f9-4058-88f6-29aeac95b2f4
# Compress (only store every 100 timepts) and save the full solution for graphing
begin
if overwrite_storage # Only store the solution if we really want to overwrite
sol_compressed = cat(sol3[:, :, 1], sol3[:, :, compression], dims=3)
for i ∈ 2:70
sol_compressed = cat(sol_compressed, sol3[:, :, i*compression], dims=3)
end
save_object("sol.dat", sol_compressed)
save_object("xgrid.dat", xgrid3)
end
end
# ╔═╡ 40e01b91-bbc1-41ba-a7b3-72c859fc8ae5
md"""
# Bibliography
"""
# ╔═╡ 00000000-0000-0000-0000-000000000001
PLUTO_PROJECT_TOML_CONTENTS = """
[deps]
ExtendableGrids = "cfc395e8-590f-11e8-1f13-43a2532b2fa8"
GridVisualize = "5eed8a63-0fb0-45eb-886d-8d5a387d12b8"
HypertextLiteral = "ac1192a8-f4b3-4bfe-ba22-af5b92cd3ab2"
JLD2 = "033835bb-8acc-5ee8-8aae-3f567f8a3819"
LinearAlgebra = "37e2e46d-f89d-539d-b4ee-838fcccc9c8e"
Plots = "91a5bcdd-55d7-5caf-9e0b-520d859cae80"
PlutoUI = "7f904dfe-b85e-4ff6-b463-dae2292396a8"
PlutoVista = "646e1f28-b900-46d7-9d87-d554eb38a413"
PyPlot = "d330b81b-6aea-500a-939a-2ce795aea3ee"
Roots = "f2b01f46-fcfa-551c-844a-d8ac1e96c665"
ShortCodes = "f62ebe17-55c5-4640-972f-b59c0dd11ccf"
VoronoiFVM = "82b139dc-5afc-11e9-35da-9b9bdfd336f3"
[compat]
ExtendableGrids = "~0.9.5"
GridVisualize = "~0.5.1"
HypertextLiteral = "~0.9.3"
JLD2 = "~0.4.22"
Plots = "~1.27.5"
PlutoUI = "~0.7.38"
PlutoVista = "~0.8.13"
PyPlot = "~2.10.0"
Roots = "~2.0.0"
ShortCodes = "~0.3.3"
VoronoiFVM = "~0.16.3"
"""
# ╔═╡ 00000000-0000-0000-0000-000000000002
PLUTO_MANIFEST_TOML_CONTENTS = """
# This file is machine-generated - editing it directly is not advised
[[AbstractPlutoDingetjes]]
deps = ["Pkg"]
git-tree-sha1 = "8eaf9f1b4921132a4cff3f36a1d9ba923b14a481"
uuid = "6e696c72-6542-2067-7265-42206c756150"
version = "1.1.4"
[[AbstractTrees]]
git-tree-sha1 = "03e0550477d86222521d254b741d470ba17ea0b5"
uuid = "1520ce14-60c1-5f80-bbc7-55ef81b5835c"
version = "0.3.4"
[[Adapt]]
deps = ["LinearAlgebra"]
git-tree-sha1 = "af92965fb30777147966f58acb05da51c5616b5f"
uuid = "79e6a3ab-5dfb-504d-930d-738a2a938a0e"
version = "3.3.3"
[[ArgCheck]]
git-tree-sha1 = "a3a402a35a2f7e0b87828ccabbd5ebfbebe356b4"
uuid = "dce04be8-c92d-5529-be00-80e4d2c0e197"
version = "2.3.0"
[[ArgTools]]
uuid = "0dad84c5-d112-42e6-8d28-ef12dabb789f"
[[ArnoldiMethod]]
deps = ["LinearAlgebra", "Random", "StaticArrays"]
git-tree-sha1 = "62e51b39331de8911e4a7ff6f5aaf38a5f4cc0ae"
uuid = "ec485272-7323-5ecc-a04f-4719b315124d"
version = "0.2.0"
[[ArrayInterface]]
deps = ["Compat", "IfElse", "LinearAlgebra", "Requires", "SparseArrays", "Static"]
git-tree-sha1 = "c933ce606f6535a7c7b98e1d86d5d1014f730596"
uuid = "4fba245c-0d91-5ea0-9b3e-6abc04ee57a9"
version = "5.0.7"
[[Artifacts]]
uuid = "56f22d72-fd6d-98f1-02f0-08ddc0907c33"
[[AutoHashEquals]]
git-tree-sha1 = "45bb6705d93be619b81451bb2006b7ee5d4e4453"
uuid = "15f4f7f2-30c1-5605-9d31-71845cf9641f"
version = "0.2.0"
[[BangBang]]
deps = ["Compat", "ConstructionBase", "Future", "InitialValues", "LinearAlgebra", "Requires", "Setfield", "Tables", "ZygoteRules"]
git-tree-sha1 = "b15a6bc52594f5e4a3b825858d1089618871bf9d"
uuid = "198e06fe-97b7-11e9-32a5-e1d131e6ad66"
version = "0.3.36"
[[Base64]]
uuid = "2a0f44e3-6c83-55bd-87e4-b1978d98bd5f"
[[Baselet]]
git-tree-sha1 = "aebf55e6d7795e02ca500a689d326ac979aaf89e"
uuid = "9718e550-a3fa-408a-8086-8db961cd8217"
version = "0.1.1"
[[Bijections]]
git-tree-sha1 = "705e7822597b432ebe152baa844b49f8026df090"
uuid = "e2ed5e7c-b2de-5872-ae92-c73ca462fb04"
version = "0.1.3"
[[Bzip2_jll]]
deps = ["Artifacts", "JLLWrappers", "Libdl", "Pkg"]
git-tree-sha1 = "19a35467a82e236ff51bc17a3a44b69ef35185a2"
uuid = "6e34b625-4abd-537c-b88f-471c36dfa7a0"
version = "1.0.8+0"
[[Cairo_jll]]
deps = ["Artifacts", "Bzip2_jll", "Fontconfig_jll", "FreeType2_jll", "Glib_jll", "JLLWrappers", "LZO_jll", "Libdl", "Pixman_jll", "Pkg", "Xorg_libXext_jll", "Xorg_libXrender_jll", "Zlib_jll", "libpng_jll"]
git-tree-sha1 = "4b859a208b2397a7a623a03449e4636bdb17bcf2"
uuid = "83423d85-b0ee-5818-9007-b63ccbeb887a"
version = "1.16.1+1"
[[ChainRulesCore]]
deps = ["Compat", "LinearAlgebra", "SparseArrays"]
git-tree-sha1 = "9950387274246d08af38f6eef8cb5480862a435f"
uuid = "d360d2e6-b24c-11e9-a2a3-2a2ae2dbcce4"
version = "1.14.0"
[[ChangesOfVariables]]
deps = ["ChainRulesCore", "LinearAlgebra", "Test"]
git-tree-sha1 = "bf98fa45a0a4cee295de98d4c1462be26345b9a1"
uuid = "9e997f8a-9a97-42d5-a9f1-ce6bfc15e2c0"
version = "0.1.2"
[[CodecZlib]]
deps = ["TranscodingStreams", "Zlib_jll"]
git-tree-sha1 = "ded953804d019afa9a3f98981d99b33e3db7b6da"
uuid = "944b1d66-785c-5afd-91f1-9de20f533193"
version = "0.7.0"
[[ColorSchemes]]
deps = ["ColorTypes", "Colors", "FixedPointNumbers", "Random"]
git-tree-sha1 = "12fc73e5e0af68ad3137b886e3f7c1eacfca2640"
uuid = "35d6a980-a343-548e-a6ea-1d62b119f2f4"
version = "3.17.1"
[[ColorTypes]]
deps = ["FixedPointNumbers", "Random"]
git-tree-sha1 = "024fe24d83e4a5bf5fc80501a314ce0d1aa35597"
uuid = "3da002f7-5984-5a60-b8a6-cbb66c0b333f"
version = "0.11.0"
[[Colors]]
deps = ["ColorTypes", "FixedPointNumbers", "Reexport"]
git-tree-sha1 = "417b0ed7b8b838aa6ca0a87aadf1bb9eb111ce40"
uuid = "5ae59095-9a9b-59fe-a467-6f913c188581"
version = "0.12.8"
[[Combinatorics]]
git-tree-sha1 = "08c8b6831dc00bfea825826be0bc8336fc369860"
uuid = "861a8166-3701-5b0c-9a16-15d98fcdc6aa"
version = "1.0.2"
[[CommonSolve]]
git-tree-sha1 = "68a0743f578349ada8bc911a5cbd5a2ef6ed6d1f"
uuid = "38540f10-b2f7-11e9-35d8-d573e4eb0ff2"
version = "0.2.0"
[[CommonSubexpressions]]
deps = ["MacroTools", "Test"]
git-tree-sha1 = "7b8a93dba8af7e3b42fecabf646260105ac373f7"
uuid = "bbf7d656-a473-5ed7-a52c-81e309532950"
version = "0.3.0"
[[Compat]]
deps = ["Base64", "Dates", "DelimitedFiles", "Distributed", "InteractiveUtils", "LibGit2", "Libdl", "LinearAlgebra", "Markdown", "Mmap", "Pkg", "Printf", "REPL", "Random", "SHA", "Serialization", "SharedArrays", "Sockets", "SparseArrays", "Statistics", "Test", "UUIDs", "Unicode"]
git-tree-sha1 = "b153278a25dd42c65abbf4e62344f9d22e59191b"
uuid = "34da2185-b29b-5c13-b0c7-acf172513d20"
version = "3.43.0"
[[CompilerSupportLibraries_jll]]
deps = ["Artifacts", "Libdl"]
uuid = "e66e0078-7015-5450-92f7-15fbd957f2ae"
[[CompositeTypes]]
git-tree-sha1 = "d5b014b216dc891e81fea299638e4c10c657b582"
uuid = "b152e2b5-7a66-4b01-a709-34e65c35f657"
version = "0.1.2"
[[CompositionsBase]]
git-tree-sha1 = "455419f7e328a1a2493cabc6428d79e951349769"
uuid = "a33af91c-f02d-484b-be07-31d278c5ca2b"
version = "0.1.1"
[[Conda]]
deps = ["Downloads", "JSON", "VersionParsing"]
git-tree-sha1 = "6e47d11ea2776bc5627421d59cdcc1296c058071"
uuid = "8f4d0f93-b110-5947-807f-2305c1781a2d"
version = "1.7.0"
[[ConstructionBase]]
deps = ["LinearAlgebra"]
git-tree-sha1 = "f74e9d5388b8620b4cee35d4c5a618dd4dc547f4"
uuid = "187b0558-2788-49d3-abe0-74a17ed4e7c9"
version = "1.3.0"
[[Contour]]
deps = ["StaticArrays"]
git-tree-sha1 = "9f02045d934dc030edad45944ea80dbd1f0ebea7"
uuid = "d38c429a-6771-53c6-b99e-75d170b6e991"
version = "0.5.7"
[[DataAPI]]
git-tree-sha1 = "cc70b17275652eb47bc9e5f81635981f13cea5c8"
uuid = "9a962f9c-6df0-11e9-0e5d-c546b8b5ee8a"
version = "1.9.0"
[[DataStructures]]
deps = ["Compat", "InteractiveUtils", "OrderedCollections"]
git-tree-sha1 = "3daef5523dd2e769dad2365274f760ff5f282c7d"
uuid = "864edb3b-99cc-5e75-8d2d-829cb0a9cfe8"
version = "0.18.11"
[[DataValueInterfaces]]
git-tree-sha1 = "bfc1187b79289637fa0ef6d4436ebdfe6905cbd6"
uuid = "e2d170a0-9d28-54be-80f0-106bbe20a464"
version = "1.0.0"
[[Dates]]
deps = ["Printf"]
uuid = "ade2ca70-3891-5945-98fb-dc099432e06a"
[[DefineSingletons]]
git-tree-sha1 = "0fba8b706d0178b4dc7fd44a96a92382c9065c2c"
uuid = "244e2a9f-e319-4986-a169-4d1fe445cd52"
version = "0.1.2"
[[DelimitedFiles]]
deps = ["Mmap"]
uuid = "8bb1440f-4735-579b-a4ab-409b98df4dab"
[[DensityInterface]]
deps = ["InverseFunctions", "Test"]
git-tree-sha1 = "80c3e8639e3353e5d2912fb3a1916b8455e2494b"
uuid = "b429d917-457f-4dbc-8f4c-0cc954292b1d"
version = "0.4.0"
[[DiffResults]]
deps = ["StaticArrays"]
git-tree-sha1 = "c18e98cba888c6c25d1c3b048e4b3380ca956805"
uuid = "163ba53b-c6d8-5494-b064-1a9d43ac40c5"
version = "1.0.3"
[[DiffRules]]
deps = ["IrrationalConstants", "LogExpFunctions", "NaNMath", "Random", "SpecialFunctions"]
git-tree-sha1 = "dd933c4ef7b4c270aacd4eb88fa64c147492acf0"
uuid = "b552c78f-8df3-52c6-915a-8e097449b14b"
version = "1.10.0"
[[Distributed]]
deps = ["Random", "Serialization", "Sockets"]
uuid = "8ba89e20-285c-5b6f-9357-94700520ee1b"
[[Distributions]]
deps = ["ChainRulesCore", "DensityInterface", "FillArrays", "LinearAlgebra", "PDMats", "Printf", "QuadGK", "Random", "SparseArrays", "SpecialFunctions", "Statistics", "StatsBase", "StatsFuns", "Test"]
git-tree-sha1 = "5a4168170ede913a2cd679e53c2123cb4b889795"
uuid = "31c24e10-a181-5473-b8eb-7969acd0382f"
version = "0.25.53"
[[DocStringExtensions]]
deps = ["LibGit2"]
git-tree-sha1 = "b19534d1895d702889b219c382a6e18010797f0b"
uuid = "ffbed154-4ef7-542d-bbb7-c09d3a79fcae"
version = "0.8.6"
[[DomainSets]]
deps = ["CompositeTypes", "IntervalSets", "LinearAlgebra", "StaticArrays", "Statistics"]
git-tree-sha1 = "5f5f0b750ac576bcf2ab1d7782959894b304923e"
uuid = "5b8099bc-c8ec-5219-889f-1d9e522a28bf"
version = "0.5.9"
[[Downloads]]
deps = ["ArgTools", "LibCURL", "NetworkOptions"]
uuid = "f43a241f-c20a-4ad4-852c-f6b1247861c6"
[[DynamicPolynomials]]
deps = ["DataStructures", "Future", "LinearAlgebra", "MultivariatePolynomials", "MutableArithmetics", "Pkg", "Reexport", "Test"]
git-tree-sha1 = "d0fa82f39c2a5cdb3ee385ad52bc05c42cb4b9f0"
uuid = "7c1d4256-1411-5781-91ec-d7bc3513ac07"
version = "0.4.5"
[[EarCut_jll]]
deps = ["Artifacts", "JLLWrappers", "Libdl", "Pkg"]
git-tree-sha1 = "3f3a2501fa7236e9b911e0f7a588c657e822bb6d"
uuid = "5ae413db-bbd1-5e63-b57d-d24a61df00f5"
version = "2.2.3+0"
[[ElasticArrays]]
deps = ["Adapt"]
git-tree-sha1 = "a0fcc1bb3c9ceaf07e1d0529c9806ce94be6adf9"
uuid = "fdbdab4c-e67f-52f5-8c3f-e7b388dad3d4"