00001
00003
00019 module mod_linear_equation
00020 #include "fml_constants.h"
00021 use mod_maths
00022 use mod_vector
00023 use mod_matrix
00024 implicit none
00025
00034 type t_soleq
00035 type(vector) :: v_sol
00036 type(vector) :: v_err
00037 integer :: iter
00038 end type t_soleq
00039
00046 interface linsolve
00047 module procedure leq_solve
00048 end interface
00049
00055 interface ilinsolve
00056 module procedure leq_solve_iter
00057 end interface
00058
00064 interface leq_getBlock
00065 module procedure gsb_getBlock_m, leq_getBlock_v
00066 end interface
00067
00073 interface destruct
00074 module procedure m_destruct_soleq
00075 end interface
00076
00077 CONTAINS
00078
00079
00080
00081
00082
00083
00089 subroutine m_destruct_soleq(soleq_)
00090 implicit none
00091 type(t_soleq) :: soleq_;
00092 call destruct(soleq_%v_sol)
00093 call destruct(soleq_%v_err)
00094 end subroutine m_destruct_soleq
00095
00096
00097
00098
00099
00100
00101
00114 function leq_gaussJourdan(m,v) result(res)
00115 implicit none
00116 type(matrix), intent(in) :: m
00117 type(vector), intent(in) :: v
00118
00119 type(type_exception) :: err_exception;
00120 type(vector) :: res
00121 type(matrix) :: m_tmp
00122 integer :: i,j,k,l,maxi;
00123 type_precision, dimension(:), allocatable :: tmp_vec;
00124 type_precision :: val_tmp, machine_eps
00125 #ifdef DEBUG_EXCEPTION
00126 if(.not.(m%is_allocate)) then
00127 m_what_exception='leq_gaussJourdan::matrix not allocated'
00128 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
00129 endif
00130 if(.not.(v%is_allocate)) then
00131 v_what_exception='leq_gaussJourdan::vector not allocated'
00132 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
00133 endif
00134 if (m%rows.ne.v%size) then
00135 v_what_exception='leq_gaussJourdan::vector size incompatible with the matrix'
00136 err_exception=e_error(stop_array_compatible,v_what_exception,stop_array_compatible)
00137 endif
00138 if(m%rows.ne.m%cols) then
00139 m_what_exception='leq_gaussJourdan::square matrix expected'
00140 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
00141 endif
00142 #endif
00143 machine_eps=math_machine_eps()
00144 res = v
00145 m_tmp= m
00146 allocate(tmp_vec(m%cols))
00147
00148 i=1
00149 j=1
00150 do while(i.le.m_tmp%rows .and. j.le.m_tmp%cols)
00151
00152 maxi=i
00153 do k=i+1,m_tmp%rows
00154 if(abs(m_tmp%ptr_container(k,j)).gt.abs(m_tmp%ptr_container(maxi,j))) then
00155 maxi=k
00156 end if
00157 end do
00158
00159 if(abs(m_tmp%ptr_container(maxi,j)).le.machine_eps) then
00160 m_what_exception='leq_gaussJourdan::matrix is close to singular or badly scaled'
00161 err_exception=e_error(stop_array_singular,m_what_exception,stop_array_singular)
00162 else
00163
00164
00165 tmp_vec(:)=m_tmp%ptr_container(i,:)
00166 m_tmp%ptr_container(i,:)=m_tmp%ptr_container(maxi,:)
00167 m_tmp%ptr_container(maxi,:)=tmp_vec(:)
00168
00169 val_tmp=res%ptr_container(i)
00170 res%ptr_container(i)=res%ptr_container(maxi)
00171 res%ptr_container(maxi)=val_tmp
00172
00173 val_tmp=(1/m_tmp%ptr_container(i,j))
00174 m_tmp%ptr_container(i,:)=val_tmp*m_tmp%ptr_container(i,:);
00175 res%ptr_container(i)=val_tmp*res%ptr_container(i);
00176
00177 do l=1,m_tmp%rows
00178 if(l.ne.j) then
00179 val_tmp=m_tmp%ptr_container(l,j);
00180 m_tmp%ptr_container(l,:)=m_tmp%ptr_container(l,:)-val_tmp*m_tmp%ptr_container(i,:)
00181 res%ptr_container(l)=res%ptr_container(l)-val_tmp*res%ptr_container(i)
00182 end if
00183 end do
00184 end if
00185 i = i + 1
00186 j = j + 1
00187 end do
00188
00189 deallocate(tmp_vec)
00190 call destruct(m_tmp)
00191 end function leq_gaussJourdan
00192
00193
00194
00195
00196
00197
00208 function leq_leastSquare(m,v) result(res)
00209 implicit none
00210 type(matrix), intent(in) :: m
00211 type(vector), intent(in) :: v
00212
00213 type(type_exception) :: err_exception;
00214 type(vector) :: res
00215 #ifdef DEBUG_EXCEPTION
00216 if(.not.(m%is_allocate)) then
00217 m_what_exception='leq_leastSquare::matrix not allocated'
00218 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
00219 endif
00220 if(.not.(v%is_allocate)) then
00221 v_what_exception='leq_leastSquare::vector not allocated'
00222 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
00223 endif
00224 if (m%rows.ne.v%size) then
00225 v_what_exception='leq_leastSquare::vector size incompatible with the matrix'
00226 err_exception=e_error(stop_array_compatible,v_what_exception,stop_array_compatible)
00227 endif
00228 #endif
00229 res=(.inv.(.tr.m*m))*.tr.m*v
00230 end function leq_leastSquare
00231
00232
00246 function leq_leastSquare_QR(m,v,meth_qr,eps_gsortho,is_permuted) result(res)
00247 implicit none
00248 type(matrix), intent(in) :: m
00249 type(vector), intent(in) :: v
00250 character*(*), intent(in), optional :: meth_qr
00251 type_precision, intent(in), optional :: eps_gsortho
00252 logical, intent(in), optional :: is_permuted
00253
00254 type(type_exception) :: err_exception;
00255 type(vector) :: res
00256 type(t_qr) :: qr_decomp
00257 #ifdef DEBUG_EXCEPTION
00258 if(.not.(m%is_allocate)) then
00259 m_what_exception='leq_leastSquare_QR::matrix not allocated'
00260 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
00261 endif
00262 if(.not.(v%is_allocate)) then
00263 v_what_exception='leq_leastSquare_QR::vector not allocated'
00264 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
00265 endif
00266 if (m%rows.ne.v%size) then
00267 v_what_exception='leq_leastSquare_QR::vector size incompatible with the matrix'
00268 err_exception=e_error(stop_array_compatible,v_what_exception,stop_array_compatible)
00269 endif
00270 #endif
00271 qr_decomp=qr(m,meth_qr,eps_gsortho,is_permuted);
00272 if(present(is_permuted) .and. is_permuted) then
00273 res=.inv.(qr_decomp%R)*.tr.(qr_decomp%Q)*qr_decomp%P*v;
00274 else
00275 res=.inv.(qr_decomp%R)*.tr.(qr_decomp%Q)*v;
00276 end if
00277
00278 call destruct(qr_decomp)
00279 end function leq_leastSquare_QR
00280
00281
00282
00283
00284
00285
00303 function leq_pseudoinverse(m,v,eps_svd,eps_chol,iter_max,meth_qr,eps_gsortho,meth_pinv) result(res)
00304 implicit none
00305 type(matrix), intent(in) :: m
00306 type(vector), intent(in) :: v
00307 type_precision, intent(in), optional :: eps_svd
00308 type_precision, intent(in), optional :: eps_chol
00309 integer, intent(in), optional :: iter_max
00310 character*(*), intent(in), optional :: meth_qr
00311 type_precision, intent(in), optional :: eps_gsortho
00312 character*(*), intent(in), optional :: meth_pinv
00313
00314 type(type_exception) :: err_exception;
00315 type(vector) :: res
00316 #ifdef DEBUG_EXCEPTION
00317 if(.not.(m%is_allocate)) then
00318 m_what_exception='leq_pseudoinverse::matrix not allocated'
00319 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
00320 endif
00321 if(.not.(v%is_allocate)) then
00322 v_what_exception='leq_pseudoinverse::vector not allocated'
00323 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
00324 endif
00325 if (m%rows.ne.v%size) then
00326 v_what_exception='leq_pseudoinverse::vector size incompatible with the matrix'
00327 err_exception=e_error(stop_array_compatible,v_what_exception,stop_array_compatible)
00328 endif
00329 #endif
00330 res=pinv(m,eps_svd,eps_chol,iter_max,meth_qr,eps_gsortho,meth_pinv)*v
00331
00332 end function leq_pseudoinverse
00333
00334
00335
00336
00337
00338
00339
00340
00341
00342
00343
00359 function leq_jacobi_elt(m,v,v0,eps,iter_max,is_precond) result(res)
00360 implicit none
00361 type(matrix), intent(in) :: m
00362 type(vector), intent(in) :: v
00363 type(vector), intent(in), optional :: v0
00364 type_precision, intent(in), optional :: eps
00365 integer, intent(in), optional :: iter_max
00366 logical, intent(in), optional :: is_precond
00367
00368 type(type_exception) :: err_exception;
00369 type(vector) :: res_k
00370 type(vector) :: in_v0
00371 type(t_soleq) :: res
00372 integer :: iter, i, j, in_iter_max
00373 type_precision :: norm_tmp, sum1, sum2, in_eps
00374 type_precision, allocatable, dimension(:) :: err_tmp
00375 logical :: is_converg
00376 #ifdef DEBUG_EXCEPTION
00377 if(.not.(m%is_allocate)) then
00378 m_what_exception='leq_jacobi_elt::matrix not allocated'
00379 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
00380 endif
00381 if(.not.(v%is_allocate)) then
00382 v_what_exception='leq_jacobi_elt::vector not allocated'
00383 err_exception=e_error(stop_aloc,v_what_exception,stop_aloc)
00384 endif
00385 if(present(v0)) then
00386 if(.not.(v0%is_allocate)) then
00387 v_what_exception='leq_jacobi_elt::vector initial not allocated'
00388 err_exception=e_error(stop_aloc,v_what_exception,stop_aloc)
00389 endif
00390 if(v0%size.ne.v%size) then
00391 v_what_exception='leq_jacobi_elt::vector initial is incompatible'
00392 err_exception=e_error(stop_array_compatible,v_what_exception,stop_array_compatible)
00393 endif
00394 end if
00395 if (m%rows.ne.v%size) then
00396 v_what_exception='leq_jacobi_elt::vector size incompatible with the matrix'
00397 err_exception=e_error(stop_array_compatible,v_what_exception,stop_array_compatible)
00398 endif
00399 if(m%rows.ne.m%cols) then
00400 m_what_exception='leq_jacobi_elt::square matrix expected'
00401 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
00402 endif
00403 #endif
00404
00405 if(present(eps)) then
00406 in_eps=eps
00407 else
00408 in_eps=p_iter_eps
00409 end if
00410
00411 if(present(iter_max)) then
00412 in_iter_max=iter_max
00413 else
00414 in_iter_max=100*max(m%rows,m%cols)
00415 end if
00416 if(.not.present(v0) .and. present(is_precond) .and. is_precond) then
00417 in_v0=leq_precond_diag(m,v);
00418 else if(present(v0)) then
00419 in_v0=v0
00420 else
00421 in_v0=v
00422 end if
00423
00424 norm_tmp=0.0
00425 is_converg = .false.
00426 allocate(err_tmp(in_iter_max))
00427 res%v_sol=in_v0
00428 res_k=in_v0
00429
00430 iter=1
00431 do while (iter.lt.in_iter_max .and. .not.(is_converg))
00432 do i=1,v%size
00433 sum1=0
00434 sum2=0
00435 do j=1,i-1
00436 sum1 = sum1 + m%ptr_container(i,j)*res%v_sol%ptr_container(j)
00437 end do
00438 do j=i+1,v%size
00439 sum2 = sum2 + m%ptr_container(i,j)*res_k%ptr_container(j)
00440 end do
00441 res%v_sol%ptr_container(i)=(1/m%ptr_container(i,i))*(v%ptr_container(i)-sum1-sum2)
00442 end do
00443
00444 norm_tmp=.norm.(m*res%v_sol-v)/.norm.v
00445
00446 err_tmp(iter+1)=norm_tmp
00447 if(norm_tmp.gt.in_eps) then
00448 res_k%ptr_container(:)=res%v_sol%ptr_container(:)
00449 else
00450 is_converg=.true.
00451 endif
00452
00453 iter = iter + 1
00454 end do
00455 res%iter=iter-1
00456 call v_init(res%v_err,res%iter)
00457 res%v_err%ptr_container=err_tmp(1:res%iter)
00458 deallocate(err_tmp)
00459 call destruct(res_k)
00460 call destruct(in_v0)
00461 if(.not.(is_converg)) then
00462 v_what_exception='leq_jacobi_elt::Jacobi elt diverge'
00463 err_exception=e_error(stop_array_diverge,v_what_exception,stop_array_diverge)
00464 end if
00465 end function leq_jacobi_elt
00466
00467
00468
00469
00470
00471
00488 function leq_sor_elt(m,v,v0,eps,iter_max,w_relax,is_precond) result(res)
00489 implicit none
00490 type(matrix), intent(in) :: m
00491 type(vector), intent(in) :: v
00492 type(vector), intent(in), optional :: v0
00493 type_precision, intent(in), optional :: eps
00494 integer, intent(in), optional :: iter_max
00495 logical, intent(in), optional :: is_precond
00496 type_precision, intent(in), optional :: w_relax
00497
00498 type(type_exception) :: err_exception;
00499 type(vector) :: res_k, res_GS
00500 type(vector) :: in_v0
00501 type(t_soleq) :: res
00502 integer :: iter, i, j, in_iter_max
00503 type_precision :: sum1, sum2, norm_tmp, in_eps, in_wrelax
00504 type_precision, allocatable, dimension(:) :: err_tmp
00505 logical :: is_converg
00506
00507 if(.not.(m%is_allocate)) then
00508 m_what_exception='leq_sor_elt::matrix not allocated'
00509 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
00510 endif
00511 if(.not.(v%is_allocate)) then
00512 v_what_exception='leq_sor_elt::vector not allocated'
00513 err_exception=e_error(stop_aloc,v_what_exception,stop_aloc)
00514 endif
00515 #ifdef DEBUG_EXCEPTION
00516 if(present(v0)) then
00517 if(.not.(v0%is_allocate)) then
00518 v_what_exception='leq_sor_eltt::vector initial not allocated'
00519 err_exception=e_error(stop_aloc,v_what_exception,stop_aloc)
00520 endif
00521 if(v0%size.ne.v%size) then
00522 v_what_exception='leq_sor_elt::vector initial is incompatible'
00523 err_exception=e_error(stop_array_compatible,v_what_exception,stop_array_compatible)
00524 endif
00525 end if
00526 if (m%rows.ne.v%size) then
00527 v_what_exception='leq_sor_elt::vector size incompatible with the matrix'
00528 err_exception=e_error(stop_array_compatible,v_what_exception,stop_array_compatible)
00529 endif
00530 if(m%rows.ne.m%cols) then
00531 m_what_exception='leq_sor_elt::square matrix expected'
00532 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
00533 endif
00534 #endif
00535
00536 if(present(w_relax)) then
00537 in_wrelax=w_relax
00538 else
00539 in_wrelax=1
00540 end if
00541
00542
00543 if(present(eps)) then
00544 in_eps=eps
00545 else
00546 in_eps=p_iter_eps
00547 end if
00548
00549 if(present(iter_max)) then
00550 in_iter_max=iter_max
00551 else
00552 in_iter_max=100*max(m%rows,m%cols)
00553 end if
00554 if(.not.present(v0) .and. present(is_precond) .and. is_precond) then
00555 in_v0=leq_precond_diag(m,v);
00556 else if(present(v0)) then
00557 in_v0=v0
00558 else
00559 in_v0=v
00560 end if
00561
00562 norm_tmp=0.0
00563 call v_init(res_GS,v%size)
00564 is_converg = .false.
00565 allocate(err_tmp(in_iter_max))
00566 res%v_sol=in_v0
00567 res_k=in_v0
00568
00569 iter=1
00570 do while (iter.lt.in_iter_max .and. .not.(is_converg))
00571 do i=1,v%size
00572 sum1=0
00573 sum2=0
00574 do j=1,i-1
00575 sum1 = sum1 + m%ptr_container(i,j)*res%v_sol%ptr_container(j)
00576 end do
00577 do j=i+1,v%size
00578 sum2 = sum2 + m%ptr_container(i,j)*res_k%ptr_container(j)
00579 end do
00580
00581 res_GS%ptr_container(i)=(1/m%ptr_container(i,i))*(v%ptr_container(i)-sum1-sum2)
00582
00583 res%v_sol%ptr_container(i)=in_wrelax*res_GS%ptr_container(i)+(1-in_wrelax)*res_k%ptr_container(i)
00584 end do
00585
00586 norm_tmp=.norm.(m*res%v_sol-v)/.norm.v
00587 err_tmp(iter+1)=norm_tmp
00588 if(norm_tmp.gt.in_eps) then
00589 res_k%ptr_container(:)=res%v_sol%ptr_container(:)
00590 else
00591 is_converg=.true.
00592 endif
00593
00594 iter = iter + 1
00595 end do
00596 res%iter=iter-1
00597 call v_init(res%v_err,res%iter)
00598 res%v_err%ptr_container=err_tmp(1:res%iter)
00599 deallocate(err_tmp)
00600 call destruct(res_k)
00601 call destruct(res_GS)
00602 call destruct(in_v0)
00603
00604 if(.not.(is_converg)) then
00605 v_what_exception='sor elt diverge'
00606 err_exception=e_error(stop_array_diverge,v_what_exception,stop_array_diverge)
00607 end if
00608 end function leq_sor_elt
00609
00610
00611
00612
00613
00614
00615
00616
00617
00618
00634 function leq_gaussSeidel_elt(m,v,v0,eps,iter_max,is_precond) result(res)
00635 implicit none
00636 type(matrix), intent(in) :: m
00637 type(vector), intent(in) :: v
00638 type(vector), intent(in), optional :: v0
00639 type_precision, intent(in), optional :: eps
00640 integer, intent(in), optional :: iter_max
00641 logical, intent(in), optional :: is_precond
00642
00643 type(type_exception) :: err_exception;
00644 type(vector) :: res_k
00645 type(vector) :: in_v0
00646 type(t_soleq) :: res
00647
00648 integer :: iter, i, j, in_iter_max
00649 type_precision :: sum1, sum2, norm_tmp, in_eps
00650 type_precision, allocatable, dimension(:) :: err_tmp
00651 logical :: is_converg
00652 #ifdef DEBUG_EXCEPTION
00653 if(.not.(m%is_allocate)) then
00654 m_what_exception='leq_gaussSeidel_elt::matrix not allocated'
00655 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
00656 endif
00657 if(.not.(v%is_allocate)) then
00658 v_what_exception='leq_gaussSeidel_elt::vector not allocated'
00659 err_exception=e_error(stop_aloc,v_what_exception,stop_aloc)
00660 endif
00661 if(present(v0)) then
00662 if(.not.(v0%is_allocate)) then
00663 v_what_exception='leq_gaussSeidel_elt::vector initial not allocated'
00664 err_exception=e_error(stop_aloc,v_what_exception,stop_aloc)
00665 endif
00666 if(v0%size.ne.v%size) then
00667 v_what_exception='leq_gaussSeidel_elt::vector initial is incompatible'
00668 err_exception=e_error(stop_array_compatible,v_what_exception,stop_array_compatible)
00669 endif
00670 end if
00671 if (m%rows.ne.v%size) then
00672 v_what_exception='leq_gaussSeidel_elt::vector size incompatible with the matrix'
00673 err_exception=e_error(stop_array_compatible,v_what_exception,stop_array_compatible)
00674 endif
00675 if(m%rows.ne.m%cols) then
00676 m_what_exception='leq_gaussSeidel_elt::square matrix expected'
00677 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
00678 endif
00679 #endif
00680
00681
00682 if(present(eps)) then
00683 in_eps=eps
00684 else
00685 in_eps=p_iter_eps
00686 end if
00687
00688 if(present(iter_max)) then
00689 in_iter_max=iter_max
00690 else
00691 in_iter_max=100*max(m%rows,m%cols)
00692 end if
00693 if(.not.present(v0) .and. present(is_precond) .and. is_precond) then
00694 in_v0=leq_precond_diag(m,v);
00695 else if(present(v0)) then
00696 in_v0=v0
00697 else
00698 in_v0=v
00699 end if
00700
00701 is_converg = .false.
00702 allocate(err_tmp(in_iter_max))
00703 res%v_sol=in_v0
00704 res_k=in_v0
00705 norm_tmp=0.0
00706 iter=1
00707
00708 do while (iter.lt.in_iter_max .and. .not.(is_converg))
00709 do i=1,v%size
00710 sum1=0
00711 sum2=0
00712 do j=1,i-1
00713 sum1 = sum1 + m%ptr_container(i,j)*res%v_sol%ptr_container(j)
00714 end do
00715 do j=i+1,v%size
00716 sum2 = sum2 + m%ptr_container(i,j)*res_k%ptr_container(j)
00717 end do
00718 res%v_sol%ptr_container(i)=(1/m%ptr_container(i,i))*(v%ptr_container(i)-sum1-sum2)
00719 end do
00720 norm_tmp=.norm.(m*res%v_sol-v)/.norm.v
00721 err_tmp(iter+1)=norm_tmp
00722 if(norm_tmp.gt.in_eps) then
00723 res_k%ptr_container(:)=res%v_sol%ptr_container(:)
00724 else
00725 is_converg=.true.
00726 endif
00727
00728 iter = iter + 1
00729 end do
00730
00731 res%iter=iter-1
00732 call v_init(res%v_err,res%iter)
00733 res%v_err%ptr_container=err_tmp(1:res%iter)
00734 deallocate(err_tmp)
00735 call destruct(res_k)
00736 call destruct(in_v0)
00737 if(.not.(is_converg)) then
00738 v_what_exception='leq_gaussSeidel_elt::Gauss-Seidel elt diverge'
00739 err_exception=e_error(stop_array_diverge,v_what_exception,stop_array_diverge)
00740 end if
00741 end function leq_gaussSeidel_elt
00742
00743
00744
00745
00746
00747
00755 function gsb_index_blocks(tab_block, size_block) result (res)
00756 integer, intent(in) :: size_block
00757 integer, intent(in) :: tab_block(size_block)
00758
00759 type(type_exception) :: err_exception;
00760 type_precision :: tmp
00761 integer :: res(size_block)
00762 integer :: i
00763
00764 tmp=0
00765 do i=1,size(res)
00766 res(i) = tmp + tab_block(i)
00767 tmp=res(i)
00768 end do
00769 end function gsb_index_blocks
00770
00771
00788 function gsb_getBlock_m(m,r_s,r_e,c_s,c_e) result (res)
00789 type(matrix), intent(in) :: m
00790 integer, intent(in) :: r_s, r_e
00791 integer, intent(in) :: c_s, c_e
00792
00793 type(type_exception) :: err_exception;
00794 type(matrix) :: res
00795 integer :: i
00796 #ifdef DEBUG_EXCEPTION
00797 if(.not.(m%is_allocate)) then
00798 m_what_exception='gsb_getBlock_m::matrix not allocated'
00799 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
00800 endif
00801 if (r_s.lt.1 .or. r_e.lt.1 .or. r_e.gt.m%rows .or. r_s.gt.m%rows) then
00802 v_what_exception='gsb_getBlock_m::index row incompatible for getting block'
00803 err_exception=e_error(stop_array_indice_exceed,v_what_exception,stop_array_indice_exceed)
00804 endif
00805 if (c_s.lt.1 .or. c_e.lt.1 .or. c_e.gt.m%cols .or. c_s.gt.m%cols) then
00806 v_what_exception='gsb_getBlock_m::index col incompatible for getting matrix block'
00807 err_exception=e_error(stop_array_indice_exceed,v_what_exception,stop_array_indice_exceed)
00808 endif
00809 #endif
00810 call m_init(res,r_e-r_s+1,c_e-c_s+1)
00811 if(r_s.eq.r_e) then
00812 if(c_s.eq.c_e) then
00813 res%ptr_container(1,1)=m%ptr_container(r_s,c_s)
00814 else
00815 res%ptr_container(1,:)=m%ptr_container(r_s,c_s:c_e)
00816 end if
00817 else if(c_s.eq.c_e) then
00818 res%ptr_container(:,1)=m%ptr_container(r_s:r_e,c_s)
00819 else
00820 res%ptr_container(:,:)=m%ptr_container(r_s:r_e,c_s:c_e)
00821 end if
00822
00823 end function gsb_getBlock_m
00824
00825
00839 function leq_getBlock_v(v,r_s,r_e) result (res)
00840 type(vector), intent(in) :: v
00841 integer, intent(in) :: r_s, r_e
00842
00843 type(type_exception) :: err_exception;
00844 type(vector) :: res
00845 integer :: i
00846 type_precision :: res_pb(v%size)
00847 #ifdef DEBUG_EXCEPTION
00848 if(.not.(v%is_allocate)) then
00849 v_what_exception='leq_getBlock_v::vector not allocated'
00850 err_exception=e_error(stop_aloc,v_what_exception,stop_aloc)
00851 endif
00852 if (r_s.lt.1 .or. r_e.lt.1 .or. r_e.gt.v%size .or. r_s.gt.v%size) then
00853 v_what_exception='leq_getBlock_v::index row incompatible for getting vector block'
00854 err_exception=e_error(stop_array_indice_exceed,v_what_exception,stop_array_indice_exceed)
00855 endif
00856 #endif
00857 call v_init(res,r_e-r_s+1)
00858 res_pb=v%ptr_container
00859 res%ptr_container(:)=res_pb(r_s:r_e)
00860 end function leq_getBlock_v
00861
00862
00880 function leq_gaussSeidel_block(m,v,v0,tab_block_i,size_block_i,eps,iter_max,is_precond) result(res)
00881 implicit none
00882 type(matrix), intent(in) :: m
00883 type(vector), intent(in) :: v
00884 type(vector), intent(in), optional :: v0
00885 integer, intent(in), optional :: tab_block_i(:)
00886 integer, intent(in), optional :: size_block_i
00887 type_precision, intent(in), optional :: eps
00888 integer, intent(in), optional :: iter_max
00889 logical, intent(in), optional :: is_precond
00890
00891 type(type_exception) :: err_exception;
00892 type(t_soleq) :: res
00893 type(vector) :: res_k
00894 type(vector) :: in_v0
00895 type(vector) :: v_tmp, v_sum1, v_sum2
00896 type(matrix) :: block_M, block_D
00897 type(vector) :: block_Res, block_B
00898 integer, allocatable, dimension(:) :: tab_index, tab_block_index
00899 integer :: iter, i, j, k, nb_block, in_iter_max, mod_size_block
00900 integer :: b_inf, b_sup
00901 type_precision :: norm_tmp, in_eps
00902 type_precision, allocatable, dimension(:) :: err_tmp
00903 logical :: is_converg
00904 #ifdef DEBUG_EXCEPTION
00905 if(.not.(m%is_allocate)) then
00906 m_what_exception='leq_gaussSeidel_block::matrix not allocated'
00907 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
00908 endif
00909 if(.not.(v%is_allocate)) then
00910 v_what_exception='leq_gaussSeidel_block::vector not allocated'
00911 err_exception=e_error(stop_aloc,v_what_exception,stop_aloc)
00912 endif
00913 if(present(v0)) then
00914 if(.not.(v0%is_allocate)) then
00915 v_what_exception='leq_gaussSeidel_block::vector initial not allocated'
00916 err_exception=e_error(stop_aloc,v_what_exception,stop_aloc)
00917 endif
00918 if(v0%size.ne.v%size) then
00919 v_what_exception='leq_gaussSeidel_block::vector initial is incompatible'
00920 err_exception=e_error(stop_array_compatible,v_what_exception,stop_array_compatible)
00921 endif
00922 end if
00923 if (m%rows.ne.v%size) then
00924 v_what_exception='leq_gaussSeidel_block::vector size incompatible with the matrix'
00925 err_exception=e_error(stop_array_compatible,v_what_exception,stop_array_compatible)
00926 endif
00927 if(m%rows.ne.m%cols) then
00928 m_what_exception='leq_gaussSeidel_bloc::square matrix expected'
00929 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
00930 endif
00931 #endif
00932
00933
00934 if(present(tab_block_i)) then
00935 nb_block=size_block_i
00936 allocate(tab_block_index(nb_block))
00937 tab_block_index(:)=tab_block_i(:)
00938 else
00939 nb_block=2
00940 mod_size_block=modulo(m%rows,nb_block)
00941 if(mod_size_block.eq.0) then
00942 allocate(tab_block_index(nb_block))
00943 tab_block_index(:)=m%rows/nb_block
00944 else
00945 allocate(tab_block_index(nb_block))
00946 tab_block_index(1:nb_block-1)=int(m%rows/nb_block)
00947 tab_block_index(nb_block)=int(m%rows/nb_block) + mod_size_block
00948 end if
00949
00950 end if
00951
00952 #ifdef DEBUG_EXCEPTION
00953 if(sum(tab_block_index).ne.v%size) then
00954 v_what_exception='leq_gaussSeidel_block::block index incompatible'
00955 err_exception=e_error(stop_array_compatible,v_what_exception,stop_array_compatible)
00956 endif
00957 #endif
00958
00959
00960 if(present(iter_max)) then
00961 in_iter_max=iter_max
00962 else
00963 in_iter_max=100*max(m%rows,m%cols)
00964 end if
00965 if(.not.present(v0) .and. present(is_precond) .and. is_precond) then
00966 in_v0=leq_precond_diag(m,v);
00967 else if(present(v0)) then
00968 in_v0=v0
00969 else
00970 in_v0=v
00971 end if
00972 if(present(eps)) then
00973 in_eps=eps
00974 else
00975 in_eps=p_iter_eps
00976 end if
00977
00978
00979 allocate(tab_index(nb_block))
00980 tab_index=gsb_index_blocks(tab_block_index,nb_block)
00981 deallocate(tab_block_index);
00982
00983 allocate(err_tmp(in_iter_max))
00984 norm_tmp=0.0
00985 is_converg = .false.
00986 res%v_sol=in_v0
00987 res_k=in_v0
00988 iter=1
00989
00990 do while (iter.lt.in_iter_max .and. .not.(is_converg))
00991
00992 do i=1,nb_block
00993 if(i.eq.1) then
00994 b_inf=1;
00995 else
00996 b_inf=tab_index(i-1)+1
00997 end if
00998
00999 b_sup=tab_index(i)
01000
01001 call v_resize(v_sum1,b_sup-b_inf+1)
01002 call v_resize(v_sum2,b_sup-b_inf+1)
01003
01004 do j=1,i-1
01005 if((i-1).eq.0) then
01006 if((j-1).eq.0) then
01007 block_M=leq_getBlock(m,1,tab_index(i),1,tab_index(j))
01008 block_Res=leq_getBlock(res%v_sol,1,tab_index(j))
01009
01010 else
01011 block_M=leq_getBlock(m,1,tab_index(i),tab_index(j-1)+1,tab_index(j))
01012
01013 block_Res=leq_getBlock(res%v_sol,tab_index(j-1)+1,tab_index(j))
01014 end if
01015 else
01016 if((j-1).eq.0) then
01017 block_M=leq_getBlock(m,tab_index(i-1)+1,tab_index(i),1,tab_index(j))
01018 block_Res=leq_getBlock(res%v_sol,1,tab_index(j))
01019 else
01020 block_M=leq_getBlock(m,tab_index(i-1)+1,tab_index(i),tab_index(j-1)+1,tab_index(j))
01021 block_Res=leq_getBlock(res%v_sol,tab_index(j-1)+1,tab_index(j))
01022 end if
01023 end if
01024
01025 v_sum1 = v_sum1 + block_M*block_Res
01026 call destruct(block_M)
01027 call destruct(block_Res)
01028 end do
01029
01030
01031 do j=i+1,nb_block
01032 if((i-1).eq.0) then
01033 block_M=leq_getBlock(m,1,tab_index(i),tab_index(j-1)+1,tab_index(j))
01034 block_Res=leq_getBlock(res_k,tab_index(j-1)+1,tab_index(j))
01035 else
01036 block_M=leq_getBlock(m,tab_index(i-1)+1,tab_index(i),tab_index(j-1)+1,tab_index(j))
01037 block_Res=leq_getBlock(res_k,tab_index(j-1)+1,tab_index(j))
01038 end if
01039
01040 v_sum2= v_sum2 + block_M*block_Res
01041
01042 call destruct(block_M)
01043 call destruct(block_Res)
01044 end do
01045
01046
01047 if((i-1).eq.0) then
01048 block_D=leq_getBlock(m,1,tab_index(i),1,tab_index(i))
01049 block_B=leq_getBlock(v,1,tab_index(i))
01050 else
01051 block_D=leq_getBlock(m,tab_index(i-1)+1,tab_index(i),tab_index(i-1)+1,tab_index(i))
01052 block_B=leq_getBlock(v,tab_index(i-1)+1,tab_index(i))
01053 end if
01054
01055
01056 v_tmp=.inv.block_D*(block_B-v_sum1-v_sum2)
01057 res%v_sol%ptr_container(b_inf:b_sup)=v_tmp%ptr_container
01058 call destruct(block_D)
01059 call destruct(block_B)
01060 call destruct(v_tmp)
01061 end do
01062
01063
01064 norm_tmp=.norm.(m*res%v_sol-v)/.norm.v;
01065 err_tmp(iter+1)=norm_tmp
01066 if(norm_tmp.gt.in_eps) then
01067 res_k%ptr_container(:)=res%v_sol%ptr_container(:)
01068 else
01069 is_converg=.true.
01070 endif
01071 iter = iter + 1
01072 end do
01073 res%iter=iter-1
01074 deallocate(tab_index)
01075 call v_init(res%v_err,res%iter)
01076 res%v_err%ptr_container=err_tmp(1:res%iter)
01077 deallocate(err_tmp)
01078
01079 call destruct(v_sum1);
01080 call destruct(v_sum2);
01081 call destruct(block_M)
01082 call destruct(block_Res)
01083 call destruct(block_D)
01084 call destruct(block_B)
01085 call destruct(v_tmp)
01086 call destruct(in_v0)
01087
01088 if(.not.(is_converg)) then
01089 v_what_exception='leq_gaussSeidel_block::Gauss-Seidel block diverge'
01090 err_exception=e_error(stop_array_diverge,v_what_exception,stop_array_diverge)
01091 end if
01092
01093 end function leq_gaussSeidel_block
01094
01095
01096
01097
01098
01099
01100
01101
01117 function leq_cg_elt(m,v,v0,eps,iter_max,is_precond) result(res)
01118 implicit none
01119 type(matrix), intent(in) :: m
01120 type(vector), intent(in) :: v
01121 type(vector), intent(in), optional :: v0
01122 type_precision, intent(in), optional :: eps
01123 integer, intent(in), optional :: iter_max
01124 logical, intent(in), optional :: is_precond
01125
01126 type(type_exception) :: err_exception;
01127 type(vector) ::p, r, y
01128 type(vector) ::res_k, p_k, r_k
01129 type(vector) :: in_v0
01130 type(t_soleq) :: res
01131 integer :: iter, in_iter_max
01132 type_precision :: alpha, beta
01133 type_precision :: norm_tmp, in_eps
01134 type_precision, allocatable, dimension(:) :: err_tmp
01135 logical :: is_converg
01136 #ifdef DEBUG_EXCEPTION
01137 if(.not.(m%is_allocate)) then
01138 m_what_exception='leq_cg_elt::matrix not allocated'
01139 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01140 endif
01141 if(.not.(v%is_allocate)) then
01142 v_what_exception='leq_cg_elt:vector not allocated'
01143 err_exception=e_error(stop_aloc,v_what_exception,stop_aloc)
01144 endif
01145 if(present(v0)) then
01146 if(.not.(v0%is_allocate)) then
01147 v_what_exception='leq_cg_elt::vector initial not allocated'
01148 err_exception=e_error(stop_aloc,v_what_exception,stop_aloc)
01149 endif
01150 if(v0%size.ne.v%size) then
01151 v_what_exception='leq_cg_elt::vector initial is incompatible'
01152 err_exception=e_error(stop_array_compatible,v_what_exception,stop_array_compatible)
01153 endif
01154 end if
01155 if (m%rows.ne.v%size) then
01156 v_what_exception='leq_cg_elt:vector size incompatible with the matrix'
01157 err_exception=e_error(stop_array_compatible,v_what_exception,stop_array_compatible)
01158 endif
01159 if(m%rows.ne.m%cols) then
01160 m_what_exception='leq_cg_elt::square matrix expected'
01161 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
01162 endif
01163 if (.not.m_isSymmetric(m)) then
01164 v_what_exception='leq_cg_elt:symmetric matrix expected'
01165 err_exception=e_error(stop_array_symmetric,v_what_exception,stop_array_symmetric)
01166 endif
01167 #endif
01168
01169
01170 if(present(eps)) then
01171 in_eps=eps
01172 else
01173 in_eps=p_iter_eps
01174 end if
01175
01176 if(present(iter_max)) then
01177 in_iter_max=iter_max
01178 else
01179 in_iter_max=100*max(m%rows,m%cols)
01180 end if
01181 if(.not.present(v0) .and. present(is_precond) .and. is_precond) then
01182 in_v0=leq_precond_diag(m,v);
01183 else if(present(v0)) then
01184 in_v0=v0
01185 else
01186 in_v0=v
01187 end if
01188
01189 is_converg = .false.
01190 norm_tmp=0.0
01191 allocate(err_tmp(in_iter_max))
01192 p=v-m*in_v0
01193 r=p
01194 y=m*p
01195 alpha=(.norm.r)**2/(y.dot.p)
01196 res%v_sol=in_v0+p_notcast(alpha)*p
01197
01198
01199 iter=1
01200 do while (iter.lt.in_iter_max .and. .not.(is_converg))
01201 res_k=res%v_sol
01202 p_k=p
01203 r_k=r
01204 r=r_k-p_notcast(alpha)*y
01205 beta=(.norm.r)**2/(.norm.r_k)**2
01206
01207 p=r+p_notcast(beta)*p_k
01208 y=m*p
01209 alpha=(.norm.r)**2/(y.dot.p)
01210 res%v_sol=res_k+p_notcast(alpha)*p
01211 norm_tmp=.norm.r_k/.norm.v;
01212 err_tmp(iter+1)=norm_tmp
01213 if(norm_tmp.gt.in_eps) then
01214 res_k%ptr_container(:)=res%v_sol%ptr_container(:)
01215 else
01216 is_converg=.true.
01217 endif
01218 iter = iter + 1
01219 end do
01220 res%iter=iter-1
01221 call v_init(res%v_err,res%iter)
01222 res%v_err%ptr_container=err_tmp(1:res%iter)
01223 deallocate(err_tmp)
01224 call destruct(res_k)
01225 call destruct(y)
01226 call destruct(p)
01227 call destruct(r)
01228 call destruct(in_v0)
01229 if(.not.(is_converg)) then
01230 v_what_exception='leq_cg_elt:Conjugate gradient diverge'
01231 err_exception=e_error(stop_array_diverge,v_what_exception,stop_array_diverge)
01232 end if
01233 end function leq_cg_elt
01234
01235
01236
01237
01238
01239
01240
01256 function leq_gradient_elt(m,v,v0,eps,iter_max,is_precond) result(res)
01257 implicit none
01258 type(matrix), intent(in) :: m
01259 type(vector), intent(in) :: v
01260 type(vector), intent(in), optional :: v0
01261 type_precision, intent(in), optional :: eps
01262 integer, intent(in), optional :: iter_max
01263 logical, intent(in), optional :: is_precond
01264
01265 type(type_exception) :: err_exception;
01266 type(vector) ::p, r, y
01267 type(vector) ::res_k, p_k, r_k
01268 type(vector) :: in_v0
01269 type(t_soleq) :: res
01270 integer :: iter, in_iter_max
01271 type_precision :: alpha, beta
01272 type_precision :: norm_tmp, in_eps
01273 type_precision, allocatable, dimension(:) :: err_tmp
01274 logical :: is_converg
01275 #ifdef DEBUG_EXCEPTION
01276 if(.not.(m%is_allocate)) then
01277 m_what_exception='leq_gradient_elt::matrix not allocated'
01278 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01279 endif
01280 if(.not.(v%is_allocate)) then
01281 v_what_exception='leq_gradient_elt::vector not allocated'
01282 err_exception=e_error(stop_aloc,v_what_exception,stop_aloc)
01283 endif
01284 if(present(v0)) then
01285 if(.not.(v0%is_allocate)) then
01286 v_what_exception='leq_gradient_elt::vector initial not allocated'
01287 err_exception=e_error(stop_aloc,v_what_exception,stop_aloc)
01288 endif
01289 if(v0%size.ne.v%size) then
01290 v_what_exception='leq_gradient_elt::vector initial is incompatible'
01291 err_exception=e_error(stop_array_compatible,v_what_exception,stop_array_compatible)
01292 endif
01293 end if
01294 if (m%rows.ne.v%size) then
01295 v_what_exception='leq_gradient_elt::vector size incompatible with the matrix'
01296 err_exception=e_error(stop_array_compatible,v_what_exception,stop_array_compatible)
01297 endif
01298 if(m%rows.ne.m%cols) then
01299 m_what_exception='leq_gradient_elt::square matrix expected'
01300 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
01301 endif
01302 if (.not.m_isSymmetric(m)) then
01303 v_what_exception='leq_gradient_elt::symmetric matrix expected'
01304 err_exception=e_error(stop_array_symmetric,v_what_exception,stop_array_symmetric)
01305 endif
01306 #endif
01307
01308 if(present(eps)) then
01309 in_eps=eps
01310 else
01311 in_eps=p_iter_eps
01312 end if
01313
01314 if(present(iter_max)) then
01315 in_iter_max=iter_max
01316 else
01317 in_iter_max=100*max(m%rows,m%cols)
01318 end if
01319 if(.not.present(v0) .and. present(is_precond) .and. is_precond) then
01320 in_v0=leq_precond_diag(m,v);
01321 else if(present(v0)) then
01322 in_v0=v0
01323 else
01324 in_v0=v
01325 end if
01326
01327
01328 allocate(err_tmp(in_iter_max))
01329 norm_tmp=0.0
01330 is_converg = .false.
01331 p=v-m*in_v0
01332 r=p
01333 y=m*p
01334 alpha=(.norm.r)**2/(y.dot.p)
01335 res%v_sol=in_v0+p_notcast(alpha)*p
01336
01337
01338 iter=1
01339 do while (iter.lt.in_iter_max .and. .not.(is_converg))
01340 res_k=res%v_sol
01341 p_k=p
01342 r_k=r
01343 r=r_k-p_notcast(alpha)*y
01344
01345 p=r
01346 y=m*p
01347 alpha=(.norm.r)**2/(y.dot.p)
01348 res%v_sol=res_k+p_notcast(alpha)*p
01349 norm_tmp=.norm.r_k/.norm.v;
01350 err_tmp(iter+1)=norm_tmp
01351 if(norm_tmp.gt.in_eps) then
01352 res_k%ptr_container(:)=res%v_sol%ptr_container(:)
01353
01354 else
01355 is_converg=.true.
01356 endif
01357
01358 iter = iter + 1
01359 end do
01360 res%iter=iter-1
01361 call v_init(res%v_err,res%iter)
01362 res%v_err%ptr_container=err_tmp(1:res%iter)
01363 deallocate(err_tmp)
01364 call destruct(res_k)
01365 call destruct(y)
01366 call destruct(p)
01367 call destruct(r)
01368 call destruct(in_v0)
01369 if(.not.(is_converg)) then
01370 v_what_exception='leq_gradient_elt::gradient method diverge'
01371 err_exception=e_error(stop_array_diverge,v_what_exception,stop_array_diverge)
01372 end if
01373 end function leq_gradient_elt
01374
01375
01376
01377
01378
01379
01380
01381
01393 function leq_thomas(m,v) result(res)
01394 implicit none
01395 type(matrix), intent(in) :: m
01396 type(vector), intent(in) :: v
01397
01398 type(type_exception) :: err_exception;
01399 type(vector) :: res, Betaa, Gammaa
01400 type(vector) :: diagCentral, diagU, diagL
01401 integer :: i, k
01402 #ifdef DEBUG_EXCEPTION
01403 if(.not.(m%is_allocate)) then
01404 m_what_exception='leq_thomas::matrix not allocated'
01405 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01406 endif
01407 if(.not.(v%is_allocate)) then
01408 v_what_exception='leq_thomas::vector not allocated'
01409 err_exception=e_error(stop_aloc,v_what_exception,stop_aloc)
01410 endif
01411 if (m%rows.ne.v%size) then
01412 v_what_exception='leq_thomas::vector size incompatible with the matrix'
01413 err_exception=e_error(stop_array_compatible,v_what_exception,stop_array_compatible)
01414 endif
01415 if(m%rows.ne.m%cols) then
01416 m_what_exception='leq_thomas::square matrix expected'
01417 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
01418 endif
01419 #endif
01420
01421 diagCentral=diag(m)
01422 diagU=diag(m,1)
01423 diagL=diag(m,-1)
01424
01425 call v_init(res,v%size)
01426 call v_init(Betaa,v%size)
01427 call v_init(Gammaa,v%size)
01428 Betaa%ptr_container(1)=m%ptr_container(1,1)
01429 Gammaa%ptr_container(1)=v%ptr_container(1)/m%ptr_container(1,1)
01430 do i=2,v%size
01431 Betaa%ptr_container(i)=diagCentral%ptr_container(i) &
01432 - ((diagU%ptr_container(i-1)/Betaa%ptr_container(i-1))*diagL%ptr_container(i-1))
01433 Gammaa%ptr_container(i)= (1/Betaa%ptr_container(i))* &
01434 (v%ptr_container(i)-(diagL%ptr_container(i-1)*Gammaa%ptr_container(i-1)))
01435 end do
01436
01437 res%ptr_container(v%size)=Gammaa%ptr_container(v%size)
01438
01439 do k=v%size-1,1,-1
01440 res%ptr_container(k)=Gammaa%ptr_container(k) - &
01441 ((diagU%ptr_container(k)/Betaa%ptr_container(k))*res%ptr_container(k+1))
01442 end do
01443 call destruct(Betaa)
01444 call destruct(Gammaa)
01445 call destruct(diagCentral)
01446 call destruct(diagL)
01447 call destruct(diagU)
01448 end function leq_thomas
01449
01450
01451
01452
01453
01454
01466 function leq_tril(m,v) result(res)
01467 implicit none
01468 type(matrix), intent(in) :: m
01469 type(vector), intent(in) :: v
01470
01471 type(type_exception) :: err_exception;
01472 type(vector) :: res
01473 integer :: i, j, N
01474 #ifdef DEBUG_EXCEPTION
01475 if(.not.(m%is_allocate)) then
01476 m_what_exception='leq_tril::matrix not allocated'
01477 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01478 endif
01479 if(.not.(v%is_allocate)) then
01480 v_what_exception='leq_tril::vector not allocated'
01481 err_exception=e_error(stop_aloc,v_what_exception,stop_aloc)
01482 endif
01483 if (m%rows.ne.v%size) then
01484 v_what_exception='leq_tril::vector size incompatible with the matrix'
01485 err_exception=e_error(stop_array_compatible,v_what_exception,stop_array_compatible)
01486 endif
01487 if(m%rows.ne.m%cols) then
01488 m_what_exception='leq_tril::square matrix expected'
01489 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
01490 endif
01491 #endif
01492 res=v;
01493 N=v%size
01494 do i=1,N
01495 do j=1,i-1
01496 res%ptr_container(i)=res%ptr_container(i)-m%ptr_container(i,j)*res%ptr_container(j)
01497 end do
01498 res%ptr_container(i)=res%ptr_container(i)/m%ptr_container(i,i)
01499 end do
01500
01501 end function leq_tril
01502
01503
01504
01505
01506
01507
01519 function leq_triu(m,v) result(res)
01520 implicit none
01521 type(matrix), intent(in) :: m
01522 type(vector), intent(in) :: v
01523
01524 type(type_exception) :: err_exception;
01525 type(vector) :: res
01526 integer :: i, j, N
01527 #ifdef DEBUG_EXCEPTION
01528 if(.not.(m%is_allocate)) then
01529 m_what_exception='leq_tril::matrix not allocated'
01530 err_exception=e_error(stop_aloc,m_what_exception,stop_aloc)
01531 endif
01532 if(.not.(v%is_allocate)) then
01533 v_what_exception='leq_tril::vector not allocated'
01534 err_exception=e_error(stop_aloc,v_what_exception,stop_aloc)
01535 endif
01536 if (m%rows.ne.v%size) then
01537 v_what_exception='leq_triu::vector size incompatible with the matrix'
01538 err_exception=e_error(stop_array_compatible,v_what_exception,stop_array_compatible)
01539 endif
01540 if(m%rows.ne.m%cols) then
01541 m_what_exception='leq_triu::square matrix expected'
01542 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
01543 endif
01544 #endif
01545 res=v;
01546 N=v%size
01547 do i=N,1,-1
01548 do j=i+1,N
01549 res%ptr_container(i)=res%ptr_container(i)-m%ptr_container(i,j)*res%ptr_container(j)
01550 end do
01551 res%ptr_container(i)=res%ptr_container(i)/m%ptr_container(i,i)
01552 end do
01553
01554 end function leq_triu
01555
01556
01557
01558
01559
01560
01574 function leq_lu(m,v,is_permuted) result(res)
01575 implicit none
01576 type(matrix), intent(in) :: m
01577 type(vector), intent(in) :: v
01578 logical, intent(in), optional :: is_permuted
01579
01580 type(type_exception) :: err_exception;
01581 type(vector) :: res
01582 type(t_m_and_p) :: lu_decomp;
01583 integer :: i, j
01584 type_precision :: tp_s, det
01585 #ifdef DEBUG_EXCEPTION
01586 if (m%rows.ne.v%size) then
01587 v_what_exception='leq_lu::vector size incompatible with the matrix'
01588 err_exception=e_error(stop_array_compatible,v_what_exception,stop_array_compatible)
01589 endif
01590 if(m%rows.ne.m%cols) then
01591 m_what_exception='leq_lu::square matrix expected'
01592 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
01593 endif
01594 #endif
01595 lu_decomp=m_lu(m,is_permuted);
01596
01597 det=v_prod(diag(lu_decomp%M))
01598 if(present(is_permuted) .and. is_permuted) then
01599
01600 det=((-1)**(lu_decomp%swops))*det
01601 end if
01602 #ifdef DEBUG_EXCEPTION
01603 if(abs(det).le.math_machine_eps()) then
01604 m_what_exception='leq_lu::regular matrix expected'
01605 err_exception=e_error(stop_array_singular,m_what_exception,stop_array_singular)
01606 end if
01607 #endif
01608
01609 if(present(is_permuted) .and. is_permuted) then
01610 res=lu_decomp%P*v;
01611 else
01612 res=v;
01613 end if
01614
01615 do i=1,lu_decomp%M%rows
01616 tp_s=res%ptr_container(i);
01617 do j=1,i-1
01618 tp_s=tp_s-lu_decomp%M%ptr_container(i,j)*res%ptr_container(j)
01619 end do
01620 res%ptr_container(i)=tp_s
01621 end do
01622
01623 do i=lu_decomp%M%rows,1,-1
01624 tp_s=res%ptr_container(i);
01625 do j=i+1,lu_decomp%M%rows
01626 tp_s=tp_s-lu_decomp%M%ptr_container(i,j)*res%ptr_container(j)
01627 end do
01628 res%ptr_container(i)=tp_s/lu_decomp%M%ptr_container(i,i)
01629 end do
01630 call destruct(lu_decomp)
01631 end function leq_lu
01632
01633
01634
01635
01636
01637
01652 function leq_Cholesky(m,v,is_permuted) result(res)
01653 implicit none
01654 type(matrix), intent(in) :: m
01655 type(vector), intent(in) :: v
01656 logical, intent(in), optional :: is_permuted
01657
01658 type(type_exception) :: err_exception;
01659 type(vector) :: res
01660 type(t_m_and_p) :: m_chol
01661 type(matrix) :: m_chol_tr
01662 integer :: i, j, N
01663 #ifdef DEBUG_EXCEPTION
01664 if (m%rows.ne.v%size) then
01665 v_what_exception='leq_Cholesky::vector size incompatible with the matrix'
01666 err_exception=e_error(stop_array_compatible,v_what_exception,stop_array_compatible)
01667 endif
01668 if(m%rows.ne.m%cols) then
01669 m_what_exception='leq_Cholesky::square matrix expected'
01670 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
01671 endif
01672
01673 if(abs(.det.m).eq.math_machine_eps()) then
01674 m_what_exception='leq_Cholesky::regular matrix expected'
01675 err_exception=e_error(stop_array_singular,m_what_exception,stop_array_singular)
01676 end if
01677 #endif
01678 m_chol=chol(m,is_permuted);
01679 m_chol_tr=.tr.m_chol%M
01680
01681 if(present(is_permuted) .and. is_permuted) then
01682 res=m_chol%P*v;
01683 else
01684 res=v;
01685 end if
01686 N=res%size
01687 do i=1,N
01688 do j=1,i-1
01689 res%ptr_container(i)=res%ptr_container(i)-m_chol%M%ptr_container(i,j)*res%ptr_container(j)
01690 end do
01691 res%ptr_container(i)=res%ptr_container(i)/m_chol%M%ptr_container(i,i)
01692 end do
01693
01694 do i=N,1,-1
01695 do j=i+1,N
01696 res%ptr_container(i)=res%ptr_container(i)-m_chol_tr%ptr_container(i,j)*res%ptr_container(j)
01697 end do
01698 res%ptr_container(i)=res%ptr_container(i)/m_chol_tr%ptr_container(i,i)
01699 end do
01700
01701 call destruct(m_chol)
01702 call destruct(m_chol_tr)
01703 end function leq_Cholesky
01704
01705
01706
01707
01708
01709
01725 function leq_qr(m,v,meth_qr,eps_gsortho,is_permuted) result(res)
01726 implicit none
01727 type(matrix), intent(in) :: m
01728 type(vector), intent(in) :: v
01729 character*(*), intent(in), optional :: meth_qr
01730 type_precision, intent(in), optional :: eps_gsortho
01731 logical, intent(in), optional :: is_permuted
01732
01733 type(type_exception) :: err_exception;
01734 type(vector) :: res
01735 type(t_qr) :: qr_decomp
01736 integer :: i, j, N
01737 #ifdef DEBUG_EXCEPTION
01738 if (m%rows.ne.v%size) then
01739 v_what_exception='leq_qr::vector size incompatible with the matrix'
01740 err_exception=e_error(stop_array_compatible,v_what_exception,stop_array_compatible)
01741 endif
01742 if (m%rows.lt.m%cols) then
01743 m_what_exception='m_decompQR_GramSchimdt:: more unknown than equations'
01744 err_exception=e_error(stop_array_compatible,m_what_exception,stop_array_compatible)
01745 endif
01746 #endif
01747 qr_decomp=qr(m,meth_qr,eps_gsortho,is_permuted);
01748
01749 if(present(is_permuted) .and. is_permuted) then
01750 res=.tr.qr_decomp%Q*qr_decomp%P*v;
01751 else
01752 res=.tr.qr_decomp%Q*v;
01753 end if
01754 N=res%size
01755 do i=N,1,-1
01756 do j=i+1,N
01757 res%ptr_container(i)=res%ptr_container(i)-qr_decomp%R%ptr_container(i,j)*res%ptr_container(j)
01758 end do
01759 res%ptr_container(i)=res%ptr_container(i)/qr_decomp%R%ptr_container(i,i)
01760 end do
01761
01762 call destruct(qr_decomp)
01763 end function leq_qr
01764
01765
01766
01767
01768
01769
01770
01771
01772
01773
01774
01782 function leq_precond_diag(m,v) result (res)
01783 type(matrix), intent(in) :: m
01784 type(vector), intent(in) :: v
01785
01786 type(vector) :: res
01787 res=.inv.(diag(m))*v
01788 end function leq_precond_diag
01789
01790
01791
01792
01793
01794
01804 function leq_precond_sor(m,v,v0,w_relax) result (res)
01805 type(matrix), intent(in) :: m
01806 type(vector), intent(in) :: v
01807 type(vector), intent(in) :: v0
01808 type_precision, intent(in) :: w_relax
01809
01810 type(vector) :: res
01811 res=(1-w_relax)*v0 + w_relax*v
01812 end function leq_precond_sor
01813
01814
01815
01816
01817
01818
01836 function leq_solve(m,v,meth_solve,v0,eps_svd,eps_chol,iter_max,meth_qr,meth_pinv,eps_gsortho,is_permuted) result(res)
01837 implicit none
01838 type(matrix), intent(in) :: m
01839 type(vector), intent(in) :: v
01840 character*(*), intent(in), optional :: meth_solve
01841 type(vector), intent(in), optional :: v0
01842 type_precision, intent(in), optional :: eps_svd
01843 type_precision, intent(in), optional :: eps_chol
01844 integer, intent(in), optional :: iter_max
01845 character*(*), intent(in), optional :: meth_qr
01846 character*(*), intent(in), optional :: meth_pinv
01847 type_precision, intent(in), optional :: eps_gsortho
01848 logical, intent(in), optional :: is_permuted
01849
01850 type(type_exception) :: err_exception;
01851 type(vector) :: res
01852
01853 if(present(meth_solve)) then
01854 select case(meth_solve)
01855 case ('gj')
01856 res=leq_gaussJourdan(m,v);
01857 case ('ls')
01858 res=leq_leastSquare(m,v);
01859 case ('lsqr')
01860 res=leq_leastSquare_QR(m,v,meth_qr,eps_gsortho,is_permuted);
01861 case ('pinvsvd')
01862 res=leq_pseudoinverse(m,v,eps_svd,eps_chol,iter_max,meth_qr,eps_gsortho,meth_pinv='svd');
01863 case ('pinvchol')
01864 res=leq_pseudoinverse(m,v,eps_svd,eps_chol,iter_max,meth_qr,eps_gsortho,meth_pinv='chol');
01865 case ('thomaselt')
01866 res=leq_thomas(m,v);
01867 case ('tril')
01868 res=leq_tril(m,v);
01869 case ('triu')
01870 res=leq_triu(m,v);
01871 case ('lu')
01872 res=leq_lu(m,v,is_permuted);
01873 case ('chol')
01874 res=leq_Cholesky(m,v,is_permuted);
01875 case ('qr')
01876 res=leq_qr(m,v,meth_qr,eps_gsortho,is_permuted);
01877 case default
01878 #ifdef DEBUG_EXCEPTION
01879 print*, "This method is not exist, your linear equation will be solved by gauss-jourdan element method"
01880 #endif
01881 res=leq_gaussJourdan(m,v);
01882 end select
01883 else
01884 res=leq_gaussJourdan(m,v)
01885 end if
01886
01887 end function leq_solve
01888
01889
01905 function leq_solve_iter(m,v,meth_solve,v0,eps,iter_max,is_precond,w_relax,tab_block_i,size_block_i) result(res)
01906 implicit none
01907 type(matrix), intent(in) :: m
01908 type(vector), intent(in) :: v
01909 character*(*), intent(in), optional :: meth_solve
01910 type(vector), intent(in),optional :: v0
01911 type_precision, intent(in), optional :: eps
01912 integer, intent(in), optional :: iter_max
01913 logical, intent(in), optional :: is_precond
01914 type_precision, intent(in), optional :: w_relax
01915 integer, intent(in), optional :: tab_block_i(:)
01916 integer, intent(in), optional :: size_block_i
01917
01918 type(type_exception) :: err_exception;
01919 type(t_soleq) :: res;
01920
01921 if(present(meth_solve)) then
01922 select case(meth_solve)
01923 case ('jacelt')
01924 res=leq_jacobi_elt(m,v,v0,eps,iter_max,is_precond);
01925 case ('sorelt')
01926 res=leq_sor_elt(m,v,v0,eps,iter_max,w_relax,is_precond);
01927 case ('gselt')
01928 res=leq_gaussSeidel_elt(m,v,v0,eps,iter_max,is_precond);
01929 case ('gsblock')
01930 res=leq_gaussSeidel_block(m,v,v0,tab_block_i,size_block_i,eps,iter_max,is_precond);
01931 case ('cgelt')
01932 res=leq_cg_elt(m,v,v0,eps,iter_max,is_precond);
01933 case ('gradelt')
01934 res=leq_gradient_elt(m,v,v0,eps,iter_max,is_precond);
01935 case default
01936 #ifdef DEBUG_EXCEPTION
01937 print*, "This method is not exist, your linear equation will be solved gauss-seidel element method"
01938 #endif
01939 res=leq_gaussSeidel_elt(m,v,v0,eps,iter_max,is_precond);
01940 end select
01941 else
01942 res=leq_gaussSeidel_elt(m,v,v0,eps,iter_max,is_precond);
01943 end if
01944
01945 end function leq_solve_iter
01946 end module mod_linear_equation
01947
01948
01949
01950
01956
01962
01966
01970
01971