Mercurial > fem-fenics-eugenio
view 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 |
line wrap: on
line source
## Copyright (C) 2013 Marco Vassallo <gedeone-octave@users.sourceforge.net> ## ## This program is free software; you can redistribute it and/or modify it under ## the terms of the GNU General Public License as published by the Free Software ## Foundation; either version 3 of the License, or (at your option) any later ## version. ## ## This program is distributed in the hope that it will be useful, but WITHOUT ## ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or ## FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more ## details. ## ## You should have received a copy of the GNU General Public License along with ## this program; if not, see <http://www.gnu.org/licenses/>. pkg load fem-fenics msh ufl start Poisson ufl 'element = FiniteElement("Lagrange", triangle, 1)' ufl ufl u = TrialFunction(element) ufl v = TestFunction(element) ufl f = Coefficient(element) ufl g = Coefficient(element) ufl ufl "a = inner(grad(u), grad(v))*dx" ufl L = f*v*dx + g*v*ds ufl end # Create 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 bc = DirichletBC(V, @(x, y) 0.0, [2;4]); f = Expression ('f', @(x,y) 10*exp(-((x - 0.5)^2 + (y - 0.5)^2) / 0.02)); g = Expression ('g', @(x,y) sin (5.0 * x)); a = BilinearForm ('Poisson', V, V); L = LinearForm ('Poisson', V, f, g); # Compute solution [A, b] = assemble_system (a, L, bc); sol = A \ b; u = Function ('u', V, sol); # Save solution in VTK format save (u, 'poisson'); # Plot solution plot (u);