rpqr79.f 2.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103
  1. *DECK RPQR79
  2. SUBROUTINE RPQR79 (NDEG, COEFF, ROOT, IERR, WORK)
  3. C***BEGIN PROLOGUE RPQR79
  4. C***PURPOSE Find the zeros of a polynomial with real coefficients.
  5. C***LIBRARY SLATEC
  6. C***CATEGORY F1A1A
  7. C***TYPE SINGLE PRECISION (RPQR79-S, CPQR79-C)
  8. C***KEYWORDS COMPLEX POLYNOMIAL, POLYNOMIAL ROOTS, POLYNOMIAL ZEROS
  9. C***AUTHOR Vandevender, W. H., (SNLA)
  10. C***DESCRIPTION
  11. C
  12. C Abstract
  13. C This routine computes all zeros of a polynomial of degree NDEG
  14. C with real coefficients by computing the eigenvalues of the
  15. C companion matrix.
  16. C
  17. C Description of Parameters
  18. C The user must dimension all arrays appearing in the call list
  19. C COEFF(NDEG+1), ROOT(NDEG), WORK(NDEG*(NDEG+2))
  20. C
  21. C --Input--
  22. C NDEG degree of polynomial
  23. C
  24. C COEFF REAL coefficients in descending order. i.e.,
  25. C P(Z)= COEFF(1)*(Z**NDEG) + COEFF(NDEG)*Z + COEFF(NDEG+1)
  26. C
  27. C WORK REAL work array of dimension at least NDEG*(NDEG+2)
  28. C
  29. C --Output--
  30. C ROOT COMPLEX vector of roots
  31. C
  32. C IERR Output Error Code
  33. C - Normal Code
  34. C 0 means the roots were computed.
  35. C - Abnormal Codes
  36. C 1 more than 30 QR iterations on some eigenvalue of the
  37. C companion matrix
  38. C 2 COEFF(1)=0.0
  39. C 3 NDEG is invalid (less than or equal to 0)
  40. C
  41. C***REFERENCES (NONE)
  42. C***ROUTINES CALLED HQR, XERMSG
  43. C***REVISION HISTORY (YYMMDD)
  44. C 800601 DATE WRITTEN
  45. C 890505 REVISION DATE from Version 3.2
  46. C 891214 Prologue converted to Version 4.0 format. (BAB)
  47. C 900315 CALLs to XERROR changed to CALLs to XERMSG. (THJ)
  48. C 911010 Code reworked and simplified. (RWC and WRB)
  49. C***END PROLOGUE RPQR79
  50. REAL COEFF(*), WORK(*), SCALE
  51. COMPLEX ROOT(*)
  52. INTEGER NDEG, IERR, K, KH, KWR, KWI, KCOL
  53. C***FIRST EXECUTABLE STATEMENT RPQR79
  54. IERR = 0
  55. IF (ABS(COEFF(1)) .EQ. 0.0) THEN
  56. IERR = 2
  57. CALL XERMSG ('SLATEC', 'RPQR79',
  58. + 'LEADING COEFFICIENT IS ZERO.', 2, 1)
  59. RETURN
  60. ENDIF
  61. C
  62. IF (NDEG .LE. 0) THEN
  63. IERR = 3
  64. CALL XERMSG ('SLATEC', 'RPQR79', 'DEGREE INVALID.', 3, 1)
  65. RETURN
  66. ENDIF
  67. C
  68. IF (NDEG .EQ. 1) THEN
  69. ROOT(1) = CMPLX(-COEFF(2)/COEFF(1),0.0)
  70. RETURN
  71. ENDIF
  72. C
  73. SCALE = 1.0E0/COEFF(1)
  74. KH = 1
  75. KWR = KH+NDEG*NDEG
  76. KWI = KWR+NDEG
  77. KWEND = KWI+NDEG-1
  78. C
  79. DO 10 K=1,KWEND
  80. WORK(K) = 0.0E0
  81. 10 CONTINUE
  82. C
  83. DO 20 K=1,NDEG
  84. KCOL = (K-1)*NDEG+1
  85. WORK(KCOL) = -COEFF(K+1)*SCALE
  86. IF (K .NE. NDEG) WORK(KCOL+K) = 1.0E0
  87. 20 CONTINUE
  88. C
  89. CALL HQR (NDEG,NDEG,1,NDEG,WORK(KH),WORK(KWR),WORK(KWI),IERR)
  90. C
  91. IF (IERR .NE. 0) THEN
  92. IERR = 1
  93. CALL XERMSG ('SLATEC', 'CPQR79',
  94. + 'NO CONVERGENCE IN 30 QR ITERATIONS.', 1, 1)
  95. RETURN
  96. ENDIF
  97. C
  98. DO 30 K=1,NDEG
  99. KM1 = K-1
  100. ROOT(K) = CMPLX(WORK(KWR+KM1),WORK(KWI+KM1))
  101. 30 CONTINUE
  102. RETURN
  103. END