s_spline_matrix_banded__init Subroutine

private subroutine s_spline_matrix_banded__init(self, n, kl, ku)

Type Bound

sll_t_spline_matrix_banded

Arguments

Type IntentOptional Attributes Name
class(sll_t_spline_matrix_banded), intent(out) :: self
integer, intent(in) :: n
integer, intent(in) :: kl
integer, intent(in) :: ku

Calls

proc~~s_spline_matrix_banded__init~~CallsGraph proc~s_spline_matrix_banded__init sll_t_spline_matrix_banded%s_spline_matrix_banded__init sll_assert sll_assert proc~s_spline_matrix_banded__init->sll_assert

Source Code

  subroutine s_spline_matrix_banded__init( self, n, kl, ku )
    class(sll_t_spline_matrix_banded), intent(  out) :: self
    integer                          , intent(in   ) :: n
    integer                          , intent(in   ) :: kl
    integer                          , intent(in   ) :: ku

    SLL_ASSERT( n  >  0 )
    SLL_ASSERT( kl >= 0 )
    SLL_ASSERT( ku >= 0 )
    SLL_ASSERT( kl < n  )
    SLL_ASSERT( ku < n  )

    self%n  = n
    self%kl = kl
    self%ku = ku

    ! Given the linear system A*x=b, we assume that A is a square (n by n)
    ! with ku super-diagonals and kl sub-diagonals.
    ! All non-zero elements of A are stored in the rectangular matrix q, using
    ! the format required by DGBTRF (LAPACK): diagonals of A are rows of q.
    ! q has 2*kl rows for the subdiagonals, 1 row for the diagonal, and ku rows
    ! for the superdiagonals. (The kl additional rows are needed for pivoting.)
    ! The term A(i,j) of the full matrix is stored in q(i-j+2*kl+1,j).

    allocate( self%ipiv(n) )
    allocate( self%q(2*kl+ku+1,n) )
    self%q(:,:) = 0.0_wp
    self%ipiv(:) = 0
    self%factorized=.FALSE.

  end subroutine s_spline_matrix_banded__init