Mercurial > fem-fenics-eugenio
annotate inst/example/Poisson/Poisson.m @ 238:b96f6b12f8ca
Move ufl to m-file in Poisson example
author | Eugenio Gianniti <eugenio.gianniti@mail.polimi.it> |
---|---|
date | Tue, 24 Jun 2014 16:43:13 +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 | |
145
93a4ee13aa75
Update the examples.
gedeone-octave <marcovass89@hotmail.it>
parents:
94
diff
changeset
|
16 pkg load fem-fenics msh |
238
b96f6b12f8ca
Move ufl to m-file in Poisson example
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
174
diff
changeset
|
17 |
b96f6b12f8ca
Move ufl to m-file in Poisson example
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
174
diff
changeset
|
18 ufl start Poisson |
b96f6b12f8ca
Move ufl to m-file in Poisson example
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
174
diff
changeset
|
19 ufl 'element = FiniteElement("Lagrange", triangle, 1)' |
b96f6b12f8ca
Move ufl to m-file in Poisson example
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
174
diff
changeset
|
20 ufl |
b96f6b12f8ca
Move ufl to m-file in Poisson example
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
174
diff
changeset
|
21 ufl u = TrialFunction(element) |
b96f6b12f8ca
Move ufl to m-file in Poisson example
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
174
diff
changeset
|
22 ufl v = TestFunction(element) |
b96f6b12f8ca
Move ufl to m-file in Poisson example
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
174
diff
changeset
|
23 ufl f = Coefficient(element) |
b96f6b12f8ca
Move ufl to m-file in Poisson example
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
174
diff
changeset
|
24 ufl g = Coefficient(element) |
b96f6b12f8ca
Move ufl to m-file in Poisson example
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
174
diff
changeset
|
25 ufl |
b96f6b12f8ca
Move ufl to m-file in Poisson example
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
174
diff
changeset
|
26 ufl "a = inner(grad(u), grad(v))*dx" |
b96f6b12f8ca
Move ufl to m-file in Poisson example
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
174
diff
changeset
|
27 ufl L = f*v*dx + g*v*ds |
b96f6b12f8ca
Move ufl to m-file in Poisson example
Eugenio Gianniti <eugenio.gianniti@mail.polimi.it>
parents:
174
diff
changeset
|
28 ufl end |
73 | 29 |
145
93a4ee13aa75
Update the examples.
gedeone-octave <marcovass89@hotmail.it>
parents:
94
diff
changeset
|
30 # Create mesh and define function space |
93a4ee13aa75
Update the examples.
gedeone-octave <marcovass89@hotmail.it>
parents:
94
diff
changeset
|
31 x = y = linspace (0, 1, 33); |
93a4ee13aa75
Update the examples.
gedeone-octave <marcovass89@hotmail.it>
parents:
94
diff
changeset
|
32 mesh = Mesh(msh2m_structured_mesh (x, y, 1, 1:4)); |
73 | 33 |
145
93a4ee13aa75
Update the examples.
gedeone-octave <marcovass89@hotmail.it>
parents:
94
diff
changeset
|
34 V = FunctionSpace('Poisson', mesh); |
73 | 35 |
145
93a4ee13aa75
Update the examples.
gedeone-octave <marcovass89@hotmail.it>
parents:
94
diff
changeset
|
36 # Define boundary condition |
93a4ee13aa75
Update the examples.
gedeone-octave <marcovass89@hotmail.it>
parents:
94
diff
changeset
|
37 bc = DirichletBC(V, @(x, y) 0.0, [2;4]); |
73 | 38 |
89
d06b423169fa
The examples now follow the new naming convention.
gedeone-octave <marcovass89@hotmail.it>
parents:
73
diff
changeset
|
39 f = Expression ('f', @(x,y) 10*exp(-((x - 0.5)^2 + (y - 0.5)^2) / 0.02)); |
d06b423169fa
The examples now follow the new naming convention.
gedeone-octave <marcovass89@hotmail.it>
parents:
73
diff
changeset
|
40 g = Expression ('g', @(x,y) sin (5.0 * x)); |
73 | 41 |
164
7190852cf57f
Update the examples giving two input parameters to the BilinearForm.
gedeone-octave <marcovass89@hotmail.it>
parents:
159
diff
changeset
|
42 a = BilinearForm ('Poisson', V, V); |
145
93a4ee13aa75
Update the examples.
gedeone-octave <marcovass89@hotmail.it>
parents:
94
diff
changeset
|
43 L = LinearForm ('Poisson', V, f, g); |
73 | 44 |
145
93a4ee13aa75
Update the examples.
gedeone-octave <marcovass89@hotmail.it>
parents:
94
diff
changeset
|
45 # Compute solution |
93a4ee13aa75
Update the examples.
gedeone-octave <marcovass89@hotmail.it>
parents:
94
diff
changeset
|
46 [A, b] = assemble_system (a, L, bc); |
93a4ee13aa75
Update the examples.
gedeone-octave <marcovass89@hotmail.it>
parents:
94
diff
changeset
|
47 sol = A \ b; |
93a4ee13aa75
Update the examples.
gedeone-octave <marcovass89@hotmail.it>
parents:
94
diff
changeset
|
48 u = Function ('u', V, sol); |
73 | 49 |
145
93a4ee13aa75
Update the examples.
gedeone-octave <marcovass89@hotmail.it>
parents:
94
diff
changeset
|
50 # Save solution in VTK format |
174 | 51 save (u, 'poisson'); |
73 | 52 |
145
93a4ee13aa75
Update the examples.
gedeone-octave <marcovass89@hotmail.it>
parents:
94
diff
changeset
|
53 # Plot solution |
93a4ee13aa75
Update the examples.
gedeone-octave <marcovass89@hotmail.it>
parents:
94
diff
changeset
|
54 plot (u); |