leq_isolve.f90

This example show how use linear_equation module for solve linear equation. We will try to solve Ax=b by use iteratives solves methods More details about this example.

00001 
00016 program leq_isolve
00017 #include "fml_constants.h" 
00018  use mod_linear_equation     ! use linear_equation module
00019   implicit none
00020 
00021   !********************************************* declaration  
00022   integer, parameter :: size_ab=30 
00023   type(matrix) :: A  !main square matrix
00024   type(vector) :: b  !second member
00025   type(t_soleq) :: x  !contains  solution (v_sol), iterations (iter) and errors (v_err)
00026     
00027   !********************************************* body 
00028   call init(A,size_ab,size_ab); !init:=m_init
00029   call init(b,size_ab); !init:=v_init 
00030   !initialize A by random values between 1.0 and 10.0 (symmetric dominant diagonal)
00031   call mc_diagDominantSymmetric(A,size_ab,low=p_notcast(1.0),high=p_notcast(6.5))
00032   !initialize b by random values between 1.0 and 5.0
00033   call random(b,low=p_notcast(1.0),high=p_notcast(6.5)) 
00034  
00035   !if a method diverge, try to change the option (tolerance error, number of iterations)
00036   !see function prototype for more information
00037   
00038   print*; !newline  
00039   print*, "********************************************* display initial data"
00040   print*,"* A="; call print(A);   !print A  (print:=m_print)
00041   write(*,'(A5)',advance='no')"* b="; call print(b);   !print b  (print:=v_print)
00042   print*; !newline  
00043   print*, "********************************************* Solve linear equation"  
00044   print*; !newline 
00045   !let's see the function prototype for more options 
00046   print*, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> solve :: Jacobi element"
00047   x=ilinsolve(A,b,meth_solve='jacelt')  
00048   write(*,'(A5)',advance='no')"* x="; call print(x%v_sol); 
00049   write(*,'(A24,I3)')"# number of iterations=", x%iter
00050   write(*,'(A9)',advance='no')"# A*x-b="; call print(A*x%v_sol-b);
00051   print*; !newline
00052   print*, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> solve :: SOR element"
00053   x=ilinsolve(A,b,meth_solve='sorelt')  
00054   write(*,'(A5)',advance='no')"* x="; call print(x%v_sol); 
00055   write(*,'(A24,I3)')"# number of iterations=", x%iter
00056   write(*,'(A9)',advance='no')"# A*x-b="; call print(A*x%v_sol-b);
00057   print*; !newline   
00058   print*, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> solve :: Gauss-Seidel element"
00059   x=ilinsolve(A,b,meth_solve='gselt')  
00060   write(*,'(A5)',advance='no')"* x="; call print(x%v_sol); 
00061   write(*,'(A24,I3)')"# number of iterations=", x%iter
00062   write(*,'(A9)',advance='no')"# A*x-b="; call print(A*x%v_sol-b);
00063   print*; !newline   
00064   print*, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> solve :: Gradient element "
00065   x=ilinsolve(A,b,meth_solve='gradelt')  
00066   write(*,'(A5)',advance='no')"* x="; call print(x%v_sol); 
00067   write(*,'(A24,I3)')"# number of iterations=", x%iter
00068   write(*,'(A9)',advance='no')"# A*x-b="; call print(A*x%v_sol-b);
00069   print*; !newline     
00070   print*, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> solve :: Conjugate Gradient element  "
00071   x=ilinsolve(A,b,meth_solve='cgelt')  
00072   write(*,'(A5)',advance='no')"* x="; call print(x%v_sol); 
00073   write(*,'(A24,I3)')"# number of iterations=", x%iter
00074   write(*,'(A9)',advance='no')"# A*x-b="; call print(A*x%v_sol-b); 
00075   print*; !newline  
00076   print*, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> solve :: Gauss-Seidel block  "
00077   !see function prototype for more information (example block composition) 
00078   x=ilinsolve(A,b,meth_solve='gsblock',eps=p_notcast(1.e-5))  
00079   write(*,'(A5)',advance='no')"* x="; call print(x%v_sol); 
00080   write(*,'(A24,I3)')"# number of iterations=", x%iter
00081   write(*,'(A9)',advance='no')"# A*x-b="; call print(A*x%v_sol-b);
00082   print*; !newline    
00083  
00084   print*, "********************************************* Solve linear equation preconditionning" 
00085   print*; !newline     
00086   print*, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> solve :: Conjugate Gradient element preconditionning"
00087   x=ilinsolve(A,b,meth_solve='cgelt',is_precond=.true.)  
00088   write(*,'(A5)',advance='no')"* x="; call print(x%v_sol); 
00089   write(*,'(A24,I3)')"# number of iterations=", x%iter
00090   write(*,'(A9)',advance='no')"# A*x-b="; call print(A*x%v_sol-b); 
00091   print*; !newline  
00092   
00093   call destruct(A)  ! destruct the matrix A  (don't forget to destruct the matrix)
00094   call destruct(b)
00095   call destruct(x)
00096 end program leq_isolve
 All Classes Namespaces Files Functions Variables Defines