Mercurial > fem-fenics-eugenio
annotate inst/femfenics_examples.m @ 204:0faaf956bd3a
Fix the HyperElasticity example
author | gedeone-octave <marcovass89@hotmail.it> |
---|---|
date | Thu, 30 Jan 2014 08:55:40 +0100 |
parents | 3ae8a7398f43 |
children |
rev | line source |
---|---|
162
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
1 ## Copyright (C) 2013 Marco Vassallo <gedeone-octave@users.sourceforge.net> |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
2 ## |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
3 ## This program is free software; you can redistribute it and/or modify it under |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
4 ## the terms of the GNU General Public License as published by the Free Software |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
5 ## Foundation; either version 3 of the License, or (at your option) any later |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
6 ## version. |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
7 ## |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
8 ## This program is distributed in the hope that it will be useful, but WITHOUT |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
9 ## ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
10 ## FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
11 ## details. |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
12 ## |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
13 ## You should have received a copy of the GNU General Public License along with |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
14 ## this program; if not, see <http://www.gnu.org/licenses/>. |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
15 |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
16 %# -*- texinfo -*- |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
17 %# @deftypefn {Function File} {[@var{}] =} femfenics_examples (@var{}) |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
18 %# Open the PDE examples menu and allow the user to select a demo that will be evaluated. |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
19 %# @end deftypefn |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
20 |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
21 function [] = femfenics_examples () |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
22 |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
23 vode = 1; while (vode > 0) |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
24 clc; |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
25 fprintf (1, ... |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
26 ['FEM-FENICS examples menu:\n', ... |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
27 '==================\n', ... |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
28 '\n', ... |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
29 ' (1) Solve the Poisson equation\n', ... |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
30 ' (2) Solve the Poisson equation with a Mixed formulation\n', ... |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
31 ' (3) Solve the Elasticity problem for an HyperElastic material\n', ... |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
32 ' (4) Solve the Navier-Stokes equation with Chorin-Teman algorithm\n', ... |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
33 ' (5) Solve the Biharmonic equation\n', ... |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
34 '\n', ... |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
35 ' Note: There are further FEM-FeNICS examples available with the pkg\n', ... |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
36 ' testsuite functions.\n', ... |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
37 '\n', ... |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
38 ' If you have another interesting example that you would like\n', ... |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
39 ' to share then please modify this file, create a patch and push\n', ... |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
40 ' your patch with your added example to the Fem-Fenics repository.\n', ... |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
41 '\n' ]); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
42 vode = input ('Please choose a number from above or press <Enter> to return: '); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
43 clc; if (vode > 0 && vode < 6) |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
44 %# We can't use the function 'demo' directly here because it does |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
45 %# not allow to run other functions within a demo. |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
46 vexa = example (mfilename (), vode); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
47 disp (vexa); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
48 eval (vexa); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
49 input ('Press <Enter> to continue: '); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
50 end %# if (vode > 0) |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
51 end %# while (vode > 0) |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
52 |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
53 %!demo |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
54 %! # In this example the Poisson equation is solved using P1 element |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
55 %! # on a uniform squared mesh. |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
56 %! # Read about the Poisson equation at |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
57 %! # http://en.wikipedia.org/wiki/Poisson%27s_equation. |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
58 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
59 %! n = length (mfilename ("fullpath")) - length (mfilename()); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
60 %! path = strtrunc(mfilename ("fullpath"), n); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
61 %! private = fullfile (path, "example/Poisson/Poisson.ufl"); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
62 %! [status, msg, msgid] = copyfile (private, pwd, 'f'); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
63 %! if (status != 1) |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
64 %! error ("couldn't copy file: %s", msg); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
65 %! endif |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
66 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
67 %! pkg load fem-fenics msh |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
68 %! import_ufl_Problem ('Poisson') |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
69 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
70 %! # Create mesh and define function space |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
71 %! x = y = linspace (0, 1, 33); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
72 %! mesh = Mesh(msh2m_structured_mesh (x, y, 1, 1:4)); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
73 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
74 %! V = FunctionSpace('Poisson', mesh); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
75 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
76 %! # Define boundary condition |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
77 %! bc = DirichletBC(V, @(x, y) 0.0, [2;4]); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
78 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
79 %! f = Expression ('f', @(x,y) 10*exp(-((x - 0.5)^2 + (y - 0.5)^2) / 0.02)); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
80 %! g = Expression ('g', @(x,y) sin (5.0 * x)); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
81 %! |
165
3ae8a7398f43
Update the examples for the tutorial.
gedeone-octave <marcovass89@hotmail.it>
parents:
162
diff
changeset
|
82 %! a = BilinearForm ('Poisson', V, V); |
162
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
83 %! L = LinearForm ('Poisson', V, f, g); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
84 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
85 %! # Compute solution |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
86 %! [A, b] = assemble_system (a, L, bc); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
87 %! sol = A \ b; |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
88 %! u = Function ('u', V, sol); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
89 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
90 %! # Save solution in VTK format |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
91 %! save (u, 'poisson') |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
92 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
93 %! # Plot solution |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
94 %! plot (u); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
95 |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
96 %!demo |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
97 %! # In this example the Poisson equation is solved using a Mixed Approach. |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
98 %! # We use the stable FE space obtained using Brezzi-Douglas-Marini polynomial |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
99 %! # of order 1 and Dicontinuos element of order 0. |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
100 %! # Read about the Poisson equation at |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
101 %! # http://en.wikipedia.org/wiki/Poisson%27s_equation. |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
102 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
103 %! n = length (mfilename ("fullpath")) - length (mfilename()); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
104 %! path = strtrunc(mfilename ("fullpath"), n); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
105 %! private = fullfile (path, "example/Mixed-Poisson/MixedPoisson.ufl"); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
106 %! [status, msg, msgid] = copyfile (private, pwd, 'f'); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
107 %! if (status != 1) |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
108 %! error ("couldn't copy file: %s", msg); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
109 %! endif |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
110 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
111 %! pkg load fem-fenics msh |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
112 %! import_ufl_Problem ('MixedPoisson') |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
113 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
114 %! # Create mesh |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
115 %! x = y = linspace (0, 1, 33); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
116 %! mesh = Mesh(msh2m_structured_mesh (x, y, 1, 1:4)); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
117 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
118 %! V = FunctionSpace('MixedPoisson', mesh); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
119 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
120 %! f = Expression ('f', @(x,y) 10*exp(-((x - 0.5)^2 + (y - 0.5)^2) / 0.02)); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
121 %! |
165
3ae8a7398f43
Update the examples for the tutorial.
gedeone-octave <marcovass89@hotmail.it>
parents:
162
diff
changeset
|
122 %! a = BilinearForm ('MixedPoisson', V, V); |
162
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
123 %! L = LinearForm ('MixedPoisson', V, f); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
124 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
125 %! # Define essential boundary |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
126 %! bc1 = DirichletBC (SubSpace (V, 1), @(x,y) [0; -sin(5.0*x)], 1); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
127 %! bc2 = DirichletBC (SubSpace (V, 1), @(x,y) [0; sin(5.0*x)], 3); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
128 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
129 %! # Compute solution |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
130 %! [A, b] = assemble_system (a, L, bc1, bc2); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
131 %! sol = A \ b; |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
132 %! func = Function ('func', V, sol); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
133 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
134 %! sigma = Function ('sigma', func, 1); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
135 %! u = Function ('u', func, 2); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
136 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
137 %! # Plot solution |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
138 %! plot (sigma); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
139 %! plot (u); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
140 |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
141 %!demo |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
142 %! # In this example we solve the equation for an Hyperelastic material |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
143 %! # minimizing the potential energy of the material. |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
144 %! # Read about the Hyperelastic materials at |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
145 %! # http://en.wikipedia.org/wiki/Hyperelastic_material |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
146 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
147 %! n = length (mfilename ("fullpath")) - length (mfilename()); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
148 %! path = strtrunc(mfilename ("fullpath"), n); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
149 %! private = fullfile (path, "example/Elasticity/HyperElasticity.ufl"); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
150 %! [status, msg, msgid] = copyfile (private, pwd, 'f'); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
151 %! if (status != 1) |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
152 %! error ("couldn't copy file: %s", msg); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
153 %! endif |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
154 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
155 %! pkg load fem-fenics msh |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
156 %! problem = 'HyperElasticity'; |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
157 %! import_ufl_Problem (problem); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
158 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
159 %! Rotation = @(x,y,z) ... |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
160 %! [0; ... |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
161 %! 0.5*(0.5 + (y - 0.5)*cos(pi/3) - (z-0.5)*sin(pi/3) - y);... |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
162 %! 0.5*(0.5 + (y - 0.5)*sin(pi/3) + (z-0.5)*cos(pi/3) - z)]; |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
163 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
164 %! #Create mesh and define function space |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
165 %! x = y = z = linspace (0, 1, 17); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
166 %! mshd = Mesh (msh3m_structured_mesh (x, y, z, 1, 1:6)); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
167 %! V = FunctionSpace (problem, mshd); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
168 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
169 %! # Create Dirichlet boundary conditions |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
170 %! bcl = DirichletBC (V, @(x,y,z) [0; 0; 0], 1); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
171 %! bcr = DirichletBC (V, Rotation, 2); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
172 %! bcs = {bcl, bcr}; |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
173 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
174 %! # Define source and boundary traction functions |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
175 %! B = Constant ('B', [0.0; -0.5; 0.0]); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
176 %! T = Constant ('T', [0.1; 0.0; 0.0]); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
177 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
178 %! # Set material parameters |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
179 %! E = 10.0; |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
180 %! nu = 0.3; |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
181 %! mu = Constant ('mu', E./(2*(1+nu))); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
182 %! lmbda = Constant ('lmbda', E*nu./((1+nu)*(1-2*nu))); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
183 %! u = Expression ('u', @(x,y,z) [0; 0; 0]); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
184 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
185 %! # Create (linear) form defining (nonlinear) variational problem |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
186 %! L = ResidualForm (problem, V, mu, lmbda, B, T, u); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
187 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
188 %! # Solve nonlinear variational problem F(u; v) = 0 |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
189 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
190 %! # Create function for the resolution of the NL problem |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
191 %! function [y, jac] = f (problem, xx, V, bc1, bc2, B, T, mu, lmbda) |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
192 %! u = Function ('u', V, xx); |
204
0faaf956bd3a
Fix the HyperElasticity example
gedeone-octave <marcovass89@hotmail.it>
parents:
165
diff
changeset
|
193 %! a = JacobianForm (problem, V, V, mu, lmbda, u); |
162
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
194 %! L = ResidualForm (problem, V, mu, lmbda, B, T, u); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
195 %! if (nargout == 1) |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
196 %! [y, xx] = assemble (L, xx, bc1, bc2); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
197 %! elseif (nargout == 2) |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
198 %! [jac, y, xx] = assemble_system (a, L, xx, bc1, bc2); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
199 %! endif |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
200 %! endfunction |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
201 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
202 %! u0 = assemble (L, bcs{:}); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
203 %! fs = @(xx) f (problem, xx, V, bcl, bcr, B, T, mu, lmbda); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
204 %! [x, fval, info] = fsolve (fs, u0, optimset ("jacobian", "on")); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
205 %! func = Function ('u', V, x); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
206 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
207 %! # Save solution in VTK format |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
208 %! save (func, 'displacement'); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
209 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
210 %! # Plot solution |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
211 %! plot (func); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
212 |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
213 %!demo |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
214 %! # In this example we solve the Navier-Stokes equation in a L-shape domain. |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
215 %! # We impose a sinusoidal pressure at the inflow and a constant pressure at |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
216 %! # the output. On the other sides, no-slip BC are applied. |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
217 %! # We solve the problem with the Chorin-Temam algorithm. |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
218 %! # The pressure is discretized with P1 element while the velocity with P2. |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
219 %! # Read about the NS equation at |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
220 %! # http://en.wikipedia.org/wiki/Navier%E2%80%93Stokes_equations |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
221 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
222 %! n = length (mfilename ("fullpath")) - length (mfilename()); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
223 %! path = strtrunc(mfilename ("fullpath"), n); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
224 %! private = fullfile (path, "example/NavierStokes/*.ufl"); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
225 %! [status, msg, msgid] = copyfile (private, pwd, 'f'); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
226 %! if (status != 1) |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
227 %! error ("couldn't copy file: %s", msg); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
228 %! endif |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
229 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
230 %! pkg load fem-fenics msh |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
231 %! import_ufl_Problem ("TentativeVelocity"); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
232 %! import_ufl_Problem ("VelocityUpdate"); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
233 %! import_ufl_Problem ("PressureUpdate"); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
234 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
235 %! # We use the msh pkg to generate the L-shape domain |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
236 %! name = [tmpnam ".geo"]; |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
237 %! fid = fopen (name, "w"); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
238 %! fputs (fid,"Point (1) = {0, 0, 0, 0.1};\n"); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
239 %! fputs (fid,"Point (2) = {1, 0, 0, 0.1};\n"); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
240 %! fputs (fid,"Point (3) = {1, 0.5, 0, 0.1};\n"); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
241 %! fputs (fid,"Point (4) = {0.5, 0.5, 0, 0.1};\n"); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
242 %! fputs (fid,"Point (5) = {0.5, 1, 0, 0.1};\n"); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
243 %! fputs (fid,"Point (6) = {0, 1, 0,0.1};\n"); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
244 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
245 %! fputs (fid,"Line (1) = {5, 6};\n"); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
246 %! fputs (fid,"Line (2) = {2, 3};\n"); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
247 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
248 %! fputs (fid,"Line(3) = {6,1,2};\n"); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
249 %! fputs (fid,"Line(4) = {5,4,3};\n"); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
250 %! fputs (fid,"Line Loop(7) = {3,2,-4,1};\n"); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
251 %! fputs (fid,"Plane Surface(8) = {7};\n"); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
252 %! fclose (fid); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
253 %! msho = msh2m_gmsh (canonicalize_file_name (name)(1:end-4),"scale", 1,"clscale", .4); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
254 %! unlink (canonicalize_file_name (name)); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
255 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
256 %! mesh = Mesh (msho); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
257 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
258 %! # Define function spaces (P2-P1). From ufl file |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
259 %! V = FunctionSpace ('VelocityUpdate', mesh); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
260 %! Q = FunctionSpace ('PressureUpdate', mesh); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
261 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
262 %! # Set parameter values |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
263 %! dt = 0.01; |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
264 %! T = 3.; |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
265 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
266 %! # Define boundary conditions |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
267 %! noslip = DirichletBC (V, @(x,y) [0; 0], [3, 4]); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
268 %! outflow = DirichletBC (Q, @(x,y) 0, 2); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
269 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
270 %! # Create functions |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
271 %! u0 = Expression ('u0', @(x,y) [0; 0]); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
272 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
273 %! # Define coefficients |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
274 %! k = Constant ('k', dt); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
275 %! f = Constant ('f', [0; 0]); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
276 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
277 %! # Tentative velocity step. |
165
3ae8a7398f43
Update the examples for the tutorial.
gedeone-octave <marcovass89@hotmail.it>
parents:
162
diff
changeset
|
278 %! a1 = BilinearForm ('TentativeVelocity', V, V, k); |
162
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
279 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
280 %! # Pressure update. |
165
3ae8a7398f43
Update the examples for the tutorial.
gedeone-octave <marcovass89@hotmail.it>
parents:
162
diff
changeset
|
281 %! a2 = BilinearForm ('PressureUpdate', Q, Q); |
162
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
282 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
283 %! # Velocity update |
165
3ae8a7398f43
Update the examples for the tutorial.
gedeone-octave <marcovass89@hotmail.it>
parents:
162
diff
changeset
|
284 %! a3 = BilinearForm ('VelocityUpdate', V, V); |
162
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
285 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
286 %! # Assemble matrices |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
287 %! A1 = assemble (a1, noslip); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
288 %! A3 = assemble (a3, noslip); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
289 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
290 %! # Time-stepping |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
291 %! t = dt; i = 0; |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
292 %! while t < T |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
293 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
294 %! # Update pressure boundary condition |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
295 %! inflow = DirichletBC (Q, @(x,y) sin(3.0*t), 1); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
296 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
297 %! # Compute tentative velocity step |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
298 %! L1 = LinearForm ('TentativeVelocity', V, k, u0, f); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
299 %! b1 = assemble (L1, noslip); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
300 %! utmp = A1 \ b1; |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
301 %! u1 = Function ('u1', V, utmp); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
302 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
303 %! # Pressure correction |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
304 %! L2 = LinearForm ('PressureUpdate', Q, u1, k); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
305 %! [A2, b2] = assemble_system (a2, L2, inflow, outflow); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
306 %! ptmp = A2 \ b2; |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
307 %! p1 = Function ('p1', Q, ptmp); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
308 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
309 %! # Velocity correction |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
310 %! L3 = LinearForm ('VelocityUpdate', V, k, u1, p1); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
311 %! b3 = assemble (L3, noslip); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
312 %! ut = A3 \ b3; |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
313 %! u1 = Function ('u0', V, ut); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
314 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
315 %! # Plot solution |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
316 %! plot (p1); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
317 %! plot (u1); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
318 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
319 %! # Save to file |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
320 %! save (p1, sprintf ("p_%3.3d", ++i)); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
321 %! delete (sprintf ("p_%3.3d.pvd", i)); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
322 %! save (u1, sprintf ("u_%3.3d", i)); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
323 %! delete (sprintf ("u_%3.3d.pvd", i)); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
324 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
325 %! # Move to next time step |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
326 %! u0 = u1; |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
327 %! t += dt |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
328 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
329 %! end |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
330 |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
331 %!demo |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
332 %! # In this example the Biharmonic equation is solved using P2 element |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
333 %! # on a uniform squared mesh. |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
334 %! # Read more about the Biharmonic equation at |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
335 %! # http://en.wikipedia.org/wiki/Biharmonic |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
336 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
337 %! n = length (mfilename ("fullpath")) - length (mfilename()); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
338 %! path = strtrunc(mfilename ("fullpath"), n); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
339 %! private = fullfile (path, "example/Biharmonic/Biharmonic.ufl"); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
340 %! [status, msg, msgid] = copyfile (private, pwd, 'f'); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
341 %! if (status != 1) |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
342 %! error ("couldn't copy file: %s", msg); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
343 %! endif |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
344 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
345 %! pkg load fem-fenics msh |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
346 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
347 %! problem = 'Biharmonic'; |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
348 %! import_ufl_Problem (problem); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
349 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
350 %! # Create mesh and define function space |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
351 %! x = y = linspace (0, 1, 33); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
352 %! mesh = Mesh(msh2m_structured_mesh (x, y, 1, 1:4)); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
353 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
354 %! V = FunctionSpace(problem, mesh); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
355 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
356 %! bc = DirichletBC (V, @(x,y) 0, 1:4); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
357 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
358 %! f = Expression ('f', @(x,y) 4.0*pi^4.*sin(pi.*x).*sin(pi.*y)); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
359 %! g = Expression ('alpha', @(x,y) 8); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
360 %! |
165
3ae8a7398f43
Update the examples for the tutorial.
gedeone-octave <marcovass89@hotmail.it>
parents:
162
diff
changeset
|
361 %! a = BilinearForm (problem, V, V, g); |
162
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
362 %! L = LinearForm (problem, V, f); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
363 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
364 %! [A, b] = assemble_system (a, L, bc); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
365 %! u = A \ b; |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
366 %! |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
367 %! func = Function ('u', V, u); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
368 %! plot (func); |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
369 |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
370 %# Local Variables: *** |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
371 %# mode: octave *** |
344e9bc30b03
New file for the ecxamples associated with the pkg.
gedeone-octave <marcovass89@hotmail.it>
parents:
diff
changeset
|
372 %# End: *** |