f1_abscissas_ab Subroutine

subroutine f1_abscissas_ab(a, b, n, x)

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

! F1_ABSCISSAS_AB computes Fejer type 1 abscissas for the interval [A,B].

Licensing:

This code is distributed under the GNU LGPL license.

Modified:

29 December 2007

Author:

John Burkardt

Reference:

Philip Davis, Philip Rabinowitz,
Methods of Numerical Integration,
Second Edition,
Dover, 2007,
ISBN: 0486453391,
LC: QA299.3.D28.

Walter Gautschi,
Numerical Quadrature in the Presence of a Singularity,
SIAM Journal on Numerical Analysis,
Volume 4, Number 3, 1967, pages 357-362.

Joerg Waldvogel,
Fast Construction of the Fejer and Clenshaw-Curtis Quadrature Rules,
BIT Numerical Mathematics,
Volume 43, Number 1, 2003, pages 1-18.

Parameters:

Input, real ( kind = 8 ) A, B, the endpoints of the interval.

Input, integer ( kind = 4 ) N, the order of the rule.

Output, real ( kind = 8 ) X(N), the abscissas.

Arguments

Type IntentOptional Attributes Name
real(kind=8) :: a
real(kind=8) :: b
integer(kind=4) :: n
real(kind=8) :: x(n)

Source Code

subroutine f1_abscissas_ab ( a, b, n, x )

  !*****************************************************************************80
  !
  !! F1_ABSCISSAS_AB computes Fejer type 1 abscissas for the interval [A,B].
  !
  !  Licensing:
  !
  !    This code is distributed under the GNU LGPL license. 
  !
  !  Modified:
  !
  !    29 December 2007
  !
  !  Author:
  !
  !    John Burkardt
  !
  !  Reference:
  !
  !    Philip Davis, Philip Rabinowitz,
  !    Methods of Numerical Integration,
  !    Second Edition,
  !    Dover, 2007,
  !    ISBN: 0486453391,
  !    LC: QA299.3.D28.
  !
  !    Walter Gautschi,
  !    Numerical Quadrature in the Presence of a Singularity,
  !    SIAM Journal on Numerical Analysis,
  !    Volume 4, Number 3, 1967, pages 357-362.
  !
  !    Joerg Waldvogel,
  !    Fast Construction of the Fejer and Clenshaw-Curtis Quadrature Rules,
  !    BIT Numerical Mathematics,
  !    Volume 43, Number 1, 2003, pages 1-18.
  !
  !  Parameters:
  !
  !    Input, real ( kind = 8 ) A, B, the endpoints of the interval.
  !
  !    Input, integer ( kind = 4 ) N, the order of the rule.
  !
  !    Output, real ( kind = 8 ) X(N), the abscissas.
  !
  implicit none

  integer ( kind = 4 ) n

  real    ( kind = 8 ) a
  real    ( kind = 8 ) b
  integer ( kind = 4 ) i
  real    ( kind = 8 ) :: pi = 3.141592653589793D+00
  real    ( kind = 8 ) theta(n)
  real    ( kind = 8 ) x(n)

  if ( n < 1 ) then
     write ( *, '(a)' ) ' '
     write ( *, '(a)' ) 'F1_ABSCISSAS_AB - Fatal error!'
     write ( *, '(a)' ) '  N < 1.'
     stop
  end if

  if ( n == 1 ) then
     x(1) = 0.5D+00 * ( b + a )
     return
  end if

  do i = 1, n
     theta(i) = real ( 2 * n - 2 * i + 1, kind = 8 ) * pi &
          / real ( 2 * n,             kind = 8 )
  end do

  x(1:n) = 0.5D+00 * ( ( b + a ) + ( b - a ) * cos ( theta(1:n) ) )

  return
end subroutine f1_abscissas_ab