[6b82ac]: doc / doc.tex Maximize Restore History

Download this file

doc.tex    1355 lines (1133 with data), 54.1 kB

   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
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
\documentclass[a4paper,10pt]{book}
\usepackage[utf8x]{inputenc}
\usepackage{graphicx}
\usepackage{amsmath}
\usepackage{amssymb}
\usepackage{caption}
\usepackage{subcaption}
\usepackage{cclicenses}
\usepackage{url}
\usepackage{listings}
\usepackage{hyperref}
\usepackage{listings}
\usepackage[cmyk]{xcolor}
\usepackage{import}
\usepackage[english]{babel}
\usepackage[T1]{fontenc}
\usepackage[framemethod=TikZ]{mdframed}
\definecolor{dkgreen}{rgb}{0,0.6,0}
\definecolor{gray}{rgb}{0.5,0.5,0.5}
\definecolor{mauve}{rgb}{0.58,0,0.82}
\definecolor{BlueLUH}{cmyk}{1.0,0.7,0,0}
\colorlet{LightBlue}{BlueLUH!7!white}
\lstset{
backgroundcolor=\color{LightBlue},
basicstyle={\small\ttfamily},
language=C++,
aboveskip=3mm,
belowskip=3mm,
showstringspaces=false,
columns=flexible,
numbers=left,
numberstyle=\tiny\color{gray},
keywordstyle=\color{blue},
commentstyle=\color{dkgreen},
stringstyle=\color{mauve},
breaklines=true,
breakatwhitespace=true
tabsize=2
}
\author{Marco Vassallo}
\title{ \textbf{Fem-fenics} \\ \bigskip \textit{General purpose Finite Element library \\ for GNU-Octave}\\
\bigskip
\textsc{ Work in progress \\(help and remarks are welcome)}}
\begin{document}
\maketitle
\tableofcontents
\chapter{Introduction}
Fem-fenics is an open source package (pkg) for the resolution of partial differential equations with Octave.
The project has been developed during the Google Summer of Code 2013 with the help and the sustain of the GNU-Octave
community under the supervision of prof. De Falco.
The report is structured as follows:
\begin{itemize}
\item in chapter \ref{intr} we provide a simple reference guide for beginners
\item in chapter \ref{impl} is presented a detailed explanation of the relevant parts of the program. In this way, the
interested reader can see what there is ``behind'' and expecially anyone interested in it can learn quickly how
it is possible to extend the code and contribute to the project.
\item in chapter \ref{exem} more examples are provided. For a lot of them, we present the octave script
alongside with the code for Fenics (in C++ and/or Python) in order to provide the user with a quick reference
guide.
\end{itemize}
If you think that going inside the report could be boring, it is available a wiki at
\begin{center}
\url{http://wiki.octave.org/Fem-fenics}
\end{center}
while if you want to see how the project has grown during the time you can give a look at
\begin{center}
\url{http://gedeone-gsoc.blogspot.com/}
\end{center}
Finally, the API is available as Appendix but also at the following address
\begin{center}
\url{http://octave.sourceforge.net/fem-fenics/overview.html}
\end{center}
\chapter{Introduction to Fem-fenics}\label{intr}
\section{Installation}
Fem-fenics is an external package for Octave, which means that it can be installed only once that Octave has been
successfully installed on the PC. Furthermore, as Fem-fenics is based on Fenics,
it is also needed a running version of the latter. They can be easily installed following the guidelines provided
on the official Octave \cite{instoctave} and Fenics \cite{instfenics} websites.
Once that Octave and Fenics are correctly installed, to install Fem-fenics open Octave (which now is provided with a new
amazing GUI) and type
\begin{verbatim}
>> pkg install fem-fenics -forge
\end{verbatim}
That's all! For any problem during the installation don't hesitate to contact us.
To be sure that everything is working fine, load the fem-fenics pkg and run
one of the examples provided within the package:
\begin{verbatim}
>> pkg load fem-fenics
>> femfenics_examples()
\end{verbatim}
For a description of the examples, look at chapter \ref{exem}.
\paragraph*{NOTE} For completing the installation process successfully,
the form compiler FFC and the header file dolfin.h should also be available on the machine.
They are managed automatically by Fenics if it is installed as a binary package or with Dorsal.
If it has been done manually, please be sure that they are available before starting the
installation of Fem-fenics.
\section{General layout and first example}\label{genlayout}
A generic problem has to be solved in two steps:
\begin{enumerate}
\item a \textbf{.ufl file} where the abstract problem is described: this file has to be written in Unified Form Language (UFL),
which is a domain specific language for defining discrete variational forms and functionals in a notation
close to pen-and-paper formulation. UFL is easy to learn, and the User manual provides explanations
and examples \cite{ufl}.
\item a script file \textbf{.m} where the abstract problem is imported and a specific problem is implemented and solved:
this is the script file where the fem-fenics functions described in the following chapters are used.
\end{enumerate}
We provide immediately a simple example in order to familiarize the user with the code.
\paragraph{The Poisson equation}
In this example, we show how it is possible to solve the Poisson equation with mixed Boundary Conditions.
If we indicate with $\Omega$ the domain and with $\Gamma = \Gamma_{N} \cup \Gamma_{D}$ the
boundaries, the problem can be expressed as
\begin{align*}
\Delta u &= f \qquad \text{on } \Omega \\
u &= 0 \qquad \text{on } \Gamma_{D} \\
\nabla u \cdot n &= g \qquad \text{on } \Gamma_{N}
\end{align*}
where $f, \, g$ are data which represent the source and the flux
of the scalar variable $u$.
A possible variational formulation of the problem is: \\
find $u \in H_{0, \Gamma_{D}}^{1} :$
\begin{align*}
a(u, v) &= L(v) \qquad \forall v \in H_{0, \Gamma_{D}}^{1} \\
a(u, v) &= \int_{\Omega} \nabla u \cdot \nabla v \\
L(v) &= \int_{\Omega} f v + \int_{\Gamma_{N}} g v \\
\end{align*}
The abstract problem can thus be written in the \verb|Poisson.ufl| file immediately.
The only thing that has to be specified at this stage is the space of Finite Elements
used for the discretization of $H_{0, \Gamma_{D}}^{1}$. In this example,
we choose the space of continuous lagrangian polynomial of degree one
\begin{lstlisting}[numbers=none]
FiniteElement("Lagrange", triangle, 1)
\end{lstlisting}
but many more possibilities are available.
\subparagraph{Poisson.ufl}
\begin{lstlisting}
element = FiniteElement("Lagrange", triangle, 1)
u = TrialFunction(element)
v = TestFunction(element)
f = Coefficient(element)
g = Coefficient(element)
a = inner(grad(u), grad(v))*dx
L = f*v*dx + g*v*ds
\end{lstlisting}
It is always a good idea to check if the ufl code is correctly written before importing it into Octave. Typing
\begin{lstlisting}[numbers=none, language = Octave]
>> ffc -l dolfin Poisson.ufl
\end{lstlisting}
in the shell shouldn't produce any error.
We can now implement and solve a specific instance of the Poisson problem with Octave.
The parameters are setted as follow
\begin{itemize}
\item $\Omega = [0, 1]\times[0, 1]$
\item $\Gamma_{D} = {(0, y) \cup(1, y)} \ \subset \partial\Omega$
\item $\Gamma_{N} = {(x, 0) \cup(x, 1)} \ \subset \partial\Omega$
\item $f = 10 \exp \dfrac{(x-0.5)^{2} + (y-0.5)^{2}}{0.02}$
\item $g = \sin(5x)$
\end{itemize}
As a first thing we need to load into Octave the pkgs previously installed
\begin{lstlisting}[numbers=none, language = Octave]
pkg load fem-fenics msh
\end{lstlisting}
The ufl file can thus be imported inside Octave. For every specific element defined inside the ufl file
there is a specific function which stores it for later use
\begin{itemize}
\item \verb$ufl_import_FunctionSpace ('Poisson')$ is a function which looks for the finite element
space defined inside the file called Poisson.ufl; if everything is ok, it generates a function
which we will use later
\item \verb$ufl_import_BilinearForm ('Poisson')$ is a function which looks for the rhs of the
equation, i.e. for the bilinear form defined inside Poisson.ufl
\item \verb$ufl_import_LinearForm ('Poisson')$ is a function which looks for the linear
form.
\end{itemize}
In some cases one could be interested in using these functions separately but if,
as in our example, all the three elements are defined in the same ufl file (and only in this case),
the \verb$import_ufl_Problem ('Poisson')$ can be used, which generates at once all
the three functions described above
\begin{lstlisting}[numbers=none, language = Octave]
ufl_import_Problem ('Poisson');
\end{lstlisting}
To set the concrete elements which define the problem,
the first things to do is to create a mesh.
It can be managed easily using the msh pkg. For a structured squared mesh
\begin{lstlisting}[numbers=none, language = Octave]
x = y = linspace (0, 1, 33);
msho = msh2m_structured_mesh (x, y, 1, 1:4);
\end{lstlisting}
Once that the mesh is available, we can thus initialize the
Fem-fenics mesh using the function \verb$Mesh ()$:
\begin{lstlisting}[numbers=none, language = Octave]
mesh = Mesh (msho);
\end{lstlisting}
To initialize the functional space, we have to specify as argument only the fem-fenics mesh,
because the finite element type and the polynomial degree have yet been specified in the ufl file:
\begin{lstlisting}[numbers=none, language = Octave]
V = FunctionSpace('Poisson', mesh);
\end{lstlisting}
Essential BC can now be applied using \verb$DirichletBC ()$; this function receives as argument the functional space,
a function handle which specifies the value to set, and the label of the sides where the BC applies.
In this case, homogenous boundary conditions hold on the left and right side of the square
\begin{lstlisting}[numbers=none, language = Octave]
bc = DirichletBC(V, @(x, y) 0.0, [2; 4]);
\end{lstlisting}
The last thing to do before solving the problem, is to set the coefficients specified
in the ufl file.
To set them, the function \verb$Expression ()$ can be used passing as argument a string
which specifies the name of the coefficient
(it is important that they are called in the same way as in the ufl file:
the source term 'f' and the normal flux 'g'),
and a function handle with the value prescribed:
\begin{lstlisting}[numbers=none, language = Octave]
ff = Expression ('f',
@(x,y) 10*exp(-((x - 0.5)^2 + (y - 0.5)^2) / 0.02));
gg = Expression ('g', @(x,y) sin (5.0 * x));
\end{lstlisting}
Another possibility for dealing with the coefficients defined in the ufl file would have been to use
the function \verb$Constant ()$ or \verb$Function ()$.
The coefficients can thus be used together with the FunctionSpace to set
the Bilinear and the Linear form
\begin{lstlisting}[numbers=none, language = Octave]
a = BilinearForm ('Poisson', V, V);
L = LinearForm ('Poisson', V, ff, gg);
\end{lstlisting}
The discretized representation of our operator is obtained using the
functions \verb$assemble ()$ or \verb$assemble_system ()$, which also allow
to specify the BC(s) to apply
\begin{lstlisting}[numbers=none, language = Octave]
[A, b] = assemble_system (a, L, bc);
\end{lstlisting}
Here A is a sparse matrix and b is a column vector. All the
functionalities available within Octave can now be exploited to solve the linear system.
The easisest possibility is the backslash command:
\begin{lstlisting}[numbers=none, language = Octave]
u = A \ b;
\end{lstlisting}
Once that the solution has been obtained, the \verb$u$ vector is converted into a
Fem-fenics function and plotted \verb$plot ()$ or saved \verb$save ()$ in the vtu
format
\begin{lstlisting}[numbers=none, language = Octave]
u = Function ('u', V, sol);
save (u, 'poisson')
plot (u);
\end{lstlisting}
The complete code for the Poisson problem is reported below, while
in figure \ref{Poissonfig} the output is presented.
\begin{figure}
\begin{center}
\includegraphics[height=7 cm,keepaspectratio=true]{./Fem-fenics_poisson.png}
\caption{The result for the Poisson equation}
\label{Poissonfig}
\end{center}
\end{figure}
\subparagraph{Poisson.m}
\begin{lstlisting}[language=Octave]
#load the pkg and import the ufl problem
pkg load fem-fenics msh
import_ufl_Problem ('Poisson')
# Create the mesh and define function space
x = y = linspace (0, 1, 33);
mesh = Mesh(msh2m_structured_mesh (x, y, 1, 1:4));
V = FunctionSpace('Poisson', mesh);
# Define boundary condition and source term
bc = DirichletBC(V, @(x, y) 0.0, [2;4]);
ff = Expression ('f', @(x,y) 10*exp(-((x - 0.5)^2 + (y - 0.5)^2) / 0.02));
gg = Expression ('g', @(x,y) sin (5.0 * x));
#Create the Bilinear and the Linear form
a = BilinearForm ('Poisson', V, V);
L = LinearForm ('Poisson', V, ff, gg);
#Extract the matrix and compute the solution
[A, b] = assemble_system (a, L, bc);
sol = A \ b;
u = Function ('u', V, sol);
# Save solution in VTK format and plot it
save (u, 'poisson')
plot (u);
\end{lstlisting}
\chapter{Implementation}\label{impl}
Fem-fenics aims to fill a gap in Octave: even if there are packages for the creation of mesh \cite{msh},
for the postprocessing of data \cite{fpl} and for the resolution of some specific pde \cite{secs1d} \cite{bim},
no general purpose finite element library is available.
The goal of the project is thus to provide a package which can be used to solve user defined problems
and which is able to exploit the functionality provided with Octave.
\begin{figure}
\begin{center}
\includegraphics[height=10 cm,keepaspectratio=true]{./code_layout.png}
\caption{General layout of the package}
\label{Codelayout}
\end{center}
\end{figure}
Instead of writing a library from scratch, an interface to one of the finite element library
which are already available has been created.
Among the many libraries taken into account, the one whch was best suited for our
purposes seemed to be the FEniCS project. It ``is a collection of free, open source, software
components with the common goal to enable automated solution of pde.''
In particular, Dolfin is the C++/Python interface of FEniCS, providing a consistent Problem
Solving Environment for ODE and PDE. The idea has been to create wrappers in Octave for C++ Dolfin,
in a similar way to what it has been done for Python.
This is a very natural choice, because Octave is mainly written in script language
and in C++. It is in fact possible to implement an Octave interpreter function in C++ through the
native oct-file interface or, conversely, to use Octave's Matrix/Array Classes in a C++ application
\cite{whatoctave}.
The works can be summirized as follows (fig. \ref{Codelayout}):
the elements already available in Octave for the resolution of PDE (Mesh and Linear
Algebra) have been exploited, and wrappers to the other FEniCS functions added.
To allow exchanges between this programs, the necessary functions
for converting an Octave mesh/matrix into a FEniCS one and viceversa have been written.
Two main ideas have guided us throughout the realization of the pkg:
\begin{itemize}
\item keep the syntax as close as possible to the original one in Fenics (Python)
\item make the interface as simple as possible.
\end{itemize}
\section{General layout of a class}
Seven new classes are implemented for dealing with FEniCS objects and for using them inside Octave:
\begin{itemize}
\item \textbf{boundarycondition} stores and builds a dolfin::DirichletBC
\item \textbf{coefficient} stores an expression object which is used for the
evaluation of user defined values
\item \textbf{expression} is needed for internal use only as explained below
\item \textbf{form} stores a general dolfin::Form and can be used either for
a dolfn::BilinearForm as well as for a dolfin::LinearForm
\item \textbf{function} for the dolfin::Function objects
\item \textbf{functionspace} stores the user defined FunctionSpace
\item \textbf{mesh} converts a PDE-tool like mesh structure in a dolfin::Mesh
\end{itemize}
The classes are written with the ``usual'' C++ style, but they need to be derived publicly
from octave\_base\_value and to be added to the Octave interpreter \cite{whatoctave}.
When a type is used for the first time during a session, it is also temporarily
registered in the interpreter after all the other basic types (int, double, ...).
The general layout of a class can thus be kept simple and with the main purpose of storing
the associated FEniCS objects, which is done throughout
boost::shared\_ptr< > to the corresponding FEniCS type.
All the classes also implement at least two constructor:
a default constructor which is necessary to register a type in the Octave interpreter,
and a constructor which takes as argument the corresponding dolfin type.
As an example, the form class implementation follows, while classes which differs from the general
layout are presented below in more details.
\begin{lstlisting}
#ifndef _FORM_OCTAVE_
#define _FORM_OCTAVE_
#include <memory>
#include <vector>
#include <dolfin.h>
#include <octave/oct.h>
class form : public octave_base_value
{
public:
form () : octave_base_value () {}
form (const dolfin::Form _frm)
: octave_base_value (), frm (new dolfin::Form (_frm)) {}
form (boost::shared_ptr <const dolfin::Form> _frm)
: octave_base_value (), frm (_frm) {}
void
print (std::ostream& os, bool pr_as_read_syntax = false) const
{
os << "Form " << ": is a form of rank " << frm->rank ()
<< " with " << frm->num_coefficients ()
<< " coefficients" << std::endl;
}
~form(void) {}
bool is_defined (void) const { return true; }
const dolfin::Form & get_form (void) const { return (*frm); }
const boost::shared_ptr <const dolfin::Form> &
get_pform (void) const { return frm; }
private:
boost::shared_ptr <const dolfin::Form> frm;
DECLARE_OCTAVE_ALLOCATOR;
DECLARE_OV_TYPEID_FUNCTIONS_AND_DATA;
};
static bool form_type_loaded = false;
DEFINE_OCTAVE_ALLOCATOR (form);
DEFINE_OV_TYPEID_FUNCTIONS_AND_DATA (form, "form", "form");
#endif
\end{lstlisting}
\subsection{Shared pointer}
In all the classes presented above, the private members are stored using
a boost::shared\_ptr< > to the corresponding FEniCS type.
This is done because we have to refer in several places to resources which are built dynamically
and we want that it is destroyed only
when the last references is destroyed \cite{Formaggia_smart}.
For example, if we have two different functional spaces in the same problem,
like with Navier-Stokes for the velocity and the pressure, the mesh is shared between
them and no one has its own copy.
Furthermore, they are widely supported inside DOLFIN, and it can thus be avoided to have a copy of the
same object for FEniCS and another one for DOLFIN: there is just one copy which is shared among DOLFIN
and FEniCS.
\subsection{The mesh class}
In addition to usual methods, the mesh class implemens functionalities which allow to deal with meshes
as they are currently available with the msh pkg, i.e. in the (p, e, t) format, and in Fenics, i.e.
in the xml Dolfin format.
It is therefore necessary to have two different constructors
\begin{lstlisting}[numbers=none]
mesh (Array<double>& p, Array<octave_idx_type>& e,
Array<octave_idx_type>& t);
mesh (std::string _filename)
: octave_base_value (), pmsh (new dolfin::Mesh(_filename)) {}
\end{lstlisting}
where the first one accepts as input a mesh in (p, e, t) format and convert it into a xml one, while
the latter load the mesh stored in the \_filename.xml file.
The constructor are used within the Mesh () function, which therefore accepts as argument
either a mesh generated with the msh pkg or a string with the name of the file
where the dolfin mesh is stored.
Furthermore, if a mesh is stored in another different format,
the program dolfin-convert can try to convert it to the dolfin xml format.
For example, for a mesh generated with Metis:
\begin{lstlisting}[numbers=none, language=bash]
Shell:
>> dolfin-convert msh.gra msh.xml
\end{lstlisting}
and then inside the Octave script:
\begin{lstlisting}[numbers=none, language=Octave]
mesh = Mesh ('msh.xml');
\end{lstlisting}
Before exploring the code in more details, the main differences between the two storing formats are presented using
the very simple, but rather instructive, example of a unit square mesh with just two elements, fig. \ref{mesh}.
\paragraph{pet}
\begin{figure}
\begin{center}
\includegraphics[height=5 cm,keepaspectratio=true]{./mesh_1.png}
\caption{The (very) simple mesh for our example}
\label{mesh}
\end{center}
\end{figure}
A mesh is represented using the three matrices $p$, $e$, $t$, and, using msh,
we can easily obtain the mesh for our example typing
\begin{lstlisting}[numbers=none, language=octave]
mesh = msh2m_structured_mesh ([0 1], [0 1], 1, [11 12 12 13])
\end{lstlisting}
The matrix $p$ stores information about the coordinates of the vertices
\begin{mdframed}[backgroundcolor=LightBlue, outerlinewidth=0.25pt,linecolor=LightBlue]
\begin{lstlisting}[numbers=none, language=octave]
>> mesh.p
\end{lstlisting}
$
\begin{array}{rrrrl}
\; 0 & 0 & 1 & 1 & \quad \text{x-coordinates} \\
\; 0 & 1 & 0 & 1 & \quad \text{y-coordinates} \\
\end{array}
$
\end{mdframed}
Thus the vertex in the $n^{th}$ column is labelled as the vertex number $n$, and so on.
The matrix $t$ stores information about the connectivity
\begin{mdframed}[backgroundcolor=LightBlue, outerlinewidth=0.25pt,linecolor=LightBlue]
\begin{lstlisting}[numbers=none, language=octave]
>> mesh.t
\end{lstlisting}
$
\begin{array}{rrl}
\; 1 & 1 & \quad \text{number of the first vertex of the element} \\
\; 3 & 4 & \quad \text{number of the second vertex of the element} \\
\; 4 & 2 & \quad \text{number of the third vertex of the element} \\
\; 0 & 0 & \\
\end{array}
$
\end{mdframed}
The first element is thus the one obtained connecting vertices 1-3-4 and so on.
The matrix $e$ stores information related to every side edge,
like the number of the vertices of the boundary elements,
and the number of the geometrical border containing the edge,
which is a convinient way to trate boundary conditions in a problem.
\begin{mdframed}[backgroundcolor=LightBlue, outerlinewidth=0.25pt,linecolor=LightBlue]
\begin{lstlisting}[numbers=none, language=octave]
>> mesh.e
\end{lstlisting}
$
\begin{array}{rrrrl}
\; 1 & 3 & 2 & 1 & \quad \text{first vertex of the side edge} \\
\; 3 & 4 & 4 & 2 & \quad \text{second vertex of the side edge} \\
\; 0 & 0 & 0 & 0 & \text{} \\
\; 0 & 0 & 0 & 0 & \text{} \\
\; 11 & 12 & 12 & 13 & \quad \text{label of the geometrical border containing the edge} \\
\; 0 & 0 & 0 & 0 & \text{}\\
\; 1 & 1 & 1 & 1 & \text{} \\
\end{array}
$
\end{mdframed}
The side edge between vertex 1-3 is labelled 11, between 3-4 is 12...
\paragraph{dolfin xml} A mesh is an object of the dolfin::Mesh class which stores information only about
the coordinates of the vertices (like $p$) and the information about the connectivity (like $t$).
A mesh can thus be manipulated using the functions and the methods of the class, which are presented below.
Instead, the information about boundaries is not directly stored in the mesh.
The mesh used in the example is stored as
\begin{lstlisting}[numbers=none, language=xml]
<?xml version="1.0"?>
<dolfin xmlns:dolfin="http://fenicsproject.org">
<mesh celltype="triangle" dim="2">
<vertices size="4">
<vertex index="0" x="0.000e+00" y="0.000e+00" />
<vertex index="1" x="0.000e+00" y="1.000e+00" />
<vertex index="2" x="1.000e+00" y="0.000e+00" />
<vertex index="3" x="1.000e+00" y="1.000e+00" />
</vertices>
<cells size="2">
<triangle index="0" v0="0" v1="2" v2="3" />
<triangle index="1" v0="0" v1="1" v2="3" />
</cells>
</mesh>
</dolfin>
\end{lstlisting}
\paragraph{Conversion between the formats}
The first necessary step in our way to a package which links Octave and FEniCS
is to convert a mesh from the $(p, e, t)$ format
into the dolfin xml one.
Furthermore, as dolfin provides methods and functions which
allow to manipulate a mesh and which don't have a conterpart in the msh pkg,
we have also created wrappers for them (specifically for mesh::refine).
As it has been showed above, the main difference between $(p, e, t)$ and $dolfin xml$
is the way in which the boundaries are distinguished.
The former stores all the information in the $e$ matrix, while the latter uses
the functions and the methods of the dolfin::mesh class to set/get informations about a mesh.
The most useful classes available in dolfin are recalled
\begin{itemize}
\item \textbf{MeshIterator}
To know whether an edge belongs or not to the boundary, we can iterate over all the
edges of our mesh using the classes provided by dolfin:
\begin{lstlisting}[numbers=none]
for (dolfin::FacetIterator f (mesh); ! f.end (); ++f)
{
if ((*f).exterior () == true)
{
//do something with the boundary cells
}
}
\end{lstlisting}
\item \textbf{MeshFunction}
To store data related to a mesh, dolfin provides the template class MeshFunctions.
"A MeshFunction is a function that can be evaluated at a set of mesh entities.
A MeshFunction is discrete and is only defined at the set of mesh entities of a fixed topological dimension.
A MeshFunction may for example be used to store a global numbering scheme for the entities of a (parallel) mesh,
marking sub domains or boolean markers for mesh refinement." \cite{meshfunction}
For example, in the function \verb$mshm_refine$ of the msh package, the list of cells to be refined
is stored as a MeshFunction, which for every cell says whether or not it has to be refined:
\begin{lstlisting}[numbers=none]
dolfin::CellFunction<bool> cell_markers (mesh);
cell_markers.set_all (false);
for (octave_idx_type i = 0;
i < cells_to_refine.length (); ++i)
cell_markers.set_value (cells_to_refine (i) , true);
\end{lstlisting}
\item \textbf{MeshValueCollection} "It differs from the MeshFunction class in two ways.
First, data does not need to be associated with all entities (only a subset).
Second, data is associated with entities through the corresponding cell index and local entity number
(relative to the cell), not by global entity index, which means that data may be stored robustly to file."\cite{meshvalue}
It is thus obvious that it is better to use the MeshValueCollection whenever saving or writing a mesh.
\end{itemize}
The containers classes presented above can be used by their own,
but to set/get data from a mesh it is better to use the methods provided by the classes:
\begin{itemize}
\item \textbf{MeshDomains} "The class MeshDomains stores the division of a Mesh into subdomains.
For each topological dimension 0 <= d <= D, where D is the topological dimension of the Mesh,
a set of integer markers are stored for a subset of the entities of dimension d,
indicating for each entity in the subset the number of the subdomain.
It should be noted that the subset does not need to contain all entities of any given dimension;
entities not contained in the subset are ���unmarked���." \cite{meshdomain}
\item \textbf{MeshData} "The class MeshData is a container for auxiliary mesh data,
represented either as MeshFunction over topological mesh entities, arrays or maps.
Each dataset is identified by a unique user-specified string." \cite{meshdata}
\end{itemize}
\subparagraph{Geometry from (p, e, t) to xml dolfin}
Converting the vertices and cells from (p, e, t) in the xml format can be done using the
dolfin editor, while caution has to be taken for storing information associated with boundaries
and subdomains, as presented in the next paragraph.
\begin{lstlisting}[numbers=none]
dolfin::MeshEditor editor;
boost::shared_ptr<dolfin::Mesh> msh (new dolfin::Mesh ());
editor.open (*msh, D, D);
editor.init_vertices (p.cols ());
editor.init_cells (t.cols ());
if (D == 2)
{
for (uint i = 0; i < p.cols (); ++i)
editor.add_vertex (i,
p.xelem (0, i),
p.xelem (1, i));
for (uint i = 0; i < t.cols (); ++i)
editor.add_cell (i,
t.xelem (0, i) - 1,
t.xelem (1, i) - 1,
t.xelem (2, i) - 1);
}
if (D == 3)
{
...
}
editor.close ();
\end{lstlisting}
\subparagraph{Subdomain markers: from (p, e, t) to dolfin xml}
There are no fundamental differences between the 2D and 3D case,
and they are thus treated together referring to the
general dimension D.
The subdomain information is contained in the t matrix,
and is temporarily copyed to a MeshValueCollection.
For every column of the t matrix, i.e. for every element of the mesh,
we have to look for the corresponding element in the DOLFIN mesh.
We use the class MeshIterator for moving around on the DOLFIN mesh:
\begin{lstlisting}[numbers=none]
dolfin::MeshValueCollection<uint> my_cell_marker (D);
for (uint i = 0; i < num_cells; ++i)
dolfin::Vertex v (mesh, t(0, i));
for (dolfin::CellIterator f (v); ! f.end (); ++f)
{
if ((*f) == all_vertices_in_the_ith_column)
{
my_cell_marker.set_value
((*f).index (), t(last_row, i), mesh);
break;
}
}
\end{lstlisting}
The \verb$all_vertices_in_the_ith_column$ is just like a pseudo code:
we have to be sure that the Cell pointed by f is the one corresponding
to the $i^{th}$ column of the matrix, checking one-by-one the vertices:
in $2D$ the cell is a triangle, and we thus have to check $3$ vertices.
As we don't know the order in which vertices are visited, we have to check all the $3! = 6$
different combinations:
\begin{lstlisting}[numbers=none]
...
if ((*f).entities(0)[0] == t(0, i)
&& (*f).entities(0)[1] == t(1, i)
&& (*f).entities(0)[2] == t(2, i)
|| ... check the other 5 possibilities... )
....
\end{lstlisting}
where the \verb$entities(std::size_t dim)$ method returns an array with the indexes of the elements
of dimension dim. Thus we use $dim = 0$ as we are looking for vertices.
In the $3D$ case, our cell is a tetrahedron, and we have to check all the $4! = 24$ possibilities,
each of which is composed by $4$ assertions; in total we have almost one hundred conditions!
Now that the information is stored in our function, it can be associated to the mesh
\begin{lstlisting}[numbers=none]
*(mesh.domains ().markers (D)) = my_marked_cell;
\end{lstlisting}
\subparagraph{Subdomain markers: from dolfin xml to (p, e, t)}
In the DOLFIN .xml file, the information is stored like:
\begin{lstlisting}[numbers=none]
...
<mesh_value_collection name="m" type="uint" dim="2" size="2">
<value cell_index="0" local_entity="0" value="1"/>
<value cell_index="1" local_entity="0" value="2"/>
...
\end{lstlisting}
When the file is read using DOLFIN, the information is automatically associated with the mesh
as a MeshValueCollection named \verb$cell_domains$, which can be accessed to extract the information
using the MeshDomains class.
Obviously we have to be sure that the information is available within
the file that we are reading, and that it is related to Cell, i.e. to elements of dimension D,
before it is associated with the last row of the t matrix:
\begin{lstlisting}[numbers=none]
dolfin::MeshFunction<uint> my_cell_marker;
if (! mesh.domains ().is_empty ())
if (mesh.domains ().num_marked (D) != 0)
my_cell_marker = *(mesh.domains ().cell_domains ());
for (j = 0; j < t.cols (); ++j)
t(D + 1, j) = my_cell_marker[j];
\end{lstlisting}
\subparagraph{Boundary Markers}
For boundary markers, things work in a similar way, as long as we remember that we are working with objects of
dimension D - 1.
In this case, the main difference is in the .xml file: it is no longer enough to say
to what cell element the label is referred to, but we have to specify to which $D - 1$
entity (a side or a face) the label is referred.
For example:
\begin{lstlisting}[numbers=none]
....
mesh_value_collection name="m" type="uint" dim="1" size="4">
<value cell_index="0" local_entity="0" value="12"/>
<value cell_index="0" local_entity="2" value="11"/>
<value cell_index="1" local_entity="0" value="12"/>
<value cell_index="1" local_entity="2" value="13"/>
...
\end{lstlisting}
The cell number $"0"$ is a triangle,
and to the \verb$local_entity$ number $"0"$, i.e. to the side number $"0"$,
is associated the label $"12"$, while to the side number $"2"$ is associated the label $"11"$.
To the side number $"1"$, there are no labels associated.
The number of the \verb$local_entity$ refers to the enumeration of the reference element.
In any case, it is DOLFIN which takes care of the conversion of indeces from this format to the usual one,
and we can thus use methods and functions as explained for the subdomain markers.
\subparagraph{Mesh refine}
Now that it is possible to convert meshes between Octave and DOLFIN,
the functions availables in the dolfin::mesh class can be used to improve the
functionality of the msh package.
For the moment, it has been added the possibility of refining a mesh,
either uniformly or specifying the list of the vertices we want to be refined.
The function is now part of the msh pkg\cite{msh}, and a more detailed desciption has been
provided previously \cite{refine}.
\subsection{The functionspace class}
A dolfin::FunctionSpace is defined by specifying a mesh and the type of the finite element which we want to use.
The mesh is handled as presented above, while the FE are specified inside the .ufl file. Possible choices are
\cite{logg2012automated}:
\begin{center}
\begin{tabular}{ l | l }
\hline
\textbf{Finite Element Space} & \textbf{Symbol} \\ \hline \hline
Argyris & ARG * \\ \hline
Arnold���Winther & AW * \\ \hline
Brezzi���Douglas���Marini & BDM\\ \hline
Crouzeix���Raviart & CR\\ \hline
Discontinuous Lagrange & DG\\ \hline
Hermite & HER*\\ \hline
Lagrange & CG\\ \hline
Mardal���Tai���Winther & MTW *\\ \hline
Morley & MOR*\\ \hline
N��d��lec 1st kind H (curl) & N1curl\\ \hline
N��d��lec 2nd kind H (curl) & N2curl\\ \hline
Raviart���Thomas & RT\\
\end{tabular}
\end{center}
where the Finite Elements denoted with * are not yet fully supported inside FEniCS.
\section{General layout of a function}
There are three general kind of functions in the code: functions which create an abstract problem
(wrappers to UFL),
functions which create the specific instance of a problem (wrapper to FEniCS) and functions
which discretize the problem and genertates the matrices.
\section{Wrappers to UFL}
As stated in section \ref{genlayout}, a problem is divided in two files:a \textit{.ufl} file
where the abstract problem is described in Unified Form Language (UFL),
and a script file \textit{.m} where a specific problem is implemented and solved.
We suppose that they are called Poisson.ufl and Poisson.m .
In order to use the information stored in the UFL file, i.e. the bilinear and the linear form,
they have to be imported inside Octave.
When the UFL file is compiled using the ffc compiler, a header file \textit{Poisson.h} is generated.
In this header file, it is defined the Poisson class, which derives from dolfin::Form,
and the constructor for the bilinear and linear form are setted.
This file is thus available only at compilation time, but it has to be included somehow
in the wrapper function for the Bilinear and the Linear form.
An easy solution would have been to write a set of pre established problems where the user could only
change the values of the coefficient for a specific problem;
but, as we want to let the user free to write its own
variational problem, a different approach has been adopted.
The ufl file is compiled at run time and generates an header file.
Then, a Poisson.cc file is written from a template which take as input the name
of the header file and is compiled including the Poisson.h file;
now the corresponding octave functions for the specific problem is available and will be used from
BilinearForm, LinearForm, FunctionSpace, ... .
As an example it is presented the import\_ufl\_BilinearForm function.
\begin{lstlisting}[language=Octave]
function import_ufl_BilinearForm (var_prob)
...
%the function which writes the var-prob.cc file
generate_rhs (var_prob);
%the function which writes the makefile
generate_makefile (var_prob, private);
% the makefile is executed in a terminal:
% 1) generate the header file from ufl
% ffc -l dolfin var_prob.ufl
% 2) compile the var_prob.cc
% mkoctfile var_prob.cc -I.
system (sprintf ("make -f Makefile_%s rhs", var_prob));
...
endfunction
\end{lstlisting}
\begin{lstlisting}[language=Octave]
function output = generate_rhs (ufl_name)
STRING ="
#include "@@UFL_NAME@@.h"
...
DEFUN_DLD (@@UFL_NAME@@_BilinearForm, args, , ""A = fem_rhs_@@UFL_NAME@@ (FUNCTIONAL SPACE, COEFF)"")
{
...
const functionspace & fspo1
= static_cast<const functionspace&> (args(0).get_rep ());
const functionspace & fspo2
= static_cast<const functionspace&> (args(1).get_rep ());
const dolfin::FunctionSpace & U = fspo1.get_fsp ();
const dolfin::FunctionSpace & V = fspo2.get_fsp ();
@@UFL_NAME@@::BilinearForm a (U, V);
...
}";
STRING = strrep (STRING, "@@UFL_NAME@@", ufl_name);
fid = fopen (sprintf ("%s_BilinearForm.cc", ufl_name), 'w');
fputs (fid, STRING);
output = fclose (fid);
endfunction
\end{lstlisting}
\section{Wrappers to DOLFIN}
The general layout of a function is very simple and is composed of 4 steps which we describe using an example:
\begin{lstlisting}
DEFUN_DLD (fem_fs, args, , "initialize a fs from a mesh")
{
// 1 read data
const mesh & msho = static_cast<const mesh&> (args(0).get_rep ());
// 2 convert the data from octave to dolfin
const dolfin::Mesh & mshd = msho.get_msh ();
// 3 build the new object using dolfin
boost::shared_ptr <const dolfin::FunctionSpace> g (new Laplace::FunctionSpace (mshd));
// 4 convert the new object from dolfin to Octave and return it
octave_value retval = new functionspace(g);
return retval;
}
\end{lstlisting}
All the functions presented above follow this general structure, and thus here we present
in detail only functions which present some differences.
\subsection{Sparse Matrices}
\subsection{Polymorphism}
\subsection{DirichletBC and Coefficient}
These two functions take as input a function handle which cannot be directly evaluated by
a dolfin function to set, respectively, the value on the boundary or the value of the coefficient.
It has thus been derived from dolfin::Expression a class "expression" which has as private member
an octave function handle and which overloads the function eval(). In this way, an object of
the class expression can be initialized throughout a function handle and can be used inside dolfin because
"it is" a dolfin::Expression
\begin{lstlisting}
class expression : public dolfin::Expression
{
...
void
eval (dolfin::Array<double>& values,
const dolfin::Array<double>& x) const
{
octave_value_list b;
b.resize (x.size ());
for (std::size_t i = 0; i < x.size (); ++i)
b(i) = x[i];
octave_value_list tmp = feval (f->function_value (), b);
Array<double> res = tmp(0).array_value ();
for (std::size_t i = 0; i < values.size (); ++i)
values[i] = res(i);
}
private:
octave_fcn_handle * f;
};
\end{lstlisting}
\paragraph{DirichletBC}
The BC are imposed directly to the mesh setting to zero all the off diagonal elements
in the corresponding line. This means that we could loose the symmetry of the matrix, if any.
To avoid this problem, instead of the method apply() it is possible to use the
function \verb$assemble_system()$ , which preserves the symmetry of the system but which needs to build
together the lhs and the rhs.
\paragraph{Coefficient}
The coefficient of the variational problem can be specified using either a Coefficient
or a Function. They are different objects which behave in different ways: a Coefficient, as exlained above,
overloads the eval() method of the dolfin::Expression class and it is evaluated at
run time using the octave function feval (). A Function instead doesn't need to be evaluated
because it is assembled copying element-by-element the values contained in the input vector.
\iffalse
\paragraph{other function}
SubSpace allows to extract a subspace from a vectorial one.
For example, if our space is P2 x P0 we can extract the one or
the other and then apply BC only where it is necessary.
\verb$fem_eval$ takes as input a Function and a coordinate and returns a
vector representing the value of the function at this point.
for dealing with form of rank 0, i.e. with functional, we have now
added the functions \verb$fem_create_functional$ to create it from a .ufl file.
We have thus extended the function assemble which returns the corresponding double value.
\verb$plot_2d$ and \verb$plot_3d$: these functions allow us to plot a function specifying
a mesh and the value of the function at every node of the mesh.
This is something which could be useful also outside of fem-fenics.
\section{Implementation Details}
The relevant implementation details which the user should know are:
We have split the construction of the form into two steps:
We set all the coefficients of the form using the function which we create on the fly.
They will be named \verb$ProblemName_BilinearForm$ or \verb$ProblemName_LinearForm$.
Then we apply specific BC to the form using the assemble() function and we get back the matrix.
If we are assembling the whole system and we want to keep the symmetry of the matrix (if any),
we can instead use the command \verb$assemble_system$ (). Finally, if we are solving a non-linear problem
and we need to apply essential BC, we should provide to the function also the vector with the
tentative solution in order to modify the entries corresponding to the boundary values.
This will be illustrated below in the HyperElasticity example.
\fi
\section{Wrapper to FEniCS}
\subsection{Code on the fly}
\iffalse
For the creation on the fly of the code from the header file,
Juan Pablo provided me a python code which makes a great job. I have spent some
time adapting it to my problem, but when I finally got a working code, we realized
that it was probably enough to use the functions available inside Octave because
the problem was rather simple. The pyton code is however available here, while the
final solution adopted can be found there.
\fi
\chapter{More Advanced Examples}\label{exem}
\iffalse
In the following examples we can see directly in action the classes and the functions presented in the
chapters before. A comparison with DOLFIN is given only for the first example, while more extensive case can
be found online. We do not report the code for all the examples but only the relevant parts.
With the following examples, we can see directly in action the new features and understand how they work.
Navier-Stokes: we learn how to deal with a vector-field problem and how we can save the solution using the
\verb$fem_save$ () function. We also use the fem pkg to generate a mesh using gmesh.
Mixed-Poisson: we solve the Poisson problem presented in the previous posts using a mixed formulation,
and we see how we can extract a scalar field from a vector one.
HyperElasticity: we exploit the fsolve () command to solve a non-linear problem. In particular,
we see how to use the assemble() function to apply BC also in this situation.
Advection-Diffusion: we solve a time dependent problem using the lsode () command and save
the solution using the pkg flp.
For each problem, we refer the reader to the complete desciption on FEniCS or bim web-page,
while here we highlight only the implementation detail relevant for our pkg.
\fi
\section{Navier-Stokes equation with Chorin-Temam projection algorithm}
Navier-Stokes: we learn how to deal with a vector-field problem and how we can save the solution using the
\verb$fem_save$ () function. We also use the msh pkg to generate a mesh using gmesh.
\paragraph{TentativeVelocity.ufl}
\begin{lstlisting}[language=Octave]
# Copyright (C) 2010 Anders Logg
# Define function spaces (P2-P1)
V = VectorElement("CG", triangle, 2)
Q = FiniteElement("CG", triangle, 1)
# Define trial and test functions
u = TrialFunction(V)
v = TestFunction(V)
# Define coefficients
k = Constant(triangle)
u0 = Coefficient(V)
f = Coefficient(V)
nu = 0.01
# Define bilinear and linear forms
eq = (1/k)*inner(u - u0, v)*dx + inner(grad(u0)*u0, v)*dx + \
nu*inner(grad(u), grad(v))*dx - inner(f, v)*dx
a = lhs(eq)
L = rhs(eq)
\end{lstlisting}
\paragraph{PressureUpdate.ufl}
\begin{lstlisting}[language=Octave]
# Copyright (C) 2010 Anders Logg
# Define function spaces (P2-P1)
V = VectorElement("CG", triangle, 2)
Q = FiniteElement("CG", triangle, 1)
# Define trial and test functions
p = TrialFunction(Q)
q = TestFunction(Q)
# Define coefficients
k = Constant(triangle)
u1 = Coefficient(V)
# Define bilinear and linear forms
a = inner(grad(p), grad(q))*dx
L = -(1/k)*div(u1)*q*dx
\end{lstlisting}
\paragraph{VelocityUpdate.ufl}
\begin{lstlisting}[language=Octave]
# Copyright (C) 2010 Anders Logg
# Define function spaces (P2-P1)
V = VectorElement("CG", triangle, 2)
Q = FiniteElement("CG", triangle, 1)
# Define trial and test functions
u = TrialFunction(V)
v = TestFunction(V)
# Define coefficients
k = Constant(triangle)
u1 = Coefficient(V)
p1 = Coefficient(Q)
# Define bilinear and linear forms
a = inner(u, v)*dx
L = inner(u1, v)*dx - k*inner(grad(p1), v)*dx
\end{lstlisting}
\paragraph{NS.m}
\begin{lstlisting}[language=Octave]
pkg load fem-fenics msh
import_ufl_Problem ("TentativeVelocity");
import_ufl_Problem ("VelocityUpdate");
import_ufl_Problem ("PressureUpdate");
# We can either load the mesh from the file as in Dolfin but
# we can also use the msh pkg to generate the L-shape domain
L-shape-domain;
mesh = Mesh (msho);
# Define function spaces (P2-P1).
V = FunctionSpace ('VelocityUpdate', mesh);
Q = FunctionSpace ('PressureUpdate', mesh);
# Set parameter values and define coefficients
dt = 0.01;
T = 3.;
k = Constant ('k', dt);
f = Constant ('f', [0; 0]);
u0 = Expression ('u0', @(x,y) [0; 0]);
# Define boundary conditions
noslip = DirichletBC (V, @(x,y) [0; 0], [3, 4]);
outflow = DirichletBC (Q, @(x,y) 0, 2);
# Assemble matrices
a1 = BilinearForm ('TentativeVelocity', V, V, k);
a2 = BilinearForm ('PressureUpdate', Q, Q);
a3 = BilinearForm ('VelocityUpdate', V, V);
A1 = assemble (a1, noslip);
A3 = assemble (a3, noslip);
# Time-stepping
t = dt; i = 0;
while t < T
# Update pressure boundary condition
inflow = DirichletBC (Q, @(x,y) sin(3.0*t), 1);
# Compute tentative velocity step
L1 = LinearForm ('TentativeVelocity', V, k, u0, f);
b1 = assemble (L1, noslip);
utmp = A1 \ b1;
u1 = Function ('u1', V, utmp);
# Pressure correction
L2 = LinearForm ('PressureUpdate', Q, u1, k);
[A2, b2] = assemble_system (a2, L2, inflow, outflow);
ptmp = A2 \ b2;
p1 = Function ('p1', Q, ptmp);
# Velocity correction
L3 = LinearForm ('VelocityUpdate', V, k, u1, p1);
b3 = assemble (L3, noslip);
ut = A3 \ b3;
u1 = Function ('u0', V, ut);
# Save to file
save (p1, sprintf ("p_%3.3d", ++i));
save (u1, sprintf ("u_%3.3d", i));
# Move to next time step
u0 = u1;
t += dt
end
\end{lstlisting}
\paragraph{L-shape-domain.m}
\begin{lstlisting}[language=Octave]
name = [tmpnam ".geo"];
fid = fopen (name, "w");
fputs (fid,"Point (1) = {0, 0, 0, 0.1};\n");
fputs (fid,"Point (2) = {1, 0, 0, 0.1};\n");
fputs (fid,"Point (3) = {1, 0.5, 0, 0.1};\n");
fputs (fid,"Point (4) = {0.5, 0.5, 0, 0.1};\n");
fputs (fid,"Point (5) = {0.5, 1, 0, 0.1};\n");
fputs (fid,"Point (6) = {0, 1, 0,0.1};\n");
fputs (fid,"Line (1) = {5, 6};\n");
fputs (fid,"Line (2) = {2, 3};\n");
fputs (fid,"Line(3) = {6,1,2};\n");
fputs (fid,"Line(4) = {5,4,3};\n");
fputs (fid,"Line Loop(7) = {3,2,-4,1};\n");
fputs (fid,"Plane Surface(8) = {7};\n");
fclose (fid);
msho = msh2m_gmsh (canonicalize_file_name (name)(1:end-4),...
"scale", 1,"clscale", .2);
unlink (canonicalize_file_name (name));
\end{lstlisting}
\section{A penalization method to take into account obstacles in incompressible viscous flows}
\newpage
\appendix
\chapter{API reference}
\section{Import problem defined with ufl}
\subsection*{import\_ufl\_BilinearForm}
\subimport{latex/}{API/import_ufl_BilinearForm.tex}
\subsection*{import\_ufl\_LinearForm}
\subimport{latex/}{API/import_ufl_LinearForm.tex}
\subsection*{ import\_ufl\_Functional}
\subimport{latex/}{API/import_ufl_Functional.tex}
\subsection*{ import\_ufl\_FunctionSpace}
\subimport{latex/}{API/import_ufl_FunctionSpace.tex}
\subsection*{ import\_ufl\_Problem}
\subimport{latex/}{API/import_ufl_Problem.tex}
\section{Problem geometry and FE space}
\subsection*{ Mesh}
\subimport{latex/}{API/Mesh.tex}
\subsection*{ FunctionSpace}
\subimport{latex/}{API/FunctionSpace.tex}
\subsection*{ SubSpace}
\subimport{latex/}{API/SubSpace.tex}
\section{Problem variables}
\subsection*{ Constant}
\subimport{latex/}{API/Constant.tex}
\subsection*{ Expression}
\subimport{latex/}{API/Expression.tex}
\subsection*{ Function}
\subimport{latex/}{API/Function.tex}
\subsection*{ DirichletBC}
\subimport{latex/}{API/DirichletBC.tex}
\section{Definition of the abstract Variational problem}
\subsection*{ BilinearForm}
\subimport{latex/}{API/BilinearForm.tex}
\subsection*{ LinearForm}
\subimport{latex/}{API/LinearForm.tex}
\subsection*{ ResidualForm}
\subimport{latex/}{API/ResidualForm.tex}
\subsection*{ JacobianForm}
\subimport{latex/}{API/JacobianForm.tex}
\subsection*{ Functional}
\subimport{latex/}{API/Functional.tex}
\section{Creation of the discretized problem}
\subsection*{ assemble}
\subimport{latex/}{API/assemble.tex}
\subsection*{ assemble\_system}
\subimport{latex/}{API/assemble_system.tex}
\section{Post processing}
\subsection*{ @function/save}
\subimport{latex/}{API/save.tex}
\subsection*{ @function/plot}
\subimport{latex/}{API/plot.tex}
\subsection*{ @mesh/plot}
\subimport{latex/}{API/plot_m.tex}
\subsection*{ @function/feval}
\subimport{latex/}{API/feval.tex}
\chapter{Autoconf and Automake}
In this section we want to discuss how we can write a config.ac and a Makefile.in files which:
\begin{itemize}
\item check if a program is available and stop if it is not
\item check if a header file is available and issue a warning if not, but go ahead with the compilation
\end{itemize}
To reach this goal, we need two components:
\paragraph{configure.ac} Is a file which checks whether the program/header is available or not
and sets consequently the values of some variables.
\begin{lstlisting}[language=make]
# Checks if the program mkoctfile is available and sets the variable HAVE_MKOCTFILE consequently
AC_CHECK_PROG([HAVE_MKOCTFILE], [mkoctfile], [yes], [no])
# if mkoctfile is not available, it issues an error and stops the compilation
if [test $HAVE_MKOCTFILE = "no"]; then
AC_MSG_ERROR([mkoctfile required to install $PACKAGE_NAME])
fi
#Checks if the header dolfin.h is available; if it is available, the value of the ac_dolfin_cpp_flags is substituted with -DHAVE_DOLFIN_H, otherwise it is left empty and a warning message is printed
AC_CHECK_HEADER([dolfin.h],
[AC_SUBST(ac_dolfin_cpp_flags,-DHAVE_DOLFIN_H) AC_SUBST(ac_dolfin_ld_flags,-ldolfin)],
[AC_MSG_WARN([dolfin headers could not be found, some functionalities will be disabled, don't worry your package will still be working, though.])] ).
# It generates the Makefile, using the template described below
AC_CONFIG_FILES([Makefile])
\end{lstlisting}
\paragraph{Makefile.ac} This file is a template for the Makefile, which will be automatically generated when the configure.ac
file is executed. The values of the variable \verb$ac_dolfin_cpp_flags$ and \verb$ac_dolfin_ld_flags$ are substituted with the
results obtained above:
\begin{lstlisting}[language=make]
CPPFLAGS += @ac_dolfin_cpp_flags@
LDFLAGS += @ac_dolfin_ld_flags@
\end{lstlisting}
In this way, if dolfin.h is available, CPPFLAGS contains also the flag -DHAVE\_DOLFIN\_H.
\paragraph {program.cc} Our .cc program, should thus include the header dolfin.h only if
\verb$-DHAVE_DOLFIN_H$ is defined at compilation time.
For example
\begin{lstlisting}
#ifdef HAVE_DOLFIN_H
#include <dolfin.h>
#endif
int main ()
{
#ifndef HAVE_DOLFIN_H
error("program: the program was built without support for dolfin");
#else
/* Body of your function */
#endif
return 0;
}
\end{lstlisting}
\paragraph {Warning} If in the Makefile.in you write something like
\begin{lstlisting}[language=make]
HAVE_DOLFIN_H = @HAVE_DOLFIN_H@
ifdef HAVE_DOLFIN_H
CPPFLAGS += -DHAVE_DOLFIN_H
LIBS += -ldolfin
endif
\end{lstlisting}
it doesn't work because the variable \verb$HAVE_DOLFIN_H$ seems to be always defined, even if the header is not available.
\bibliographystyle{unsrt}
\bibliography{doc}
\end{document}