This example show how use linear_equation module for solve linear equation. We will try to solve Ax=b by classical solves methods (not iteratives method) More details about this example.
00001 00016 program leq_solve 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=4 00023 type(matrix) :: A !main square matrix 00024 type(vector) :: b !second member 00025 type(vector) :: x !solution 00026 !others 00027 type(matrix) :: LowA, UpA, TriDiag 00028 00029 !********************************************* body 00030 call init(A,size_ab,size_ab); !init:=m_init 00031 call init(b,size_ab); !init:=v_init 00032 !initialize A by random values between 1.0 and 10.0 (dominant diagonal) 00033 call mc_diagDominant(A,size_ab,low=p_notcast(1.0),high=p_notcast(6.5)) 00034 !initialize b by random values between 1.0 and 5.0 00035 call random(b,low=p_notcast(1.0),high=p_notcast(6.5)) 00036 !create tridiagonal matrix 00037 TriDiag=mc_tridiag(d=p_notcast(3.0),u=p_notcast(-5.0),l=p_notcast(-0.5),size_n=size_ab) 00038 !lower matrix of A 00039 LowA=tril(A) 00040 !upper matrix of A 00041 UpA=triu(A) 00042 00043 print*; !newline 00044 print*, "********************************************* display initial data" 00045 print*,"* A="; call print(A); !print A (print:=m_print) 00046 write(*,'(A5)',advance='no')"* b="; call print(b); !print b (print:=v_print) 00047 print*; !newline 00048 00049 print*, "********************************************* Solve linear equation" 00050 print*; !newline 00051 !let's see the function prototype for more options 00052 print*, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> solve :: Gauss-Jourdan" 00053 x=linsolve(A,b,meth_solve='gj') 00054 write(*,'(A5)',advance='no')"* x="; call print(x); 00055 write(*,'(A9)',advance='no')"# A*x-b="; call print(A*x-b); 00056 00057 print*; !newline 00058 print*, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> solve :: least-square" 00059 x=linsolve(A,b,meth_solve='ls') 00060 write(*,'(A5)',advance='no')"* x="; call print(x); 00061 write(*,'(A9)',advance='no')"# A*x-b="; call print(A*x-b); 00062 print*; !newline 00063 print*, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> solve :: least-square (QR factorization)" 00064 x=linsolve(A,b,meth_solve='lsqr') 00065 write(*,'(A5)',advance='no')"* x="; call print(x); 00066 write(*,'(A9)',advance='no')"# A*x-b="; call print(A*x-b); 00067 print*; !newline 00068 print*, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> solve :: pseudo-inverse by svd decomposition " 00069 x=linsolve(A,b,meth_solve='pinvsvd') 00070 write(*,'(A5)',advance='no')"* x="; call print(x); 00071 write(*,'(A9)',advance='no')"# A*x-b="; call print(A*x-b); 00072 print*; !newline 00073 print*, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> solve :: pseudo-inverse by cholesky decomposition " 00074 x=linsolve(A,b,meth_solve='pinvchol') 00075 write(*,'(A5)',advance='no')"* x="; call print(x); 00076 write(*,'(A9)',advance='no')"# A*x-b="; call print(A*x-b); 00077 print*; !newline 00078 print*, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> solve :: LU decomposition " 00079 x=linsolve(A,b,meth_solve='lu') 00080 write(*,'(A5)',advance='no')"* x="; call print(x); 00081 write(*,'(A9)',advance='no')"# A*x-b="; call print(A*x-b); 00082 print*; !newline 00083 print*, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> solve :: QR decomposition " 00084 x=linsolve(A,b,meth_solve='qr') 00085 write(*,'(A5)',advance='no')"* x="; call print(x); 00086 write(*,'(A9)',advance='no')"# A*x-b="; call print(A*x-b); 00087 print*; !newline 00088 print*, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> solve :: cholesky decomposition " 00089 x=linsolve(A,b,meth_solve='chol') 00090 write(*,'(A5)',advance='no')"* x="; call print(x); 00091 write(*,'(A9)',advance='no')"# A*x-b="; call print(A*x-b); 00092 print*; !newline 00093 00094 print*, "********************************************* Solve particular equation" 00095 print*; !newline 00096 print*, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> solve :: lower triangular matrix " 00097 print*,"* LowA="; call print(LowA); 00098 write(*,'(A5)',advance='no')"* b="; call print(b); 00099 print*; !newline 00100 x=linsolve(LowA,b,meth_solve='tril') 00101 write(*,'(A5)',advance='no')"* x="; call print(x); 00102 write(*,'(A12)',advance='no')"# LowA*x-b="; call print(LowA*x-b); 00103 print*; !newline 00104 print*, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> solve :: upper triangular matrix " 00105 print*,"* UpA="; call print(UpA); 00106 write(*,'(A5)',advance='no')"* b="; call print(b); 00107 x=linsolve(UpA,b,meth_solve='triu') 00108 write(*,'(A5)',advance='no')"* x="; call print(x); 00109 write(*,'(A11)',advance='no')"# UpA*x-b="; call print(UpA*x-b); 00110 print*; !newline 00111 print*, ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> solve :: tridiagonal matrix by Thomas algorithm " 00112 print*,"* TriDiag="; call print(TriDiag); 00113 write(*,'(A5)',advance='no')"* b="; call print(b); 00114 x=linsolve(TriDiag,b,meth_solve='thomaselt') 00115 write(*,'(A5)',advance='no')"* x="; call print(x); 00116 write(*,'(A15)',advance='no')"# TriDiag*x-b="; call print(TriDiag*x-b); 00117 print*; !newline 00118 00119 call destruct(A) ! destruct the matrix A (don't forget to destruct the matrix) 00120 call destruct(b) 00121 call destruct(x) 00122 call destruct(LowA) 00123 call destruct(UpA) 00124 end program leq_solve