JacobiGL Subroutine

subroutine JacobiGL(alpha, beta, n, r)

Compute the n'th order Gauss Lobatto quadrature points, , associated with the Jacobi polynomial of type .

Arguments

Type IntentOptional AttributesName
real(kind=wp), intent(in) :: alpha

The value of .

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

The value of .

integer(kind=ip), intent(in) :: n

The order of the Jacobi polynomial.

real(kind=wp), intent(out), dimension(n+1):: r

On output the location of the Gauss Lobatto quadrature points.


Called by

proc~~jacobigl~~CalledByGraph proc~jacobigl JacobiGL proc~init_ref_element init_ref_element proc~init_ref_element->proc~jacobigl interface~init_ref_element init_ref_element interface~init_ref_element->proc~init_ref_element interface~ref_element ref_element interface~ref_element->interface~init_ref_element proc~scal_schw_init scal_schw_init proc~scal_schw_init->interface~ref_element program~test test program~test->interface~ref_element interface~scal_schw_init scal_schw_init interface~scal_schw_init->proc~scal_schw_init

Contents

Source Code


Source Code

    subroutine JacobiGL(alpha, beta, n, r)
    !! Compute the n'th order Gauss Lobatto quadrature points, \(x\),
    !! associated with the Jacobi polynomial of type
    !! \((\alpha,\beta)>-1\).
      real(wp), intent(in) :: alpha
      !! The value of \(\alpha\).
      real(wp), intent(in) :: beta
      !! The value of \(\beta\).
      integer(ip), intent(in) :: n
      !! The order of the Jacobi polynomial.
      real(wp), dimension(n+1), intent(out) :: r
      !! On output the location of the Gauss Lobatto quadrature points.
      real(wp), dimension(n-1) :: w

      if ( n<=0 ) then
        print*,'JacobiGL called with n<=0. Aborting'
        stop
      end if
      r = 0.0_wp
      if ( n==1 ) then
        r(1) = -1.0_wp
        r(2) = 1.0_wp
        return
      end if

      call JacobiGQ ( alpha+1.0_wp, beta+1.0_wp, n-2, r(2:n), w )
      r(1) = -1.0_wp
      r(n+1) = 1.0_wp
      return
    end subroutine JacobiGL