Mercurial > fem-fenics-eugenio
view inst/example/Mixed-Poisson/MixedPoisson.m @ 159:d03d998e64af
Move examples in the ./inst directory.
author | gedeone-octave <marcovass89@hotmail.it> |
---|---|
date | Thu, 12 Sep 2013 15:22:19 +0200 |
parents | example/Mixed-Poisson/MixedPoisson.m@93a4ee13aa75 |
children | 7190852cf57f |
line wrap: on
line source
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); 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);