Mercurial > fem-fenics-eugenio
changeset 95:77aabcc087dc
Update the examples to the new naming convention.
author | gedeone-octave <marcovass89@hotmail.it> |
---|---|
date | Wed, 07 Aug 2013 20:45:53 +0200 |
parents | 4f688918f76a |
children | 38b31ad6ab4e |
files | example/Evolution.m example/MixedPoisson.m example/NavierStokes/NS.m |
diffstat | 3 files changed, 18 insertions(+), 10 deletions(-) [+] |
line wrap: on
line diff
--- a/example/Evolution.m Wed Aug 07 18:26:26 2013 +0200 +++ b/example/Evolution.m Wed Aug 07 20:45:53 2013 +0200 @@ -22,7 +22,7 @@ u0 = Expression ('u0', @(x,y) 10*exp(-((x - 0.5)^2 + (y - 0.5)^2) / 0.02)); a = Evolution_BilinearForm (V, k); -b = Evolution_LinearForm (V, k, u0); +L = Evolution_LinearForm (V, k, u0); A = assemble (a, bc); @@ -33,6 +33,7 @@ # we can pass u0 to the lhs indifferently as a fem_coeff or # as a fem_func + L = Evolution_LinearForm (V, k, u0); b = assemble (L, bc); u = A \ b;
--- a/example/MixedPoisson.m Wed Aug 07 18:26:26 2013 +0200 +++ b/example/MixedPoisson.m Wed Aug 07 20:45:53 2013 +0200 @@ -22,9 +22,9 @@ u = A \ b; -func = fem_func ('u', V, u); -uu = fem_func ('u', func, 0); -sigma = fem_func ('sigma', func, 1); +func = Function ('u', V, u); +uu = Function ('u', func, 0); +sigma = Function ('sigma', func, 1); fem_plot (sigma); fem_plot (uu);
--- a/example/NavierStokes/NS.m Wed Aug 07 18:26:26 2013 +0200 +++ b/example/NavierStokes/NS.m Wed Aug 07 20:45:53 2013 +0200 @@ -24,27 +24,34 @@ f = Expression ('f', @(x,y) [0; 0]); u0 = Expression ('u0', @(x,y) [0; 0]); p0 = Expression ('p1', @(x,y) 0); -A1 = TentativeVelocity_BilinearForm (V, k, bcnoslip); -A3 = VelocityUpdate_BilinearForm (V, bcnoslip); + +a1 = TentativeVelocity_BilinearForm (V, k); +A1 = assemble (a1, bcnoslip); +a3 = VelocityUpdate_BilinearForm (V); +A3 = assemble (a3, bcnoslip) + # solve the problem for each time step # We need to update only the lhs while t < T t += dt; bcin = DirichletBC (Q, @(x,y) sin(3.0*t), 1); - b1 = TentativeVelocity_LinearForm (V, k, u0, f, bcnoslip); + L1 = TentativeVelocity_LinearForm (V, k, u0, f); + b1 = assemble (L1, bcnoslip); utmp = A1 \ b1; u1 = Function ('u1', V, utmp); - A2 = PressureUpdate_BilinearForm (Q, bcin, bcout); - b2 = PressureUpdate_LinearForm (Q, bcin, bcout, u1, k); + a2 = PressureUpdate_BilinearForm (Q); + L2 = PressureUpdate_LinearForm (Q,u1, k); + [A2, b2] = assemble_system (a2, L2, bcin, bcout); ptmp = A2 \ b2; p1 = Function ('p1', Q, ptmp); - b3 = VelocityUpdate_LinearForm (V, k, u1, p1, bcnoslip); + L3 = VelocityUpdate_LinearForm (V, k, u1, p1); + b3 = assemble (L3, bcnoslip); ut = A3 \ b3; u0 = Function ('u0', V, ut);