r8vec_expand_linear2 Subroutine

subroutine r8vec_expand_linear2(n, x, before, fat, after, xfat)

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

! R8VEC_EXPAND_LINEAR2 linearly interpolates new data into an R8VEC.

Discussion:

This routine starts with a vector of data.

The intent is to "fatten" the data, that is, to insert more points
between successive values of the original data.

There will also be extra points placed BEFORE the first original
value and AFTER that last original value.

The "fattened" data is equally spaced between the original points.

The BEFORE data uses the spacing of the first original interval,
and the AFTER data uses the spacing of the last original interval.

Example:

N = 3
BEFORE = 3
FAT = 2
AFTER = 1

X(1:N)        = (/                   0.0,           6.0,             7.0       /)
XFAT(1:2*3+1) = (/ -6.0, -4.0, -2.0, 0.0, 2.0, 4.0, 6.0, 6.33, 6.66, 7.0, 7.66 /)
                   3 "BEFORE's"      Old  2 "FATS"  Old    2 "FATS"  Old  1 "AFTER"

Licensing:

This code is distributed under the GNU LGPL license.

Modified:

03 December 2007

Author:

John Burkardt

Parameters:

Input, integer ( kind = 4 ) N, the number of input data values.
N must be at least 2.

Input, real ( kind = 8 ) X(N), the original data.

Input, integer ( kind = 4 ) BEFORE, the number of "before" values.

Input, integer ( kind = 4 ) FAT, the number of data values to interpolate
between each pair of original data values.

Input, integer ( kind = 4 ) AFTER, the number of "after" values.

Output, real ( kind = 8 ) XFAT(BEFORE+(N-1)*(FAT+1)+1+AFTER), the "fattened" data.

Arguments

Type IntentOptional Attributes Name
integer(kind=4) :: n
real(kind=8) :: x(n)
integer(kind=4) :: before
integer(kind=4) :: fat
integer(kind=4) :: after
real(kind=8) :: xfat(before+(n-1)*(fat+1)+1+after)

Source Code

subroutine r8vec_expand_linear2 ( n, x, before, fat, after, xfat )

  !*****************************************************************************80
  !
  !! R8VEC_EXPAND_LINEAR2 linearly interpolates new data into an R8VEC.
  !
  !  Discussion:
  !
  !    This routine starts with a vector of data.
  !
  !    The intent is to "fatten" the data, that is, to insert more points
  !    between successive values of the original data.
  !
  !    There will also be extra points placed BEFORE the first original
  !    value and AFTER that last original value.
  !
  !    The "fattened" data is equally spaced between the original points.
  !
  !    The BEFORE data uses the spacing of the first original interval,
  !    and the AFTER data uses the spacing of the last original interval.
  !
  !  Example:
  !
  !    N = 3
  !    BEFORE = 3
  !    FAT = 2
  !    AFTER = 1
  !
  !    X(1:N)        = (/                   0.0,           6.0,             7.0       /)
  !    XFAT(1:2*3+1) = (/ -6.0, -4.0, -2.0, 0.0, 2.0, 4.0, 6.0, 6.33, 6.66, 7.0, 7.66 /)
  !                       3 "BEFORE's"      Old  2 "FATS"  Old    2 "FATS"  Old  1 "AFTER"
  !
  !  Licensing:
  !
  !    This code is distributed under the GNU LGPL license. 
  !
  !  Modified:
  !
  !    03 December 2007
  !
  !  Author:
  !
  !    John Burkardt
  !
  !  Parameters:
  !
  !    Input, integer ( kind = 4 ) N, the number of input data values.
  !    N must be at least 2.
  !
  !    Input, real ( kind = 8 ) X(N), the original data.
  !
  !    Input, integer ( kind = 4 ) BEFORE, the number of "before" values.
  !
  !    Input, integer ( kind = 4 ) FAT, the number of data values to interpolate
  !    between each pair of original data values.
  !
  !    Input, integer ( kind = 4 ) AFTER, the number of "after" values.
  !
  !    Output, real ( kind = 8 ) XFAT(BEFORE+(N-1)*(FAT+1)+1+AFTER), the "fattened" data.
  !
  implicit none

  integer ( kind = 4 ) after
  integer ( kind = 4 ) before
  integer ( kind = 4 ) fat
  integer ( kind = 4 ) n

  integer ( kind = 4 ) i
  integer ( kind = 4 ) j
  integer ( kind = 4 ) k
  real    ( kind = 8 ) x(n)
  real    ( kind = 8 ) xfat(before+(n-1)*(fat+1)+1+after)

  k = 0
  !
  !  Points BEFORE.
  !
  do j = 1 - before + fat, fat
     k = k + 1
     xfat(k) = ( real ( fat - j + 1, kind = 8 ) * ( x(1) - ( x(2) - x(1) ) ) &
          + real (       j,     kind = 8 ) *   x(1)          ) &
          / real ( fat     + 1, kind = 8 )
  end do
  !
  !  Original points and FAT points.
  !
  do i = 1, n - 1

     k = k + 1
     xfat(k) = x(i)

     do j = 1, fat
        k = k + 1
        xfat(k) = ( real ( fat - j + 1, kind = 8 ) * x(i)     &
             + real (       j,     kind = 8 ) * x(i+1) ) &
             / real ( fat     + 1, kind = 8 )
     end do

  end do

  k = k + 1
  xfat(k) = x(n)
  !
  !  Points AFTER.
  !
  do j = 1, after
     k = k + 1
     xfat(k) = ( real ( fat - j + 1, kind = 8 ) * x(n)     &
          + real (       j,     kind = 8 ) * ( x(n) + ( x(n) - x(n-1) ) ) ) &
          / real ( fat     + 1, kind = 8 )
  end do

  return
end subroutine r8vec_expand_linear2