Mercurial > fem-fenics-eugenio
view inst/example/Mixed-Poisson/MixedPoisson.m @ 174:e85390ed620d
Add copyright.
author | gedeone-octave <marcovass89@hotmail.it> |
---|---|
date | Sat, 26 Oct 2013 17:24:56 +0100 |
parents | 7190852cf57f |
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 import_ufl_Problem ('MixedPoisson') # Create mesh x = y = linspace (0, 1, 33); mesh = Mesh(msh2m_structured_mesh (x, y, 1, 1:4)); V = FunctionSpace('MixedPoisson', mesh); f = Expression ('f', @(x,y) 10*exp(-((x - 0.5)^2 + (y - 0.5)^2) / 0.02)); a = BilinearForm ('MixedPoisson', V, V); L = LinearForm ('MixedPoisson', V, f); # Define essential boundary bc1 = DirichletBC (SubSpace (V, 1), @(x,y) [0; -sin(5.0*x)], 1); bc2 = DirichletBC (SubSpace (V, 1), @(x,y) [0; sin(5.0*x)], 3); # Compute solution [A, b] = assemble_system (a, L, bc1, bc2); sol = A \ b; func = Function ('func', V, sol); sigma = Function ('sigma', func, 1); u = Function ('u', func, 2); # Plot solution plot (sigma); plot (u);