changeset 94:4f688918f76a

Used new available functions assemble and assemble-systems.
author gedeone-octave <marcovass89@hotmail.it>
date Wed, 07 Aug 2013 18:26:26 +0200
parents 9e7035e0494b
children 77aabcc087dc
files example/Biharmonic.m example/Evolution.m example/MixedPoisson.m example/Poisson.m
diffstat 4 files changed, 18 insertions(+), 10 deletions(-) [+]
line wrap: on
line diff
--- a/example/Biharmonic.m	Wed Aug 07 18:25:27 2013 +0200
+++ b/example/Biharmonic.m	Wed Aug 07 18:26:26 2013 +0200
@@ -18,8 +18,10 @@
 
 g = Expression ('alpha', @(x,y) 8);
 
-A = Biharmonic_BilinearForm (V, bc, g);
-b = Biharmonic_LinearForm (V, bc, f );
+a = Biharmonic_BilinearForm (V, g);
+L = Biharmonic_LinearForm (V, f);
+
+[A, b] = assemble_system (a, L, bc);
 
 u = A \ b;
 
--- a/example/Evolution.m	Wed Aug 07 18:25:27 2013 +0200
+++ b/example/Evolution.m	Wed Aug 07 18:26:26 2013 +0200
@@ -21,7 +21,10 @@
 
 u0 = Expression ('u0', @(x,y) 10*exp(-((x - 0.5)^2 + (y - 0.5)^2) / 0.02));
 
-A = Evolution_BilinearForm (V, bc, k);
+a = Evolution_BilinearForm (V, k);
+b = Evolution_LinearForm (V, k, u0);
+
+A = assemble (a, bc);
 
 # solve the problem for each time step
 # We need to update only the lhs
@@ -30,7 +33,7 @@
 
   # we can pass u0 to the lhs indifferently as a fem_coeff or
   # as a fem_func
-  b = Evolution_LinearForm (V, k, u0, bc);
+  b = assemble (L, bc);
 
   u = A \ b;
  
--- a/example/MixedPoisson.m	Wed Aug 07 18:25:27 2013 +0200
+++ b/example/MixedPoisson.m	Wed Aug 07 18:26:26 2013 +0200
@@ -15,8 +15,10 @@
 
 f = Expression ('f', @(x,y) 10*exp(-((x - 0.5)^2 + (y - 0.5)^2) / 0.02));
 
-A = MixedPoisson_BilinearForm (V, bc1, bc2);
-b = MixedPoisson_LinearForm (V, f, bc1, bc2);
+a = MixedPoisson_BilinearForm (V);
+L = MixedPoisson_LinearForm (V, f);
+
+[A, b] = assemble_system (a, L, bc1, bc2);
 
 u = A \ b;
 
--- a/example/Poisson.m	Wed Aug 07 18:25:27 2013 +0200
+++ b/example/Poisson.m	Wed Aug 07 18:26:26 2013 +0200
@@ -16,12 +16,13 @@
 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 = Poisson_BilinearForm (V, bc);
-b = Poisson_LinearForm (V, f, g, bc);
+a = Poisson_BilinearForm (V);
+L = Poisson_LinearForm (V, f, g);
+
+[A, b] = assemble_system (a, L, bc);
 
 u = A \ b;
 
-
-func = fem_func ('u', V, u);
+func = Function ('u', V, u);
 
 fem_plot (func);