s_spline_matrix_dense__solve_inplace Subroutine

private subroutine s_spline_matrix_dense__solve_inplace(self, nrhs, bx)

Type Bound

sll_t_spline_matrix_dense

Arguments

Type IntentOptional Attributes Name
class(sll_t_spline_matrix_dense), intent(in) :: self
integer, intent(in) :: nrhs
real(kind=wp), intent(inout), dimension(*) :: bx

Calls

proc~~s_spline_matrix_dense__solve_inplace~~CallsGraph proc~s_spline_matrix_dense__solve_inplace sll_t_spline_matrix_dense%s_spline_matrix_dense__solve_inplace dgetrs dgetrs proc~s_spline_matrix_dense__solve_inplace->dgetrs proc~sll_s_error_handler sll_s_error_handler proc~s_spline_matrix_dense__solve_inplace->proc~sll_s_error_handler interface~c_abort~2 c_abort proc~sll_s_error_handler->interface~c_abort~2 proc~errout errout proc~sll_s_error_handler->proc~errout

Source Code

  subroutine s_spline_matrix_dense__solve_inplace( self, nrhs,bx )
    class(sll_t_spline_matrix_dense), intent(in   ) :: self
    integer                         , intent(in   ) :: nrhs
    real(wp),dimension(*)           , intent(inout) :: bx

    integer :: info

    character(len=*), parameter :: &
         this_sub_name = "sll_t_spline_matrix_dense % solve_inplace"
    character(len=256) :: err_msg

    SLL_ASSERT( size(self%a,1) == self%n )
    SLL_ASSERT( size(self%a,2) == self%n )
!    SLL_ASSERT( size(bx)  == self%n*nrhs )
    SLL_ASSERT( self%factorized )

    ! Solve linear system PLU*x=b using Lapack
    call dgetrs( 'N', self%n, nrhs, self%a, self%n, self%ipiv, bx, self%n, info )

    if ( info < 0 ) then
      write( err_msg, '(i0,a)' ) abs(info), "-th argument had an illegal value"
      SLL_ERROR(this_sub_name,err_msg)
    end if

    end subroutine s_spline_matrix_dense__solve_inplace