s_bsplines_non_uniform__eval_basis Subroutine

private pure subroutine s_bsplines_non_uniform__eval_basis(self, x, values, jmin)

Evaluate value at x of all basis functions with support in local cell values[j] = B_j(x) for jmin <= j <= jmin+degree

Type Bound

sll_t_bsplines_non_uniform

Arguments

Type IntentOptional Attributes Name
class(sll_t_bsplines_non_uniform), intent(in) :: self

non-uniform B-splines

real(kind=wp), intent(in) :: x

evaluation point

real(kind=wp), intent(out) :: values(0:)
integer, intent(out) :: jmin

index of first non-zero B-spline


Calls

proc~~s_bsplines_non_uniform__eval_basis~~CallsGraph proc~s_bsplines_non_uniform__eval_basis sll_t_bsplines_non_uniform%s_bsplines_non_uniform__eval_basis proc~f_bsplines_non_uniform__find_cell sll_t_bsplines_non_uniform%f_bsplines_non_uniform__find_cell proc~s_bsplines_non_uniform__eval_basis->proc~f_bsplines_non_uniform__find_cell sll_assert sll_assert proc~s_bsplines_non_uniform__eval_basis->sll_assert

Source Code

  SLL_PURE subroutine s_bsplines_non_uniform__eval_basis( self, x, values, jmin )
    class(sll_t_bsplines_non_uniform), intent(in   ) :: self
    real(wp)                         , intent(in   ) :: x
    real(wp)                         , intent(  out) :: values(0:)
    integer                          , intent(  out) :: jmin

    integer  :: icell

    integer  :: j, r
    real(wp) :: temp
    real(wp) :: saved

    ! GFortran: to allocate on stack use -fstack-arrays
    real(wp) :: left (1:self%degree)
    real(wp) :: right(1:self%degree)

    ! Check on inputs
    SLL_ASSERT( x > self%xmin - 1.0d-14 )
    SLL_ASSERT( x < self%xmax + 1.0d-14 )
    SLL_ASSERT( size(values) == 1+self%degree )

    ! 1. Compute cell index 'icell'
!    icell = f_bsplines_non_uniform__find_cell( self, x )
    icell = self % find_cell( x )
    SLL_ASSERT( icell >= 1               )
    SLL_ASSERT( icell <= self%ncells     )
    SLL_ASSERT( self%knots(icell) <= x   )
    SLL_ASSERT( x <= self%knots(icell+1) )

    ! 2. Compute index range of B-splines with support over cell 'icell'
    jmin = icell - self%offset

    ! 3. Compute values of aforementioned B-splines
    values(0) = 1.0_wp
    do j = 1, self%degree
       left (j) = x - self%knots(icell+1-j)
       right(j) = self%knots(icell+j) - x
       saved    = 0.0_wp
       do r = 0, j-1
          temp      = values(r) / (right(r+1) + left(j-r))
          values(r) = saved + right(r+1) * temp
          saved     = left(j-r) * temp
       end do
      values(j) = saved
    end do

  end subroutine s_bsplines_non_uniform__eval_basis