Computes all non-zero basis function of all degrees from up to .
N(j,i)
is the value of the th-degree basis function, , where and . It is a simple modification of Algorithm A2.2.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(basis), | intent(in) | :: | me | |||
real(kind=wp), | intent(in) | :: | u | Given value |
||
real(kind=wp), | intent(out) | :: | N(0:me%p,0:me%p) |
|
||
integer, | intent(in), | optional | :: | span | Knot span where lies on |
pure subroutine AllBasisFuns(me,u,N,span)
class(basis), intent(in) :: me
real(wp), intent(in) :: u !< Given \( u \) value
integer, intent(in), optional :: span !< Knot span where \( u \) lies on
real(wp), intent(out) :: N(0:me%p, 0:me%p) !<
!< $$ = \begin{bmatrix} {N_{i - p,0}(u)} & \cdots &{N_{i - p,p}(u)}\\
!< \vdots & \ddots & \vdots \\ {N_{i,0}(u)} & \cdots &{N_{i,p}(u)} \end{bmatrix} $$
!locals
integer :: i !< Knot span index
integer :: j, k, r
real(wp) :: left(me%p), right(me%p), saved, temp
associate( p => me%p, Ui => me%kv ) !get rid of % signs
! set the knot span
if(present(span))then
i = span
else
i = me%FindSpan(u)
end if
N = 0.0_wp; right=0.0_wp; left=0.0_wp
N(0,0) = 1.0_wp
do k = 1, p
do j = 1, k
left(j) = u - Ui(i+1-j)
right(j) = Ui(i+j) - u
saved = 0.0_wp
do r = 0, j-1
temp = N(r,k-1) / (right(r+1) + left(j-r))
N(r,k) = saved + right(r+1) * temp
saved = left(j-r) * temp
end do
N(j,k) = saved
end do
end do
end associate
end subroutine AllBasisFuns