NLSL
lmnls.f90
Go to the documentation of this file.
1 c NLSL Version 1.5 beta 11/25/95
2 c----------------------------------------------------------------------
3 c =========================
4 c subroutine LMNLS
5 c =========================
6 c> @brief (L)evenberg-(M)arquardt (N)onlinear (L)east (S)quares
7 c> ============================================================
8 c> This is a modification of the original lmder subroutine from the
9 c> MINPACK subroutine library. It uses a Levenberg-Marquardt nonlinear
10 c> least-squares algorithm modified to carry out a local optimization
11 c> constrained to lie within a "trust region" defined by a step bound
12 c> delta using scaling of the variables.
13 c> @details
14 c> For a description of the trust region approach for least squares
15 c> problems, see J.E. Dennis and R.B. Schnabel, Numerical Methods for
16 c> Unconstrained Optimization and Nonlinear Equations, Prentice-Hall,
17 c> Englewood Cliffs, NJ (1983), sections 6.4, 7.1, and 10.2.
18 c----------------------------------------------------------------------
19  subroutine lmnls(fcn,m,n,x,fvec,fjac,ldfjac,ftol,xtol,gtol,
20  * maxfev,maxitr,diag,scale,factor,nprint,
21  * info,nfev,njev,ipvt,qtf,gnvec,gradf,
22  * wa1,wa2,wa3,wa4)
23 c
24 c ----added for EPR NLS -----
25 c (also note that fnorm and iter have been moved to common /iterat/)
26  use rnddbl
27  use nlsdim
28  use parcom
29  use iterat
30  use stdio
31 c
32  integer m,n,ldfjac,maxfev,maxitr,nprint,info,istep,nfev,njev
33  integer ipvt(n)
34  double precision ftol,xtol,gtol,factor
35  double precision x(n),fvec(m),fjac(ldfjac,njcol),diag(njcol),
36  * scale(njcol),qtf(njcol),gnvec(njcol),gradf(njcol),
37  * wa1(njcol),wa2(njcol),wa3(njcol),wa4(m)
38  external fcn
39 c
40 c----------------------------------------------------------------------
41 c
42 c> @details
43 c> The purpose of LMDER is to minimize the sum of the squares of
44 c> m nonlinear functions in n variables by a modification of
45 c> the Levenberg-Marquardt algorithm. The user must provide a
46 c> subroutine which calculates the functions and the Jacobian.
47 c>
48 c> The subroutine statement is
49 c>
50 c> subroutine lmnls(fcn,m,n,x,fvec,fjac,ldfjac,ftol,xtol,gtol,
51 c> maxfev,maxitr,diag,scale,factor,nprint,info,
52 c> nfev,njev,ipvt,qtf,wa1,wa2,wa3,wa4)
53 c
54 c where
55 c
56 c>@param FCN is the name of the user-supplied subroutine which
57 c> calculates the functions and the Jacobian. FCN must
58 c> be declared in an external statement in the user
59 c> calling program, and should be written as follows:
60 c>>
61 c>> subroutine fcn(m,n,x,fvec,fjac,ldfjac,iflag)
62 c>> integer m,n,ldfjac,iflag
63 c>> double precision x(n),fvec(m),fjac(ldfjac,n)
64 c>> ----------
65 c>> If iflag=1 calculate the functions at x and
66 c>> return this vector in fvec. Do not alter fjac.
67 c>> If iflag=2 calculate the Jacobian at x and
68 c>> return this matrix in fjac. Do not alter fvec.
69 c>> ----------
70 c>> return
71 c>> end
72 c>>
73 c> The value of IFLAG should not be changed by FCN unless
74 c> the user wants to terminate execution of LMDER.
75 c> In this case set iflag to a negative integer.
76 c>
77 c>@param M is a positive integer input variable set to the number
78 c> of functions.
79 c>
80 c>@param N is a positive integer input variable set to the number
81 c> of variables. N must not exceed M.
82 c>
83 c>@param X is an array of length N. On input X must contain
84 c> an initial estimate of the solution vector. On output X
85 c> contains the final estimate of the solution vector.
86 c>
87 c>@param FVEC is an output array of length M which contains
88 c> the functions evaluated at the output X.
89 c>
90 c>@param FJAC is an output M by N array. the upper N by N submatrix
91 c> of FJAC contains an upper triangular matrix R with
92 c> diagonal elements of nonincreasing magnitude such that
93 c>>
94 c>> T T T
95 c>> P *(JAC *JAC)*P = R *R,
96 c>>
97 c> where P is a permutation matrix and JAC is the final
98 c> calculated Jacobian. column j of P is column IPVT(j)
99 c> (see below) of the identity matrix. The lower trapezoidal
100 c> part of FJAC contains information generated during
101 c> the computation of R.
102 c>
103 c>@param LDFJAC is a positive integer input variable not less than M
104 c> which specifies the leading dimension of the array FJAC.
105 c>
106 c>@param FTOL is a nonnegative input variable. Termination
107 c> occurs when both the actual and predicted relative
108 c> reductions in the sum of squares are at most FTOL.
109 c> Therefore, FTOL measures the relative error desired
110 c> in the sum of squares.
111 c>
112 c>@param XTOL is a nonnegative input variable. Termination
113 c> occurs when the relative error between two consecutive
114 c> iterates is at most XTOL. Therefore, XTOL measures the
115 c> relative error desired in the approximate solution.
116 c>
117 c>@param GTOL is a nonnegative input variable. Termination
118 c> occurs when the cosine of the angle between FVEC and
119 c> any column of the Jacobian is at most GTOL in absolute
120 c> value. therefore, GTOL measures the orthogonality
121 c> desired between the function vector and the columns
122 c> of the Jacobian.
123 c>
124 c>@param MAXFEV is a positive integer input variable. Termination
125 c> occurs when the number of calls to FCN with IFLAG=1
126 c> has reached MAXFEV.
127 c>
128 c>@param SCALE is an array of length N containing multiplicative scale
129 c> factors for each of the variables in X. If an element of SCALE
130 c> is non-positive, it will be reset internally to unity.
131 c> Positive entries in the SCALE array will be retained as
132 c> user-specified scaling factors for the trust-region search
133 c> of the algorithm. The step for the Ith parameter will be scaled
134 c> using the norm of the Ith column of the Jacobian *divided* by
135 c> SCALE(I). This produces larger steps along the ith dimension,
136 c> at least initialy. The default value for all parameters is unity
137 c> (i.e., the column norms of the Jacobian will be used).
138 c>> NB: This convention differs from the original
139 c>> specifications of LMDER in MINPACK.
140 c>
141 c>@param FACTOR is a positive input variable used in determining the
142 c> initial trust region bound. This bound is set to the product of
143 c> FACTOR and the Euclidean norm of DIAG*X if nonzero, or else
144 c> to FACTOR itself. In most cases FACTOR should lie in the
145 c> interval (.1,100.). 100 is a generally recommended value.
146 c>
147 c>@param NPRINT is an integer input variable that enables controlled
148 c> printing of iterates if it is positive. In this case,
149 c> fcn is called with IFLAG=0 at the beginning of the first
150 c> iteration and every NPRINT iterations thereafter and
151 c> immediately prior to return, with X, FVEC, and FJAC
152 c> available for printing. FVEC and FJAC should not be
153 c> altered. If NPRINT is not positive, no special calls
154 c> of FCN with IFLAG=0 are made.
155 c>
156 c>@param INFO @parblock is an integer output variable. If the user has
157 c> terminated execution, INFFO is set to the (negative)
158 c> value of IFLAG. See description of FCN. Otherwise,
159 c> INFO is set as follows:
160 c>
161 c> INFO=0 Improper input parameters.
162 c>
163 c> INFO=1 Both actual and predicted relative reductions
164 c> in the sum of squares are at most FTOL.
165 c>
166 c> INFO=2 Relative error between two consecutive iterates
167 c> is at most XTOL.
168 c>
169 c> INFO=3 conditions for INFO=1 and INFO=2 both hold.
170 c>
171 c> INFO=4 The cosine of the angle between FVEC and any
172 c> column of the Jacobian is at most GTOL in
173 c> absolute value.
174 c>
175 c> INFO=5 number of calls to FCN with IFLAG=1 has
176 c> reached MAXFEV.
177 c>
178 c> INFO=6 FTOL is too small. No further reduction in
179 c> the sum of squares is possible.
180 c>
181 c> INFO=7 XTOL is too small. No further improvement in
182 c> the approximate solution X is possible.
183 c>
184 c> INFO=8 GTOL is too small. FVEC is orthogonal to the
185 c> columns of the Jacobian to machine precision.
186 c> @endparblock
187 c>@param NFEV is an integer output variable set to the number of
188 c> calls to FCN with IFLAG=1.
189 c>
190 c>@param NJEV is an integer output variable set to the number of
191 c> calls to NJEV with IFLAG=2.
192 c>
193 c>@param IPVT is an integer output array of length N. IPVT
194 c> defines a permutation matrix p such that JAC*P=Q*R,
195 c> where JAC is the final calculated Jacobian, Q is
196 c> orthogonal (not stored), and R is upper triangular
197 c> with diagonal elements of nonincreasing magnitude.
198 c> column j of P is column IPVT(j) of the identity matrix.
199 c>
200 c>@param QTF is an output array of length N which contains
201 c> the first N elements of the vector (Q transpose)*FVEC.
202 c>
203 c>@param WA1, WA2, and WA3 are work arrays of length N.
204 c>
205 c>@param WA4 is a work array of length M.
206 c
207 c Subprograms called
208 c
209 c User-supplied ...... FCN
210 c
211 c MINPACK-supplied ... DPMPAR,ENORM,LMPAR,QRFAC
212 c
213 c FORTRAN-supplied ... DABS,DMAX1,DMIN1,DSQRT,MOD
214 c
215 c Argonne National Laboratory. MINPACK Project. March 1980.
216 c Burton S. Garbow, Kenneth E. Hillstrom, Jorge J. More
217 c
218 c----------------------------------------------------------------------
219 cc
220  integer i,iflag,j,l
221  double precision actred,delta,dirder,epsmch,fnorm1,gnorm,par,
222  * pnorm,prered,ratio,sum,temp,temp1,temp2,xnorm
223  parameter(epsmch=rndoff)
224  double precision enorm
225 c
226  double precision ONE,P1,P5,P25,P75,P0001,ZERO
227 c
228  character*1 trstr
229  logical grdclc
230 c
231  integer itrim
232  external itrim
233 c ---------------------------
234  data one,p1,p5,p25,p75,p0001,zero
235  * /1.0d0,1.0d-1,5.0d-1,2.5d-1,7.5d-1,1.0d-4,0.0d0/
236 c
237 c######################################################################
238 c
239 c epsmch is the machine precision.
240 c
241 c
242  info=0
243  iflag=0
244  nfev=0
245  njev=0
246 c
247 c----------------------------------------------------------------------
248 c Check the input parameters for errors.
249 c----------------------------------------------------------------------
250  if (n.le.0 .or. m.lt.n .or. ldfjac.lt.m
251  * .or. ftol.lt.zero .or. xtol.lt.zero .or. gtol.lt.zero
252  * .or. maxfev.le.0 .or.maxitr.le.0 .or. factor.le.zero)
253  * go to 300
254 c
255 c**********************************************************************
256 c
257 c----------------------------------------------------------------------
258 c Initialize Levenberg-Marquardt parameter and iteration counter
259 c----------------------------------------------------------------------
260  par=zero
261  iter=1
262 c
263 c----------------------------------------------------------------------
264 c Evaluate the function at the starting point
265 c and calculate its norm.
266 c----------------------------------------------------------------------
267  iflag=1
268  call fcn(m,n,x,fvec,fjac,ldfjac,iflag)
269  nfev=1
270  if (iflag.lt.0) go to 300
271  fnorm=enorm(m,fvec)
272 c
273 c----------------------------------------------------------------------
274 c ********* Beginning of the outer loop *******************************
275 c----------------------------------------------------------------------
276 c
277  30 continue
278 c
279 c----------------------------------------------------------------------
280 c Calculate the Jacobian matrix (iflag=2)
281 c----------------------------------------------------------------------
282  iflag=2
283  call fcn(m,n,x,fvec,fjac,ldfjac,iflag)
284  njev=njev+1
285  nfev=nfev+n*njev
286  if (iflag.lt.0) go to 300
287 c
288 c----------------------------------------------------------------------
289 c If requested, call fcn to enable printing of iterates
290 c----------------------------------------------------------------------
291  if (nprint.gt.0) then
292  iflag=0
293  if (mod(iter-1,nprint).eq.0)
294  * call fcn(m,n,x,fvec,fjac,ldfjac,iflag)
295  if (iflag.lt.0) go to 300
296  end if
297 c
298 c----------------------------------------------------------------------
299 c Compute the QR factorization of the Jacobian
300 c----------------------------------------------------------------------
301 
302  call qrfac(m,njcol,fjac,ldfjac,.true.,ipvt,n,wa1,wa2,wa3)
303 c
304 c----------------------------------------------------------------------
305 c On the first iteration, set each non-positive element of the
306 c SCALE scaling array according to the norms of the columns of
307 c the initial Jacobian
308 c----------------------------------------------------------------------
309  if (iter.eq.1) then
310  do j=1, n
311  if (scale(j).le.zero) then
312  diag(j)=wa2(j)
313  else
314  diag(j)=wa2(j)/scale(j)
315  end if
316  if (diag(j).eq.zero) diag(j)=one
317  end do
318 c
319  if (itrace.ne.0) then
320  write(itrace,1000) (tag(j)(:itrim(tag(j))),j=1,n)
321  write(itrace,1001) (wa2(j),j=1,n)
322  write(itrace,1002) (diag(j),j=1,n)
323  end if
324 c
325 c----------------------------------------------------------------------
326 c On the first iteration, calculate the norm of the scaled x
327 c and initialize the trust region bound delta
328 c----------------------------------------------------------------------
329  do j=1, n
330  wa3(j)=diag(j)*x(j)
331  end do
332  xnorm=enorm(n,wa3)
333  delta=factor*xnorm
334  if (delta.eq.zero) delta=factor
335  if (itrace.ne.0) write(itrace,1003) xnorm,delta,factor
336  end if
337 c
338 c----------------------------------------------------------------------
339 c Form (Q transpose)*fvec and store the first njcol components in
340 c QtF.
341 c----------------------------------------------------------------------
342  do i=1, m
343  wa4(i)=fvec(i)
344  end do
345 c
346  do j=1, njcol
347  if (fjac(j,j).ne.zero) then
348  sum=zero
349  do i=j, m
350  sum=sum + fjac(i,j)*wa4(i)
351  end do
352  temp=-sum/fjac(j,j)
353  do i=j, m
354  wa4(i)=wa4(i) + fjac(i,j)*temp
355  end do
356  end if
357 c
358  fjac(j,j)=wa1(j)
359  qtf(j)=wa4(j)
360  end do
361 c
362 c----------------------------------------------------------------------
363 c Compute the norm of the scaled gradient.
364 c----------------------------------------------------------------------
365  gnorm=zero
366  if (fnorm.ne.zero) then
367  do j=1, n
368  l=ipvt(j)
369  if (wa2(l).ne.zero) then
370  sum=zero
371  do i=1, j
372  sum=sum + fjac(i,j)*(qtf(i)/fnorm)
373  end do
374  gnorm=dmax1(gnorm,dabs(sum/wa2(l)))
375  end if
376  end do
377  end if
378 c
379 c----------------------------------------------------------------------
380 c Test for convergence of the gradient norm
381 c----------------------------------------------------------------------
382  if (gnorm.le.gtol) info=4
383  if (info.ne.0) go to 300
384 c
385 c----------------------------------------------------------------------
386 c Rescale diag array
387 c----------------------------------------------------------------------
388  do j=1,n
389  if (scale(j).gt.zero) then
390  temp=wa2(j)/scale(j)
391  diag(j)=dmax1(diag(j),temp)
392  end if
393  end do
394 c
395 c----------------------------------------------------------------------
396 c ******** Beginning of the inner loop ******************************
397 c----------------------------------------------------------------------
398  istep=0
399  grdclc=.false.
400  200 continue
401 c
402 c----------------------------------------------------------------------
403 c Determine the Levenberg-Marquardt parameter.
404 c----------------------------------------------------------------------
405  call lmpar(n,fjac,ldfjac,ipvt,diag,qtf,delta,par,wa1,wa2,
406  * wa3,wa4,gnvec,gradf)
407 c
408  grdclc=grdclc.or.(par.ne.zero)
409 c----------------------------------------------------------------------
410 c Store the direction p and X + p. Calculate the norm of p.
411 c----------------------------------------------------------------------
412  do j=1, n
413  wa1(j)=-wa1(j)
414  wa2(j)=x(j) + wa1(j)
415  wa3(j)=diag(j)*wa1(j)
416  end do
417  pnorm=enorm(n,wa3)
418 c
419 c----------------------------------------------------------------------
420 c On the first iteration, adjust the initial trust region bound
421 c to the size of the initial step.
422 c----------------------------------------------------------------------
423  trstr=' '
424  if (iter.eq.1) then
425  if (delta.gt.pnorm) then
426  trstr='*'
427  else
428  trstr='s'
429  end if
430  delta=dmin1(delta,pnorm)
431  end if
432 c
433  if (istep.eq.0 .and. itrace.ne.0) then
434  write (itrace,1012) iter,(tag(j),j=1,n)
435  write (itrace,1013) fnorm,(x(j),j=1,n)
436  end if
437 c
438 c----------------------------------------------------------------------
439 c Evaluate the function at x + p and calculate its norm.
440 c----------------------------------------------------------------------
441  iflag=1
442  call fcn(m,n,wa2,wa4,fjac,ldfjac,iflag)
443  nfev=nfev+1
444  if (iflag.lt.0) go to 300
445  fnorm1=enorm(m,wa4)
446  istep=istep+1
447 c
448 c----------------------------------------------------------------------
449 c Compute the scaled actual reduction.
450 c----------------------------------------------------------------------
451  actred=-one
452  if (p1*fnorm1.lt.fnorm) actred=one-(fnorm1/fnorm)**2
453 c
454 c----------------------------------------------------------------------
455 c Compute the scaled predicted reduction and
456 c the scaled directional derivative.
457 c----------------------------------------------------------------------
458  do j=1,n
459  wa3(j)=zero
460  l=ipvt(j)
461  temp=wa1(l)
462  do i=1, j
463  wa3(i)=wa3(i) + fjac(i,j)*temp
464  end do
465  end do
466  temp1=enorm(n,wa3)/fnorm
467  temp2=(dsqrt(par)*pnorm)/fnorm
468  prered=temp1**2 + temp2**2/p5
469  dirder=-(temp1**2 + temp2**2)
470 c
471 c----------------------------------------------------------------------
472 c Compute the ratio of the actual to the predicted
473 c reduction.
474 c----------------------------------------------------------------------
475  ratio=zero
476  if (prered.ne.zero) ratio=actred/prered
477 c
478 c----------------------------------------------------------------------
479 c Update the step bound.
480 c----------------------------------------------------------------------
481 c
482 c----------------------------------------------------------------------
483 c If actual reduction is too much smaller than the predicted
484 c reduction (i.e. actred/prered ratio is too small)
485 c the function is not well-approximated by a quadratic
486 c equation. Reduce the size of the trust region by a
487 c factor of 0.1 to 0.5 and increase the L-M parameter.
488 c----------------------------------------------------------------------
489  if (ratio.le.p25) then
490  if (actred.ge.zero) temp=p5
491  if (actred.lt.zero) temp=p5*dirder/(dirder + p5*actred)
492  if (p1*fnorm1.ge.fnorm .or. temp.lt.p1) temp=p1
493  if (delta.gt.pnorm/p1) then
494  trstr='x'
495  delta=pnorm/p1
496  endif
497  delta=delta*temp
498  par=par/temp
499 c
500  else
501 c
502 c----------------------------------------------------------------------
503 c If ratio of actual to predicted reduction is close to 1,
504 c the quadratic model is a good approximation to the function,
505 c and we can try increasing the trust region to twice the
506 c last step size in order to check whether a better solution
507 c is available. Otherwise, the size of the trust region
508 c is left unchanged.
509 c----------------------------------------------------------------------
510  if (par.eq.zero .or. ratio.ge.p75) then
511  delta=pnorm/p5
512  par=p5*par
513  temp=one/p5
514  else
515  temp=one
516  end if
517  end if
518 c
519  if (itrace.ne.0) write (itrace,1014) istep,par,ratio,
520  # one/temp,trstr,fnorm1,(wa2(j)-x(j),j=1,n)
521 c
522 c----------------------------------------------------------------------
523 c Test for successful iteration.
524 c----------------------------------------------------------------------
525  if (ratio.ge.p0001) then
526 c
527 c----------------------------------------------------------------------
528 c Successful iteration. Update X, FVEC, and their norms.
529 c----------------------------------------------------------------------
530  do j=1, n
531  x(j)=wa2(j)
532  wa2(j)=diag(j)*x(j)
533  end do
534 c
535  do i=1, m
536  fvec(i)=wa4(i)
537  end do
538  xnorm=enorm(n,wa2)
539  fnorm=fnorm1
540  iter=iter+1
541  if (itrace.ne.0) then
542  write(itrace,1006) (diag(j),j=1,n)
543  if (grdclc) write(itrace,1007) (gradf(j),j=1,n)
544  write(itrace,1008) (gnvec(j),j=1,n)
545  write (itrace,1015) delta
546  end if
547 c
548  end if
549 c----------------------------------------------------------------------
550 c Tests for convergence.
551 c----------------------------------------------------------------------
552  info=0
553  if (dabs(actred).le.ftol .and. prered.le.ftol
554  * .and. p5*ratio.le.one) info=1
555  if (delta.le.xtol*xnorm) info=info + 2
556  if (info.ne.0) go to 300
557 c
558 c----------------------------------------------------------------------
559 c Tests for termination and stringent tolerances.
560 c----------------------------------------------------------------------
561  if (nfev.ge.maxfev) info=5
562  if (iter.ge.maxitr) info=6
563  if (dabs(actred).le.epsmch .and. prered.le.epsmch
564  * .and. p5*ratio.le.one) info=7
565  if (delta.le.epsmch*xnorm) info=8
566  if (gnorm.le.epsmch) info=9
567  if (info.ne.0) go to 300
568 c----------------------------------------------------------------------
569 c End of the inner loop. Repeat if iteration unsuccessful.
570 c----------------------------------------------------------------------
571  if (ratio.lt.p0001) go to 200
572 c----------------------------------------------------------------------
573 c End of the outer loop.
574 c----------------------------------------------------------------------
575  go to 30
576  300 continue
577 c----------------------------------------------------------------------
578 c Termination, either normal or user-imposed.
579 c----------------------------------------------------------------------
580  if (iflag.lt.0) info=10
581  iflag=3
582  if (info.gt.4) iflag=-iflag
583  if (info.ne.0.and.info.ne.10.and.nprint.gt.0)
584  # call fcn(m,n,x,fvec,fjac,ldfjac,iflag)
585  return
586 c
587 c ##### format statements for trace printout ########################
588 c
589  1000 format(/10x,41('=')/10x,
590  # 'TRACE OF LEVENBERG-MARQUARDT MINIMIZATION'/
591  # 10x,41('=')//'INITIAL SCALING:'/13x,10(1x,a9,2x))
592  1001 format('Col norms of J:',10(2x,g10.4)/)
593  1002 format(9x,'Scale:',10(2x,g10.4))
594  1003 format(/10x,'Scaled X norm: ',g11.5/6x,'Trust region (TR):',
595  # g11.5,' =(Xnorm*',g9.3,')')
596  1006 format(79('-')/t26,'Scale:',10(1x,g10.4))
597  1007 format(t23,'Gradient:',10(1x,g10.4))
598  1008 format(t21,'G-N vector:',10(1x,g10.4))
599  1012 format(/'##### Iteration',i3,1x,59('#')/
600  #'Stp LMpar Ratio Trscl Fnorm ',12(2x,a9))
601  1013 format(79('-')/' 0',t22,g11.5,12g11.4)
602  1014 format(i3,2f6.2,f5.2,a1,g11.5,sp,12g11.4)
603  1015 format(t23,'Final TR: ',g10.4)
604  end
char info[81]
Definition: genio.c:45
subroutine qrfac(m, n, a, lda, pivot, ipvt, lipvt, rdiag, acnorm, wa)
Definition: qrfac.f90:2
integer, save iter
Definition: iterat.f90:14
Definition: stdio.f90:26
double precision, save fnorm
Definition: iterat.f90:15
subroutine lmpar(n, r, ldr, ipvt, diag, qtf, delta, par, x, sdiag, wa1, wa2, gnvec, gradf)
Definition: lmpar.f90:3
double factor
Definition: genio.c:44
subroutine lmnls(fcn, m, n, x, fvec, fjac, ldfjac, ftol, xtol, gtol, maxfev, maxitr, diag, scale, factor, nprint, info, nfev, njev, ipvt, qtf, gnvec, gradf, wa1, wa2, wa3, wa4)
(L)evenberg-(M)arquardt (N)onlinear (L)east (S)quares This is a modification of the original lmder su...
Definition: lmnls.f90:23
integer, save itrace
Definition: stdio.f90:32
character *9, dimension(mxjcol), save tag
Definition: parcom.f90:67
double precision, parameter rndoff
Definition: rnddbl.f90:86
integer, save njcol
Definition: parcom.f90:62