00001
00003
00020 module mod_matrix
00021 #include "fml_constants.h"
00022 use mod_exception
00023 use mod_utility
00024 use mod_maths
00025 use mod_vector
00026 implicit none
00027
00028
00029
00030
00031
00041 type matrix
00042 integer :: rows=2;
00043 integer :: cols=2;
00044 type_precision, dimension(:,:), allocatable :: ptr_container(:,:)
00045 logical :: is_allocate = .false.
00046 end type matrix
00047
00056 type t_lu
00057 type(matrix) :: L
00058 type(matrix) :: U
00059 type(matrix) :: P
00060 integer :: swops=0
00061 end type t_lu
00062
00071 type t_qr
00072 type(matrix) :: Q
00073 type(matrix) :: R
00074 type(matrix) :: P
00075 integer :: swops=0
00076 end type t_qr
00077
00086 type t_svd
00087 type(matrix) :: S
00088 type(matrix) :: V
00089 type(matrix) :: U
00090 end type t_svd
00091
00092
00100 type t_m_and_p
00101 type(matrix) :: M
00102 type(matrix) :: P
00103 integer :: swops=0
00104 end type t_m_and_p
00105
00115 type t_poweig
00116 type_precision :: lambda
00117 type(vector) :: v_lambda
00118 type(vector) :: v_err
00119 integer :: iter
00120 end type t_poweig
00121
00129 type t_eig
00130 type(vector) :: v_eigvalues
00131 type(matrix) :: m_eigvectors
00132 end type t_eig
00133
00134
00140 interface operator(*)
00141 module procedure m_prod_mat, m_prod_vec_c, m_prod_vec2, m_prod_scalar1, m_prod_scalar2
00142 end interface
00143
00149 interface operator(/)
00150 module procedure m_div_scalar
00151 end interface
00152
00158 interface operator(==)
00159 module procedure m_isEqual, m_isEqual_scalar
00160 end interface
00161
00167 interface operator(+)
00168 module procedure m_add
00169 end interface
00170
00176 interface operator(-)
00177 module procedure m_minus
00178 end interface
00179
00185 interface assignment(=)
00186 module procedure m_affect
00187 end interface
00188
00189
00195 interface operator(.tr.)
00196 module procedure m_trans
00197 end interface
00198
00204 interface operator(.det.)
00205 module procedure m_det_gaussj
00206 end interface
00207
00208
00214 interface operator(.rank.)
00215 module procedure m_rank_gaussj
00216 end interface
00217
00223 interface operator(.inv.)
00224 module procedure m_inverse_gaussj
00225 end interface
00226
00232 interface operator(.cond.)
00233 module procedure m_cond
00234 end interface
00235
00241 interface destruct
00242 module procedure m_destruct
00243 module procedure m_destruct_qr, m_destruct_lu
00244 module procedure m_destruct_m_and_p, m_destruct_t_poweig
00245 module procedure m_destruct_t_eig, m_destruct_t_svd
00246 end interface
00247
00253 interface sum
00254 module procedure m_sum
00255 end interface
00256
00262 interface min
00263 module procedure m_min
00264 end interface
00265
00271 interface max
00272 module procedure m_max
00273 end interface
00274
00280 interface inv
00281 module procedure m_inverse_gaussj
00282 end interface
00283
00289 interface pinv
00290 module procedure m_pinv
00291 end interface
00292
00293
00299 interface norm
00300 module procedure m_norm
00301 end interface
00302
00308 interface triu
00309 module procedure m_triu
00310 end interface
00311
00317 interface tril
00318 module procedure m_tril
00319 end interface
00320
00326 interface diag
00327 module procedure m_diag
00328 end interface
00329
00335 interface get
00336 module procedure m_get , m_get_m, m_getRow
00337 end interface
00338
00344 interface set
00345 module procedure m_set
00346 end interface
00347
00353 interface spec
00354 module procedure m_eig_qr
00355 end interface
00356
00362 interface rank
00363 module procedure m_rank
00364 end interface
00365
00371 interface det
00372 module procedure m_det
00373 end interface
00374
00380 interface chol
00381 module procedure m_decompCholesky
00382 end interface
00383
00389 interface lu
00390 module procedure m_decompLU
00391 end interface
00392
00398 interface m_lu
00399 module procedure m_decompLU_m
00400 end interface
00401
00407 interface qr
00408 module procedure m_decompQR
00409 end interface
00410
00416 interface svd
00417 module procedure m_decompsvd
00418 end interface
00419
00425 interface init
00426 module procedure m_init, m_init_fromfile
00427 end interface
00428
00434 interface random
00435 module procedure mc_random
00436 end interface
00437
00443 interface print
00444 module procedure m_print, m_print_tofile, m_print_lu_tofile
00445 end interface
00446
00449 character(len=len_what_exception) :: m_what_exception
00450
00451
00452 integer, parameter :: frobenius = -1
00453
00454 CONTAINS
00455
00456
00457
00458
00459
00460
00469 subroutine m_init(m,rows_,cols_)
00470 implicit none
00471 type(matrix) :: m;
00472 integer, intent(in), optional :: rows_;
00473 integer, intent(in), optional :: cols_;
00474
00475 integer(kind=2) :: ierr;
00476 type(type_exception) :: err_exception;
00477 integer :: in_rows=2, in_cols=2
00478
00479 if(present(rows_)) in_rows=rows_;
00480 if(present(cols_)) in_cols=cols_;
00481
00482 #ifdef DEBUG_EXCEPTION
00483 if (in_rows .lt. 1) then
00484 m_what_exception='m_init::row is nil or negative'
00485 err_exception=e_error(stop_array_indice_exceed,m_what_exception,stop_array_indice_exceed)
00486 endif
00487 if (in_cols .lt. 1) then
00488 m_what_exception='m_init::column is nil or negative'
00489 err_exception=e_error(stop_array_indice_exceed,m_what_exception,stop_array_indice_exceed)
00490 endif
00491 #endif
00492 call destruct(m)
00493
00494 m%rows=in_rows
00495 m%cols=in_cols
00496 allocate(m%ptr_container(m%rows,m%cols),stat=ierr)
00497
00498 #ifdef DEBUG_EXCEPTION
00499 if (ierr .ne. 0) then
00500 m_what_exception='m_init::allocate error'
00501 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
00502 endif
00503 #endif
00504 m%ptr_container(:,:)=0
00505 m%is_allocate = allocated(m%ptr_container)
00506 end subroutine m_init
00507
00508
00516 subroutine m_init_fromfile(m,filename,unit)
00517 implicit none
00518 type(matrix) :: m;
00519 character*(*), intent(in) :: filename
00520 integer, intent(in), optional :: unit
00521
00522 integer(kind=2) :: ierr;
00523 type(type_exception) :: err_exception;
00524 integer :: in_unit
00525 integer :: i
00526 integer :: size_mat(2)
00527
00528 if(u_is_exist_file(filename)) then
00529 call destruct(m)
00530
00531 if(present(unit)) then
00532 in_unit=unit;
00533 else
00534 in_unit=22
00535 end if
00536
00537 if(.not.u_is_exist_file(filename)) then
00538 m_what_exception='m_init_fromfile::'//filename//' not exist'
00539 err_exception=e_error(stop_open_file,m_what_exception,stop_open_file)
00540 end if
00541
00542 open (unit=in_unit,file=filename, form="formatted", action="read",status="unknown")
00543 read (unit=in_unit,fmt=*,end=1000)size_mat
00544 #ifdef DEBUG_EXCEPTION
00545 if (size_mat(1) .lt. 1) then
00546 m_what_exception='m_init_fromfile::row is nil or negative'
00547 err_exception=e_error(stop_array_indice_exceed,m_what_exception,stop_array_indice_exceed)
00548 endif
00549 if (size_mat(2) .lt. 1) then
00550 m_what_exception='m_init_fromfile::column is nil or negative'
00551 err_exception=e_error(stop_array_indice_exceed,m_what_exception,stop_array_indice_exceed)
00552 endif
00553 #endif
00554 m%rows=size_mat(1)
00555 m%cols=size_mat(2)
00556 allocate(m%ptr_container(m%rows,m%cols),stat=ierr)
00557 #ifdef DEBUG_EXCEPTION
00558 if (ierr .ne. 0) then
00559 m_what_exception='m_init_fromfile::matrix allocate error'
00560 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
00561 endif
00562 #endif
00563 do i=1,m%rows
00564 read (unit=in_unit,fmt=*,end=1000) m%ptr_container(i,:)
00565 end do
00566 1000 close(unit=in_unit)
00567 m%is_allocate = allocated(m%ptr_container)
00568 end if
00569 end subroutine m_init_fromfile
00570
00571
00580 subroutine mc_diagDominant(m,n,low,high)
00581 implicit none
00582 type(matrix) :: m;
00583 integer, intent(in) :: n
00584 type_precision, intent(in) :: low, high
00585
00586 integer(kind=2) :: ierr;
00587 type(type_exception) :: err_exception;
00588 integer :: i,j
00589 #ifdef DEBUG_EXCEPTION
00590 if (n .lt. 1) then
00591 m_what_exception='mc_diagDominant::matrix index is nil or negative'
00592 err_exception=e_error(stop_array_indice_exceed,m_what_exception,stop_array_indice_exceed)
00593 endif
00594 #endif
00595 call destruct(m)
00596 m%rows=n
00597 m%cols=n
00598 allocate(m%ptr_container(m%rows,m%cols),stat=ierr)
00599 #ifdef DEBUG_EXCEPTION
00600 if (ierr .ne. 0) then
00601 m_what_exception='mc_diagDominant::allocate error'
00602 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
00603 endif
00604 #endif
00605 do i=1,m%rows
00606 do j=1,m%cols
00607 m%ptr_container(i,j)=rand(0) * (high - low) + low
00608 end do
00609 end do
00610
00611 forall(i=1:m%rows) m%ptr_container(i,i)=1+sum(m%ptr_container(i,:))
00612
00613 m%is_allocate = allocated(m%ptr_container)
00614 end subroutine mc_diagDominant
00615
00616
00625 subroutine mc_diagDominantSymmetric(m,n,low,high)
00626 implicit none
00627 type(matrix) :: m;
00628 integer, intent(in) :: n
00629 type_precision, intent(in) :: low, high
00630
00631 integer(kind=2) :: ierr;
00632 type(type_exception) :: err_exception;
00633 integer :: i,j
00634 #ifdef DEBUG_EXCEPTION
00635 if (n .lt. 1) then
00636 m_what_exception='mc_diagDominantSymmetric::matrix index is nil or negative'
00637 err_exception=e_error(stop_array_indice_exceed,m_what_exception,stop_array_indice_exceed)
00638 endif
00639 #endif
00640 call destruct(m)
00641 m%rows=n
00642 m%cols=n
00643 allocate(m%ptr_container(m%rows,m%cols),stat=ierr)
00644 #ifdef DEBUG_EXCEPTION
00645 if (ierr .ne. 0) then
00646 m_what_exception='mc_diagDominantSymmetric::allocate error'
00647 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
00648 endif
00649 #endif
00650 do i=1,m%rows
00651 do j=i+1,m%cols
00652 m%ptr_container(i,j)=rand(0) * (high - low) + low
00653 end do
00654 do j=1,i-1
00655 m%ptr_container(i,j)=m%ptr_container(j,i)
00656 end do
00657 end do
00658 forall(i=1:m%rows) m%ptr_container(i,i)=1+sum(m%ptr_container(i,:))
00659
00660 m%is_allocate = allocated(m%ptr_container)
00661 end subroutine mc_diagDominantSymmetric
00662
00663
00671 subroutine m_resize(m,rows_,cols_)
00672 implicit none
00673 type(matrix), intent(inout) :: m;
00674 integer, intent(in) :: rows_;
00675 integer, intent(in) :: cols_;
00676
00677 integer(kind=2) :: ierr;
00678 type(type_exception) :: err_exception;
00679 #ifdef DEBUG_EXCEPTION
00680 if (rows_ .lt. 1) then
00681 m_what_exception='m_resize::row is nil or negative'
00682 err_exception=e_error(stop_array_indice_exceed,m_what_exception,stop_array_indice_exceed)
00683 endif
00684 if (cols_ .lt. 1) then
00685 m_what_exception='m_resize::column is nil or negative'
00686 err_exception=e_error(stop_array_indice_exceed,m_what_exception,stop_array_indice_exceed)
00687 endif
00688 #endif
00689 call destruct(m)
00690
00691 m%rows=rows_
00692 m%cols=cols_
00693 allocate(m%ptr_container(m%rows,m%cols),stat=ierr)
00694 #ifdef DEBUG_EXCEPTION
00695 if (ierr .ne. 0) then
00696 m_what_exception='m_resize::allocate error'
00697 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
00698 endif
00699 #endif
00700 m%ptr_container(:,:)=0
00701 m%is_allocate = allocated(m%ptr_container)
00702 end subroutine m_resize
00703
00704
00710 subroutine m_destruct(m)
00711 implicit none
00712 type(matrix) :: m;
00713
00714 integer(kind=2) :: ierr;
00715 type(type_exception) :: err_exception;
00716
00717 if(m%is_allocate) then
00718 deallocate(m%ptr_container,stat=ierr)
00719 #ifdef DEBUG_EXCEPTION
00720 if (ierr .ne. 0) then
00721 m_what_exception='m_destruct::matrix deallocate error'
00722 err_exception=e_error(stop_dealoc,m_what_exception,stop_dealoc)
00723 endif
00724 #endif
00725 m%is_allocate = allocated(m%ptr_container)
00726 end if
00727 end subroutine m_destruct
00728
00729
00730
00731
00737 subroutine m_destruct_lu(lu_)
00738 implicit none
00739 type(t_lu) :: lu_;
00740 call destruct(lu_%L)
00741 call destruct(lu_%U)
00742 call destruct(lu_%P)
00743 end subroutine m_destruct_lu
00744
00745
00751 subroutine m_destruct_m_and_p(m_and_p_)
00752 implicit none
00753 type(t_m_and_p) :: m_and_p_;
00754 call destruct(m_and_p_%M)
00755 call destruct(m_and_p_%P)
00756 end subroutine m_destruct_m_and_p
00757
00758
00764 subroutine m_destruct_qr(qr_)
00765 implicit none
00766 type(t_qr) :: qr_;
00767 call destruct(qr_%Q)
00768 call destruct(qr_%R)
00769 call destruct(qr_%P)
00770 end subroutine m_destruct_qr
00771
00772
00773
00779 subroutine m_destruct_t_poweig(t_poweig_)
00780 implicit none
00781 type(t_poweig) :: t_poweig_;
00782 call destruct(t_poweig_%v_err)
00783 call destruct(t_poweig_%v_lambda)
00784 end subroutine m_destruct_t_poweig
00785
00786
00792 subroutine m_destruct_t_eig(t_eig_)
00793 implicit none
00794 type(t_eig) :: t_eig_;
00795 call destruct(t_eig_%v_eigvalues)
00796 call destruct(t_eig_%m_eigvectors)
00797 end subroutine m_destruct_t_eig
00798
00799
00805 subroutine m_destruct_t_svd(t_svd_)
00806 implicit none
00807 type(t_svd) :: t_svd_;
00808 call destruct(t_svd_%S)
00809 call destruct(t_svd_%V)
00810 call destruct(t_svd_%U)
00811 end subroutine m_destruct_t_svd
00812
00813
00822 subroutine mc_random(m,low,high)
00823 implicit none
00824 type(matrix), intent(inout) :: m;
00825 type_precision, intent(in), optional :: low
00826 type_precision, intent(in), optional :: high
00827
00828 type_precision :: in_low, in_high
00829 type(type_exception) :: err_exception;
00830 integer :: i,j
00831 real :: x
00832 #ifdef DEBUG_EXCEPTION
00833 if (.not. m%is_allocate) then
00834 m_what_exception='mc_random::matrix not allocate yet'
00835 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
00836 endif
00837 #endif
00838 if(present(low)) then
00839 in_low=low;
00840 else
00841 in_low=1.0;
00842 end if
00843 if(present(high)) then
00844 in_high=high;
00845 else
00846 in_high=max(m%rows,m%cols);
00847 end if
00848 do i=1,m%rows
00849 do j=1,m%cols
00850 m%ptr_container(i,j)=rand((i+j)**100*int(secnds(x)))*(in_high-in_low)+in_low
00851 end do
00852 end do
00853 end subroutine mc_random
00854
00855
00856
00863 subroutine m_minit_value(m,value)
00864 implicit none
00865 type_precision,intent(in) :: value
00866 type(matrix) :: m;
00867
00868 type(type_exception) :: err_exception;
00869 #ifdef DEBUG_EXCEPTION
00870 if (.NOT. m%is_allocate) then
00871 m_what_exception='m_minit_value::matrix not allocate yet'
00872 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
00873 endif
00874 #endif
00875 m%ptr_container=value
00876 end subroutine m_minit_value
00877
00878
00886 function m_zeros(rows_,cols_) result(res)
00887 implicit none
00888 integer, intent(in) :: rows_;
00889 integer, intent(in) :: cols_;
00890
00891 type(matrix) ::res;
00892 type(type_exception) :: err_exception;
00893 integer(kind=2) :: ierr;
00894
00895 res%rows=rows_
00896 res%cols=cols_
00897 allocate(res%ptr_container(res%rows,res%cols),stat=ierr)
00898 #ifdef DEBUG_EXCEPTION
00899 if (ierr .ne. 0) then
00900 m_what_exception='m_zeros::allocate error'
00901 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
00902 endif
00903 #endif
00904 res%ptr_container(:,:)=0
00905 res%is_allocate = allocated(res%ptr_container)
00906 end function m_zeros
00907
00908
00924
00925 function m_extract(m,r_lbound,r_ubound,c_lbound,c_ubound) result (res)
00926 type(matrix), intent(in) :: m
00927 integer, intent(in), optional :: r_lbound, r_ubound
00928 integer, intent(in), optional :: c_lbound, c_ubound
00929
00930 type(type_exception) :: err_exception;
00931 type(matrix) :: res
00932 integer :: in_r_lbound=1, in_r_ubound=1
00933 integer :: in_c_lbound=1, in_c_ubound=1
00934
00935 if(present(r_lbound)) in_r_lbound=r_lbound;
00936 if(present(r_ubound)) in_r_ubound=r_ubound;
00937 if(present(c_lbound)) in_c_lbound=c_lbound;
00938 if(present(c_ubound)) in_c_ubound=c_ubound;
00939
00940 #ifdef DEBUG_EXCEPTION
00941 if (in_r_lbound.le.0 .or. in_r_ubound.le.0 .or. in_c_lbound.le.0 .or. in_c_ubound.le.0) then
00942 m_what_exception='m_extract::positive boundary expected'
00943 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
00944 endif
00945 if (in_r_lbound.gt.m%rows .or. in_r_ubound.gt.m%rows .or. in_c_lbound.gt.m%cols .or. in_c_ubound.gt.m%cols) then
00946 m_what_exception='m_extract::inconsistent boundary'
00947 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
00948 endif
00949 #endif
00950
00951 call m_init(res,in_r_ubound-in_r_lbound+1,in_c_ubound-in_c_lbound+1);
00952 res%ptr_container=m%ptr_container(in_r_lbound:in_r_ubound,in_c_lbound:in_c_ubound)
00953 end function m_extract
00954
00955
00973 subroutine m_setsub(m,r_lbound,r_ubound,c_lbound,c_ubound,m_sub)
00974 type(matrix), intent(inout) :: m
00975 integer, intent(in), optional :: r_lbound, r_ubound
00976 integer, intent(in), optional :: c_lbound, c_ubound
00977 type(matrix), intent(in) :: m_sub
00978
00979 type(type_exception) :: err_exception;
00980 type(matrix) :: res
00981 integer :: i, j
00982 integer :: in_r_lbound, in_r_ubound
00983 integer :: in_c_lbound, in_c_ubound
00984
00985
00986 in_r_lbound=1; in_r_ubound=m%rows
00987 in_c_lbound=1; in_c_ubound=m%cols
00988 if(present(r_lbound)) in_r_lbound=r_lbound;
00989 if(present(r_ubound)) in_r_ubound=r_ubound;
00990 if(present(c_lbound)) in_c_lbound=c_lbound;
00991 if(present(c_ubound)) in_c_ubound=c_ubound;
00992
00993 #ifdef DEBUG_EXCEPTION
00994 if(.not.(m%is_allocate) .or. .not.(m_sub%is_allocate)) then
00995 m_what_exception='m_setsub::matrix not allocated'
00996 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
00997 endif
00998 if (in_r_lbound.le.0 .or. in_r_ubound.le.0 .or. in_c_lbound.le.0 .or. in_c_ubound.le.0) then
00999 m_what_exception='m_setsub::positive boundary expected'
01000 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01001 endif
01002 if (in_r_lbound.gt.m%rows .or. in_r_ubound.gt.m%rows .or. in_c_lbound.gt.m%cols .or. in_c_ubound.gt.m%cols) then
01003 m_what_exception='m_setsub::inconsistent boundary'
01004 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01005 endif
01006 if ((in_r_ubound-in_r_lbound+1).ne.m_sub%rows .or. (in_c_ubound-in_c_lbound+1).ne.m_sub%cols) then
01007 m_what_exception='m_setsub::inconsistent boundary between m and m_sub'
01008 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01009 endif
01010 #endif
01011
01012 forall(i=r_lbound:r_ubound,j=c_lbound:c_ubound)
01013 m%ptr_container(i,j)=m_sub%ptr_container(i-r_lbound+1,j-c_lbound+1)
01014 end forall
01015
01016 end subroutine m_setsub
01017
01018
01026 function m_matrixTOvector(m) result (res)
01027 type(matrix), intent(in) :: m
01028
01029 integer :: i, s, e
01030 type(type_exception) :: err_exception;
01031 type(vector) :: res
01032 #ifdef DEBUG_EXCEPTION
01033 if(.not.(m%is_allocate)) then
01034 m_what_exception='m_matrixTOvector::matrix not allocated'
01035 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01036 endif
01037 #endif
01038 call v_init(res,m%rows*m%cols);
01039 do i=1,m%rows
01040 s=(i-1)*m%cols+1
01041 e=i*m%cols
01042 res%ptr_container(s:e)=m%ptr_container(i,:)
01043 end do
01044 end function m_matrixTOvector
01045
01046
01055 function m_get(m,i,j) result (res)
01056 type(matrix), intent(in) :: m
01057 integer, intent(in) :: i, j
01058
01059 type(type_exception) :: err_exception;
01060 type_precision :: res
01061 #ifdef DEBUG_EXCEPTION
01062 if(.not.(m%is_allocate)) then
01063 m_what_exception='m_get::matrix not allocated'
01064 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01065 endif
01066 if ((i .lt. 1) .or. (i .GT. m%rows) .or. (j .lt. 1) .or. (j .GT. m%cols)) then
01067 m_what_exception='m_get::matrix index nil or negative'
01068 err_exception=e_error(stop_array_indice_exceed,m_what_exception,stop_array_indice_exceed)
01069 endif
01070 #endif
01071 res=m%ptr_container(i,j)
01072 end function m_get
01073
01074
01081 function m_get_m(m) result (res)
01082 type(matrix), intent(in) :: m
01083
01084 type(type_exception) :: err_exception;
01085 type_precision :: res(m%rows,m%cols)
01086 #ifdef DEBUG_EXCEPTION
01087 if(.not.(m%is_allocate)) then
01088 m_what_exception='m_get_m::matrix not allocated'
01089 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01090 endif
01091 #endif
01092 res=m%ptr_container
01093 end function m_get_m
01094
01095
01103 function m_getRow(m,i) result (res)
01104 type(matrix), intent(in) :: m
01105 integer, intent(in) :: i
01106
01107 type(type_exception) :: err_exception;
01108 type_precision :: res(m%cols)
01109 #ifdef DEBUG_EXCEPTION
01110 if(.not.(m%is_allocate)) then
01111 m_what_exception='m_getRow::matrix not allocated'
01112 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01113 endif
01114 if ((i .lt. 1) .or. (i .GT. m%rows)) then
01115 m_what_exception='m_getRow::matrix index nil or negative'
01116 err_exception=e_error(stop_array_indice_exceed,m_what_exception,stop_array_indice_exceed)
01117 endif
01118 #endif
01119 res(:)=m%ptr_container(i,:)
01120 end function m_getRow
01121
01122
01130 function m_getCol(m,j) result (res)
01131 type(matrix), intent(in) :: m
01132 integer, intent(in) :: j
01133
01134 type(type_exception) :: err_exception;
01135 type_precision :: res(m%rows)
01136 #ifdef DEBUG_EXCEPTION
01137 if(.not.(m%is_allocate)) then
01138 m_what_exception='m_getCol::matrix not allocated'
01139 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01140 endif
01141 if ((j .lt. 1) .or. (j .GT. m%cols)) then
01142 m_what_exception='m_getCol::matrix index nil or negative'
01143 err_exception=e_error(stop_array_indice_exceed,m_what_exception,stop_array_indice_exceed)
01144 endif
01145 #endif
01146 res=m%ptr_container(:,j)
01147 end function m_getCol
01148
01149
01158 subroutine m_set(m,i,j,value)
01159 type(matrix), intent(inout) :: m
01160 integer, intent(in) :: i, j
01161 type_precision, intent(in) :: value
01162
01163 type(type_exception) :: err_exception;
01164 #ifdef DEBUG_EXCEPTION
01165 if(.not.(m%is_allocate)) then
01166 m_what_exception='m_set::matrix not allocated'
01167 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01168 endif
01169 if ((i .lt. 1) .or. (i .GT. m%rows) .or. (j .lt. 1) .or. (j .GT. m%cols)) then
01170 m_what_exception='m_set::matrix size nil or negative'
01171 err_exception=e_error(stop_array_indice_exceed,m_what_exception,stop_array_indice_exceed)
01172 endif
01173 #endif
01174 m%ptr_container(i,j)=value
01175 end subroutine m_set
01176
01177
01184 subroutine m_affect(m,value)
01185 type(matrix), intent(inout) :: m
01186 type_precision, intent(in) :: value
01187
01188 type(type_exception) :: err_exception;
01189 #ifdef DEBUG_EXCEPTION
01190 if(.not.(m%is_allocate)) then
01191 m_what_exception='m_affect::matrix not allocated'
01192 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01193 endif
01194 #endif
01195 m%ptr_container(:,:)=value
01196 end subroutine m_affect
01197
01198
01199
01200
01201
01202
01209 function m_getSize(m) result (nb)
01210 type(matrix), intent(in) :: m
01211
01212 integer :: nb
01213 type(type_exception) :: err_exception;
01214 #ifdef DEBUG_EXCEPTION
01215 if(.not.(m%is_allocate)) then
01216 m_what_exception='m_getSize::matrix not allocated'
01217 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01218 endif
01219 #endif
01220 nb=m%cols*m%rows
01221 end function m_getSize
01222
01223
01230 function m_getSizeCols(m) result (nb)
01231 type(matrix), intent(in) :: m
01232
01233 integer :: nb
01234 type(type_exception) :: err_exception;
01235 #ifdef DEBUG_EXCEPTION
01236 if(.not.(m%is_allocate)) then
01237 m_what_exception='m_getSizeCols::matrix not allocated'
01238 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01239 endif
01240 #endif
01241 nb=m%cols
01242 end function m_getSizeCols
01243
01244
01251 function m_getSizeRows(m) result (nb)
01252 type(matrix), intent(in) :: m
01253
01254 integer :: nb
01255 type(type_exception) :: err_exception;
01256 #ifdef DEBUG_EXCEPTION
01257 if(.not.(m%is_allocate)) then
01258 m_what_exception='m_getSizeRows::matrix not allocated'
01259 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01260 endif
01261 #endif
01262 nb=m%rows
01263 end function m_getSizeRows
01264
01265
01266
01273 function m_nbnegative(m) result (nb)
01274 type(matrix), intent(in) :: m
01275
01276 integer :: nb
01277 type(type_exception) :: err_exception;
01278 #ifdef DEBUG_EXCEPTION
01279 if(.not.(m%is_allocate)) then
01280 m_what_exception='m_nbnegative::matrix not allocated'
01281 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01282 endif
01283 #endif
01284 nb = count(m%ptr_container<0);
01285 end function m_nbnegative
01286
01287
01295 function m_nbnegativeRow(m,i) result (nb)
01296 type(matrix), intent(in) :: m
01297 integer, intent(in) :: i
01298
01299 integer :: nb
01300 type(type_exception) :: err_exception;
01301 #ifdef DEBUG_EXCEPTION
01302 if(.not.(m%is_allocate)) then
01303 m_what_exception='m_nbnegativeRow::matrix not allocated'
01304 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01305 endif
01306 if ((i .lt. 1) .or. (i .GT. m%rows)) then
01307 m_what_exception='m_nbnegativeRow::row is nil or negative'
01308 err_exception=e_error(stop_array_indice_exceed,m_what_exception,stop_array_indice_exceed)
01309 endif
01310 #endif
01311 nb = count(m%ptr_container(i,:)<0);
01312 end function m_nbnegativeRow
01313
01314
01322 function m_nbnegativeCol(m,j) result (nb)
01323 type(matrix), intent(in) :: m
01324 integer, intent(in) :: j
01325
01326 integer :: nb
01327 type(type_exception) :: err_exception;
01328 #ifdef DEBUG_EXCEPTION
01329 if(.not.(m%is_allocate)) then
01330 m_what_exception='m_nbnegativeCol::matrix not allocated'
01331 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01332 endif
01333 if ((j .lt. 1) .or. (j .GT. m%cols)) then
01334 m_what_exception='m_nbnegativeCol::column nil or negative'
01335 err_exception=e_error(stop_array_indice_exceed,m_what_exception,stop_array_indice_exceed)
01336 endif
01337 #endif
01338 nb = count(m%ptr_container(:,j)<0);
01339 end function m_nbnegativeCol
01340
01341
01348 function m_nbpositive(m) result (nb)
01349 type(matrix), intent(in) :: m
01350
01351 integer :: nb
01352 type(type_exception) :: err_exception;
01353 #ifdef DEBUG_EXCEPTION
01354 if(.not.(m%is_allocate)) then
01355 m_what_exception='m_nbpositive::matrix not allocated'
01356 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01357 endif
01358 #endif
01359 nb= count(m%ptr_container>0);
01360 end function m_nbpositive
01361
01362
01370 function m_nbpositiveRow(m,i) result (nb)
01371 type(matrix), intent(in) :: m
01372 integer, intent(in) :: i
01373
01374 integer :: nb
01375 type(type_exception) :: err_exception;
01376 #ifdef DEBUG_EXCEPTION
01377 if(.not.(m%is_allocate)) then
01378 m_what_exception='m_nbpositiveRow::matrix not allocated'
01379 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01380 endif
01381 if ((i .lt. 1) .or. (i .GT. m%rows)) then
01382 m_what_exception='m_nbpositiveRow::row is nil or negative'
01383 err_exception=e_error(stop_array_indice_exceed,m_what_exception,stop_array_indice_exceed)
01384 endif
01385 #endif
01386 nb = count(m%ptr_container(i,:)>0);
01387 end function m_nbpositiveRow
01388
01389
01390
01398 function m_nbpositiveCol(m,j) result (nb)
01399 type(matrix), intent(in) :: m
01400 integer, intent(in) :: j
01401
01402 integer :: nb
01403 type(type_exception) :: err_exception;
01404 #ifdef DEBUG_EXCEPTION
01405 if(.not.(m%is_allocate)) then
01406 m_what_exception='m_nbpositiveCol::matrix not allocated'
01407 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01408 endif
01409 if ((j .lt. 1) .or. (j .GT. m%cols)) then
01410 m_what_exception='m_nbpositiveCol::column nil or negative'
01411 err_exception=e_error(stop_array_indice_exceed,m_what_exception,stop_array_indice_exceed)
01412 endif
01413 #endif
01414 nb = count(m%ptr_container(:,j)>0);
01415 end function m_nbpositiveCol
01416
01417
01424 function m_nbzeros(m) result (nb)
01425 type(matrix), intent(in) :: m
01426
01427 integer :: nb
01428 type(type_exception) :: err_exception;
01429 #ifdef DEBUG_EXCEPTION
01430 if(.not.(m%is_allocate)) then
01431 m_what_exception='m_nbzeros::matrix not allocated'
01432 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01433 endif
01434 #endif
01435 nb = count(m%ptr_container==0);
01436 end function m_nbzeros
01437
01438
01446 function m_nbzerosCol(m,j) result (nb)
01447 type(matrix), intent(in) :: m
01448 integer, intent(in) :: j
01449
01450 integer :: nb
01451 type(type_exception) :: err_exception;
01452 #ifdef DEBUG_EXCEPTION
01453 if(.not.(m%is_allocate)) then
01454 m_what_exception='m_nbzerosCol::matrix not allocated'
01455 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01456 endif
01457 if ((j .lt. 1) .or. (j .GT. m%cols)) then
01458 m_what_exception='m_nbzerosCol::column nil or negative'
01459 err_exception=e_error(stop_array_indice_exceed,m_what_exception,stop_array_indice_exceed)
01460 endif
01461 #endif
01462 nb = count(m%ptr_container(:,j)==0);
01463 end function m_nbzerosCol
01464
01465
01473 function m_nbzerosRow(m,i) result (nb)
01474 type(matrix), intent(in) :: m
01475 integer, intent(in) :: i
01476
01477 integer :: nb
01478 type(type_exception) :: err_exception;
01479 #ifdef DEBUG_EXCEPTION
01480 if(.not.(m%is_allocate)) then
01481 m_what_exception='m_nbzerosRow::matrix not allocated'
01482 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01483 endif
01484 if ((i .lt. 1) .or. (i .GT. m%rows)) then
01485 m_what_exception='m_nbzerosRow::row nil or negative'
01486 err_exception=e_error(stop_array_indice_exceed,m_what_exception,stop_array_indice_exceed)
01487 endif
01488 #endif
01489 nb = count(m%ptr_container(i,:)==0);
01490 end function m_nbzerosRow
01491
01492
01499 function m_max(m) result (val)
01500 type(matrix), intent(in) :: m
01501
01502 type_precision :: val
01503 type(type_exception) :: err_exception;
01504 #ifdef DEBUG_EXCEPTION
01505 if(.not.(m%is_allocate)) then
01506 m_what_exception='m_max::matrix not allocated'
01507 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01508 endif
01509 #endif
01510 val= maxval(m%ptr_container);
01511 end function m_max
01512
01513
01521 function m_maxRow(m,i) result (val)
01522 type(matrix), intent(in) :: m
01523 integer, intent(in) :: i
01524
01525 type_precision :: val
01526 type(type_exception) :: err_exception;
01527 #ifdef DEBUG_EXCEPTION
01528 if(.not.(m%is_allocate)) then
01529 m_what_exception='m_maxRow::matrix not allocated'
01530 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01531 endif
01532 if ((i .lt. 1) .or. (i .GT. m%rows)) then
01533 m_what_exception='m_maxRow::row is nil or negative'
01534 err_exception=e_error(stop_array_indice_exceed,m_what_exception,stop_array_indice_exceed)
01535 endif
01536 #endif
01537 val = maxval(m%ptr_container(i,:));
01538 end function m_maxRow
01539
01540
01541
01549 function m_maxCol(m,j) result (val)
01550 type(matrix), intent(in) :: m
01551 integer, intent(in) :: j
01552
01553 type_precision :: val
01554 type(type_exception) :: err_exception;
01555 #ifdef DEBUG_EXCEPTION
01556 if(.not.(m%is_allocate)) then
01557 m_what_exception='m_maxCol::matrix not allocated'
01558 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01559 endif
01560 if ((j .lt. 1) .or. (j .GT. m%cols)) then
01561 m_what_exception='m_maxCol::column nil or negative'
01562 err_exception=e_error(stop_array_indice_exceed,m_what_exception,stop_array_indice_exceed)
01563 endif
01564 #endif
01565 val = maxval(m%ptr_container(:,j));
01566 end function m_maxCol
01567
01568
01569
01576 function m_min(m) result (val)
01577 type(matrix), intent(in) :: m
01578
01579 type_precision :: val
01580 type(type_exception) :: err_exception;
01581 #ifdef DEBUG_EXCEPTION
01582 if(.not.(m%is_allocate)) then
01583 m_what_exception='m_min::matrix not allocated'
01584 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01585 endif
01586 #endif
01587 val= minval(m%ptr_container);
01588 end function m_min
01589
01590
01598 function m_minRow(m,i) result (val)
01599 type(matrix), intent(in) :: m
01600 integer, intent(in) :: i
01601
01602 type_precision :: val
01603 type(type_exception) :: err_exception;
01604 #ifdef DEBUG_EXCEPTION
01605 if(.not.(m%is_allocate)) then
01606 m_what_exception='m_minRow::matrix not allocated'
01607 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01608 endif
01609 if ((i .lt. 1) .or. (i .GT. m%rows)) then
01610 m_what_exception='m_minRow::row is nil or negative'
01611 err_exception=e_error(stop_array_indice_exceed,m_what_exception,stop_array_indice_exceed)
01612 endif
01613 #endif
01614 val = minval(m%ptr_container(i,:));
01615 end function m_minRow
01616
01617
01618
01619
01627 function m_minCol(m,j) result (val)
01628 type(matrix), intent(in) :: m
01629 integer, intent(in) :: j
01630
01631 type_precision :: val
01632 type(type_exception) :: err_exception;
01633 #ifdef DEBUG_EXCEPTION
01634 if(.not.(m%is_allocate)) then
01635 m_what_exception='m_minCol::matrix not allocated'
01636 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01637 endif
01638 if ((j .lt. 1) .or. (j .GT. m%cols)) then
01639 m_what_exception='m_minCol::column nil or negative'
01640 err_exception=e_error(stop_array_indice_exceed,m_what_exception,stop_array_indice_exceed)
01641 endif
01642 #endif
01643 val = minval(m%ptr_container(:,j));
01644 end function m_minCol
01645
01646
01653 function m_sum(m) result (val)
01654 type(matrix), intent(in) :: m
01655
01656 type_precision :: val
01657 type(type_exception) :: err_exception;
01658 #ifdef DEBUG_EXCEPTION
01659 if(.not.(m%is_allocate)) then
01660 m_what_exception='m_sum::matrix not allocated'
01661 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01662 endif
01663 #endif
01664 val= sum(m%ptr_container);
01665 end function m_sum
01666
01667
01675 function m_sumRow(m,i) result (val)
01676 type(matrix), intent(in) :: m
01677 integer, intent(in) :: i
01678
01679 type_precision :: val
01680 type(type_exception) :: err_exception;
01681 #ifdef DEBUG_EXCEPTION
01682 if(.not.(m%is_allocate)) then
01683 m_what_exception='m_sumRow::matrix not allocated'
01684 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01685 endif
01686 if ((i .lt. 1) .or. (i .GT. m%rows)) then
01687 m_what_exception='m_sumRow::row is nil or negative'
01688 err_exception=e_error(stop_array_indice_exceed,m_what_exception,stop_array_indice_exceed)
01689 endif
01690 #endif
01691 val = sum(m%ptr_container(i,:));
01692 end function m_sumRow
01693
01694
01702 function m_sumCol(m,j) result (val)
01703 type(matrix), intent(in) :: m
01704 integer, intent(in) :: j
01705
01706 type_precision :: val
01707 type(type_exception) :: err_exception;
01708 #ifdef DEBUG_EXCEPTION
01709 if(.not.(m%is_allocate)) then
01710 m_what_exception='m_sumCol::matrix not allocated'
01711 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01712 endif
01713 if ((j .lt. 1) .or. (j .GT. m%cols)) then
01714 m_what_exception='m_sumCol::column nil or negative'
01715 err_exception=e_error(stop_array_indice_exceed,m_what_exception,stop_array_indice_exceed)
01716 endif
01717 #endif
01718 val = sum(m%ptr_container(:,j));
01719 end function m_sumCol
01720
01721
01722
01723
01724
01725
01726
01727
01728
01729
01730
01731
01732
01740 function m_prod_mat(m1,m2) result (res)
01741 implicit none
01742 type(matrix),intent(in) :: m1,m2
01743 type(matrix) :: res
01744
01745 type(type_exception) :: err_exception;
01746 type_precision :: sum_tmp = 0
01747 integer :: i,j,k
01748 #ifdef DEBUG_EXCEPTION
01749 if(.not.(m1%is_allocate) .or. .not.(m2%is_allocate)) then
01750 m_what_exception='m_prod_mat::matrix not allocated'
01751 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01752 endif
01753 if (m1%cols.ne.m2%rows) then
01754 m_what_exception='m_prod_mat::The number cols of the first matrix must be equal to the number rows of the second matrix'
01755 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
01756 endif
01757 #endif
01758 call m_init(res,m1%rows,m2%cols)
01759
01760 do i=1,m1%rows
01761 do j=1,m2%cols
01762 sum_tmp=0
01763 do k=1,m1%cols
01764 sum_tmp = sum_tmp + m1%ptr_container(i,k)*m2%ptr_container(k,j)
01765 end do
01766 res%ptr_container(i,j)=sum_tmp
01767 end do
01768 end do
01769
01770 end function m_prod_mat
01771
01772
01780 function m_prod_vec1(m,v) result (res)
01781 implicit none
01782 type(matrix),intent(in) :: m
01783 type(vector),intent(in) :: v
01784 type(vector) :: res
01785
01786 type(type_exception) :: err_exception;
01787 integer :: i,j
01788 type_precision :: sum_tmp = 0
01789 #ifdef DEBUG_EXCEPTION
01790 if(.not.(m%is_allocate) .or. .not.(v%is_allocate)) then
01791 m_what_exception='m_prod_vec1::matrix not allocated'
01792 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01793 endif
01794 if (m%cols.ne.v%size) then
01795 m_what_exception='m_prod_vec1::The number cols of the first matrix must be equal to vector size'
01796 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
01797 endif
01798 #endif
01799 call v_init(res,m%rows)
01800
01801 do i=1,m%rows
01802 sum_tmp=0
01803 do j=1,v%size
01804 sum_tmp = sum_tmp + m%ptr_container(i,j)*v%ptr_container(j)
01805 end do
01806 res%ptr_container(i)=sum_tmp
01807 end do
01808
01809
01810 end function m_prod_vec1
01811
01812
01821 function m_prod_vec2(v,m) result (res)
01822 implicit none
01823 type(vector),intent(in) :: v
01824 type(matrix),intent(in) :: m
01825 type(vector) :: res
01826
01827 type(type_exception) :: err_exception;
01828 integer :: i,j
01829 type_precision :: sum_tmp = 0
01830 #ifdef DEBUG_EXCEPTION
01831 if(.not.(m%is_allocate) .or. .not.(v%is_allocate)) then
01832 m_what_exception='m_prod_vec2::matrix not allocated'
01833 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01834 endif
01835 if (v%size.ne.m%rows) then
01836 m_what_exception='m_prod_vec2::matrix size must be equal to matrix row'
01837 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
01838 endif
01839 #endif
01840 call v_init(res,m%cols)
01841
01842 do j=1,m%cols
01843 sum_tmp=0
01844 do i=1,v%size
01845 sum_tmp = sum_tmp + v%ptr_container(j)*m%ptr_container(i,j)
01846 end do
01847 res%ptr_container(j) = sum_tmp
01848 end do
01849
01850 end function m_prod_vec2
01851
01852
01861 function m_prod_vec_c(m,v) result (res)
01862 implicit none
01863 type(matrix),intent(in) :: m
01864 type(vector),intent(in) :: v
01865 type(vector) :: res
01866
01867 type(type_exception) :: err_exception;
01868 integer :: k
01869 #ifdef DEBUG_EXCEPTION
01870 if(.not.(m%is_allocate) .or. .not.(v%is_allocate)) then
01871 m_what_exception='m_prod_vec_c::matrix not allocated'
01872 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01873 endif
01874 if (m%cols.ne.v%size) then
01875 m_what_exception='m_prod_vec_c::The number cols of the first matrix must be equal to vector size'
01876 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
01877 endif
01878 #endif
01879 call v_init(res,m%rows)
01880
01881 do k=1,v%size
01882 res%ptr_container(:) = res%ptr_container(:) + m%ptr_container(:,k)*v%ptr_container(k)
01883 end do
01884
01885 end function m_prod_vec_c
01886
01887
01895 function m_prod_scalar1(alpha,m) result (res)
01896 implicit none
01897 type_precision,intent(in) :: alpha
01898 type(matrix),intent(in) :: m
01899 type(matrix) :: res
01900
01901 type(type_exception) :: err_exception;
01902 #ifdef DEBUG_EXCEPTION
01903 if(.not.(m%is_allocate)) then
01904 m_what_exception='m_prod_scalar1::matrix not allocated'
01905 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01906 endif
01907 #endif
01908 call m_init(res,m%rows,m%cols)
01909 res%ptr_container=alpha*m%ptr_container
01910 end function m_prod_scalar1
01911
01912
01920 function m_prod_scalar2(m,alpha) result (res)
01921 implicit none
01922 type(matrix),intent(in) :: m
01923 type_precision,intent(in) :: alpha
01924 type(matrix) :: res
01925
01926 type(type_exception) :: err_exception;
01927 #ifdef DEBUG_EXCEPTION
01928 if(.not.(m%is_allocate)) then
01929 m_what_exception='m_prod_scalar2::matrix not allocated'
01930 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01931 endif
01932 #endif
01933 call m_init(res,m%rows,m%cols)
01934 res%ptr_container=m%ptr_container*alpha
01935 end function m_prod_scalar2
01936
01937
01945 function m_div_scalar(m,alpha) result (res)
01946 implicit none
01947 type_precision,intent(in) :: alpha
01948 type(matrix),intent(in) :: m
01949 type(matrix) :: res
01950
01951 type(type_exception) :: err_exception;
01952 #ifdef DEBUG_EXCEPTION
01953 if(.not.(m%is_allocate)) then
01954 m_what_exception='m_div_scalar::matrix not allocated'
01955 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01956 endif
01957 if(alpha.eq.0) then
01958 m_what_exception='m_div_scalar::can not divide by zero'
01959 err_exception=e_error(stop_div0,m_what_exception,stop_div0)
01960 end if
01961 #endif
01962 call m_init(res,m%rows,m%cols)
01963 res%ptr_container=m%ptr_container/alpha
01964 end function m_div_scalar
01965
01966
01967
01968
01969
01970
01978 function m_add(m1,m2) result (res)
01979 implicit none
01980 type(matrix),intent(in) :: m1,m2
01981 type(matrix) :: res
01982
01983 type(type_exception) :: err_exception;
01984 #ifdef DEBUG_EXCEPTION
01985 if(.not.(m1%is_allocate) .or. .not.(m2%is_allocate)) then
01986 m_what_exception='m_add::matrix not allocated'
01987 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01988 endif
01989 if ((m1%rows.ne.m2%rows) .or. (m1%cols.ne.m2%cols)) then
01990 m_what_exception='m_add::matrices are not the same size'
01991 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
01992 endif
01993 #endif
01994 call m_init(res,m1%rows,m1%cols)
01995 res%ptr_container=m1%ptr_container+m2%ptr_container
01996 end function m_add
01997
01998
01999
02000
02001
02002
02010 function m_minus(m1,m2) result (res)
02011 implicit none
02012 type(matrix),intent(in) :: m1,m2
02013 type(matrix) :: res
02014
02015 type(type_exception) :: err_exception;
02016 #ifdef DEBUG_EXCEPTION
02017 if(.not.(m1%is_allocate) .or. .not.(m2%is_allocate)) then
02018 m_what_exception='m_minus::matrix not allocated'
02019 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
02020 endif
02021 if ((m1%rows.ne.m2%rows) .or. (m1%cols.ne.m2%cols)) then
02022 m_what_exception='m_minus::matrices are not the same size'
02023 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
02024 endif
02025 #endif
02026 call m_init(res,m1%rows,m1%cols)
02027 res%ptr_container=m1%ptr_container-m2%ptr_container
02028 end function m_minus
02029
02030
02031
02032
02033
02034
02035
02042 function m_trans(m) result (res)
02043
02044 type(matrix), intent(in) :: m
02045
02046 type(type_exception) :: err_exception;
02047 type(matrix) :: res
02048 #ifdef DEBUG_EXCEPTION
02049 if(.not.(m%is_allocate)) then
02050 m_what_exception='m_trans::matrix not allocated'
02051 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
02052 endif
02053 #endif
02054 call m_init(res,m%cols,m%rows)
02055 res%ptr_container = transpose(m%ptr_container)
02056 end function m_trans
02057
02058
02059
02060
02061
02062
02070 function m_isEqual(m1,m2) result (res)
02071 type(matrix), intent(in) :: m1, m2
02072
02073 type(type_exception) :: err_exception;
02074 logical :: res
02075 #ifdef DEBUG_EXCEPTION
02076 if(.not.(m1%is_allocate) .or. .not.(m2%is_allocate)) then
02077 m_what_exception='m_isEqual::matrix not allocated'
02078 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
02079 endif
02080 if ((m1%rows.ne.m2%rows) .or. (m1%cols.ne.m2%cols)) then
02081 m_what_exception='m_isEqual::matrices are not the same size'
02082 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
02083 endif
02084 #endif
02085 res=all(m1%ptr_container.eq.m2%ptr_container);
02086 end function m_isEqual
02087
02088
02096 function m_isEqual_scalar(m,val) result (res)
02097 type(matrix), intent(in) :: m
02098 type_precision, intent(in) :: val
02099
02100 type(type_exception) :: err_exception;
02101 logical :: res
02102 res = .false.
02103 #ifdef DEBUG_EXCEPTION
02104 if(.not.(m%is_allocate)) then
02105 m_what_exception='m_isEqual_scalar::matrix not allocated'
02106 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
02107 endif
02108 #endif
02109 if( count(m%ptr_container.ne.val)==0) res = .true.
02110 end function m_isEqual_scalar
02111
02112
02119 function m_isSymmetric(m) result (res)
02120 type(matrix), intent(in) :: m
02121
02122 type(type_exception) :: err_exception;
02123 logical :: res
02124 type(matrix) :: tmp
02125 #ifdef DEBUG_EXCEPTION
02126 if(.not.(m%is_allocate)) then
02127 m_what_exception=' m_isSymmetric::matrix not allocated'
02128 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
02129 endif
02130 #endif
02131 tmp=m_trans(m);
02132 res=all(m%ptr_container.eq.tmp%ptr_container);
02133 end function m_isSymmetric
02134
02135
02136
02137
02138
02139
02140
02147 function m_identity(n) result(res)
02148 implicit none
02149 integer, intent(in) :: n;
02150
02151 type(matrix) ::res;
02152 type(type_exception) :: err_exception;
02153 integer(kind=2) :: ierr;
02154 integer :: i
02155 #ifdef DEBUG_EXCEPTION
02156 if (n.lt.1) then
02157 m_what_exception='m_identity::matrix size is nil or negative'
02158 err_exception=e_error(stop_array_indice_exceed,m_what_exception,stop_array_indice_exceed)
02159 endif
02160 #endif
02161 res%rows=n
02162 res%cols=n
02163 allocate(res%ptr_container(n,n),stat=ierr)
02164 #ifdef DEBUG_EXCEPTION
02165 if (ierr .ne. 0) then
02166 m_what_exception='m_identity::allocate error'
02167 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
02168 endif
02169 #endif
02170 res%ptr_container=0
02171 forall(i=1:res%rows) res%ptr_container(i,i)=1
02172 res%is_allocate = allocated(res%ptr_container)
02173 end function m_identity
02174
02175
02176
02186 function mc_tridiag(d,u,l,size_n) result (res)
02187 type_precision, intent(in) :: d;
02188 type_precision, intent(in) :: u;
02189 type_precision, intent(in) :: l;
02190 integer, intent(in) :: size_n;
02191
02192 type(matrix) ::res;
02193 type(type_exception) :: err_exception;
02194 integer(kind=2) :: ierr;
02195 integer :: i,j
02196 #ifdef DEBUG_EXCEPTION
02197 if (size_n.lt.2) then
02198 m_what_exception='mc_tridiag::size vector must be superior or equal 2'
02199 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
02200 endif
02201 #endif
02202 res%rows=size_n
02203 res%cols=size_n
02204 allocate(res%ptr_container(res%rows,res%cols),stat=ierr)
02205 #ifdef DEBUG_EXCEPTION
02206 if (ierr .ne. 0) then
02207 m_what_exception='mc_tridiag::allocate error'
02208 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
02209 endif
02210 #endif
02211 res%ptr_container(:,:)=0
02212 do i=1,res%rows
02213 do j=1,res%cols
02214 if(i.eq.(j-1)) res%ptr_container(i,j) = l
02215 if(i.eq.j) res%ptr_container(i,j) = d
02216 if(i.eq.(j+1)) res%ptr_container(i,j) = u
02217 end do
02218 end do
02219 res%is_allocate = allocated(res%ptr_container)
02220 end function mc_tridiag
02221
02222
02231 function mc_bidiag_up(d,u,size_n) result (res)
02232 type_precision, intent(in) :: d;
02233 type_precision, intent(in) :: u;
02234 integer, intent(in) :: size_n;
02235
02236 type(matrix) ::res;
02237 type(type_exception) :: err_exception;
02238 integer(kind=2) :: ierr;
02239 integer :: i,j
02240 #ifdef DEBUG_EXCEPTION
02241 if (size_n.lt.2) then
02242 m_what_exception='mc_bidiag_up::size vector must be superior or equal 2'
02243 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
02244 endif
02245 #endif
02246 res%rows=size_n
02247 res%cols=size_n
02248 allocate(res%ptr_container(res%rows,res%cols),stat=ierr)
02249 #ifdef DEBUG_EXCEPTION
02250 if (ierr .ne. 0) then
02251 m_what_exception='mc_bidiag_up::allocate error'
02252 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
02253 endif
02254 #endif
02255 res%ptr_container(:,:)=0
02256 do i=1,res%rows
02257 do j=1,res%cols
02258 if(i.eq.j) res%ptr_container(i,j) = d
02259 if(i.eq.(j+1)) res%ptr_container(i,j) = u
02260 end do
02261 end do
02262 res%is_allocate = allocated(res%ptr_container)
02263 end function mc_bidiag_up
02264
02265
02274 function mc_bidiag_low(d,l,size_n) result (res)
02275 type_precision, intent(in) :: d;
02276 type_precision, intent(in) :: l;
02277 integer, intent(in) :: size_n;
02278
02279 type(matrix) ::res;
02280 type(type_exception) :: err_exception;
02281 integer(kind=2) :: ierr;
02282 integer :: i,j
02283 #ifdef DEBUG_EXCEPTION
02284 if (size_n.lt.2) then
02285 m_what_exception='mc_bidiag_low::size vector must be superior or equal 2'
02286 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
02287 endif
02288 #endif
02289 res%rows=size_n
02290 res%cols=size_n
02291 allocate(res%ptr_container(res%rows,res%cols),stat=ierr)
02292 #ifdef DEBUG_EXCEPTION
02293 if (ierr .ne. 0) then
02294 m_what_exception='mc_bidiag_low::allocate error'
02295 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
02296 endif
02297 #endif
02298 res%ptr_container(:,:)=0
02299 do i=1,res%rows
02300 do j=1,res%cols
02301 if(i.eq.j) res%ptr_container(i,j) = d
02302 if(i.eq.(j-1)) res%ptr_container(i,j) = l
02303 end do
02304 end do
02305 res%is_allocate = allocated(res%ptr_container)
02306 end function mc_bidiag_low
02307
02308
02315 function mc_diag(vect_data) result (res)
02316 type_precision :: vect_data(:);
02317
02318 type(matrix) ::res;
02319 type(type_exception) :: err_exception;
02320 integer(kind=2) :: ierr;
02321 integer :: i,j
02322 integer :: size_vect
02323
02324 size_vect=size(vect_data)
02325 #ifdef DEBUG_EXCEPTION
02326 if (size_vect.lt.1) then
02327 m_what_exception='mc_diag::size vector must be superior or equal 1'
02328 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
02329 endif
02330 #endif
02331 res%rows=size_vect
02332 res%cols=size_vect
02333 allocate(res%ptr_container(res%rows,res%cols),stat=ierr)
02334 #ifdef DEBUG_EXCEPTION
02335 if (ierr .ne. 0) then
02336 m_what_exception='allocate error'
02337 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
02338 endif
02339 #endif
02340 res%ptr_container(:,:)=0
02341 do i=1,res%rows
02342 do j=1,res%cols
02343 if(i.eq.j) res%ptr_container(i,j) = vect_data(i)
02344 end do
02345 end do
02346 res%is_allocate = allocated(res%ptr_container)
02347 end function mc_diag
02348
02349
02357 function m_tril(m,swap_diag) result (res)
02358 type(matrix), intent(in) :: m
02359 integer, intent(in), optional :: swap_diag
02360
02361 type(matrix) :: res
02362 integer :: i,j
02363 integer :: in_sd
02364 type(type_exception) :: err_exception;
02365 #ifdef DEBUG_EXCEPTION
02366 if(.not.(m%is_allocate)) then
02367 m_what_exception='m_tril::matrix not allocated'
02368 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
02369 endif
02370 #endif
02371 if(present(swap_diag)) then
02372 in_sd=swap_diag;
02373 else
02374 in_sd=0;
02375 end if
02376
02377 call m_init(res,m%rows,m%cols)
02378 do i=1,m%rows
02379 do j=1,m%cols
02380 if(j.le.(i+in_sd)) res%ptr_container(i,j) = m%ptr_container(i,j)
02381 end do
02382 end do
02383 end function m_tril
02384
02385
02393 function m_triu(m,swap_diag) result (res)
02394 type(matrix), intent(in) :: m
02395 integer, intent(in), optional :: swap_diag
02396
02397 type(matrix) :: res
02398 integer :: i,j
02399 integer :: in_sd
02400 type(type_exception) :: err_exception;
02401 #ifdef DEBUG_EXCEPTION
02402 if(.not.(m%is_allocate)) then
02403 m_what_exception='m_triu::matrix not allocated'
02404 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
02405 endif
02406 #endif
02407 if(present(swap_diag)) then
02408 in_sd=swap_diag;
02409 else
02410 in_sd=0;
02411 end if
02412 call m_init(res,m%rows,m%cols)
02413 do i=1,m%rows
02414 do j=1,m%cols
02415 if(i.le.(j-in_sd)) res%ptr_container(i,j) = m%ptr_container(i,j)
02416 end do
02417 end do
02418 end function m_triu
02419
02420
02427 function m_bidiag_low(m) result (res)
02428 type(matrix), intent(in) :: m
02429
02430 type(matrix) :: res
02431 integer :: i,j
02432 type(type_exception) :: err_exception;
02433 #ifdef DEBUG_EXCEPTION
02434 if(.not.(m%is_allocate)) then
02435 m_what_exception='m_bidiag_low::matrix not allocated'
02436 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
02437 endif
02438 #endif
02439 call m_init(res,m%rows,m%cols)
02440 do i=1,m%rows
02441 do j=1,m%cols
02442 if( i.eq.j .or. j.eq.(i-1) ) res%ptr_container(i,j) = m%ptr_container(i,j)
02443 end do
02444 end do
02445 end function m_bidiag_low
02446
02447
02454 function m_bidiag_up(m) result (res)
02455 type(matrix), intent(in) :: m
02456
02457 type(matrix) :: res
02458 integer :: i,j
02459 type(type_exception) :: err_exception;
02460 #ifdef DEBUG_EXCEPTION
02461 if(.not.(m%is_allocate)) then
02462 m_what_exception='m_bidiag_up::matrix not allocated'
02463 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
02464 endif
02465 #endif
02466 call m_init(res,m%rows,m%cols)
02467 do i=1,m%rows
02468 do j=1,m%cols
02469 if( i.eq.j .or. j.eq.(i+1) ) res%ptr_container(i,j) = m%ptr_container(i,j)
02470 end do
02471 end do
02472 end function m_bidiag_up
02473
02474
02475
02476
02485 function m_diag(m,i) result(res)
02486 implicit none
02487 type(matrix), intent(in) :: m
02488 integer, intent(in), optional :: i
02489
02490 type(type_exception) :: err_exception;
02491 type(vector) :: res
02492 integer :: size_diag
02493 integer :: k, t1=0, t2=0
02494 integer :: in_i
02495 #ifdef DEBUG_EXCEPTION
02496 if(.not.(m%is_allocate)) then
02497 m_what_exception='m_diag::matrix not allocated'
02498 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
02499 endif
02500 #endif
02501 if(.not.present(i)) then
02502 in_i=0
02503 else
02504 in_i=i
02505 end if
02506 #ifdef DEBUG_EXCEPTION
02507 if ((in_i .le. -m%rows) .or. (in_i .ge. m%cols)) then
02508 m_what_exception='m_diag::diag is not defined'
02509 err_exception=e_error(stop_array_diag_compatible,m_what_exception,stop_array_diag_compatible)
02510 endif
02511 #endif
02512 if(m%rows.ge.m%cols) then
02513 if(in_i.le.0 .and. abs(in_i).le.(m%rows-m%cols)) then
02514 size_diag=min(m%rows,m%cols)
02515 else if(in_i.gt.0) then
02516 size_diag=(min(m%rows,m%cols))-in_i
02517 else
02518 size_diag=(min(m%rows,m%cols)+(m%rows-m%cols))+in_i
02519 end if
02520 else
02521 if(in_i.ge.0 .and. abs(in_i).le.(m%cols-m%rows)) then
02522 size_diag=min(m%rows,m%cols)
02523 else if(in_i.lt.0) then
02524 size_diag=(min(m%rows,m%cols))+in_i
02525 else
02526 size_diag=(min(m%rows,m%cols)+(m%cols-m%rows))-in_i
02527 end if
02528 end if
02529
02530 call v_init(res,size_diag)
02531
02532 if(in_i.gt.0) then
02533 t1=0; t2=1
02534 else if(in_i.lt.0) then
02535 t1=1; t2=0
02536 else
02537 t1=0; t2=0
02538 end if
02539
02540 forall (k=1:size_diag) res%ptr_container(k)=m%ptr_container(t1*abs(in_i)+k,t2*abs(in_i)+k)
02541 end function m_diag
02542
02543
02544
02545
02546
02547
02558 function m_det(m,meth_det,is_permuted) result(det)
02559 implicit none
02560 type(matrix), intent(in) :: m
02561 character*(*), intent(in), optional :: meth_det
02562 logical, intent(in), optional :: is_permuted
02563 type_precision :: det
02564
02565 if(present(meth_det) .and. meth_det.eq.'lu') then
02566 det=m_det_lu(m,is_permuted)
02567 else if(present(meth_det) .and. meth_det.eq.'chol') then
02568 det=m_det_chol(m,is_permuted)
02569 else
02570 det=m_det_gaussj(m)
02571 end if
02572
02573 end function m_det
02574
02575
02576
02584 function m_det_lu(m,is_permuted) result(res)
02585 implicit none
02586 type(matrix), intent(in) :: m
02587 logical, intent(in), optional :: is_permuted
02588
02589 type(type_exception) :: err_exception;
02590 type_precision :: res
02591 type(t_m_and_p) :: decomp_lu_m
02592 #ifdef DEBUG_EXCEPTION
02593 if(.not.(m%is_allocate)) then
02594 m_what_exception='m_det_lu::matrix not allocated'
02595 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
02596 endif
02597 if (m%rows.ne.m%cols) then
02598 m_what_exception='m_det_lu::square matrix expected'
02599 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
02600 endif
02601 #endif
02602 decomp_lu_m=m_lu(m,is_permuted);
02603 res=v_prod(diag(decomp_lu_m%M))
02604
02605 if(present(is_permuted) .and. is_permuted) then
02606
02607 res=((-1)**(decomp_lu_m%swops))*res
02608 end if
02609 call destruct(decomp_lu_m)
02610 end function m_det_lu
02611
02612
02620 function m_det_lu_all(m,is_permuted) result(res)
02621 implicit none
02622 type(matrix), intent(in) :: m
02623 logical, intent(in), optional :: is_permuted
02624
02625 type(type_exception) :: err_exception;
02626 type_precision :: res
02627 type(t_lu) :: decomp_lu
02628 #ifdef DEBUG_EXCEPTION
02629 if(.not.(m%is_allocate)) then
02630 m_what_exception='m_det_lu_all::matrix not allocated'
02631 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
02632 endif
02633 #endif
02634 decomp_lu=lu(m,is_permuted);
02635 res=v_prod(diag(decomp_lu%U))
02636
02637 if(present(is_permuted) .and. is_permuted) then
02638
02639 res=((-1)**(decomp_lu%swops))*res
02640 end if
02641 call destruct(decomp_lu)
02642 end function m_det_lu_all
02643
02644
02653 function m_det_chol(m,is_permuted) result(res)
02654 implicit none
02655 type(matrix), intent(in) :: m
02656 logical, intent(in), optional :: is_permuted
02657
02658 type(type_exception) :: err_exception;
02659 type_precision :: res
02660 type(t_m_and_p) :: m_chol
02661 #ifdef DEBUG_EXCEPTION
02662 if(.not.(m%is_allocate)) then
02663 m_what_exception='m_det_chol::matrix not allocated'
02664 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
02665 endif
02666 if (m%rows.ne.m%cols) then
02667 m_what_exception='m_det_chol::square matrix expected'
02668 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
02669 endif
02670 #endif
02671 m_chol=chol(m,is_permuted);
02672 res=(v_prod(diag(m_chol%M)))**2
02673 if(present(is_permuted) .and. is_permuted) then
02674
02675 res=((-1)**(m_chol%swops))*res
02676 end if
02677 call destruct(m_chol)
02678 end function m_det_chol
02679
02680
02694
02695 function m_det_gaussj(m) result(det)
02696 implicit none
02697 type(matrix), intent(in) :: m
02698
02699 type(type_exception) :: err_exception;
02700 type(matrix) :: m_tmp, res
02701 integer :: i,j,k,l,maxi;
02702 type_precision, dimension(:), allocatable :: tmp_vec;
02703 logical :: eq_zero
02704 type_precision :: det, val_tmp, val_tmp1, machine_eps
02705 #ifdef DEBUG_EXCEPTION
02706 if(.not.(m%is_allocate)) then
02707 m_what_exception='m_det_gaussj::matrix not allocated'
02708 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
02709 endif
02710 if (m%rows.ne.m%cols) then
02711 m_what_exception='m_det_gaussj::square matrix expected'
02712 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
02713 endif
02714 #endif
02715 m_tmp=m
02716 allocate(tmp_vec(m%cols))
02717 machine_eps=math_machine_eps()
02718 det=1
02719 i=1
02720 j=1
02721 do while(i.le.m_tmp%rows .and. j.le.m_tmp%cols)
02722
02723 maxi=i
02724 do k=i+1,m_tmp%rows
02725 if(abs(m_tmp%ptr_container(k,j)).gt.abs(m_tmp%ptr_container(maxi,j))) then
02726 maxi=k
02727 eq_zero = .false.
02728 end if
02729 end do
02730
02731 if(abs(m_tmp%ptr_container(maxi,j)).le.machine_eps) then
02732 det=0
02733 i=m_tmp%rows
02734 j=m_tmp%cols
02735 else
02736
02737 if(i.ne.maxi) then
02738
02739 tmp_vec(:)=m_tmp%ptr_container(i,:)
02740 m_tmp%ptr_container(i,:)=m_tmp%ptr_container(maxi,:)
02741 m_tmp%ptr_container(maxi,:)=tmp_vec(:)
02742 det=(-1)*det
02743 end if
02744
02745 val_tmp=(1/m_tmp%ptr_container(i,j))
02746 m_tmp%ptr_container(i,:)=val_tmp*m_tmp%ptr_container(i,:);
02747 det=det/val_tmp
02748
02749 do l=1,m_tmp%rows
02750 if(l.ne.i) then
02751 val_tmp=m_tmp%ptr_container(l,j);
02752 m_tmp%ptr_container(l,:)=m_tmp%ptr_container(l,:)-val_tmp*m_tmp%ptr_container(i,:)
02753 end if
02754 end do
02755 end if
02756 i = i + 1
02757 j = j + 1
02758 end do
02759
02760 deallocate(tmp_vec)
02761 call destruct(res)
02762 call destruct(m_tmp)
02763 end function m_det_gaussj
02764
02765
02766
02767
02768
02769
02778 function m_rank_gaussj(m) result(rank)
02779 implicit none
02780 type(matrix), intent(in) :: m
02781
02782 type(type_exception) :: err_exception;
02783 type(matrix) :: m_tmp
02784 integer :: i,j,k,l,maxi;
02785 type_precision, dimension(:), allocatable :: tmp_vec;
02786 type_precision :: val_tmp, machine_eps
02787 integer :: rank
02788 #ifdef DEBUG_EXCEPTION
02789 if(.not.(m%is_allocate)) then
02790 m_what_exception='m_rank_gaussj::matrix not allocated'
02791 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
02792 endif
02793 #endif
02794 rank=min(m%rows,m%cols)
02795 m_tmp=m
02796 allocate(tmp_vec(m%cols))
02797 machine_eps=math_machine_eps()
02798 i=1
02799 j=1
02800 do while(i.le.m_tmp%rows .and. j.le.m_tmp%cols)
02801
02802
02803 maxi=i
02804 do k=i+1,m_tmp%rows
02805 if(abs(m_tmp%ptr_container(k,j)).gt.abs(m_tmp%ptr_container(maxi,j))) then
02806 maxi=k
02807 end if
02808 end do
02809
02810 if(abs(m_tmp%ptr_container(maxi,j)).gt.machine_eps) then
02811
02812
02813 if(j.eq.m_tmp%cols) then
02814 tmp_vec(:)=m_tmp%ptr_container(i,:)
02815 m_tmp%ptr_container(i,:)=m_tmp%ptr_container(maxi,:)
02816 m_tmp%ptr_container(maxi,:)=tmp_vec(:)
02817 end if
02818
02819 val_tmp=m_tmp%ptr_container(i,j)
02820 m_tmp%ptr_container(i,:)=(1/val_tmp)*m_tmp%ptr_container(i,:);
02821
02822 do l=i+1,m_tmp%rows
02823 val_tmp=m_tmp%ptr_container(l,j);
02824 m_tmp%ptr_container(l,:)=m_tmp%ptr_container(l,:)-val_tmp*m_tmp%ptr_container(i,:)
02825 end do
02826 end if
02827 i = i + 1
02828 j = j + 1
02829 end do
02830 do i=1,min(m%rows,m%cols)
02831 if(abs(m_tmp%ptr_container(i,i)).le.machine_eps) rank=rank-1;
02832 end do
02833
02834
02835 deallocate(tmp_vec)
02836 call destruct(m_tmp)
02837 end function m_rank_gaussj
02838
02839
02852 function m_rank_svd(m,tol_rank,eps_svd,iter_max,meth_qr,eps_gsortho,is_permuted) result(rank)
02853 implicit none
02854 type(matrix), intent(in) :: m
02855 type_precision, intent(in), optional :: tol_rank
02856 type_precision, intent(in), optional :: eps_svd
02857 integer, intent(in), optional :: iter_max
02858 character*(*), intent(in), optional :: meth_qr
02859 type_precision, intent(in), optional :: eps_gsortho
02860 logical, intent(in), optional :: is_permuted
02861
02862 integer :: rank
02863 type_precision :: in_tol_rank
02864 type(vector) :: v_tmp
02865
02866
02867
02868 v_tmp=m_decompsvd_s(m,eps_svd,iter_max,meth_qr,eps_gsortho,is_permuted);
02869 if(present(tol_rank)) then
02870 in_tol_rank=tol_rank;
02871 else
02872 in_tol_rank=p_notcast(1.d-8);
02873 end if
02874
02875 rank=count(v_tmp%ptr_container.gt.in_tol_rank);
02876
02877 call destruct(v_tmp)
02878 end function m_rank_svd
02879
02880
02895 function m_rank(m,tol_rank,meth_rk,eps_svd,iter_max,meth_qr,eps_gsortho,is_permuted) result(rank)
02896 implicit none
02897 type(matrix), intent(in) :: m
02898 type_precision, intent(in), optional :: tol_rank
02899 type_precision, intent(in), optional :: eps_svd
02900 character*(*), intent(in), optional :: meth_rk
02901 character*(*), intent(in), optional :: meth_qr
02902 type_precision, intent(in), optional :: eps_gsortho
02903 integer, intent(in), optional :: iter_max
02904 logical, intent(in), optional :: is_permuted
02905 integer :: rank
02906
02907 if(present(meth_rk) .and. meth_rk.eq.'svd') then
02908 rank=m_rank_svd(m,tol_rank,eps_svd,iter_max,meth_qr,eps_gsortho,is_permuted)
02909 else
02910 rank=m_rank_gaussj(m)
02911 end if
02912
02913 end function m_rank
02914
02915
02916
02917
02918
02919
02920
02928 function m_inverse_gaussj(m) result(res)
02929 implicit none
02930 type(matrix), intent(in) :: m
02931
02932 type(type_exception) :: err_exception;
02933 type(matrix) :: m_tmp, res
02934 integer :: i,j,k,l,maxi;
02935 type_precision, dimension(:), allocatable :: tmp_vec;
02936 type_precision :: val_tmp, machine_eps
02937 #ifdef DEBUG_EXCEPTION
02938 if(.not.(m%is_allocate)) then
02939 m_what_exception='m_inverse_gaussj::matrix not allocated'
02940 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
02941 endif
02942 if (m%rows.ne.m%cols) then
02943 m_what_exception='m_inverse_gaussj::square matrix expected'
02944 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
02945 endif
02946 #endif
02947 res=m_identity(m%rows)
02948 m_tmp=m
02949 machine_eps=math_machine_eps()
02950 allocate(tmp_vec(m%cols))
02951
02952 i=1
02953 j=1
02954 do while(i.le.m_tmp%rows .and. j.le.m_tmp%cols)
02955
02956 maxi=i
02957 do k=i+1,m_tmp%rows
02958 if(abs(m_tmp%ptr_container(k,j)).gt.abs(m_tmp%ptr_container(maxi,j))) then
02959 maxi=k
02960 end if
02961 end do
02962
02963 if(abs(m_tmp%ptr_container(maxi,j)).le.machine_eps) then
02964 m_what_exception='m_inverse_gaussj::matrix is close to singular or badly scaled'
02965 err_exception=e_error(stop_array_singular,m_what_exception,stop_array_singular)
02966 else
02967
02968
02969 tmp_vec(:)=m_tmp%ptr_container(i,:)
02970 m_tmp%ptr_container(i,:)=m_tmp%ptr_container(maxi,:)
02971 m_tmp%ptr_container(maxi,:)=tmp_vec(:)
02972
02973 tmp_vec(:)=res%ptr_container(i,:)
02974 res%ptr_container(i,:)=res%ptr_container(maxi,:)
02975 res%ptr_container(maxi,:)=tmp_vec(:)
02976
02977 val_tmp=(1/m_tmp%ptr_container(i,j))
02978 m_tmp%ptr_container(i,:)=val_tmp*m_tmp%ptr_container(i,:);
02979 res%ptr_container(i,:)=val_tmp*res%ptr_container(i,:);
02980
02981 do l=1,m_tmp%rows
02982 if(l.ne.i) then
02983 val_tmp=m_tmp%ptr_container(l,j);
02984 m_tmp%ptr_container(l,:)=m_tmp%ptr_container(l,:)-val_tmp*m_tmp%ptr_container(i,:)
02985 res%ptr_container(l,:)=res%ptr_container(l,:)-val_tmp*res%ptr_container(i,:)
02986 end if
02987 end do
02988 end if
02989 i = i + 1
02990 j = j + 1
02991 end do
02992
02993 deallocate(tmp_vec)
02994 call destruct(m_tmp)
02995 end function m_inverse_gaussj
02996
02997
02998
02999
03000
03001
03026 function m_pinv(m,eps_svd,eps_chol,iter_max,meth_qr,eps_gsortho,meth_pinv) result(res)
03027 implicit none
03028 type(matrix), intent(in) :: m
03029 type_precision, intent(in), optional :: eps_svd
03030 type_precision, intent(in), optional :: eps_chol
03031 integer, intent(in), optional :: iter_max
03032 character*(*), intent(in), optional :: meth_qr
03033 type_precision, intent(in), optional :: eps_gsortho
03034 character*(*), intent(in), optional :: meth_pinv
03035
03036 type(matrix) :: res
03037
03038 if(present(meth_pinv) .and. meth_pinv.eq.'svd') then
03039 res=m_pseudoinv_svd(m,eps_svd,iter_max,meth_qr,eps_gsortho)
03040 else
03041 res=m_pseudoinv_chol(m,eps_chol)
03042 end if
03043 end function m_pinv
03044
03045
03046
03047
03048
03049
03060 function m_pseudoinv_svd(m,eps_svd,iter_max,meth_qr,eps_gsortho) result(res)
03061 implicit none
03062 type(matrix), intent(in) :: m
03063 type_precision, intent(in), optional :: eps_svd
03064 integer, intent(in), optional :: iter_max
03065 character*(*), intent(in), optional :: meth_qr
03066 type_precision, intent(in), optional :: eps_gsortho
03067
03068 type(t_svd) :: decomp_svd
03069 type(matrix) :: res
03070 integer :: i
03071 type_precision :: machine_eps
03072
03073 machine_eps=math_machine_eps()
03074
03075 decomp_svd=m_decompsvd(m,eps_svd,iter_max,meth_qr,eps_gsortho);
03076
03077
03078 do i=1,min(decomp_svd%S%rows,decomp_svd%S%cols)
03079 if(decomp_svd%S%ptr_container(i,i).ge.machine_eps) then
03080 decomp_svd%S%ptr_container(i,i)=1.0/decomp_svd%S%ptr_container(i,i)
03081 else
03082 decomp_svd%S%ptr_container(i,i)=0.0
03083 end if
03084 end do
03085
03086 res=decomp_svd%V*(.tr.decomp_svd%S)*(.tr.decomp_svd%U)
03087 call destruct(decomp_svd)
03088 end function m_pseudoinv_svd
03089
03090
03091
03092
03093
03094
03103 function m_pseudoinv_chol(m,eps_chol) result(res)
03104 implicit none
03105 type(matrix), intent(in) :: m
03106 type_precision, intent(in), optional :: eps_chol
03107
03108 type(matrix) :: m_tmp, L
03109 type(vector) :: diag_m_tmp
03110 integer :: size_r, size_c, k, i, j, r
03111 type_precision :: sum_tmp, in_eps_chol
03112 logical :: is_transpose
03113 type(matrix) :: res
03114
03115
03116 size_r=m%rows;
03117 size_c=m%cols;
03118 is_transpose=.false.
03119 if(size_r.lt.size_c) then
03120 is_transpose=.true.
03121 m_tmp=m*.tr.m
03122 size_c=size_r;
03123 else
03124 m_tmp=.tr.m*m
03125 end if
03126
03127
03128 diag_m_tmp=diag(m_tmp);
03129 if(present(eps_chol)) then
03130 in_eps_chol=eps_chol;
03131 else
03132 in_eps_chol=minval(abs(diag_m_tmp%ptr_container))*1.d-9;
03133 end if
03134 call destruct(diag_m_tmp);
03135 call m_init(L,m_tmp%rows,m_tmp%cols);
03136
03137 r=0;
03138 do k=1,size_c
03139 r=r+1
03140
03141 if(r.eq.1) then
03142 call m_setsub(L,r_lbound=k,r_ubound=size_c, c_lbound=r,c_ubound=r, &
03143 m_sub=m_extract(m_tmp,r_lbound=k,r_ubound=size_c,c_lbound=k,c_ubound=k) );
03144 else
03145 call m_setsub(L,r_lbound=k,r_ubound=size_c, c_lbound=r,c_ubound=r, &
03146 m_sub=m_extract(m_tmp,r_lbound=k,r_ubound=size_c,c_lbound=k,c_ubound=k)- &
03147 m_extract(L,r_lbound=k,r_ubound=size_c,c_lbound=1,c_ubound=r-1) * &
03148 .tr.m_extract(L,r_lbound=k,r_ubound=k,c_lbound=1,c_ubound=r-1) );
03149 end if
03150 if(L%ptr_container(k,r).gt.in_eps_chol) then
03151 L%ptr_container(k,r)=sqrt(L%ptr_container(k,r))
03152 if(k.lt.size_c) then
03153 call m_setsub(L,r_lbound=k+1,r_ubound=size_c, c_lbound=r,c_ubound=r, &
03154 m_sub=m_extract(L,r_lbound=k+1,r_ubound=size_c, c_lbound=r,c_ubound=r)/L%ptr_container(k,r) );
03155 end if
03156 else
03157 r=r-1;
03158 end if
03159 end do
03160 call destruct(m_tmp);
03161 L=m_extract(L,r_lbound=1,r_ubound=L%cols,c_lbound=1,c_ubound=r)
03162
03163
03164 res=inv(.tr.L*L);
03165 res=res*res
03166 if(is_transpose) then
03167 res=(.tr.m)*L*res*(.tr.L);
03168 else
03169 res=L*res*(.tr.L)*(.tr.m);
03170 end if
03171 call destruct(L);
03172 end function m_pseudoinv_chol
03173
03174
03175
03176
03177
03178
03179
03186 function m_trace(m) result(res)
03187 implicit none
03188 type(matrix), intent(in) :: m
03189
03190 type(type_exception) :: err_exception;
03191 type_precision :: res
03192 #ifdef DEBUG_EXCEPTION
03193 if(.not.(m%is_allocate)) then
03194 m_what_exception='m_trace::matrix not allocated'
03195 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
03196 endif
03197 if (m%rows.ne.m%cols) then
03198 m_what_exception='m_trace::square matrix expected'
03199 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
03200 endif
03201 #endif
03202 res=sum(diag(m));
03203
03204 end function m_trace
03205
03206
03207
03208
03209
03210
03218 function m_permut(m,is_permuted) result(res)
03219 implicit none
03220 type(matrix), intent(in) :: m
03221 logical, intent(in), optional :: is_permuted
03222
03223 type(type_exception) :: err_exception;
03224 type(t_m_and_p) :: res
03225 integer :: i, maxi, j, p;
03226 type_precision, dimension(:), allocatable :: tmp_vec_M, tmp_vec_P ;
03227 #ifdef DEBUG_EXCEPTION
03228 if(.not.(m%is_allocate)) then
03229 m_what_exception='m_permut::matrix not allocated'
03230 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
03231 endif
03232 #endif
03233 res%M=m
03234
03235 if(present(is_permuted) .and. is_permuted) then
03236 res%P=m_identity(m%rows);
03237 allocate(tmp_vec_M(m%cols))
03238 allocate(tmp_vec_P(m%rows))
03239 do i=1,res%M%rows
03240 do j=1,res%M%cols
03241
03242 maxi=i
03243 do p=i+1,res%M%rows
03244 if(abs(res%M%ptr_container(p,j)).gt.abs(res%M%ptr_container(maxi,j))) then
03245 maxi=p
03246 end if
03247 end do
03248
03249
03250
03251 if(i.ne.maxi) then
03252 tmp_vec_M(:)=res%M%ptr_container(i,:)
03253 res%M%ptr_container(i,:)=res%M%ptr_container(maxi,:)
03254 res%M%ptr_container(maxi,:)=tmp_vec_M(:)
03255
03256 tmp_vec_P(:)=res%P%ptr_container(i,:)
03257 res%P%ptr_container(i,:)=res%P%ptr_container(maxi,:)
03258 res%P%ptr_container(maxi,:)=tmp_vec_P(:)
03259
03260 res%swops=res%swops+1
03261 end if
03262 end do
03263 end do
03264 deallocate(tmp_vec_M)
03265 deallocate(tmp_vec_P)
03266 end if
03267 end function m_permut
03268
03269
03270
03271
03272
03273
03281 function m_permut_col(m,is_permuted) result(res)
03282 implicit none
03283 type(matrix), intent(in) :: m
03284 logical, intent(in), optional :: is_permuted
03285
03286 type(type_exception) :: err_exception;
03287 type(t_m_and_p) :: res
03288 integer :: i, maxj, j, p;
03289 type_precision, dimension(:), allocatable :: tmp_vec_M, tmp_vec_P ;
03290 #ifdef DEBUG_EXCEPTION
03291 if(.not.(m%is_allocate)) then
03292 m_what_exception='m_permut_col::matrix not allocated'
03293 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
03294 endif
03295 #endif
03296 res%M=m
03297
03298 if(present(is_permuted) .and. is_permuted) then
03299 res%P=m_identity(m%cols);
03300 allocate(tmp_vec_M(m%rows))
03301 allocate(tmp_vec_P(m%cols))
03302 do i=1,res%M%rows
03303 do j=1,res%M%cols
03304
03305 maxj=j
03306 do p=j+1,res%M%cols
03307 if(abs(res%M%ptr_container(i,p)).gt.abs(res%M%ptr_container(i,maxj))) then
03308 maxj=p
03309 end if
03310 end do
03311
03312
03313
03314 if(j.ne.maxj) then
03315 tmp_vec_M(:)=res%M%ptr_container(:,j)
03316 res%M%ptr_container(:,j)=res%M%ptr_container(:,maxj)
03317 res%M%ptr_container(:,maxj)=tmp_vec_M(:)
03318
03319 tmp_vec_P(:)=res%P%ptr_container(:,j)
03320 res%P%ptr_container(:,j)=res%P%ptr_container(:,maxj)
03321 res%P%ptr_container(:,maxj)=tmp_vec_P(:)
03322
03323 res%swops=res%swops+1
03324 end if
03325 end do
03326 end do
03327 deallocate(tmp_vec_M)
03328 deallocate(tmp_vec_P)
03329 end if
03330 end function m_permut_col
03331
03332
03333
03334
03335
03336
03344 function m_decompLU(m,is_permuted) result(res)
03345 implicit none
03346 type(matrix), intent(in) :: m
03347 logical, intent(in), optional :: is_permuted
03348
03349 type(type_exception) :: err_exception;
03350 type(t_lu) :: res
03351 type(t_m_and_p) :: res_mp
03352 integer :: i, maxi, j, k, p
03353 type_precision :: tp_s
03354 type(matrix) :: m_tmp
03355 type_precision, dimension(:), allocatable :: tmp_vec;
03356 #ifdef DEBUG_EXCEPTION
03357 if(.not.(m%is_allocate)) then
03358 m_what_exception='m_decompLU::matrix not allocated'
03359 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
03360 endif
03361 #endif
03362
03363 call m_init(res%L,m%rows,min(m%rows,m%cols))
03364
03365 call m_init(res%U,min(m%rows,m%cols),m%cols)
03366
03367 if(present(is_permuted) .and. is_permuted) then
03368 res_mp=m_permut(m,is_permuted);
03369 res%P=res_mp%P
03370 res%swops=res_mp%swops
03371 m_tmp=res_mp%M
03372 call destruct(res_mp)
03373 else
03374 m_tmp=m
03375 end if
03376
03377 do i=1,m_tmp%rows
03378 do j=1,m_tmp%cols
03379 tp_s=m_tmp%ptr_container(i,j)
03380 do k=1,min(i,j)-1
03381 tp_s=tp_s-res%L%ptr_container(i,k)*res%U%ptr_container(k,j)
03382 end do
03383 if(i.gt.j) then
03384 res%L%ptr_container(i,j)=tp_s/res%U%ptr_container(j,j)
03385 else
03386 res%U%ptr_container(i,j)=tp_s
03387 end if
03388 end do
03389 end do
03390 forall(i=1:m%cols) res%L%ptr_container(i,i)=1
03391 call destruct(m_tmp)
03392 end function m_decompLU
03393
03394
03402 function m_decompLU_m(m,is_permuted) result(res)
03403 implicit none
03404 type(matrix), intent(in) :: m
03405 logical, intent(in), optional :: is_permuted
03406
03407 type(type_exception) :: err_exception;
03408 integer :: i, maxi, p, j, k
03409 type_precision :: tp_s
03410 type(t_m_and_p) :: res
03411 type_precision, dimension(:), allocatable :: tmp_vec;
03412 #ifdef DEBUG_EXCEPTION
03413 if(.not.(m%is_allocate)) then
03414 m_what_exception='m_decompLU_m::matrix not allocated'
03415 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
03416 endif
03417 #endif
03418 res=m_permut(m,is_permuted);
03419
03420 do i=1,res%M%rows
03421 do j=1,res%M%cols
03422 tp_s=res%M%ptr_container(i,j)
03423 do k=1,min(i,j)-1
03424 tp_s=tp_s-res%M%ptr_container(i,k)*res%M%ptr_container(k,j)
03425 end do
03426 if(i.gt.j) tp_s=tp_s/res%M%ptr_container(j,j)
03427
03428 res%M%ptr_container(i,j)=tp_s
03429 end do
03430 end do
03431 end function m_decompLU_m
03432
03433
03434
03435
03436
03437
03447 function m_decompCholesky(m,is_permuted) result(res)
03448 implicit none
03449 type(matrix), intent(in) :: m
03450 logical, intent(in), optional :: is_permuted
03451
03452 type(type_exception) :: err_exception;
03453 integer :: i, j, k, s, N
03454 type_precision :: tp_s, machine_eps
03455 type(matrix) :: m_tmp
03456 type(t_m_and_p) :: res
03457 #ifdef DEBUG_EXCEPTION
03458 if(.not.(m%is_allocate)) then
03459 m_what_exception='m_decompCholesky::matrix not allocated'
03460 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
03461 endif
03462 if (m%rows.ne.m%cols) then
03463 m_what_exception='m_decompCholesky:: square matrix expected'
03464 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
03465 endif
03466 #endif
03467 machine_eps=math_machine_eps()
03468 res=m_permut(m,is_permuted);
03469 N=res%M%rows;
03470 m_tmp=res%M
03471 res%M%ptr_container=0
03472 do k=1,N-1
03473 res%M%ptr_container(k,k)=sqrt(m_tmp%ptr_container(k,k))
03474 do s=k+1,N
03475 res%M%ptr_container(s,k)=m_tmp%ptr_container(s,k)/res%M%ptr_container(k,k)
03476 end do
03477 do j=k+1,N
03478 do i=j,N
03479 m_tmp%ptr_container(i,j)=m_tmp%ptr_container(i,j) - &
03480 res%M%ptr_container(i,k)*res%M%ptr_container(j,k)
03481 if(i.eq.j .and. m_tmp%ptr_container(i,j).le.machine_eps) then
03482 m_what_exception='m_decompCholesky::matrix is not positive definite'
03483 err_exception=e_error(stop_array_positive,m_what_exception,stop_array_positive)
03484 end if
03485 end do
03486 end do
03487 end do
03488 res%M%ptr_container(N,N)=sqrt(m_tmp%ptr_container(N,N))
03489
03490 call destruct(m_tmp)
03491 end function m_decompCholesky
03492
03493
03494
03495
03496
03497
03513 function m_decompQR_GramSchimdt(m,is_permuted) result(res)
03514 implicit none
03515 type(matrix), intent(in) :: m
03516 logical, intent(in), optional :: is_permuted
03517
03518 type(type_exception) :: err_exception;
03519 type(t_qr) :: res
03520 integer :: j,i
03521 logical :: is_converg
03522 type(vector) :: v_tmp, v_tmp2
03523 type(t_m_and_p) :: res_mp
03524 type(matrix) :: m_tmp
03525 type_precision :: rank_tmp, test, machine_eps
03526 #ifdef DEBUG_EXCEPTION
03527 if(.not.(m%is_allocate)) then
03528 m_what_exception='m_decompQR_GramSchimdt::matrix not allocated'
03529 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
03530 endif
03531 #endif
03532 machine_eps=math_machine_eps()
03533 res_mp=m_permut(m,is_permuted);
03534
03535
03536 call m_init(res%Q,m%rows,m%cols)
03537
03538 call m_init(res%R,m%cols,m%cols)
03539
03540 is_converg = .true.
03541 if(present(is_permuted) .and. is_permuted) then
03542 res%P=res_mp%P
03543 res%swops=res_mp%swops
03544 m_tmp=res_mp%M
03545 call destruct(res_mp);
03546 else
03547 m_tmp=m
03548 end if
03549 rank_tmp=rank(m_tmp)
03550 call v_init(v_tmp,m_tmp%rows)
03551 call v_init(v_tmp2,m_tmp%rows)
03552
03553 v_tmp%ptr_container=m_tmp%ptr_container(:,1)
03554 res%R%ptr_container(1,1)=norm(v_tmp,type_norm=2)
03555 if(res%R%ptr_container(1,1).le.machine_eps) then
03556 is_converg=.false.
03557 else
03558 res%Q%ptr_container(:,1)=(1.0/res%R%ptr_container(1,1))*v_tmp%ptr_container
03559 end if
03560
03561 j=2
03562 do while (j.le.m_tmp%cols .and. is_converg)
03563 v_tmp%ptr_container=m_tmp%ptr_container(:,j)
03564 do i=1,j-1
03565 v_tmp2%ptr_container=res%Q%ptr_container(:,i)
03566 res%R%ptr_container(i,j)=v_dot(v_tmp,v_tmp2)
03567 v_tmp=v_tmp-res%R%ptr_container(i,j)*v_tmp2
03568 end do
03569
03570 if(j.le.rank_tmp) then
03571
03572 res%R%ptr_container(j,j)=norm(v_tmp,type_norm=2)
03573
03574
03575 if(res%R%ptr_container(j,j).le.machine_eps) then
03576 is_converg=.false.
03577 else
03578 res%Q%ptr_container(:,j)=(1.0/res%R%ptr_container(j,j))*v_tmp%ptr_container
03579 end if
03580 end if
03581 j = j + 1
03582 end do
03583
03584 call destruct(m_tmp)
03585 call destruct(v_tmp)
03586 call destruct(v_tmp2)
03587 if(.not.(is_converg)) then
03588 v_what_exception='m_decompQR_GramSchimdt::QR decomposition not exist'
03589 err_exception=e_error(stop_array_factorisation ,v_what_exception,stop_array_factorisation)
03590 end if
03591 end function m_decompQR_GramSchimdt
03592
03593
03594
03595
03596
03597
03620 function m_decompQR_GramSchimdt_Reortho(m,eps,is_permuted) result(res)
03621 implicit none
03622 type(matrix), intent(in) :: m
03623 type_precision, intent(in), optional :: eps
03624 logical, intent(in), optional :: is_permuted
03625
03626 type(type_exception) :: err_exception;
03627 type(t_qr) :: res
03628 integer :: k, i
03629 logical :: is_nach
03630 type(vector) :: z
03631 type(matrix) :: m_tmp
03632 type(t_m_and_p) :: res_mp
03633 type_precision :: in_eps, s, t, tt
03634 #ifdef DEBUG_EXCEPTION
03635 if(.not.(m%is_allocate)) then
03636 m_what_exception='m_decompQR_GramSchimdt::matrix not allocated'
03637 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
03638 endif
03639 #endif
03640
03641 if(present(eps)) then
03642 in_eps=eps
03643 else
03644 in_eps=math_machine_eps()
03645 end if
03646 in_eps=1.e-25
03647
03648 if(present(is_permuted) .and. is_permuted) then
03649 res%P=res_mp%P
03650 res%swops=res_mp%swops
03651 res%Q=res_mp%M
03652 call destruct(res_mp);
03653 else
03654 res%Q=m
03655 end if
03656
03657 call m_init(res%R,m%cols,m%cols);
03658 call v_init(z,m%cols);
03659
03660 do k=1,m%cols
03661
03662 t=norm(m_extract(res%Q,r_lbound=1,r_ubound=m%rows,c_lbound=k,c_ubound=k),frobenius)
03663 is_nach= .false.
03664 do while(.not. is_nach)
03665 z%ptr_container(k) = z%ptr_container(k) + 1
03666 do i=1,k-1
03667
03668 m_tmp=.tr.m_extract(res%Q,r_lbound=1,r_ubound=m%rows,c_lbound=i,c_ubound=i) &
03669 *m_extract(res%Q,r_lbound=1,r_ubound=m%rows,c_lbound=k,c_ubound=k)
03670 s=m_tmp%ptr_container(1,1)
03671
03672 res%R%ptr_container(i,k)=res%R%ptr_container(i,k) + s
03673
03674 call m_setsub(res%Q,r_lbound=1,r_ubound=m%rows, c_lbound=k,c_ubound=k, &
03675 m_sub=m_extract(res%Q,r_lbound=1,r_ubound=m%rows, c_lbound=k,c_ubound=k) - &
03676 s*m_extract(res%Q,r_lbound=1,r_ubound=m%rows, c_lbound=i,c_ubound=i))
03677 end do
03678
03679
03680 tt=norm(m_matrixTOvector(m_extract(res%Q,r_lbound=1,r_ubound=m%rows, &
03681 c_lbound=k, c_ubound=k)))
03682 if(tt.gt.10*in_eps*t .and. tt.lt.(t/10.0)) then
03683 is_nach = .false.
03684 t=tt;
03685 else
03686 is_nach = .true.
03687 if(tt.lt.10*in_eps*t) tt=0;
03688 end if
03689 end do
03690 res%R%ptr_container(k,k)=tt;
03691 if (tt*in_eps.ne.0) then
03692 tt=1/tt;
03693 else
03694 tt=0;
03695 end if
03696
03697 call m_setsub(res%Q,r_lbound=1,r_ubound=m%rows, c_lbound=k,c_ubound=k, &
03698 m_sub=tt*m_extract(res%Q,r_lbound=1,r_ubound=m%rows, c_lbound=k,c_ubound=k))
03699 end do
03700
03701 call destruct(z)
03702 end function m_decompQR_GramSchimdt_Reortho
03703
03704
03705
03706
03707
03708
03709
03722 function m_decompQR_HouseHolder(m,is_permuted) result(res)
03723 implicit none
03724 type(matrix), intent(in) :: m
03725 logical, intent(in), optional :: is_permuted
03726
03727 type(type_exception) :: err_exception;
03728 type(t_qr) :: res
03729 type(t_m_and_p) :: res_mp
03730 type(matrix) :: a_k, a_temp, u, Qi
03731 integer :: i, j, k
03732 integer :: size_r, size_c, size_max
03733 type_precision :: alpha, machine_eps
03734 logical :: is_converg
03735 #ifdef DEBUG_EXCEPTION
03736 if(.not.(m%is_allocate)) then
03737 m_what_exception='m_decompQR_HouseHolder::matrix not allocated'
03738 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
03739 endif
03740 #endif
03741 if(present(is_permuted) .and. is_permuted) then
03742 res_mp=m_permut(m,is_permuted);
03743 res%P=res_mp%P
03744 res%swops=res_mp%swops
03745 res%R=res_mp%M
03746 call destruct(res_mp);
03747 else
03748
03749 res%R=m
03750 end if
03751 machine_eps=math_machine_eps()
03752 size_r=m%rows
03753 size_c=m%cols
03754 size_max=min(size_r-1,size_c)
03755
03756 res%Q=m_identity(size_r)
03757
03758 do k=1, size_max
03759
03760 a_k=m_extract(res%R,r_lbound=k,r_ubound=res%R%rows,c_lbound=k,c_ubound=k);
03761
03762
03763 a_temp=.tr.a_k*a_k
03764 alpha=sqrt(a_temp%ptr_container(1,1))
03765
03766
03767 if(a_k%ptr_container(1,1).lt.0) alpha = (-1.0)*alpha;
03768
03769 u=a_k
03770 u%ptr_container(1,1)=u%ptr_container(1,1)-alpha
03771
03772 a_temp=.tr.u*u
03773 if(abs(a_temp%ptr_container(1,1)).lt.machine_eps) then
03774 a_temp%ptr_container(1,1)=1
03775 end if
03776 u=u/(sqrt(a_temp%ptr_container(1,1)))
03777
03778 Qi=m_identity(size_r-k+1) - p_notcast(2)*u*.tr.u
03779
03780 call m_setsub(res%Q,r_lbound=1,r_ubound=res%Q%rows, &
03781 c_lbound=k,c_ubound=res%Q%cols, &
03782 m_sub=m_extract(res%Q,r_lbound=1,r_ubound=res%Q%rows, &
03783 c_lbound=k,c_ubound=res%Q%cols)*Qi)
03784
03785
03786 call m_setsub(res%R,r_lbound=k,r_ubound=res%R%rows, &
03787 c_lbound=k,c_ubound=res%R%cols, &
03788 m_sub=Qi*m_extract(res%R,r_lbound=k,r_ubound=res%R%rows, &
03789 c_lbound=k,c_ubound=res%R%cols))
03790 end do
03791
03792 call destruct(a_temp)
03793 call destruct(a_k)
03794 call destruct(u)
03795 call destruct(Qi)
03796 end function m_decompQR_HouseHolder
03797
03798
03799
03800
03812 function m_decompQR(m,meth_qr,eps_gsortho,is_permuted) result(res)
03813 implicit none
03814 type(matrix), intent(in) :: m
03815 character*(*), intent(in), optional :: meth_qr
03816 type_precision, intent(in), optional :: eps_gsortho
03817 logical, intent(in), optional :: is_permuted
03818
03819 type(type_exception) :: err_exception;
03820 type(t_qr) :: res
03821 integer :: i,j
03822 if(present(meth_qr) .and. meth_qr.eq.'ho') then
03823 res=m_decompQR_HouseHolder(m,is_permuted)
03824 else if(present(meth_qr) .and. meth_qr.eq.'gsro') then
03825 res=m_decompQR_GramSchimdt_Reortho(m,eps_gsortho,is_permuted)
03826 else
03827 res=m_decompQR_GramSchimdt(m,is_permuted)
03828 end if
03829
03830 end function m_decompQR
03831
03832
03833
03834
03835
03836
03853 function m_decompsvd(m,eps,iter_max,meth_qr,eps_gsortho,is_permuted) result(res)
03854 implicit none
03855 type(matrix), intent(in) :: m
03856 type_precision, intent(in), optional :: eps
03857 integer, intent(in), optional :: iter_max
03858 character*(*), intent(in), optional :: meth_qr
03859 type_precision, intent(in), optional :: eps_gsortho
03860 logical, intent(in), optional :: is_permuted
03861
03862 type(type_exception) :: err_exception;
03863 type(t_svd) :: res
03864 integer :: i, j, iter, in_iter_max
03865 type_precision :: norm_tmp=0, norm_fro=0, norm_diag=0, in_eps, machine_eps
03866 type_precision :: fro=0
03867 type_precision :: s_tmp=0
03868 type(t_qr) :: qr_decomp
03869 #ifdef DEBUG_EXCEPTION
03870 if(.not.(m%is_allocate)) then
03871 m_what_exception='m_decompsvd::matrix not allocated'
03872 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
03873 endif
03874 #endif
03875 machine_eps=math_machine_eps()
03876
03877 if(present(eps)) then
03878 in_eps=eps
03879 else
03880 in_eps=machine_eps
03881 end if
03882
03883 if(present(iter_max)) then
03884 in_iter_max=iter_max
03885 else
03886 in_iter_max=100*max(m%rows,m%cols)
03887 end if
03888
03889
03890 res%U=m_identity(m%rows);
03891 qr_decomp%R=.tr.m
03892 res%V=m_identity(m%cols);
03893
03894 norm_tmp=10000
03895 iter=0
03896 do while(norm_tmp.gt.in_eps .and. iter.le.in_iter_max)
03897 qr_decomp=qr(.tr.qr_decomp%R,meth_qr,eps_gsortho,is_permuted);
03898 res%U=res%U*qr_decomp%Q;
03899 qr_decomp=qr(.tr.qr_decomp%R,meth_qr,eps_gsortho,is_permuted);
03900 res%V=res%V*qr_decomp%Q;
03901
03902 norm_fro=norm(triu(qr_decomp%R,1),type_norm=frobenius)
03903 norm_diag=norm(diag(qr_decomp%R))
03904
03905 if(abs(norm_diag).le.machine_eps) norm_diag=1;
03906 norm_tmp=norm_fro/norm_diag;
03907 iter = iter +1
03908 end do
03909
03910
03911
03912 call m_init(res%S,m%rows,m%cols)
03913
03914 do i=1,min(m%rows,m%cols)
03915 s_tmp=qr_decomp%R%ptr_container(i,i)
03916 res%S%ptr_container(i,i)=abs(s_tmp);
03917 if(s_tmp.lt.machine_eps) then
03918 res%U%ptr_container(:,i)=-res%U%ptr_container(:,i)
03919 end if
03920 end do
03921
03922 call destruct(qr_decomp)
03923
03924 if(iter.eq.in_iter_max) then
03925 v_what_exception='m_decompsvd::svd method diverge'
03926 err_exception=e_error(stop_array_diverge,v_what_exception,stop_array_diverge)
03927 end if
03928 end function m_decompsvd
03929
03930
03943 function m_decompsvd_s(m,eps,iter_max,meth_qr,eps_gsortho,is_permuted) result(res)
03944 implicit none
03945 type(matrix), intent(in) :: m
03946 type_precision, intent(in), optional :: eps
03947 integer, intent(in), optional :: iter_max
03948 character*(*), intent(in), optional :: meth_qr
03949 type_precision, intent(in), optional :: eps_gsortho
03950 logical, intent(in), optional :: is_permuted
03951
03952 type(type_exception) :: err_exception;
03953 type(vector) :: res
03954 integer :: i, j, iter, in_iter_max
03955 type_precision :: norm_tmp=0, norm_fro=0, norm_diag=0, in_eps, machine_eps
03956 type_precision :: fro=0
03957 type_precision :: s_tmp=0
03958 type(t_qr) :: qr_decomp
03959 #ifdef DEBUG_EXCEPTION
03960 if(.not.(m%is_allocate)) then
03961 m_what_exception='m_decompsvd_s::matrix not allocated'
03962 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
03963 endif
03964 #endif
03965 machine_eps=math_machine_eps()
03966
03967 if(present(eps)) then
03968 in_eps=eps
03969 else
03970 in_eps=machine_eps
03971 end if
03972
03973
03974 if(present(iter_max)) then
03975 in_iter_max=iter_max
03976 else
03977 in_iter_max=100*max(m%rows,m%cols)
03978 end if
03979
03980 qr_decomp%R=.tr.m
03981
03982 norm_tmp=10000
03983 iter=0
03984 do while(norm_tmp.gt.in_eps .and. iter.le.in_iter_max)
03985 qr_decomp=qr(.tr.qr_decomp%R,meth_qr,eps_gsortho,is_permuted);
03986 qr_decomp=qr(.tr.qr_decomp%R,meth_qr,eps_gsortho,is_permuted);
03987
03988 norm_fro=norm(triu(qr_decomp%R,1),type_norm=frobenius)
03989 norm_diag=norm(diag(qr_decomp%R))
03990
03991 if(abs(norm_diag).le.machine_eps) norm_diag=1;
03992 norm_tmp=norm_fro/norm_diag;
03993 iter = iter +1
03994 end do
03995
03996
03997 call v_init(res,min(m%rows,m%cols))
03998
03999 do i=1,res%size
04000 res%ptr_container(i)=abs(qr_decomp%R%ptr_container(i,i));
04001 end do
04002
04003 call destruct(qr_decomp)
04004 if(iter.eq.in_iter_max) then
04005 v_what_exception='m_decompsvd_s::svd method diverge'
04006 err_exception=e_error(stop_array_diverge,v_what_exception,stop_array_diverge)
04007 end if
04008 end function m_decompsvd_s
04009
04010
04011
04012
04013
04014
04015
04027 function m_pow_eig(m,v0,eps,iter_max) result(res)
04028 implicit none
04029 type(matrix), intent(in) :: m
04030 type(vector), intent(in), optional :: v0
04031 type_precision, intent(in), optional :: eps
04032 integer, intent(in), optional :: iter_max
04033
04034 type(type_exception) :: err_exception;
04035 integer :: iter, in_iter_max
04036 type(t_poweig) :: res
04037 logical :: is_converg
04038 type(vector) :: w
04039 type_precision, allocatable, dimension(:) :: err_tmp
04040 type_precision :: norm_tmp, in_eps
04041 integer :: idx_lambda
04042 #ifdef DEBUG_EXCEPTION
04043 if(.not.(m%is_allocate)) then
04044 m_what_exception='m_pow_eig::matrix not allocated'
04045 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
04046 endif
04047 if (m%rows.ne.m%cols) then
04048 m_what_exception='m_pow_eig:: square matrix expected'
04049 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
04050 endif
04051 #endif
04052 is_converg = .false.
04053 if(present(v0)) then
04054 res%v_lambda=v0;
04055 else
04056 call v_init(res%v_lambda,m%rows);
04057 res%v_lambda=p_notcast(1.0)
04058 end if
04059
04060 if(present(eps)) then
04061 in_eps=eps
04062 else
04063 in_eps=math_machine_eps()
04064 end if
04065
04066 if(present(iter_max)) then
04067 in_iter_max=iter_max
04068 else
04069 in_iter_max=100*max(m%rows,m%cols)
04070 end if
04071
04072 allocate(err_tmp(in_iter_max))
04073 w=m*res%v_lambda
04074 idx_lambda=maxloc(abs(w%ptr_container),dim=1)
04075 res%lambda=abs(w%ptr_container(idx_lambda));
04076 res%v_lambda=(1.0/res%lambda)*w
04077 iter=1
04078 do while (iter.lt.in_iter_max .and. .not.(is_converg))
04079 w=m*res%v_lambda
04080 idx_lambda=maxloc(abs(w%ptr_container),dim=1)
04081 res%lambda=abs(w%ptr_container(idx_lambda));
04082 res%v_lambda=(1.0/res%lambda)*w
04083 norm_tmp=norm(m*res%v_lambda-w,type_norm=infty)
04084 err_tmp(iter+1)=norm_tmp
04085 if(norm_tmp>in_eps) then
04086 res%v_lambda=(1.0/res%lambda)*w
04087 else
04088 is_converg=.true.
04089 endif
04090 iter = iter + 1
04091 end do
04092 res%lambda=w%ptr_container(idx_lambda);
04093 res%iter=iter-1
04094 call v_init(res%v_err,res%iter)
04095 res%v_err%ptr_container=err_tmp(1:res%iter)
04096 deallocate(err_tmp)
04097 call destruct(w)
04098 if(.not.(is_converg)) then
04099 v_what_exception='m_pow_eig::power eig method diverge'
04100 err_exception=e_error(stop_array_diverge,v_what_exception,stop_array_diverge)
04101 end if
04102
04103 end function m_pow_eig
04104
04105
04117 function m_eig_deflation(m,v0,eps,iter_max) result(res)
04118 implicit none
04119 type(matrix), intent(in) :: m
04120 type(vector), intent(in), optional :: v0
04121 type_precision, intent(in), optional :: eps
04122 integer, intent(in), optional :: iter_max
04123
04124 type(type_exception) :: err_exception;
04125 type(t_eig) :: res
04126 type(matrix) :: m_tmp
04127 type(t_poweig) :: poweig
04128 integer :: i, j
04129 #ifdef DEBUG_EXCEPTION
04130 if (m%rows.ne.m%cols) then
04131 m_what_exception='m_eig_deflation:: square matrix expected'
04132 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
04133 endif
04134 #endif
04135 m_tmp=m
04136 call v_init(res%v_eigvalues,m_tmp%rows)
04137 call m_init(res%m_eigvectors,m_tmp%rows,m_tmp%rows)
04138 do i=1,m_tmp%rows
04139 poweig=m_pow_eig(m_tmp,v0,eps,iter_max);
04140 res%v_eigvalues%ptr_container(i)=poweig%lambda;
04141 res%m_eigvectors%ptr_container(:,i)=poweig%v_lambda%ptr_container;
04142
04143 do j=1,m_tmp%cols
04144 m_tmp%ptr_container(i,j)=m_tmp%ptr_container(i,j)- &
04145 (poweig%lambda/sqrnorm(poweig%v_lambda))* &
04146 (poweig%v_lambda%ptr_container(i)*poweig%v_lambda%ptr_container(j));
04147 end do
04148
04149
04150
04151
04152 call print((m-poweig%lambda*m_identity(m%rows))*poweig%v_lambda)
04153 end do
04154
04155 call destruct(m_tmp)
04156 call destruct(poweig)
04157 end function m_eig_deflation
04158
04159
04160
04171 function m_eig_qr(m,iter_max,meth_qr,eps_gsortho,is_permuted) result(res)
04172 implicit none
04173 type(matrix), intent(in) :: m
04174 integer, intent(in) :: iter_max
04175 character*(*), intent(in), optional :: meth_qr
04176 type_precision, intent(in), optional :: eps_gsortho
04177 logical, intent(in), optional :: is_permuted
04178
04179 type(type_exception) :: err_exception;
04180 integer :: k
04181 integer :: iter
04182 type(vector) :: res
04183 type(t_qr) :: qr_decomp
04184 type(matrix) :: m_tmp, I
04185 type_precision :: s
04186 #ifdef DEBUG_EXCEPTION
04187 if(.not.(m%is_allocate)) then
04188 m_what_exception='m_eig_qr::matrix not allocated'
04189 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
04190 endif
04191 if (m%rows.ne.m%cols) then
04192 m_what_exception='m_eig_qr:: square matrix expected'
04193 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
04194 endif
04195 #endif
04196 iter=1
04197 m_tmp=m
04198 I=m_identity(m%rows)
04199 do while(iter.le.iter_max)
04200 s=m_tmp%ptr_container(m%rows,m%cols)
04201 qr_decomp=qr(m_tmp-s*I,meth_qr,eps_gsortho,is_permuted);
04202 m_tmp=qr_decomp%R*qr_decomp%Q+s*I
04203 iter = iter + 1
04204 end do
04205 res=diag(m_tmp)
04206 call destruct(m_tmp)
04207 call destruct(qr_decomp)
04208 end function m_eig_qr
04209
04210
04211
04212
04213
04214
04215
04216
04225 function m_cond(m) result (res)
04226 type(matrix), intent(in) :: m
04227
04228 type(type_exception) :: err_exception;
04229 type_precision :: res
04230 #ifdef DEBUG_EXCEPTION
04231 if(.not.(m%is_allocate)) then
04232 m_what_exception='m_cond::matrix not allocated'
04233 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
04234 endif
04235 #endif
04236 res=norm(.inv.m,type_norm=2)*norm(m,type_norm=2)
04237
04238 end function m_cond
04239
04240
04241
04242
04243
04244
04254 function m_norm(m,type_norm) result (res)
04255 type(matrix), intent(in) :: m
04256 integer, intent(in), optional :: type_norm
04257
04258 type(type_exception) :: err_exception;
04259 integer :: i, j
04260 integer :: in_type_norm
04261 type_precision :: machine_eps, res
04262 type(t_poweig) :: poweig
04263 #ifdef DEBUG_EXCEPTION
04264 if(.not.(m%is_allocate)) then
04265 m_what_exception='m_norm::matrix not allocated'
04266 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
04267 endif
04268 #endif
04269 machine_eps=math_machine_eps()
04270 if(present(type_norm) .and. (type_norm.ge.-1 .and. type_norm.le.2) ) then
04271 in_type_norm=type_norm
04272 else
04273 in_type_norm=2
04274 end if
04275 res=0
04276 if(in_type_norm.eq.1) then
04277 res=maxval(sum(abs(m%ptr_container), dim = 1))
04278 else if(in_type_norm.eq.infty) then
04279 res=maxval(sum(abs(m%ptr_container), dim = 2))
04280 else if(in_type_norm.eq.frobenius) then
04281
04282 do j=1,m%cols
04283 do i=1,m%rows
04284 res = res + abs(m%ptr_container(i,j))**2
04285 end do
04286 end do
04287
04288 res=sqrt(res);
04289 else if(in_type_norm.eq.2) then
04290 poweig=m_pow_eig(.tr.m*m,eps=p_notcast(machine_eps),iter_max=1000)
04291 res=sqrt(abs(poweig%lambda))
04292 call destruct(poweig)
04293 end if
04294
04295 end function m_norm
04296
04297
04298
04299
04300
04301
04308 subroutine m_print(m)
04309 implicit none
04310 type(matrix),intent(in) :: m
04311
04312 integer :: i
04313 type(type_exception) :: err_exception;
04314
04315 type_precision :: res(m%cols)
04316 res=0
04317 #ifdef DEBUG_EXCEPTION
04318 if(.not.(m%is_allocate)) then
04319 m_what_exception='matrix not allocated'
04320 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
04321 endif
04322 #endif
04323 do i=1,m%rows
04324 print*, m%ptr_container(i,:)
04325 end do
04326 end subroutine m_print
04327
04328
04339 subroutine m_print_tofile(m,filename,unit,status,position)
04340 implicit none
04341 type(matrix), intent(in) :: m
04342 character*(*), intent(in) :: filename
04343 integer, intent(in), optional :: unit
04344 character*(*), intent(in), optional :: status
04345 character*(*), intent(in), optional :: position
04346
04347 type(type_exception) :: err_exception;
04348 integer :: in_unit
04349 character*10 :: in_status
04350 character*10 :: in_position
04351 integer :: i, j
04352 #ifdef DEBUG_EXCEPTION
04353 if(.not.(m%is_allocate)) then
04354 m_what_exception='m_print_tofile::matrix not allocated'
04355 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
04356 endif
04357 #endif
04358 if(present(unit)) then
04359 in_unit=unit;
04360 else
04361 in_unit=11
04362 end if
04363 if(present(status)) then
04364 in_status=status;
04365 else
04366 in_status='unknown'
04367 end if
04368 if(present(position)) then
04369 in_position=position;
04370 else
04371 in_position='rewind'
04372 end if
04373
04374 open (unit=in_unit,file=filename, form="formatted", action="write",status=in_status,position=in_position)
04375 write(in_unit,fmt='(i0)',advance='no') m%rows
04376 write(in_unit,fmt='(a1)',advance='no') " "
04377 write(in_unit,fmt='(i0)') m%cols
04378 do i=1,m%rows
04379 do j=1,m%cols-1
04380 write(in_unit,fmt=p_fmt_file,advance='no') m%ptr_container(i,j)
04381 write(in_unit,fmt='(a1)',advance='no') " "
04382 end do
04383 write(in_unit,fmt=p_fmt_file) m%ptr_container(i,m%cols)
04384 end do
04385 close(unit=in_unit)
04386 end subroutine m_print_tofile
04387
04388
04399 subroutine m_print_lu_tofile(lu_decomp,filename,unit,status,position)
04400 implicit none
04401 type(t_lu), intent(in) :: lu_decomp
04402 character*(*), intent(in) :: filename
04403 integer, intent(in), optional :: unit
04404 character*(*), intent(in), optional :: status
04405 character*(*), intent(in), optional :: position
04406
04407 integer :: in_unit
04408 character*10 :: in_status
04409 character*10 :: in_position
04410
04411 if(present(unit)) then
04412 in_unit=unit;
04413 else
04414 in_unit=11
04415 end if
04416 if(present(status)) then
04417 in_status=status;
04418 else
04419 in_status='unknown'
04420 end if
04421 if(present(position)) then
04422 in_position=position;
04423 else
04424 in_position='rewind'
04425 end if
04426
04427 open (unit=in_unit,file=filename, form="formatted", action="write",status=in_status,position=in_position)
04428 write(in_unit,fmt=*) "------> LU decomposition :"
04429 close(unit=in_unit)
04430
04431 open (unit=in_unit,file=filename, form="formatted", action="write",status=in_status,position="append")
04432 write(in_unit,fmt=*) "------> L matrix :"
04433 close(unit=in_unit)
04434 call m_print_tofile(lu_decomp%L,filename,unit,status=in_status,position="append")
04435
04436 open (unit=in_unit,file=filename, form="formatted", action="write",status=in_status,position="append")
04437 write(in_unit,fmt=*) "------> U matrix :"
04438 close(unit=in_unit)
04439 call m_print_tofile(lu_decomp%U,filename,unit,status=in_status,position="append")
04440
04441 if(lu_decomp%P%is_allocate) then
04442 open (unit=in_unit,file=filename, form="formatted", action="write",status=in_status,position="append")
04443 write(in_unit,fmt=*) "------> P (permutation) matrix :"
04444 close(unit=in_unit)
04445 call m_print_tofile(lu_decomp%P,filename,unit,status=in_status,position="append")
04446 end if
04447
04448 end subroutine m_print_lu_tofile
04449
04450
04461 subroutine m_print_qr_tofile(qr_decomp,filename,unit,status,position)
04462 implicit none
04463 type(t_qr), intent(in) :: qr_decomp
04464 character*(*), intent(in) :: filename
04465 integer, intent(in), optional :: unit
04466 character*(*), intent(in), optional :: status
04467 character*(*), intent(in), optional :: position
04468
04469 integer :: in_unit
04470 character*10 :: in_status
04471 character*10 :: in_position
04472
04473 if(present(unit)) then
04474 in_unit=unit;
04475 else
04476 in_unit=11
04477 end if
04478 if(present(status)) then
04479 in_status=status;
04480 else
04481 in_status='unknown'
04482 end if
04483 if(present(position)) then
04484 in_position=position;
04485 else
04486 in_position='rewind'
04487 end if
04488
04489 open (unit=in_unit,file=filename, form="formatted", action="write",status=in_status,position=in_position)
04490 write(in_unit,fmt=*) "------> QR decomposition :"
04491 close(unit=in_unit)
04492
04493 open (unit=in_unit,file=filename, form="formatted", action="write",status=in_status,position="append")
04494 write(in_unit,fmt=*) "------> Q matrix :"
04495 close(unit=in_unit)
04496 call m_print_tofile(qr_decomp%Q,filename,unit,status=in_status,position="append")
04497
04498 open (unit=in_unit,file=filename, form="formatted", action="write",status=in_status,position="append")
04499 write(in_unit,fmt=*) "------> R matrix :"
04500 close(unit=in_unit)
04501 call m_print_tofile(qr_decomp%R,filename,unit,status=in_status,position="append")
04502
04503 if(qr_decomp%P%is_allocate) then
04504 open (unit=in_unit,file=filename, form="formatted", action="write",status=in_status,position="append")
04505 write(in_unit,fmt=*) "------> P (permutation) matrix :"
04506 close(unit=in_unit)
04507 call m_print_tofile(qr_decomp%P,filename,unit,status=in_status,position="append")
04508 end if
04509
04510 end subroutine m_print_qr_tofile
04511
04512
04523 subroutine m_print_m_and_p_tofile(m_and_p_decomp,filename,unit,status,position)
04524 implicit none
04525 type(t_m_and_p), intent(in) :: m_and_p_decomp
04526 character*(*), intent(in) :: filename
04527 integer, intent(in), optional :: unit
04528 character*(*), intent(in), optional :: status
04529 character*(*), intent(in), optional :: position
04530
04531 integer :: in_unit
04532 character*10 :: in_status
04533 character*10 :: in_position
04534
04535
04536 if(present(unit)) then
04537 in_unit=unit;
04538 else
04539 in_unit=11
04540 end if
04541 if(present(status)) then
04542 in_status=status;
04543 else
04544 in_status='unknown'
04545 end if
04546 if(present(position)) then
04547 in_position=position;
04548 else
04549 in_position='rewind'
04550 end if
04551
04552 open (unit=in_unit,file=filename, form="formatted", action="write",status=in_status,position=in_position)
04553 write(in_unit,fmt=*) "------> M_AND_P (example:Cholesky) decomposition :"
04554 close(unit=in_unit)
04555
04556 open (unit=in_unit,file=filename, form="formatted", action="write",status=in_status,position="append")
04557 write(in_unit,fmt=*) "------> M matrix :"
04558 close(unit=in_unit)
04559 call m_print_tofile(m_and_p_decomp%M,filename,unit,status=in_status,position="append")
04560
04561 if(m_and_p_decomp%P%is_allocate) then
04562 open (unit=in_unit,file=filename, form="formatted", action="write",status=in_status,position="append")
04563 write(in_unit,fmt=*) "------> P (permutation) matrix :"
04564 close(unit=in_unit)
04565 call m_print_tofile(m_and_p_decomp%P,filename,unit,status=in_status,position="append")
04566 end if
04567
04568 end subroutine m_print_m_and_p_tofile
04569
04570 end module mod_matrix
04571
04572
04573
04577
04581
04585
04589
04596
04600
04604
04608
04612