defcmn.f 6.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236
  1. *DECK DEFCMN
  2. SUBROUTINE DEFCMN (NDATA, XDATA, YDATA, SDDATA, NORD, NBKPT,
  3. + BKPTIN, MDEIN, MDEOUT, COEFF, BF, XTEMP, PTEMP, BKPT, G, MDG,
  4. + W, MDW, LW)
  5. C***BEGIN PROLOGUE DEFCMN
  6. C***SUBSIDIARY
  7. C***PURPOSE Subsidiary to DEFC
  8. C***LIBRARY SLATEC
  9. C***TYPE DOUBLE PRECISION (EFCMN-S, DEFCMN-D)
  10. C***AUTHOR Hanson, R. J., (SNLA)
  11. C***DESCRIPTION
  12. C
  13. C This is a companion subprogram to DEFC( ).
  14. C This subprogram does weighted least squares fitting of data by
  15. C B-spline curves.
  16. C The documentation for DEFC( ) has complete usage instructions.
  17. C
  18. C***SEE ALSO DEFC
  19. C***ROUTINES CALLED DBNDAC, DBNDSL, DCOPY, DFSPVN, DSCAL, DSORT, XERMSG
  20. C***REVISION HISTORY (YYMMDD)
  21. C 800801 DATE WRITTEN
  22. C 890531 Changed all specific intrinsics to generic. (WRB)
  23. C 890618 Completely restructured and extensively revised (WRB & RWC)
  24. C 891214 Prologue converted to Version 4.0 format. (BAB)
  25. C 900315 CALLs to XERROR changed to CALLs to XERMSG. (THJ)
  26. C 900328 Added TYPE section. (WRB)
  27. C 900510 Convert XERRWV calls to XERMSG calls. (RWC)
  28. C 900604 DP version created from SP version. (RWC)
  29. C***END PROLOGUE DEFCMN
  30. INTEGER LW, MDEIN, MDEOUT, MDG, MDW, NBKPT, NDATA, NORD
  31. DOUBLE PRECISION BF(NORD,*), BKPT(*), BKPTIN(*), COEFF(*),
  32. * G(MDG,*), PTEMP(*), SDDATA(*), W(MDW,*), XDATA(*), XTEMP(*),
  33. * YDATA(*)
  34. C
  35. EXTERNAL DBNDAC, DBNDSL, DCOPY, DFSPVN, DSCAL, DSORT, XERMSG
  36. C
  37. DOUBLE PRECISION DUMMY, RNORM, XMAX, XMIN, XVAL
  38. INTEGER I, IDATA, ILEFT, INTSEQ, IP, IR, IROW, L, MT, N, NB,
  39. * NORDM1, NORDP1, NP1
  40. CHARACTER*8 XERN1, XERN2
  41. C
  42. C***FIRST EXECUTABLE STATEMENT DEFCMN
  43. C
  44. C Initialize variables and analyze input.
  45. C
  46. N = NBKPT - NORD
  47. NP1 = N + 1
  48. C
  49. C Initially set all output coefficients to zero.
  50. C
  51. CALL DCOPY (N, 0.D0, 0, COEFF, 1)
  52. MDEOUT = -1
  53. IF (NORD.LT.1 .OR. NORD.GT.20) THEN
  54. CALL XERMSG ('SLATEC', 'DEFCMN',
  55. + 'IN DEFC, THE ORDER OF THE B-SPLINE MUST BE 1 THRU 20.',
  56. + 3, 1)
  57. RETURN
  58. ENDIF
  59. C
  60. IF (NBKPT.LT.2*NORD) THEN
  61. CALL XERMSG ('SLATEC', 'DEFCMN',
  62. + 'IN DEFC, THE NUMBER OF KNOTS MUST BE AT LEAST TWICE ' //
  63. + 'THE B-SPLINE ORDER.', 4, 1)
  64. RETURN
  65. ENDIF
  66. C
  67. IF (NDATA.LT.0) THEN
  68. CALL XERMSG ('SLATEC', 'DEFCMN',
  69. + 'IN DEFC, THE NUMBER OF DATA POINTS MUST BE NONNEGATIVE.',
  70. + 5, 1)
  71. RETURN
  72. ENDIF
  73. C
  74. NB = (NBKPT-NORD+3)*(NORD+1) + (NBKPT+1)*(NORD+1) +
  75. + 2*MAX(NBKPT,NDATA) + NBKPT + NORD**2
  76. IF (LW .LT. NB) THEN
  77. WRITE (XERN1, '(I8)') NB
  78. WRITE (XERN2, '(I8)') LW
  79. CALL XERMSG ('SLATEC', 'DEFCMN',
  80. * 'IN DEFC, INSUFFICIENT STORAGE FOR W(*). CHECK FORMULA ' //
  81. * 'THAT READS LW.GE. ... . NEED = ' // XERN1 //
  82. * ' GIVEN = ' // XERN2, 6, 1)
  83. MDEOUT = -1
  84. RETURN
  85. ENDIF
  86. C
  87. IF (MDEIN.NE.1 .AND. MDEIN.NE.2) THEN
  88. CALL XERMSG ('SLATEC', 'DEFCMN',
  89. + 'IN DEFC, INPUT VALUE OF MDEIN MUST BE 1-2.', 7, 1)
  90. RETURN
  91. ENDIF
  92. C
  93. C Sort the breakpoints.
  94. C
  95. CALL DCOPY (NBKPT, BKPTIN, 1, BKPT, 1)
  96. CALL DSORT (BKPT, DUMMY, NBKPT, 1)
  97. C
  98. C Save interval containing knots.
  99. C
  100. XMIN = BKPT(NORD)
  101. XMAX = BKPT(NP1)
  102. NORDM1 = NORD - 1
  103. NORDP1 = NORD + 1
  104. C
  105. C Process least squares equations.
  106. C
  107. C Sort data and an array of pointers.
  108. C
  109. CALL DCOPY (NDATA, XDATA, 1, XTEMP, 1)
  110. DO 100 I = 1,NDATA
  111. PTEMP(I) = I
  112. 100 CONTINUE
  113. C
  114. IF (NDATA.GT.0) THEN
  115. CALL DSORT (XTEMP, PTEMP, NDATA, 2)
  116. XMIN = MIN(XMIN,XTEMP(1))
  117. XMAX = MAX(XMAX,XTEMP(NDATA))
  118. ENDIF
  119. C
  120. C Fix breakpoint array if needed. This should only involve very
  121. C minor differences with the input array of breakpoints.
  122. C
  123. DO 110 I = 1,NORD
  124. BKPT(I) = MIN(BKPT(I),XMIN)
  125. 110 CONTINUE
  126. C
  127. DO 120 I = NP1,NBKPT
  128. BKPT(I) = MAX(BKPT(I),XMAX)
  129. 120 CONTINUE
  130. C
  131. C Initialize parameters of banded matrix processor, DBNDAC( ).
  132. C
  133. MT = 0
  134. IP = 1
  135. IR = 1
  136. ILEFT = NORD
  137. INTSEQ = 1
  138. DO 150 IDATA = 1,NDATA
  139. C
  140. C Sorted indices are in PTEMP(*).
  141. C
  142. L = PTEMP(IDATA)
  143. XVAL = XDATA(L)
  144. C
  145. C When interval changes, process equations in the last block.
  146. C
  147. IF (XVAL.GE.BKPT(ILEFT+1)) THEN
  148. CALL DBNDAC (G, MDG, NORD, IP, IR, MT, ILEFT-NORDM1)
  149. MT = 0
  150. C
  151. C Move pointer up to have BKPT(ILEFT).LE.XVAL, ILEFT.LE.N.
  152. C
  153. DO 130 ILEFT = ILEFT,N
  154. IF (XVAL.LT.BKPT(ILEFT+1)) GO TO 140
  155. IF (MDEIN.EQ.2) THEN
  156. C
  157. C Data is being sequentially accumulated.
  158. C Transfer previously accumulated rows from W(*,*) to
  159. C G(*,*) and process them.
  160. C
  161. CALL DCOPY (NORDP1, W(INTSEQ,1), MDW, G(IR,1), MDG)
  162. CALL DBNDAC (G, MDG, NORD, IP, IR, 1, INTSEQ)
  163. INTSEQ = INTSEQ + 1
  164. ENDIF
  165. 130 CONTINUE
  166. ENDIF
  167. C
  168. C Obtain B-spline function value.
  169. C
  170. 140 CALL DFSPVN (BKPT, NORD, 1, XVAL, ILEFT, BF)
  171. C
  172. C Move row into place.
  173. C
  174. IROW = IR + MT
  175. MT = MT + 1
  176. CALL DCOPY (NORD, BF, 1, G(IROW,1), MDG)
  177. G(IROW,NORDP1) = YDATA(L)
  178. C
  179. C Scale data if uncertainty is nonzero.
  180. C
  181. IF (SDDATA(L).NE.0.D0) CALL DSCAL (NORDP1, 1.D0/SDDATA(L),
  182. + G(IROW,1), MDG)
  183. C
  184. C When staging work area is exhausted, process rows.
  185. C
  186. IF (IROW.EQ.MDG-1) THEN
  187. CALL DBNDAC (G, MDG, NORD, IP, IR, MT, ILEFT-NORDM1)
  188. MT = 0
  189. ENDIF
  190. 150 CONTINUE
  191. C
  192. C Process last block of equations.
  193. C
  194. CALL DBNDAC (G, MDG, NORD, IP, IR, MT, ILEFT-NORDM1)
  195. C
  196. C Finish processing any previously accumulated rows from W(*,*)
  197. C to G(*,*).
  198. C
  199. IF (MDEIN.EQ.2) THEN
  200. DO 160 I = INTSEQ,NP1
  201. CALL DCOPY (NORDP1, W(I,1), MDW, G(IR,1), MDG)
  202. CALL DBNDAC (G, MDG, NORD, IP, IR, 1, MIN(N,I))
  203. 160 CONTINUE
  204. ENDIF
  205. C
  206. C Last call to adjust block positioning.
  207. C
  208. CALL DCOPY (NORDP1, 0.D0, 0, G(IR,1), MDG)
  209. CALL DBNDAC (G, MDG, NORD, IP, IR, 1, NP1)
  210. C
  211. C Transfer accumulated rows from G(*,*) to W(*,*) for
  212. C possible later sequential accumulation.
  213. C
  214. DO 170 I = 1,NP1
  215. CALL DCOPY (NORDP1, G(I,1), MDG, W(I,1), MDW)
  216. 170 CONTINUE
  217. C
  218. C Solve for coefficients when possible.
  219. C
  220. DO 180 I = 1,N
  221. IF (G(I,1).EQ.0.D0) THEN
  222. MDEOUT = 2
  223. RETURN
  224. ENDIF
  225. 180 CONTINUE
  226. C
  227. C All the diagonal terms in the accumulated triangular
  228. C matrix are nonzero. The solution can be computed but
  229. C it may be unsuitable for further use due to poor
  230. C conditioning or the lack of constraints. No checking
  231. C for either of these is done here.
  232. C
  233. CALL DBNDSL (1, G, MDG, NORD, IP, IR, COEFF, N, RNORM)
  234. MDEOUT = 1
  235. RETURN
  236. END