Mercurial > fem-fenics-eugenio
annotate inst/example/Evolution/Evolution.m @ 239:cb6f11a09667
Move ufl to m-file in Evolution example
author | Eugenio Gianniti <eugenio.gianniti@mail.polimi.it> |
---|---|
date | Tue, 24 Jun 2014 16:55:49 +0200 |
parents | e85390ed620d |
children |
rev | line source |
---|---|
174 | 1 ## Copyright (C) 2013 Marco Vassallo <gedeone-octave@users.sourceforge.net> |
2 ## | |
3 ## This program is free software; you can redistribute it and/or modify it under | |
4 ## the terms of the GNU General Public License as published by the Free Software | |
5 ## Foundation; either version 3 of the License, or (at your option) any later | |
6 ## version. | |
7 ## | |
8 ## This program is distributed in the hope that it will be useful, but WITHOUT | |
9 ## ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or | |
10 ## FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more | |
11 ## details. | |
12 ## | |
13 ## You should have received a copy of the GNU General Public License along with | |
14 ## this program; if not, see <http://www.gnu.org/licenses/>. | |
15 | |
157
08dc5547f4d6
Update the examples to the new naming conventions.
gedeone-octave <marcovass89@hotmail.it>
parents:
95
diff
changeset
|
16 pkg load fem-fenics msh |
73 | 17 |
239
cb6f11a09667
Move ufl to m-file in Evolution example
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
174
diff
changeset
|
18 ufl start Evolution |
cb6f11a09667
Move ufl to m-file in Evolution example
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
174
diff
changeset
|
19 ufl 'element = FiniteElement("Lagrange", triangle, 1)' |
cb6f11a09667
Move ufl to m-file in Evolution example
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
174
diff
changeset
|
20 ufl |
cb6f11a09667
Move ufl to m-file in Evolution example
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
174
diff
changeset
|
21 ufl u = TrialFunction(element) |
cb6f11a09667
Move ufl to m-file in Evolution example
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
174
diff
changeset
|
22 ufl v = TestFunction(element) |
cb6f11a09667
Move ufl to m-file in Evolution example
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
174
diff
changeset
|
23 ufl |
cb6f11a09667
Move ufl to m-file in Evolution example
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
174
diff
changeset
|
24 ufl u0 = Coefficient(element) |
cb6f11a09667
Move ufl to m-file in Evolution example
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
174
diff
changeset
|
25 ufl |
cb6f11a09667
Move ufl to m-file in Evolution example
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
174
diff
changeset
|
26 ufl dt = Constant(triangle) |
cb6f11a09667
Move ufl to m-file in Evolution example
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
174
diff
changeset
|
27 ufl |
cb6f11a09667
Move ufl to m-file in Evolution example
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
174
diff
changeset
|
28 ufl eq = (1/dt)*(u-u0)*v*dx + "inner(grad(u), grad(v))*dx" |
cb6f11a09667
Move ufl to m-file in Evolution example
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
174
diff
changeset
|
29 ufl a = rhs (eq) |
cb6f11a09667
Move ufl to m-file in Evolution example
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
174
diff
changeset
|
30 ufl L = lhs (eq) |
cb6f11a09667
Move ufl to m-file in Evolution example
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
174
diff
changeset
|
31 ufl end |
cb6f11a09667
Move ufl to m-file in Evolution example
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
174
diff
changeset
|
32 |
73 | 33 problem = 'Evolution'; |
34 | |
157
08dc5547f4d6
Update the examples to the new naming conventions.
gedeone-octave <marcovass89@hotmail.it>
parents:
95
diff
changeset
|
35 mesh = Mesh (msh2m_structured_mesh (0:0.05:1, 0:0.05:1, 1, 1:4)); |
73 | 36 |
157
08dc5547f4d6
Update the examples to the new naming conventions.
gedeone-octave <marcovass89@hotmail.it>
parents:
95
diff
changeset
|
37 V = FunctionSpace(problem, mesh); |
73 | 38 |
89
d06b423169fa
The examples now follow the new naming convention.
gedeone-octave <marcovass89@hotmail.it>
parents:
73
diff
changeset
|
39 bc = DirichletBC (V, @(x,y) 1, 1:4); |
73 | 40 |
41 t = 0; | |
42 dt = 0.01; | |
43 T = 0.3; | |
44 | |
89
d06b423169fa
The examples now follow the new naming convention.
gedeone-octave <marcovass89@hotmail.it>
parents:
73
diff
changeset
|
45 k = Expression ('dt', @(x,y) dt); |
73 | 46 |
89
d06b423169fa
The examples now follow the new naming convention.
gedeone-octave <marcovass89@hotmail.it>
parents:
73
diff
changeset
|
47 u0 = Expression ('u0', @(x,y) 10*exp(-((x - 0.5)^2 + (y - 0.5)^2) / 0.02)); |
73 | 48 |
164
7190852cf57f
Update the examples giving two input parameters to the BilinearForm.
gedeone-octave <marcovass89@hotmail.it>
parents:
159
diff
changeset
|
49 a = BilinearForm (problem, V, V, k); |
157
08dc5547f4d6
Update the examples to the new naming conventions.
gedeone-octave <marcovass89@hotmail.it>
parents:
95
diff
changeset
|
50 L = LinearForm (problem, V, k, u0); |
94
4f688918f76a
Used new available functions assemble and assemble-systems.
gedeone-octave <marcovass89@hotmail.it>
parents:
89
diff
changeset
|
51 |
4f688918f76a
Used new available functions assemble and assemble-systems.
gedeone-octave <marcovass89@hotmail.it>
parents:
89
diff
changeset
|
52 A = assemble (a, bc); |
73 | 53 |
54 # solve the problem for each time step | |
55 # We need to update only the lhs | |
56 while t < T | |
57 t += dt; | |
58 | |
59 # we can pass u0 to the lhs indifferently as a fem_coeff or | |
60 # as a fem_func | |
157
08dc5547f4d6
Update the examples to the new naming conventions.
gedeone-octave <marcovass89@hotmail.it>
parents:
95
diff
changeset
|
61 L = LinearForm (problem, V, k, u0); |
94
4f688918f76a
Used new available functions assemble and assemble-systems.
gedeone-octave <marcovass89@hotmail.it>
parents:
89
diff
changeset
|
62 b = assemble (L, bc); |
73 | 63 |
64 u = A \ b; | |
65 | |
89
d06b423169fa
The examples now follow the new naming convention.
gedeone-octave <marcovass89@hotmail.it>
parents:
73
diff
changeset
|
66 u0 = Function ('u0', V, u); |
73 | 67 |
68 #press Q to make the plot continue | |
157
08dc5547f4d6
Update the examples to the new naming conventions.
gedeone-octave <marcovass89@hotmail.it>
parents:
95
diff
changeset
|
69 plot (u0); |
08dc5547f4d6
Update the examples to the new naming conventions.
gedeone-octave <marcovass89@hotmail.it>
parents:
95
diff
changeset
|
70 end |