rpzero Subroutine

public subroutine rpzero(n, a, r, iflg, s)

Find the zeros of a polynomial with real coefficients: P(X)= A(1)*X**N + A(2)*X**(N-1) +...+ A(N+1)

REVISION HISTORY (YYMMDD)

  • 810223 DATE WRITTEN. Kahaner, D. K., (NBS)
  • 890206 REVISION DATE from Version 3.2
  • 891214 Prologue converted to Version 4.0 format. (BAB)
  • Jacob Williams, 9/13/2022 : modernized this routine

Note

This is just a wrapper to cpzero

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: n

degree of P(X)

real(kind=wp), intent(in), dimension(n+1) :: a

real vector containing coefficients of P(X), A(I) = coefficient of X**(N+1-I)

complex(kind=wp), intent(inout), dimension(n) :: r

N word complex vector. On Input: containing initial estimates for zeros if these are known. On output: ith zero.

integer, intent(inout) :: iflg

On Input:

flag to indicate if initial estimates of zeros are input:

  • If IFLG == 0, no estimates are input.
  • If IFLG /= 0, the vector R contains estimates of the zeros

** WARNING ** If estimates are input, they must be separated, that is, distinct or not repeated.

On Output:

Error Diagnostics:

  • If IFLG == 0 on return, all is well
  • If IFLG == 1 on return, A(1)=0.0 or N=0 on input
  • If IFLG == 2 on return, the program failed to converge after 25*N iterations. Best current estimates of the zeros are in R(I). Error bounds are not calculated.
real(kind=wp), intent(out), dimension(n) :: s

an N word array. bound for R(I).


Calls

proc~~rpzero~~CallsGraph proc~rpzero polyroots_module::rpzero proc~cpzero polyroots_module::cpzero proc~rpzero->proc~cpzero proc~cpevl polyroots_module::cpevl proc~cpzero->proc~cpevl

Source Code

subroutine rpzero(n, a, r, iflg, s)

    implicit none

    integer, intent(in) :: n !! degree of `P(X)`
    real(wp), dimension(n+1), intent(in) :: a !! real vector containing coefficients of `P(X)`,
                                              !! `A(I)` = coefficient of `X**(N+1-I)`
    complex(wp), dimension(n), intent(inout) :: r !! `N` word complex vector. On Input: containing initial estimates for zeros
                                                  !! if these are known. On output: ith zero.
    integer, intent(inout) :: iflg !!### On Input:
                                   !!
                                   !! flag to indicate if initial estimates of zeros are input:
                                   !!
                                   !!  * If `IFLG == 0`, no estimates are input.
                                   !!  * If `IFLG /= 0`, the vector R contains estimates of the zeros
                                   !!
                                   !! ** WARNING ****** If estimates are input, they must
                                   !!                   be separated, that is, distinct or
                                   !!                   not repeated.
                                   !!### On Output:
                                   !!
                                   !! Error Diagnostics:
                                   !!
                                   !! * If `IFLG == 0` on return, all is well
                                   !! * If `IFLG == 1` on return, `A(1)=0.0` or `N=0` on input
                                   !! * If `IFLG == 2` on return, the program failed to converge
                                   !!   after `25*N` iterations.  Best current estimates of the
                                   !!   zeros are in `R(I)`.  Error bounds are not calculated.
    real(wp), dimension(n), intent(out) :: s !! an `N` word array. bound for `R(I)`.

    integer :: i
    complex(wp), dimension(:), allocatable :: p !! complex coefficients

    allocate(p(n+1))

    do i = 1, n + 1
        p(i) = cmplx(a(i), 0.0_wp, wp)
    end do
    call cpzero(n, p, r, iflg, s)

end subroutine rpzero