invert_symmetric_matrix_dp Module Subroutine

module subroutine invert_symmetric_matrix_dp(matrix_)

invert a symmetric matrix using DSYTRI method (double precision version)

Arguments

Type IntentOptional Attributes Name
real(kind=dp), intent(inout) :: matrix_(:,:)

Calls

proc~~invert_symmetric_matrix_dp~~CallsGraph proc~invert_symmetric_matrix_dp invert_symmetric_matrix_dp dsytrf dsytrf proc~invert_symmetric_matrix_dp->dsytrf proc~integer_to_character integer_to_character proc~invert_symmetric_matrix_dp->proc~integer_to_character proc~write_error write_error proc~invert_symmetric_matrix_dp->proc~write_error proc~write_message write_message proc~invert_symmetric_matrix_dp->proc~write_message dsytri dsytri proc~invert_symmetric_matrix_dp->dsytri proc~write_error->proc~write_message

Contents


Variables

Type Visibility Attributes Name Initial
integer(kind=int32), private :: ILAENV
integer(kind=int32), private :: i_size_
integer(kind=int32), private :: i_size_2_
integer(kind=int32), private :: lwork_
integer(kind=int32), private :: nb_
integer(kind=int32), private :: ok_
integer(kind=int32), private, allocatable :: pivot_(:)
integer(kind=int32), private :: size_
integer(kind=int32), private :: size_1_
integer(kind=int32), private :: size_2_
real(kind=dp), private, allocatable :: work_(:)

Source Code

      module subroutine invert_symmetric_matrix_dp(matrix_)
         !! invert a symmetric matrix using DSYTRI method (double precision version)
         real(dp), intent(inout) :: matrix_(:,:)
         !---------------------------------------------------------------------!
         integer(int32) :: size_1_, size_2_, size_, lwork_, nb_, ok_, ILAENV, i_size_, i_size_2_
         integer(int32), allocatable :: pivot_(:)
         real(dp), allocatable :: work_(:)
         !---------------------------------------------------------------------!
         size_1_ = size(matrix_, dim = 1)
         size_2_ = size(matrix_, dim = 2)
         if (size_1_ .eq. size_2_) then
            size_ = size_1_
         else
            call write_message("Error in invert_symmetric_matrix_dp: size "    &
               // "in dim = 1 ("//trim(adjustl(integer_to_character(size_1_))) &
               // ") is different than in dim = 2 (" //                        &
               trim(adjustl(integer_to_character(size_2_))) // ")")
            call write_error("Adapt this subroutine to rectangle matrices")
         endif
         !---------------------------------------------------------------------!
         if (allocated(pivot_)) deallocate(pivot_)
         allocate(pivot_(size_))
         !---------------------------------------------------------------------!
         nb_ = ILAENV(1,'DSYTRF','L',size_,-1,-1,-1)
         lwork_ = nb_ * size_
         if (allocated(work_)) deallocate(work_)
         allocate(work_(lwork_))
         work_ = 0
         !---------------------------------------------------------------------!
         call DSYTRF('L',size_,matrix_,size_,pivot_,work_,lwork_,ok_)
         if (ok_ /= 0) then
            call write_error("DSYTRF failed with status: " //                  &
               trim(adjustl(integer_to_character(ok_))) )
         endif
         !---------------------------------------------------------------------!
         call DSYTRI('L',size_,matrix_,size_,pivot_,work_,ok_)
         if (ok_ /= 0) then
            call write_error("DSYTRI failed with status: " //                  &
               trim(adjustl(integer_to_character(ok_))) )
         endif
         !---------------------------------------------------------------------!
      end subroutine invert_symmetric_matrix_dp