subfor Subroutine

subroutine subfor(w, ipivot, nrow, last, b, x)

************80

! SUBFOR carries out the forward pass of substitution for the current block.

Discussion:

The forward pass is the action on the right hand side corresponding to the
elimination carried out in FACTRB for this block.

At the end, X(1:NROW) contains the right hand side of the transformed
IPIVOT(1:NROW)-th equation in this block.

Then, since for I=1,...,NROW-LAST, B(NROW+I) is going to be used as
the right hand side of equation I in the next block (shifted over there
from this block during factorization), it is set equal to X(LAST+I) here.

Modified:

14 February 2007

Author:

Carl DeBoor

Reference:

Carl DeBoor,
A Practical Guide to Splines,
Springer, 2001,
ISBN: 0387953663,
LC: QA1.A647.v27.

Parameters:

Input, real ( kind = 8 ) W(NROW,LAST), integer IPIVOT(NROW),
integer ( kind = 4 ) NROW, integer LAST, are as on return from FACTRB.

Output, real ( kind = 8 ) B(2*NROW-LAST).  On input, B(1:NROW)
contains the right hand sides for this block.  On output,
B(NROW+1:2*NROW-LAST) contains the appropriately modified right
hand sides for the next block.

Output, real X(NROW), contains, on output, the appropriately modified
right hand sides of equations IPIVOT(1:NROW).

Arguments

Type IntentOptional Attributes Name
real(kind=8) :: w(nrow,last)
integer(kind=4) :: ipivot(nrow)
integer(kind=4) :: nrow
integer(kind=4) :: last
real(kind=8) :: b(nrow+nrow-last)
real(kind=8) :: x(nrow)

Source Code

subroutine subfor ( w, ipivot, nrow, last, b, x )

  !*****************************************************************************80
  !
  !! SUBFOR carries out the forward pass of substitution for the current block.
  !
  !  Discussion:
  !
  !    The forward pass is the action on the right hand side corresponding to the 
  !    elimination carried out in FACTRB for this block.
  !
  !    At the end, X(1:NROW) contains the right hand side of the transformed
  !    IPIVOT(1:NROW)-th equation in this block.  
  !
  !    Then, since for I=1,...,NROW-LAST, B(NROW+I) is going to be used as 
  !    the right hand side of equation I in the next block (shifted over there 
  !    from this block during factorization), it is set equal to X(LAST+I) here.
  !
  !  Modified:
  !
  !    14 February 2007
  !
  !  Author:
  !
  !    Carl DeBoor
  !
  !  Reference:
  !
  !    Carl DeBoor,
  !    A Practical Guide to Splines,
  !    Springer, 2001,
  !    ISBN: 0387953663,
  !    LC: QA1.A647.v27.
  !
  !  Parameters:
  !
  !    Input, real ( kind = 8 ) W(NROW,LAST), integer IPIVOT(NROW), 
  !    integer ( kind = 4 ) NROW, integer LAST, are as on return from FACTRB.
  !
  !    Output, real ( kind = 8 ) B(2*NROW-LAST).  On input, B(1:NROW)
  !    contains the right hand sides for this block.  On output,
  !    B(NROW+1:2*NROW-LAST) contains the appropriately modified right
  !    hand sides for the next block.
  !
  !    Output, real X(NROW), contains, on output, the appropriately modified
  !    right hand sides of equations IPIVOT(1:NROW).
  !
  implicit none

  integer ( kind = 4 ) last
  integer ( kind = 4 ) nrow

  real ( kind = 8 ) b(nrow+nrow-last)
  integer ( kind = 4 ) ip
  integer ( kind = 4 ) ipivot(nrow)
  integer ( kind = 4 ) jhi
  integer ( kind = 4 ) k
  real ( kind = 8 ) w(nrow,last)
  real ( kind = 8 ) x(nrow)

  ip = ipivot(1)
  x(1) = b(ip)

  do k = 2, nrow

     ip = ipivot(k)

     jhi = min ( k - 1, last )

     x(k) = b(ip) - dot_product ( w(ip,1:jhi), x(1:jhi) )

  end do
  !
  !  Transfer modified right hand sides of equations IPIVOT(LAST+1:NROW) 
  !  to next block.
  !
  b(nrow+1:2*nrow-last) = x(last+1:nrow)

  return
end subroutine subfor