dpchkt.f 2.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596
  1. *DECK DPCHKT
  2. SUBROUTINE DPCHKT (N, X, KNOTYP, T)
  3. C***BEGIN PROLOGUE DPCHKT
  4. C***SUBSIDIARY
  5. C***PURPOSE Compute B-spline knot sequence for DPCHBS.
  6. C***LIBRARY SLATEC (PCHIP)
  7. C***CATEGORY E3
  8. C***TYPE DOUBLE PRECISION (PCHKT-S, DPCHKT-D)
  9. C***AUTHOR Fritsch, F. N., (LLNL)
  10. C***DESCRIPTION
  11. C
  12. C Set a knot sequence for the B-spline representation of a PCH
  13. C function with breakpoints X. All knots will be at least double.
  14. C Endknots are set as:
  15. C (1) quadruple knots at endpoints if KNOTYP=0;
  16. C (2) extrapolate the length of end interval if KNOTYP=1;
  17. C (3) periodic if KNOTYP=2.
  18. C
  19. C Input arguments: N, X, KNOTYP.
  20. C Output arguments: T.
  21. C
  22. C Restrictions/assumptions:
  23. C 1. N.GE.2 . (not checked)
  24. C 2. X(i).LT.X(i+1), i=1,...,N . (not checked)
  25. C 3. 0.LE.KNOTYP.LE.2 . (Acts like KNOTYP=0 for any other value.)
  26. C
  27. C***SEE ALSO DPCHBS
  28. C***ROUTINES CALLED (NONE)
  29. C***REVISION HISTORY (YYMMDD)
  30. C 870701 DATE WRITTEN
  31. C 900405 Converted Fortran to upper case.
  32. C 900410 Converted prologue to SLATEC 4.0 format.
  33. C 900410 Minor cosmetic changes.
  34. C 900430 Produced double precision version.
  35. C 930514 Changed NKNOTS from an output to an input variable. (FNF)
  36. C 930604 Removed unused variable NKNOTS from argument list. (FNF)
  37. C***END PROLOGUE DPCHKT
  38. C
  39. C*Internal Notes:
  40. C
  41. C Since this is subsidiary to DPCHBS, which validates its input before
  42. C calling, it is unnecessary for such validation to be done here.
  43. C
  44. C**End
  45. C
  46. C Declare arguments.
  47. C
  48. INTEGER N, KNOTYP
  49. DOUBLE PRECISION X(*), T(*)
  50. C
  51. C Declare local variables.
  52. C
  53. INTEGER J, K, NDIM
  54. DOUBLE PRECISION HBEG, HEND
  55. C***FIRST EXECUTABLE STATEMENT DPCHKT
  56. C
  57. C Initialize.
  58. C
  59. NDIM = 2*N
  60. C
  61. C Set interior knots.
  62. C
  63. J = 1
  64. DO 20 K = 1, N
  65. J = J + 2
  66. T(J) = X(K)
  67. T(J+1) = T(J)
  68. 20 CONTINUE
  69. C Assertion: At this point T(3),...,T(NDIM+2) have been set and
  70. C J=NDIM+1.
  71. C
  72. C Set end knots according to KNOTYP.
  73. C
  74. HBEG = X(2) - X(1)
  75. HEND = X(N) - X(N-1)
  76. IF (KNOTYP.EQ.1 ) THEN
  77. C Extrapolate.
  78. T(2) = X(1) - HBEG
  79. T(NDIM+3) = X(N) + HEND
  80. ELSE IF ( KNOTYP.EQ.2 ) THEN
  81. C Periodic.
  82. T(2) = X(1) - HEND
  83. T(NDIM+3) = X(N) + HBEG
  84. ELSE
  85. C Quadruple end knots.
  86. T(2) = X(1)
  87. T(NDIM+3) = X(N)
  88. ENDIF
  89. T(1) = T(2)
  90. T(NDIM+4) = T(NDIM+3)
  91. C
  92. C Terminate.
  93. C
  94. RETURN
  95. C------------- LAST LINE OF DPCHKT FOLLOWS -----------------------------
  96. END