dnls1.f 37 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018
  1. *DECK DNLS1
  2. SUBROUTINE DNLS1 (FCN, IOPT, M, N, X, FVEC, FJAC, LDFJAC, FTOL,
  3. + XTOL, GTOL, MAXFEV, EPSFCN, DIAG, MODE, FACTOR, NPRINT, INFO,
  4. + NFEV, NJEV, IPVT, QTF, WA1, WA2, WA3, WA4)
  5. C***BEGIN PROLOGUE DNLS1
  6. C***PURPOSE Minimize the sum of the squares of M nonlinear functions
  7. C in N variables by a modification of the Levenberg-Marquardt
  8. C algorithm.
  9. C***LIBRARY SLATEC
  10. C***CATEGORY K1B1A1, K1B1A2
  11. C***TYPE DOUBLE PRECISION (SNLS1-S, DNLS1-D)
  12. C***KEYWORDS LEVENBERG-MARQUARDT, NONLINEAR DATA FITTING,
  13. C NONLINEAR LEAST SQUARES
  14. C***AUTHOR Hiebert, K. L., (SNLA)
  15. C***DESCRIPTION
  16. C
  17. C 1. Purpose.
  18. C
  19. C The purpose of DNLS1 is to minimize the sum of the squares of M
  20. C nonlinear functions in N variables by a modification of the
  21. C Levenberg-Marquardt algorithm. The user must provide a subrou-
  22. C tine which calculates the functions. The user has the option
  23. C of how the Jacobian will be supplied. The user can supply the
  24. C full Jacobian, or the rows of the Jacobian (to avoid storing
  25. C the full Jacobian), or let the code approximate the Jacobian by
  26. C forward-differencing. This code is the combination of the
  27. C MINPACK codes (Argonne) LMDER, LMDIF, and LMSTR.
  28. C
  29. C
  30. C 2. Subroutine and Type Statements.
  31. C
  32. C SUBROUTINE DNLS1(FCN,IOPT,M,N,X,FVEC,FJAC,LDFJAC,FTOL,XTOL,
  33. C * GTOL,MAXFEV,EPSFCN,DIAG,MODE,FACTOR,NPRINT,INFO
  34. C * ,NFEV,NJEV,IPVT,QTF,WA1,WA2,WA3,WA4)
  35. C INTEGER IOPT,M,N,LDFJAC,MAXFEV,MODE,NPRINT,INFO,NFEV,NJEV
  36. C INTEGER IPVT(N)
  37. C DOUBLE PRECISION FTOL,XTOL,GTOL,EPSFCN,FACTOR
  38. C DOUBLE PRECISION X(N),FVEC(M),FJAC(LDFJAC,N),DIAG(N),QTF(N),
  39. C * WA1(N),WA2(N),WA3(N),WA4(M)
  40. C
  41. C
  42. C 3. Parameters.
  43. C
  44. C Parameters designated as input parameters must be specified on
  45. C entry to DNLS1 and are not changed on exit, while parameters
  46. C designated as output parameters need not be specified on entry
  47. C and are set to appropriate values on exit from DNLS1.
  48. C
  49. C FCN is the name of the user-supplied subroutine which calculate
  50. C the functions. If the user wants to supply the Jacobian
  51. C (IOPT=2 or 3), then FCN must be written to calculate the
  52. C Jacobian, as well as the functions. See the explanation
  53. C of the IOPT argument below.
  54. C If the user wants the iterates printed (NPRINT positive), then
  55. C FCN must do the printing. See the explanation of NPRINT
  56. C below. FCN must be declared in an EXTERNAL statement in the
  57. C calling program and should be written as follows.
  58. C
  59. C
  60. C SUBROUTINE FCN(IFLAG,M,N,X,FVEC,FJAC,LDFJAC)
  61. C INTEGER IFLAG,LDFJAC,M,N
  62. C DOUBLE PRECISION X(N),FVEC(M)
  63. C ----------
  64. C FJAC and LDFJAC may be ignored , if IOPT=1.
  65. C DOUBLE PRECISION FJAC(LDFJAC,N) , if IOPT=2.
  66. C DOUBLE PRECISION FJAC(N) , if IOPT=3.
  67. C ----------
  68. C If IFLAG=0, the values in X and FVEC are available
  69. C for printing. See the explanation of NPRINT below.
  70. C IFLAG will never be zero unless NPRINT is positive.
  71. C The values of X and FVEC must not be changed.
  72. C RETURN
  73. C ----------
  74. C If IFLAG=1, calculate the functions at X and return
  75. C this vector in FVEC.
  76. C RETURN
  77. C ----------
  78. C If IFLAG=2, calculate the full Jacobian at X and return
  79. C this matrix in FJAC. Note that IFLAG will never be 2 unless
  80. C IOPT=2. FVEC contains the function values at X and must
  81. C not be altered. FJAC(I,J) must be set to the derivative
  82. C of FVEC(I) with respect to X(J).
  83. C RETURN
  84. C ----------
  85. C If IFLAG=3, calculate the LDFJAC-th row of the Jacobian
  86. C and return this vector in FJAC. Note that IFLAG will
  87. C never be 3 unless IOPT=3. FVEC contains the function
  88. C values at X and must not be altered. FJAC(J) must be
  89. C set to the derivative of FVEC(LDFJAC) with respect to X(J).
  90. C RETURN
  91. C ----------
  92. C END
  93. C
  94. C
  95. C The value of IFLAG should not be changed by FCN unless the
  96. C user wants to terminate execution of DNLS1. In this case, set
  97. C IFLAG to a negative integer.
  98. C
  99. C
  100. C IOPT is an input variable which specifies how the Jacobian will
  101. C be calculated. If IOPT=2 or 3, then the user must supply the
  102. C Jacobian, as well as the function values, through the
  103. C subroutine FCN. If IOPT=2, the user supplies the full
  104. C Jacobian with one call to FCN. If IOPT=3, the user supplies
  105. C one row of the Jacobian with each call. (In this manner,
  106. C storage can be saved because the full Jacobian is not stored.)
  107. C If IOPT=1, the code will approximate the Jacobian by forward
  108. C differencing.
  109. C
  110. C M is a positive integer input variable set to the number of
  111. C functions.
  112. C
  113. C N is a positive integer input variable set to the number of
  114. C variables. N must not exceed M.
  115. C
  116. C X is an array of length N. On input, X must contain an initial
  117. C estimate of the solution vector. On output, X contains the
  118. C final estimate of the solution vector.
  119. C
  120. C FVEC is an output array of length M which contains the functions
  121. C evaluated at the output X.
  122. C
  123. C FJAC is an output array. For IOPT=1 and 2, FJAC is an M by N
  124. C array. For IOPT=3, FJAC is an N by N array. The upper N by N
  125. C submatrix of FJAC contains an upper triangular matrix R with
  126. C diagonal elements of nonincreasing magnitude such that
  127. C
  128. C T T T
  129. C P *(JAC *JAC)*P = R *R,
  130. C
  131. C where P is a permutation matrix and JAC is the final calcu-
  132. C lated Jacobian. Column J of P is column IPVT(J) (see below)
  133. C of the identity matrix. The lower part of FJAC contains
  134. C information generated during the computation of R.
  135. C
  136. C LDFJAC is a positive integer input variable which specifies
  137. C the leading dimension of the array FJAC. For IOPT=1 and 2,
  138. C LDFJAC must not be less than M. For IOPT=3, LDFJAC must not
  139. C be less than N.
  140. C
  141. C FTOL is a non-negative input variable. Termination occurs when
  142. C both the actual and predicted relative reductions in the sum
  143. C of squares are at most FTOL. Therefore, FTOL measures the
  144. C relative error desired in the sum of squares. Section 4 con-
  145. C tains more details about FTOL.
  146. C
  147. C XTOL is a non-negative input variable. Termination occurs when
  148. C the relative error between two consecutive iterates is at most
  149. C XTOL. Therefore, XTOL measures the relative error desired in
  150. C the approximate solution. Section 4 contains more details
  151. C about XTOL.
  152. C
  153. C GTOL is a non-negative input variable. Termination occurs when
  154. C the cosine of the angle between FVEC and any column of the
  155. C Jacobian is at most GTOL in absolute value. Therefore, GTOL
  156. C measures the orthogonality desired between the function vector
  157. C and the columns of the Jacobian. Section 4 contains more
  158. C details about GTOL.
  159. C
  160. C MAXFEV is a positive integer input variable. Termination occurs
  161. C when the number of calls to FCN to evaluate the functions
  162. C has reached MAXFEV.
  163. C
  164. C EPSFCN is an input variable used in determining a suitable step
  165. C for the forward-difference approximation. This approximation
  166. C assumes that the relative errors in the functions are of the
  167. C order of EPSFCN. If EPSFCN is less than the machine preci-
  168. C sion, it is assumed that the relative errors in the functions
  169. C are of the order of the machine precision. If IOPT=2 or 3,
  170. C then EPSFCN can be ignored (treat it as a dummy argument).
  171. C
  172. C DIAG is an array of length N. If MODE = 1 (see below), DIAG is
  173. C internally set. If MODE = 2, DIAG must contain positive
  174. C entries that serve as implicit (multiplicative) scale factors
  175. C for the variables.
  176. C
  177. C MODE is an integer input variable. If MODE = 1, the variables
  178. C will be scaled internally. If MODE = 2, the scaling is speci-
  179. C fied by the input DIAG. Other values of MODE are equivalent
  180. C to MODE = 1.
  181. C
  182. C FACTOR is a positive input variable used in determining the ini-
  183. C tial step bound. This bound is set to the product of FACTOR
  184. C and the Euclidean norm of DIAG*X if nonzero, or else to FACTOR
  185. C itself. In most cases FACTOR should lie in the interval
  186. C (.1,100.). 100. is a generally recommended value.
  187. C
  188. C NPRINT is an integer input variable that enables controlled
  189. C printing of iterates if it is positive. In this case, FCN is
  190. C called with IFLAG = 0 at the beginning of the first iteration
  191. C and every NPRINT iterations thereafter and immediately prior
  192. C to return, with X and FVEC available for printing. Appropriate
  193. C print statements must be added to FCN (see example) and
  194. C FVEC should not be altered. If NPRINT is not positive, no
  195. C special calls to FCN with IFLAG = 0 are made.
  196. C
  197. C INFO is an integer output variable. If the user has terminated
  198. C execution, INFO is set to the (negative) value of IFLAG. See
  199. C description of FCN and JAC. Otherwise, INFO is set as follows
  200. C
  201. C INFO = 0 improper input parameters.
  202. C
  203. C INFO = 1 both actual and predicted relative reductions in the
  204. C sum of squares are at most FTOL.
  205. C
  206. C INFO = 2 relative error between two consecutive iterates is
  207. C at most XTOL.
  208. C
  209. C INFO = 3 conditions for INFO = 1 and INFO = 2 both hold.
  210. C
  211. C INFO = 4 the cosine of the angle between FVEC and any column
  212. C of the Jacobian is at most GTOL in absolute value.
  213. C
  214. C INFO = 5 number of calls to FCN for function evaluation
  215. C has reached MAXFEV.
  216. C
  217. C INFO = 6 FTOL is too small. No further reduction in the sum
  218. C of squares is possible.
  219. C
  220. C INFO = 7 XTOL is too small. No further improvement in the
  221. C approximate solution X is possible.
  222. C
  223. C INFO = 8 GTOL is too small. FVEC is orthogonal to the
  224. C columns of the Jacobian to machine precision.
  225. C
  226. C Sections 4 and 5 contain more details about INFO.
  227. C
  228. C NFEV is an integer output variable set to the number of calls to
  229. C FCN for function evaluation.
  230. C
  231. C NJEV is an integer output variable set to the number of
  232. C evaluations of the full Jacobian. If IOPT=2, only one call to
  233. C FCN is required for each evaluation of the full Jacobian.
  234. C If IOPT=3, the M calls to FCN are required.
  235. C If IOPT=1, then NJEV is set to zero.
  236. C
  237. C IPVT is an integer output array of length N. IPVT defines a
  238. C permutation matrix P such that JAC*P = Q*R, where JAC is the
  239. C final calculated Jacobian, Q is orthogonal (not stored), and R
  240. C is upper triangular with diagonal elements of nonincreasing
  241. C magnitude. Column J of P is column IPVT(J) of the identity
  242. C matrix.
  243. C
  244. C QTF is an output array of length N which contains the first N
  245. C elements of the vector (Q transpose)*FVEC.
  246. C
  247. C WA1, WA2, and WA3 are work arrays of length N.
  248. C
  249. C WA4 is a work array of length M.
  250. C
  251. C
  252. C 4. Successful Completion.
  253. C
  254. C The accuracy of DNLS1 is controlled by the convergence parame-
  255. C ters FTOL, XTOL, and GTOL. These parameters are used in tests
  256. C which make three types of comparisons between the approximation
  257. C X and a solution XSOL. DNLS1 terminates when any of the tests
  258. C is satisfied. If any of the convergence parameters is less than
  259. C the machine precision (as defined by the function R1MACH(4)),
  260. C then DNLS1 only attempts to satisfy the test defined by the
  261. C machine precision. Further progress is not usually possible.
  262. C
  263. C The tests assume that the functions are reasonably well behaved,
  264. C and, if the Jacobian is supplied by the user, that the functions
  265. C and the Jacobian are coded consistently. If these conditions
  266. C are not satisfied, then DNLS1 may incorrectly indicate conver-
  267. C gence. If the Jacobian is coded correctly or IOPT=1,
  268. C then the validity of the answer can be checked, for example, by
  269. C rerunning DNLS1 with tighter tolerances.
  270. C
  271. C First Convergence Test. If ENORM(Z) denotes the Euclidean norm
  272. C of a vector Z, then this test attempts to guarantee that
  273. C
  274. C ENORM(FVEC) .LE. (1+FTOL)*ENORM(FVECS),
  275. C
  276. C where FVECS denotes the functions evaluated at XSOL. If this
  277. C condition is satisfied with FTOL = 10**(-K), then the final
  278. C residual norm ENORM(FVEC) has K significant decimal digits and
  279. C INFO is set to 1 (or to 3 if the second test is also satis-
  280. C fied). Unless high precision solutions are required, the
  281. C recommended value for FTOL is the square root of the machine
  282. C precision.
  283. C
  284. C Second Convergence Test. If D is the diagonal matrix whose
  285. C entries are defined by the array DIAG, then this test attempts
  286. C to guarantee that
  287. C
  288. C ENORM(D*(X-XSOL)) .LE. XTOL*ENORM(D*XSOL).
  289. C
  290. C If this condition is satisfied with XTOL = 10**(-K), then the
  291. C larger components of D*X have K significant decimal digits and
  292. C INFO is set to 2 (or to 3 if the first test is also satis-
  293. C fied). There is a danger that the smaller components of D*X
  294. C may have large relative errors, but if MODE = 1, then the
  295. C accuracy of the components of X is usually related to their
  296. C sensitivity. Unless high precision solutions are required,
  297. C the recommended value for XTOL is the square root of the
  298. C machine precision.
  299. C
  300. C Third Convergence Test. This test is satisfied when the cosine
  301. C of the angle between FVEC and any column of the Jacobian at X
  302. C is at most GTOL in absolute value. There is no clear rela-
  303. C tionship between this test and the accuracy of DNLS1, and
  304. C furthermore, the test is equally well satisfied at other crit-
  305. C ical points, namely maximizers and saddle points. Therefore,
  306. C termination caused by this test (INFO = 4) should be examined
  307. C carefully. The recommended value for GTOL is zero.
  308. C
  309. C
  310. C 5. Unsuccessful Completion.
  311. C
  312. C Unsuccessful termination of DNLS1 can be due to improper input
  313. C parameters, arithmetic interrupts, or an excessive number of
  314. C function evaluations.
  315. C
  316. C Improper Input Parameters. INFO is set to 0 if IOPT .LT. 1
  317. C or IOPT .GT. 3, or N .LE. 0, or M .LT. N, or for IOPT=1 or 2
  318. C LDFJAC .LT. M, or for IOPT=3 LDFJAC .LT. N, or FTOL .LT. 0.E0,
  319. C or XTOL .LT. 0.E0, or GTOL .LT. 0.E0, or MAXFEV .LE. 0, or
  320. C FACTOR .LE. 0.E0.
  321. C
  322. C Arithmetic Interrupts. If these interrupts occur in the FCN
  323. C subroutine during an early stage of the computation, they may
  324. C be caused by an unacceptable choice of X by DNLS1. In this
  325. C case, it may be possible to remedy the situation by rerunning
  326. C DNLS1 with a smaller value of FACTOR.
  327. C
  328. C Excessive Number of Function Evaluations. A reasonable value
  329. C for MAXFEV is 100*(N+1) for IOPT=2 or 3 and 200*(N+1) for
  330. C IOPT=1. If the number of calls to FCN reaches MAXFEV, then
  331. C this indicates that the routine is converging very slowly
  332. C as measured by the progress of FVEC, and INFO is set to 5.
  333. C In this case, it may be helpful to restart DNLS1 with MODE
  334. C set to 1.
  335. C
  336. C
  337. C 6. Characteristics of the Algorithm.
  338. C
  339. C DNLS1 is a modification of the Levenberg-Marquardt algorithm.
  340. C Two of its main characteristics involve the proper use of
  341. C implicitly scaled variables (if MODE = 1) and an optimal choice
  342. C for the correction. The use of implicitly scaled variables
  343. C achieves scale invariance of DNLS1 and limits the size of the
  344. C correction in any direction where the functions are changing
  345. C rapidly. The optimal choice of the correction guarantees (under
  346. C reasonable conditions) global convergence from starting points
  347. C far from the solution and a fast rate of convergence for
  348. C problems with small residuals.
  349. C
  350. C Timing. The time required by DNLS1 to solve a given problem
  351. C depends on M and N, the behavior of the functions, the accu-
  352. C racy requested, and the starting point. The number of arith-
  353. C metic operations needed by DNLS1 is about N**3 to process each
  354. C evaluation of the functions (call to FCN) and to process each
  355. C evaluation of the Jacobian it takes M*N**2 for IOPT=2 (one
  356. C call to FCN), M*N**2 for IOPT=1 (N calls to FCN) and
  357. C 1.5*M*N**2 for IOPT=3 (M calls to FCN). Unless FCN
  358. C can be evaluated quickly, the timing of DNLS1 will be
  359. C strongly influenced by the time spent in FCN.
  360. C
  361. C Storage. DNLS1 requires (M*N + 2*M + 6*N) for IOPT=1 or 2 and
  362. C (N**2 + 2*M + 6*N) for IOPT=3 single precision storage
  363. C locations and N integer storage locations, in addition to
  364. C the storage required by the program. There are no internally
  365. C declared storage arrays.
  366. C
  367. C *Long Description:
  368. C
  369. C 7. Example.
  370. C
  371. C The problem is to determine the values of X(1), X(2), and X(3)
  372. C which provide the best fit (in the least squares sense) of
  373. C
  374. C X(1) + U(I)/(V(I)*X(2) + W(I)*X(3)), I = 1, 15
  375. C
  376. C to the data
  377. C
  378. C Y = (0.14,0.18,0.22,0.25,0.29,0.32,0.35,0.39,
  379. C 0.37,0.58,0.73,0.96,1.34,2.10,4.39),
  380. C
  381. C where U(I) = I, V(I) = 16 - I, and W(I) = MIN(U(I),V(I)). The
  382. C I-th component of FVEC is thus defined by
  383. C
  384. C Y(I) - (X(1) + U(I)/(V(I)*X(2) + W(I)*X(3))).
  385. C
  386. C **********
  387. C
  388. C PROGRAM TEST
  389. C C
  390. C C Driver for DNLS1 example.
  391. C C
  392. C INTEGER J,IOPT,M,N,LDFJAC,MAXFEV,MODE,NPRINT,INFO,NFEV,NJEV,
  393. C * NWRITE
  394. C INTEGER IPVT(3)
  395. C DOUBLE PRECISION FTOL,XTOL,GTOL,FACTOR,FNORM,EPSFCN
  396. C DOUBLE PRECISION X(3),FVEC(15),FJAC(15,3),DIAG(3),QTF(3),
  397. C * WA1(3),WA2(3),WA3(3),WA4(15)
  398. C DOUBLE PRECISION DENORM,D1MACH
  399. C EXTERNAL FCN
  400. C DATA NWRITE /6/
  401. C C
  402. C IOPT = 1
  403. C M = 15
  404. C N = 3
  405. C C
  406. C C The following starting values provide a rough fit.
  407. C C
  408. C X(1) = 1.E0
  409. C X(2) = 1.E0
  410. C X(3) = 1.E0
  411. C C
  412. C LDFJAC = 15
  413. C C
  414. C C Set FTOL and XTOL to the square root of the machine precision
  415. C C and GTOL to zero. Unless high precision solutions are
  416. C C required, these are the recommended settings.
  417. C C
  418. C FTOL = SQRT(R1MACH(4))
  419. C XTOL = SQRT(R1MACH(4))
  420. C GTOL = 0.E0
  421. C C
  422. C MAXFEV = 400
  423. C EPSFCN = 0.0
  424. C MODE = 1
  425. C FACTOR = 1.E2
  426. C NPRINT = 0
  427. C C
  428. C CALL DNLS1(FCN,IOPT,M,N,X,FVEC,FJAC,LDFJAC,FTOL,XTOL,
  429. C * GTOL,MAXFEV,EPSFCN,DIAG,MODE,FACTOR,NPRINT,
  430. C * INFO,NFEV,NJEV,IPVT,QTF,WA1,WA2,WA3,WA4)
  431. C FNORM = ENORM(M,FVEC)
  432. C WRITE (NWRITE,1000) FNORM,NFEV,NJEV,INFO,(X(J),J=1,N)
  433. C STOP
  434. C 1000 FORMAT (5X,' FINAL L2 NORM OF THE RESIDUALS',E15.7 //
  435. C * 5X,' NUMBER OF FUNCTION EVALUATIONS',I10 //
  436. C * 5X,' NUMBER OF JACOBIAN EVALUATIONS',I10 //
  437. C * 5X,' EXIT PARAMETER',16X,I10 //
  438. C * 5X,' FINAL APPROXIMATE SOLUTION' // 5X,3E15.7)
  439. C END
  440. C SUBROUTINE FCN(IFLAG,M,N,X,FVEC,DUM,IDUM)
  441. C C This is the form of the FCN routine if IOPT=1,
  442. C C that is, if the user does not calculate the Jacobian.
  443. C INTEGER I,M,N,IFLAG
  444. C DOUBLE PRECISION X(N),FVEC(M),Y(15)
  445. C DOUBLE PRECISION TMP1,TMP2,TMP3,TMP4
  446. C DATA Y(1),Y(2),Y(3),Y(4),Y(5),Y(6),Y(7),Y(8),
  447. C * Y(9),Y(10),Y(11),Y(12),Y(13),Y(14),Y(15)
  448. C * /1.4E-1,1.8E-1,2.2E-1,2.5E-1,2.9E-1,3.2E-1,3.5E-1,3.9E-1,
  449. C * 3.7E-1,5.8E-1,7.3E-1,9.6E-1,1.34E0,2.1E0,4.39E0/
  450. C C
  451. C IF (IFLAG .NE. 0) GO TO 5
  452. C C
  453. C C Insert print statements here when NPRINT is positive.
  454. C C
  455. C RETURN
  456. C 5 CONTINUE
  457. C DO 10 I = 1, M
  458. C TMP1 = I
  459. C TMP2 = 16 - I
  460. C TMP3 = TMP1
  461. C IF (I .GT. 8) TMP3 = TMP2
  462. C FVEC(I) = Y(I) - (X(1) + TMP1/(X(2)*TMP2 + X(3)*TMP3))
  463. C 10 CONTINUE
  464. C RETURN
  465. C END
  466. C
  467. C
  468. C Results obtained with different compilers or machines
  469. C may be slightly different.
  470. C
  471. C FINAL L2 NORM OF THE RESIDUALS 0.9063596E-01
  472. C
  473. C NUMBER OF FUNCTION EVALUATIONS 25
  474. C
  475. C NUMBER OF JACOBIAN EVALUATIONS 0
  476. C
  477. C EXIT PARAMETER 1
  478. C
  479. C FINAL APPROXIMATE SOLUTION
  480. C
  481. C 0.8241058E-01 0.1133037E+01 0.2343695E+01
  482. C
  483. C
  484. C For IOPT=2, FCN would be modified as follows to also
  485. C calculate the full Jacobian when IFLAG=2.
  486. C
  487. C SUBROUTINE FCN(IFLAG,M,N,X,FVEC,FJAC,LDFJAC)
  488. C C
  489. C C This is the form of the FCN routine if IOPT=2,
  490. C C that is, if the user calculates the full Jacobian.
  491. C C
  492. C INTEGER I,LDFJAC,M,N,IFLAG
  493. C DOUBLE PRECISION X(N),FVEC(M),FJAC(LDFJAC,N),Y(15)
  494. C DOUBLE PRECISION TMP1,TMP2,TMP3,TMP4
  495. C DATA Y(1),Y(2),Y(3),Y(4),Y(5),Y(6),Y(7),Y(8),
  496. C * Y(9),Y(10),Y(11),Y(12),Y(13),Y(14),Y(15)
  497. C * /1.4E-1,1.8E-1,2.2E-1,2.5E-1,2.9E-1,3.2E-1,3.5E-1,3.9E-1,
  498. C * 3.7E-1,5.8E-1,7.3E-1,9.6E-1,1.34E0,2.1E0,4.39E0/
  499. C C
  500. C IF (IFLAG .NE. 0) GO TO 5
  501. C C
  502. C C Insert print statements here when NPRINT is positive.
  503. C C
  504. C RETURN
  505. C 5 CONTINUE
  506. C IF(IFLAG.NE.1) GO TO 20
  507. C DO 10 I = 1, M
  508. C TMP1 = I
  509. C TMP2 = 16 - I
  510. C TMP3 = TMP1
  511. C IF (I .GT. 8) TMP3 = TMP2
  512. C FVEC(I) = Y(I) - (X(1) + TMP1/(X(2)*TMP2 + X(3)*TMP3))
  513. C 10 CONTINUE
  514. C RETURN
  515. C C
  516. C C Below, calculate the full Jacobian.
  517. C C
  518. C 20 CONTINUE
  519. C C
  520. C DO 30 I = 1, M
  521. C TMP1 = I
  522. C TMP2 = 16 - I
  523. C TMP3 = TMP1
  524. C IF (I .GT. 8) TMP3 = TMP2
  525. C TMP4 = (X(2)*TMP2 + X(3)*TMP3)**2
  526. C FJAC(I,1) = -1.E0
  527. C FJAC(I,2) = TMP1*TMP2/TMP4
  528. C FJAC(I,3) = TMP1*TMP3/TMP4
  529. C 30 CONTINUE
  530. C RETURN
  531. C END
  532. C
  533. C
  534. C For IOPT = 3, FJAC would be dimensioned as FJAC(3,3),
  535. C LDFJAC would be set to 3, and FCN would be written as
  536. C follows to calculate a row of the Jacobian when IFLAG=3.
  537. C
  538. C SUBROUTINE FCN(IFLAG,M,N,X,FVEC,FJAC,LDFJAC)
  539. C C This is the form of the FCN routine if IOPT=3,
  540. C C that is, if the user calculates the Jacobian row by row.
  541. C INTEGER I,M,N,IFLAG
  542. C DOUBLE PRECISION X(N),FVEC(M),FJAC(N),Y(15)
  543. C DOUBLE PRECISION TMP1,TMP2,TMP3,TMP4
  544. C DATA Y(1),Y(2),Y(3),Y(4),Y(5),Y(6),Y(7),Y(8),
  545. C * Y(9),Y(10),Y(11),Y(12),Y(13),Y(14),Y(15)
  546. C * /1.4E-1,1.8E-1,2.2E-1,2.5E-1,2.9E-1,3.2E-1,3.5E-1,3.9E-1,
  547. C * 3.7E-1,5.8E-1,7.3E-1,9.6E-1,1.34E0,2.1E0,4.39E0/
  548. C C
  549. C IF (IFLAG .NE. 0) GO TO 5
  550. C C
  551. C C Insert print statements here when NPRINT is positive.
  552. C C
  553. C RETURN
  554. C 5 CONTINUE
  555. C IF( IFLAG.NE.1) GO TO 20
  556. C DO 10 I = 1, M
  557. C TMP1 = I
  558. C TMP2 = 16 - I
  559. C TMP3 = TMP1
  560. C IF (I .GT. 8) TMP3 = TMP2
  561. C FVEC(I) = Y(I) - (X(1) + TMP1/(X(2)*TMP2 + X(3)*TMP3))
  562. C 10 CONTINUE
  563. C RETURN
  564. C C
  565. C C Below, calculate the LDFJAC-th row of the Jacobian.
  566. C C
  567. C 20 CONTINUE
  568. C
  569. C I = LDFJAC
  570. C TMP1 = I
  571. C TMP2 = 16 - I
  572. C TMP3 = TMP1
  573. C IF (I .GT. 8) TMP3 = TMP2
  574. C TMP4 = (X(2)*TMP2 + X(3)*TMP3)**2
  575. C FJAC(1) = -1.E0
  576. C FJAC(2) = TMP1*TMP2/TMP4
  577. C FJAC(3) = TMP1*TMP3/TMP4
  578. C RETURN
  579. C END
  580. C
  581. C***REFERENCES Jorge J. More, The Levenberg-Marquardt algorithm:
  582. C implementation and theory. In Numerical Analysis
  583. C Proceedings (Dundee, June 28 - July 1, 1977, G. A.
  584. C Watson, Editor), Lecture Notes in Mathematics 630,
  585. C Springer-Verlag, 1978.
  586. C***ROUTINES CALLED D1MACH, DCKDER, DENORM, DFDJC3, DMPAR, DQRFAC,
  587. C DWUPDT, XERMSG
  588. C***REVISION HISTORY (YYMMDD)
  589. C 800301 DATE WRITTEN
  590. C 890531 Changed all specific intrinsics to generic. (WRB)
  591. C 890831 Modified array declarations. (WRB)
  592. C 891006 Cosmetic changes to prologue. (WRB)
  593. C 891006 REVISION DATE from Version 3.2
  594. C 891214 Prologue converted to Version 4.0 format. (BAB)
  595. C 900315 CALLs to XERROR changed to CALLs to XERMSG. (THJ)
  596. C 900510 Convert XERRWV calls to XERMSG calls. (RWC)
  597. C 920205 Corrected XERN1 declaration. (WRB)
  598. C 920501 Reformatted the REFERENCES section. (WRB)
  599. C***END PROLOGUE DNLS1
  600. IMPLICIT DOUBLE PRECISION (A-H,O-Z)
  601. INTEGER IOPT,M,N,LDFJAC,MAXFEV,MODE,NPRINT,INFO,NFEV,NJEV
  602. INTEGER IJUNK,NROW,IPVT(*)
  603. DOUBLE PRECISION FTOL,XTOL,GTOL,FACTOR,EPSFCN
  604. DOUBLE PRECISION X(*),FVEC(*),FJAC(LDFJAC,*),DIAG(*),QTF(*),
  605. 1 WA1(*),WA2(*),WA3(*),WA4(*)
  606. LOGICAL SING
  607. EXTERNAL FCN
  608. INTEGER I,IFLAG,ITER,J,L,MODECH
  609. DOUBLE PRECISION ACTRED,DELTA,DIRDER,EPSMCH,FNORM,FNORM1,GNORM,
  610. 1 ONE,PAR,PNORM,PRERED,P1,P5,P25,P75,P0001,RATIO,SUM,TEMP,
  611. 2 TEMP1,TEMP2,XNORM,ZERO
  612. DOUBLE PRECISION D1MACH,DENORM,ERR,CHKLIM
  613. CHARACTER*8 XERN1
  614. CHARACTER*16 XERN3
  615. SAVE CHKLIM, ONE, P1, P5, P25, P75, P0001, ZERO
  616. C
  617. DATA CHKLIM/.1D0/
  618. DATA ONE,P1,P5,P25,P75,P0001,ZERO
  619. 1 /1.0D0,1.0D-1,5.0D-1,2.5D-1,7.5D-1,1.0D-4,0.0D0/
  620. C***FIRST EXECUTABLE STATEMENT DNLS1
  621. EPSMCH = D1MACH(4)
  622. C
  623. INFO = 0
  624. IFLAG = 0
  625. NFEV = 0
  626. NJEV = 0
  627. C
  628. C CHECK THE INPUT PARAMETERS FOR ERRORS.
  629. C
  630. IF (IOPT .LT. 1 .OR. IOPT .GT. 3 .OR. N .LE. 0 .OR.
  631. 1 M .LT. N .OR. LDFJAC .LT. N .OR. FTOL .LT. ZERO
  632. 2 .OR. XTOL .LT. ZERO .OR. GTOL .LT. ZERO
  633. 3 .OR. MAXFEV .LE. 0 .OR. FACTOR .LE. ZERO) GO TO 300
  634. IF (IOPT .LT. 3 .AND. LDFJAC .LT. M) GO TO 300
  635. IF (MODE .NE. 2) GO TO 20
  636. DO 10 J = 1, N
  637. IF (DIAG(J) .LE. ZERO) GO TO 300
  638. 10 CONTINUE
  639. 20 CONTINUE
  640. C
  641. C EVALUATE THE FUNCTION AT THE STARTING POINT
  642. C AND CALCULATE ITS NORM.
  643. C
  644. IFLAG = 1
  645. IJUNK = 1
  646. CALL FCN(IFLAG,M,N,X,FVEC,FJAC,IJUNK)
  647. NFEV = 1
  648. IF (IFLAG .LT. 0) GO TO 300
  649. FNORM = DENORM(M,FVEC)
  650. C
  651. C INITIALIZE LEVENBERG-MARQUARDT PARAMETER AND ITERATION COUNTER.
  652. C
  653. PAR = ZERO
  654. ITER = 1
  655. C
  656. C BEGINNING OF THE OUTER LOOP.
  657. C
  658. 30 CONTINUE
  659. C
  660. C IF REQUESTED, CALL FCN TO ENABLE PRINTING OF ITERATES.
  661. C
  662. IF (NPRINT .LE. 0) GO TO 40
  663. IFLAG = 0
  664. IF (MOD(ITER-1,NPRINT) .EQ. 0)
  665. 1 CALL FCN(IFLAG,M,N,X,FVEC,FJAC,IJUNK)
  666. IF (IFLAG .LT. 0) GO TO 300
  667. 40 CONTINUE
  668. C
  669. C CALCULATE THE JACOBIAN MATRIX.
  670. C
  671. IF (IOPT .EQ. 3) GO TO 475
  672. C
  673. C STORE THE FULL JACOBIAN USING M*N STORAGE
  674. C
  675. IF (IOPT .EQ. 1) GO TO 410
  676. C
  677. C THE USER SUPPLIES THE JACOBIAN
  678. C
  679. IFLAG = 2
  680. CALL FCN(IFLAG,M,N,X,FVEC,FJAC,LDFJAC)
  681. NJEV = NJEV + 1
  682. C
  683. C ON THE FIRST ITERATION, CHECK THE USER SUPPLIED JACOBIAN
  684. C
  685. IF (ITER .LE. 1) THEN
  686. IF (IFLAG .LT. 0) GO TO 300
  687. C
  688. C GET THE INCREMENTED X-VALUES INTO WA1(*).
  689. C
  690. MODECH = 1
  691. CALL DCKDER(M,N,X,FVEC,FJAC,LDFJAC,WA1,WA4,MODECH,ERR)
  692. C
  693. C EVALUATE FUNCTION AT INCREMENTED VALUE AND PUT IN WA4(*).
  694. C
  695. IFLAG = 1
  696. CALL FCN(IFLAG,M,N,WA1,WA4,FJAC,LDFJAC)
  697. NFEV = NFEV + 1
  698. IF(IFLAG .LT. 0) GO TO 300
  699. DO 350 I = 1, M
  700. MODECH = 2
  701. CALL DCKDER(1,N,X,FVEC(I),FJAC(I,1),LDFJAC,WA1,
  702. 1 WA4(I),MODECH,ERR)
  703. IF (ERR .LT. CHKLIM) THEN
  704. WRITE (XERN1, '(I8)') I
  705. WRITE (XERN3, '(1PE15.6)') ERR
  706. CALL XERMSG ('SLATEC', 'DNLS1', 'DERIVATIVE OF ' //
  707. * 'FUNCTION ' // XERN1 // ' MAY BE WRONG, ERR = ' //
  708. * XERN3 // ' TOO CLOSE TO 0.', 7, 0)
  709. ENDIF
  710. 350 CONTINUE
  711. ENDIF
  712. C
  713. GO TO 420
  714. C
  715. C THE CODE APPROXIMATES THE JACOBIAN
  716. C
  717. 410 IFLAG = 1
  718. CALL DFDJC3(FCN,M,N,X,FVEC,FJAC,LDFJAC,IFLAG,EPSFCN,WA4)
  719. NFEV = NFEV + N
  720. 420 IF (IFLAG .LT. 0) GO TO 300
  721. C
  722. C COMPUTE THE QR FACTORIZATION OF THE JACOBIAN.
  723. C
  724. CALL DQRFAC(M,N,FJAC,LDFJAC,.TRUE.,IPVT,N,WA1,WA2,WA3)
  725. C
  726. C FORM (Q TRANSPOSE)*FVEC AND STORE THE FIRST N COMPONENTS IN
  727. C QTF.
  728. C
  729. DO 430 I = 1, M
  730. WA4(I) = FVEC(I)
  731. 430 CONTINUE
  732. DO 470 J = 1, N
  733. IF (FJAC(J,J) .EQ. ZERO) GO TO 460
  734. SUM = ZERO
  735. DO 440 I = J, M
  736. SUM = SUM + FJAC(I,J)*WA4(I)
  737. 440 CONTINUE
  738. TEMP = -SUM/FJAC(J,J)
  739. DO 450 I = J, M
  740. WA4(I) = WA4(I) + FJAC(I,J)*TEMP
  741. 450 CONTINUE
  742. 460 CONTINUE
  743. FJAC(J,J) = WA1(J)
  744. QTF(J) = WA4(J)
  745. 470 CONTINUE
  746. GO TO 560
  747. C
  748. C ACCUMULATE THE JACOBIAN BY ROWS IN ORDER TO SAVE STORAGE.
  749. C COMPUTE THE QR FACTORIZATION OF THE JACOBIAN MATRIX
  750. C CALCULATED ONE ROW AT A TIME, WHILE SIMULTANEOUSLY
  751. C FORMING (Q TRANSPOSE)*FVEC AND STORING THE FIRST
  752. C N COMPONENTS IN QTF.
  753. C
  754. 475 DO 490 J = 1, N
  755. QTF(J) = ZERO
  756. DO 480 I = 1, N
  757. FJAC(I,J) = ZERO
  758. 480 CONTINUE
  759. 490 CONTINUE
  760. DO 500 I = 1, M
  761. NROW = I
  762. IFLAG = 3
  763. CALL FCN(IFLAG,M,N,X,FVEC,WA3,NROW)
  764. IF (IFLAG .LT. 0) GO TO 300
  765. C
  766. C ON THE FIRST ITERATION, CHECK THE USER SUPPLIED JACOBIAN.
  767. C
  768. IF(ITER .GT. 1) GO TO 498
  769. C
  770. C GET THE INCREMENTED X-VALUES INTO WA1(*).
  771. C
  772. MODECH = 1
  773. CALL DCKDER(M,N,X,FVEC,FJAC,LDFJAC,WA1,WA4,MODECH,ERR)
  774. C
  775. C EVALUATE AT INCREMENTED VALUES, IF NOT ALREADY EVALUATED.
  776. C
  777. IF(I .NE. 1) GO TO 495
  778. C
  779. C EVALUATE FUNCTION AT INCREMENTED VALUE AND PUT INTO WA4(*).
  780. C
  781. IFLAG = 1
  782. CALL FCN(IFLAG,M,N,WA1,WA4,FJAC,NROW)
  783. NFEV = NFEV + 1
  784. IF(IFLAG .LT. 0) GO TO 300
  785. 495 CONTINUE
  786. MODECH = 2
  787. CALL DCKDER(1,N,X,FVEC(I),WA3,1,WA1,WA4(I),MODECH,ERR)
  788. IF (ERR .LT. CHKLIM) THEN
  789. WRITE (XERN1, '(I8)') I
  790. WRITE (XERN3, '(1PE15.6)') ERR
  791. CALL XERMSG ('SLATEC', 'DNLS1', 'DERIVATIVE OF FUNCTION '
  792. * // XERN1 // ' MAY BE WRONG, ERR = ' // XERN3 //
  793. * ' TOO CLOSE TO 0.', 7, 0)
  794. ENDIF
  795. 498 CONTINUE
  796. C
  797. TEMP = FVEC(I)
  798. CALL DWUPDT(N,FJAC,LDFJAC,WA3,QTF,TEMP,WA1,WA2)
  799. 500 CONTINUE
  800. NJEV = NJEV + 1
  801. C
  802. C IF THE JACOBIAN IS RANK DEFICIENT, CALL DQRFAC TO
  803. C REORDER ITS COLUMNS AND UPDATE THE COMPONENTS OF QTF.
  804. C
  805. SING = .FALSE.
  806. DO 510 J = 1, N
  807. IF (FJAC(J,J) .EQ. ZERO) SING = .TRUE.
  808. IPVT(J) = J
  809. WA2(J) = DENORM(J,FJAC(1,J))
  810. 510 CONTINUE
  811. IF (.NOT.SING) GO TO 560
  812. CALL DQRFAC(N,N,FJAC,LDFJAC,.TRUE.,IPVT,N,WA1,WA2,WA3)
  813. DO 550 J = 1, N
  814. IF (FJAC(J,J) .EQ. ZERO) GO TO 540
  815. SUM = ZERO
  816. DO 520 I = J, N
  817. SUM = SUM + FJAC(I,J)*QTF(I)
  818. 520 CONTINUE
  819. TEMP = -SUM/FJAC(J,J)
  820. DO 530 I = J, N
  821. QTF(I) = QTF(I) + FJAC(I,J)*TEMP
  822. 530 CONTINUE
  823. 540 CONTINUE
  824. FJAC(J,J) = WA1(J)
  825. 550 CONTINUE
  826. 560 CONTINUE
  827. C
  828. C ON THE FIRST ITERATION AND IF MODE IS 1, SCALE ACCORDING
  829. C TO THE NORMS OF THE COLUMNS OF THE INITIAL JACOBIAN.
  830. C
  831. IF (ITER .NE. 1) GO TO 80
  832. IF (MODE .EQ. 2) GO TO 60
  833. DO 50 J = 1, N
  834. DIAG(J) = WA2(J)
  835. IF (WA2(J) .EQ. ZERO) DIAG(J) = ONE
  836. 50 CONTINUE
  837. 60 CONTINUE
  838. C
  839. C ON THE FIRST ITERATION, CALCULATE THE NORM OF THE SCALED X
  840. C AND INITIALIZE THE STEP BOUND DELTA.
  841. C
  842. DO 70 J = 1, N
  843. WA3(J) = DIAG(J)*X(J)
  844. 70 CONTINUE
  845. XNORM = DENORM(N,WA3)
  846. DELTA = FACTOR*XNORM
  847. IF (DELTA .EQ. ZERO) DELTA = FACTOR
  848. 80 CONTINUE
  849. C
  850. C COMPUTE THE NORM OF THE SCALED GRADIENT.
  851. C
  852. GNORM = ZERO
  853. IF (FNORM .EQ. ZERO) GO TO 170
  854. DO 160 J = 1, N
  855. L = IPVT(J)
  856. IF (WA2(L) .EQ. ZERO) GO TO 150
  857. SUM = ZERO
  858. DO 140 I = 1, J
  859. SUM = SUM + FJAC(I,J)*(QTF(I)/FNORM)
  860. 140 CONTINUE
  861. GNORM = MAX(GNORM,ABS(SUM/WA2(L)))
  862. 150 CONTINUE
  863. 160 CONTINUE
  864. 170 CONTINUE
  865. C
  866. C TEST FOR CONVERGENCE OF THE GRADIENT NORM.
  867. C
  868. IF (GNORM .LE. GTOL) INFO = 4
  869. IF (INFO .NE. 0) GO TO 300
  870. C
  871. C RESCALE IF NECESSARY.
  872. C
  873. IF (MODE .EQ. 2) GO TO 190
  874. DO 180 J = 1, N
  875. DIAG(J) = MAX(DIAG(J),WA2(J))
  876. 180 CONTINUE
  877. 190 CONTINUE
  878. C
  879. C BEGINNING OF THE INNER LOOP.
  880. C
  881. 200 CONTINUE
  882. C
  883. C DETERMINE THE LEVENBERG-MARQUARDT PARAMETER.
  884. C
  885. CALL DMPAR(N,FJAC,LDFJAC,IPVT,DIAG,QTF,DELTA,PAR,WA1,WA2,
  886. 1 WA3,WA4)
  887. C
  888. C STORE THE DIRECTION P AND X + P. CALCULATE THE NORM OF P.
  889. C
  890. DO 210 J = 1, N
  891. WA1(J) = -WA1(J)
  892. WA2(J) = X(J) + WA1(J)
  893. WA3(J) = DIAG(J)*WA1(J)
  894. 210 CONTINUE
  895. PNORM = DENORM(N,WA3)
  896. C
  897. C ON THE FIRST ITERATION, ADJUST THE INITIAL STEP BOUND.
  898. C
  899. IF (ITER .EQ. 1) DELTA = MIN(DELTA,PNORM)
  900. C
  901. C EVALUATE THE FUNCTION AT X + P AND CALCULATE ITS NORM.
  902. C
  903. IFLAG = 1
  904. CALL FCN(IFLAG,M,N,WA2,WA4,FJAC,IJUNK)
  905. NFEV = NFEV + 1
  906. IF (IFLAG .LT. 0) GO TO 300
  907. FNORM1 = DENORM(M,WA4)
  908. C
  909. C COMPUTE THE SCALED ACTUAL REDUCTION.
  910. C
  911. ACTRED = -ONE
  912. IF (P1*FNORM1 .LT. FNORM) ACTRED = ONE - (FNORM1/FNORM)**2
  913. C
  914. C COMPUTE THE SCALED PREDICTED REDUCTION AND
  915. C THE SCALED DIRECTIONAL DERIVATIVE.
  916. C
  917. DO 230 J = 1, N
  918. WA3(J) = ZERO
  919. L = IPVT(J)
  920. TEMP = WA1(L)
  921. DO 220 I = 1, J
  922. WA3(I) = WA3(I) + FJAC(I,J)*TEMP
  923. 220 CONTINUE
  924. 230 CONTINUE
  925. TEMP1 = DENORM(N,WA3)/FNORM
  926. TEMP2 = (SQRT(PAR)*PNORM)/FNORM
  927. PRERED = TEMP1**2 + TEMP2**2/P5
  928. DIRDER = -(TEMP1**2 + TEMP2**2)
  929. C
  930. C COMPUTE THE RATIO OF THE ACTUAL TO THE PREDICTED
  931. C REDUCTION.
  932. C
  933. RATIO = ZERO
  934. IF (PRERED .NE. ZERO) RATIO = ACTRED/PRERED
  935. C
  936. C UPDATE THE STEP BOUND.
  937. C
  938. IF (RATIO .GT. P25) GO TO 240
  939. IF (ACTRED .GE. ZERO) TEMP = P5
  940. IF (ACTRED .LT. ZERO)
  941. 1 TEMP = P5*DIRDER/(DIRDER + P5*ACTRED)
  942. IF (P1*FNORM1 .GE. FNORM .OR. TEMP .LT. P1) TEMP = P1
  943. DELTA = TEMP*MIN(DELTA,PNORM/P1)
  944. PAR = PAR/TEMP
  945. GO TO 260
  946. 240 CONTINUE
  947. IF (PAR .NE. ZERO .AND. RATIO .LT. P75) GO TO 250
  948. DELTA = PNORM/P5
  949. PAR = P5*PAR
  950. 250 CONTINUE
  951. 260 CONTINUE
  952. C
  953. C TEST FOR SUCCESSFUL ITERATION.
  954. C
  955. IF (RATIO .LT. P0001) GO TO 290
  956. C
  957. C SUCCESSFUL ITERATION. UPDATE X, FVEC, AND THEIR NORMS.
  958. C
  959. DO 270 J = 1, N
  960. X(J) = WA2(J)
  961. WA2(J) = DIAG(J)*X(J)
  962. 270 CONTINUE
  963. DO 280 I = 1, M
  964. FVEC(I) = WA4(I)
  965. 280 CONTINUE
  966. XNORM = DENORM(N,WA2)
  967. FNORM = FNORM1
  968. ITER = ITER + 1
  969. 290 CONTINUE
  970. C
  971. C TESTS FOR CONVERGENCE.
  972. C
  973. IF (ABS(ACTRED) .LE. FTOL .AND. PRERED .LE. FTOL
  974. 1 .AND. P5*RATIO .LE. ONE) INFO = 1
  975. IF (DELTA .LE. XTOL*XNORM) INFO = 2
  976. IF (ABS(ACTRED) .LE. FTOL .AND. PRERED .LE. FTOL
  977. 1 .AND. P5*RATIO .LE. ONE .AND. INFO .EQ. 2) INFO = 3
  978. IF (INFO .NE. 0) GO TO 300
  979. C
  980. C TESTS FOR TERMINATION AND STRINGENT TOLERANCES.
  981. C
  982. IF (NFEV .GE. MAXFEV) INFO = 5
  983. IF (ABS(ACTRED) .LE. EPSMCH .AND. PRERED .LE. EPSMCH
  984. 1 .AND. P5*RATIO .LE. ONE) INFO = 6
  985. IF (DELTA .LE. EPSMCH*XNORM) INFO = 7
  986. IF (GNORM .LE. EPSMCH) INFO = 8
  987. IF (INFO .NE. 0) GO TO 300
  988. C
  989. C END OF THE INNER LOOP. REPEAT IF ITERATION UNSUCCESSFUL.
  990. C
  991. IF (RATIO .LT. P0001) GO TO 200
  992. C
  993. C END OF THE OUTER LOOP.
  994. C
  995. GO TO 30
  996. 300 CONTINUE
  997. C
  998. C TERMINATION, EITHER NORMAL OR USER IMPOSED.
  999. C
  1000. IF (IFLAG .LT. 0) INFO = IFLAG
  1001. IFLAG = 0
  1002. IF (NPRINT .GT. 0) CALL FCN(IFLAG,M,N,X,FVEC,FJAC,IJUNK)
  1003. IF (INFO .LT. 0) CALL XERMSG ('SLATEC', 'DNLS1',
  1004. + 'EXECUTION TERMINATED BECAUSE USER SET IFLAG NEGATIVE.', 1, 1)
  1005. IF (INFO .EQ. 0) CALL XERMSG ('SLATEC', 'DNLS1',
  1006. + 'INVALID INPUT PARAMETER.', 2, 1)
  1007. IF (INFO .EQ. 4) CALL XERMSG ('SLATEC', 'DNLS1',
  1008. + 'THIRD CONVERGENCE CONDITION, CHECK RESULTS BEFORE ACCEPTING.',
  1009. + 1, 1)
  1010. IF (INFO .EQ. 5) CALL XERMSG ('SLATEC', 'DNLS1',
  1011. + 'TOO MANY FUNCTION EVALUATIONS.', 9, 1)
  1012. IF (INFO .GE. 6) CALL XERMSG ('SLATEC', 'DNLS1',
  1013. + 'TOLERANCES TOO SMALL, NO FURTHER IMPROVEMENT POSSIBLE.', 3, 1)
  1014. RETURN
  1015. C
  1016. C LAST CARD OF SUBROUTINE DNLS1.
  1017. C
  1018. END