NLSL
lmpar.f90
Go to the documentation of this file.
1  subroutine lmpar(n,r,ldr,ipvt,diag,qtf,delta,par,x,sdiag,wa1,
2  * wa2,gnvec,gradf)
3  use rnddbl
4  integer n,ldr
5  integer ipvt(n)
6  double precision delta,par
7  double precision r(ldr,n),diag(n),qtf(n),x(n),sdiag(n),wa1(n),
8  * wa2(n),gnvec(n), gradf(n)
9 c **********
10 c
11 c Subroutine lmpar
12 c
13 c Given an m by n matrix a, an n by n nonsingular diagonal
14 c matrix d, an m-vector b, and a positive number delta,
15 c the problem is to determine a value for the parameter
16 c par such that if x solves the system
17 c
18 c A*x = f , sqrt(par)*D*x = 0 ,
19 c
20 c in the least squares sense, and dxnorm is the Euclidean
21 c norm of D*x, then either par is zero and
22 c
23 c (dxnorm-delta) .le. 0.1*delta ,
24 c
25 c or par is positive and
26 c
27 c abs(dxnorm-delta) .le. 0.1*delta .
28 c
29 c This subroutine completes the solution of the problem
30 c if it is provided with the necessary information from the
31 c QR factorization, with column pivoting, of A. That is, if
32 c A*P = Q*R, where P is a permutation matrix, Q has orthogonal
33 c columns, and R is an upper triangular matrix with diagonal
34 c elements of nonincreasing magnitude, then lmpar expects
35 c the full upper triangle of R, the permutation matrix P,
36 c and the first n components of (Q transpose)*b. On output
37 c lmpar also provides an upper triangular matrix S such that
38 c
39 c t t t
40 c P *(A *A + par*D*D)*P = S *S .
41 c
42 c S is employed within lmpar and may be of separate interest.
43 c
44 c Only a few iterations are generally needed for convergence
45 c of the algorithm. If, however, the limit of 10 iterations
46 c is reached, then the output par will contain the best
47 c value obtained so far.
48 c
49 c The subroutine statement is
50 c
51 c subroutine lmpar(n,r,ldr,ipvt,diag,qtf,delta,par,x,sdiag,
52 c wa1,wa2)
53 c
54 c where
55 c
56 c n is a positive integer input variable set to the order of R.
57 c
58 c R is an n by n array. On input the full upper triangle
59 c must contain the full upper triangle of the matrix R.
60 c On output the full upper triangle is unaltered, and the
61 c strict lower triangle contains the strict upper triangle
62 c (transposed) of the upper triangular matrix S.
63 c
64 c ldr is a positive integer input variable not less than n
65 c which specifies the leading dimension of the array r.
66 c
67 c ipvt is an integer input array of length n which defines the
68 c permutation matrix P such that A*P = Q*R. Column j of P
69 c is column ipvt(j) of the identity matrix.
70 c
71 c diag is an input array of length n which must contain the
72 c diagonal elements of the matrix D.
73 c
74 c qtf is an input array of length n which must contain the first
75 c n elements of the vector (Q transpose)*f (f is the vector
76 c containing the function values to be minimized).
77 c
78 c delta is a positive input variable which specifies an upper
79 c bound on the Euclidean norm of D*x.
80 c
81 c par is a nonnegative variable. On input par contains an
82 c initial estimate of the Levenberg-Marquardt parameter.
83 c On output par contains the final estimate.
84 c
85 c x is an output array of length n which contains the least
86 c squares solution of the system A*x = b, sqrt(par)*D*x = 0,
87 c for the output par.
88 c
89 c sdiag is an output array of length n which contains the
90 c diagonal elements of the upper triangular matrix S.
91 c
92 c gnvec is an output array of length n which contains the
93 c Gauss-Newton vector corresponding to the input matrix R.
94 c
95 c gradf is an output array of length n which contains the
96 c gradient (steepest descent) vector. This is only calculated
97 c if the L-M parameter, par, is nonzero.
98 c
99 c wa1 and wa2 are work arrays of length n.
100 c
101 c Subprograms called
102 c
103 c MINPACK-supplied ... enorm,qrsolv
104 c
105 c Fortran-supplied ... dabs,dmax1,dmin1,dsqrt
106 c
107 c Argonne National Laboratory. MINPACK project. March 1980.
108 c Burton S. Garbow, Kenneth E. Hillstrom, Jorge J. More
109 c
110 c **********
111  integer i,iter,j,jm1,jp1,k,l,nsing
112  double precision dxnorm,dwarf,fp,gnorm,parc,parl,paru,sum,temp
113  double precision enorm
114 c
115  double precision P1,P001,ZERO
116  data p1,p001,zero /1.0d-1,1.0d-3,0.0d0/
117 c
118 c dwarf is the smallest positive magnitude.
119  parameter(dwarf=dbl_min)
120 c
121 c
122 c----------------------------------------------------------------------
123 c Compute and store in x the Gauss-Newton direction. If the
124 c Jacobian is rank-deficient, obtain a least squares solution.
125 c Store the scaled Gauss-Newton direction in gnvec.
126 c
127 c The Gauss-Newton direction is the solution of the problem given
128 c above for par=0. The solution to the above least-squares problem
129 c in this case is equivalent to finding the vector
130 c
131 c + T -1 T +
132 c x = -J f = -P R Q f (i.e. J is the pseudoinverse of the Jacobian)
133 c
134 c See exercise 24 in chapter 6 of Dennis and Schnabel.
135 c----------------------------------------------------------------------
136  nsing = n
137  do j = 1, n
138  wa1(j) = qtf(j)
139  if (r(j,j) .eq. zero .and. nsing .eq. n) nsing = j - 1
140  if (nsing .lt. n) wa1(j) = zero
141  end do
142 c
143  if (nsing .ge. 1) then
144  do k = 1, nsing
145  j = nsing - k + 1
146  wa1(j) = wa1(j)/r(j,j)
147  temp = wa1(j)
148  jm1 = j - 1
149  if (jm1 .ge. 1) then
150  do i = 1, jm1
151  wa1(i) = wa1(i) - r(i,j)*temp
152  end do
153  end if
154  end do
155  end if
156 c
157  do j = 1, n
158  l = ipvt(j)
159  x(l) = wa1(j)
160  gnvec(l)= -x(l)
161  end do
162 c
163 c----------------------------------------------------------------------
164 c Initialize the iteration counter. If the (scaled) Gauss-Newton
165 c step size is within 10% of the trust-region bound, terminate the
166 c search (par will be returned as ZERO).
167 c----------------------------------------------------------------------
168  iter = 0
169  do j = 1, n
170  wa2(j) = diag(j)*x(j)
171  end do
172  dxnorm = enorm(n,wa2)
173  fp = dxnorm - delta
174  if (fp .le. p1*delta) go to 220
175 c
176 c----------------------------------------------------------------------
177 c If the Jacobian is not rank deficient, the scaled Gauss-Newton
178 c step provides a lower bound, parl, for the zero of the function.
179 c Otherwise set this bound to zero.
180 c----------------------------------------------------------------------
181  parl = zero
182  if (nsing .ge. n) then
183  do j = 1, n
184  l = ipvt(j)
185  wa1(j) = diag(l)*(wa2(l)/dxnorm)
186  end do
187  do j = 1, n
188  sum = zero
189  jm1 = j - 1
190  if (jm1 .ge. 1) then
191  do i = 1, jm1
192  sum = sum + r(i,j)*wa1(i)
193  end do
194  end if
195  wa1(j) = (wa1(j) - sum)/r(j,j)
196  end do
197  temp = enorm(n,wa1)
198  parl = ((fp/delta)/temp)/temp
199  end if
200 c
201 c----------------------------------------------------------------------
202 c Calculate the steepest descent direction of the residuals at
203 c point where fvec has been evaluated.
204 c
205 c T T T
206 c Grad(f) = J f = P R Q f
207 c
208 c This provides an upper bound, paru, for the zero of the function.
209 c----------------------------------------------------------------------
210  do j = 1, n
211  sum = zero
212  do i = 1, j
213  sum = sum + r(i,j)*qtf(i)
214  end do
215  l = ipvt(j)
216  wa1(j) = sum/diag(l)
217  gradf(l) = wa1(j)
218  end do
219  gnorm = enorm(n,wa1)
220  paru = gnorm/delta
221  if (paru .eq. zero) paru = dwarf/dmin1(delta,p1)
222 c
223 c----------------------------------------------------------------------
224 c If the input par lies outside of the interval (parl,paru),
225 c set par to the closer endpoint.
226 c----------------------------------------------------------------------
227  par = dmax1(par,parl)
228  par = dmin1(par,paru)
229  if (par .eq. zero) par = gnorm/dxnorm
230 c
231 c----------------------------------------------------------------------
232 c *** Beginning of an iteration.
233 c----------------------------------------------------------------------
234  150 continue
235  iter = iter + 1
236 c
237 c----------------------------------------------------------------------
238 c Evaluate the function at the current value of par.
239 c----------------------------------------------------------------------
240  if (par .eq. zero) par = dmax1(dwarf,p001*paru)
241  temp = dsqrt(par)
242  do j = 1, n
243  wa1(j) = temp*diag(j)
244  end do
245  call qrsolv(n,r,ldr,ipvt,wa1,qtf,x,sdiag,wa2)
246  do j = 1, n
247  wa2(j) = diag(j)*x(j)
248  end do
249  dxnorm = enorm(n,wa2)
250  temp = fp
251  fp = dxnorm - delta
252 c
253 c----------------------------------------------------------------------
254 c If the function is small enough, accept the current value
255 c of par. Also test for the exceptional cases where parl
256 c is zero or the number of iterations has reached 10.
257 c----------------------------------------------------------------------
258  if (dabs(fp) .le. p1*delta
259  * .or. parl .eq. zero .and. fp .le. temp
260  * .and. temp .lt. zero .or. iter .eq. 10) go to 220
261 c
262 c----------------------------------------------------------------------
263 c Compute the Newton correction.
264 c----------------------------------------------------------------------
265  do j = 1, n
266  l = ipvt(j)
267  wa1(j) = diag(l)*(wa2(l)/dxnorm)
268  end do
269  do j = 1, n
270  wa1(j) = wa1(j)/sdiag(j)
271  temp = wa1(j)
272  jp1 = j + 1
273  if (n .ge. jp1) then
274  do i = jp1, n
275  wa1(i) = wa1(i) - r(i,j)*temp
276  end do
277  end if
278  end do
279  temp = enorm(n,wa1)
280  parc = ((fp/delta)/temp)/temp
281 c
282 c----------------------------------------------------------------------
283 c Depending on the sign of the function, update parl or paru.
284 c----------------------------------------------------------------------
285  if (fp .gt. zero) parl = dmax1(parl,par)
286  if (fp .lt. zero) paru = dmin1(paru,par)
287 c
288 c----------------------------------------------------------------------
289 c Compute an improved estimate for par.
290 c----------------------------------------------------------------------
291  par = dmax1(parl,par+parc)
292 c
293 c----------------------------------------------------------------------
294 c *** End of an iteration.
295 c----------------------------------------------------------------------
296  go to 150
297  220 continue
298 c
299 c----------------------------------------------------------------------
300 c Termination.
301 c----------------------------------------------------------------------
302  if (iter .eq. 0) par = zero
303 c
304  return
305  end
double precision, parameter dbl_min
Definition: rnddbl.f90:80
subroutine lmpar(n, r, ldr, ipvt, diag, qtf, delta, par, x, sdiag, wa1, wa2, gnvec, gradf)
Definition: lmpar.f90:3
subroutine qrsolv(n, r, ldr, ipvt, diag, qtb, x, sdiag, wa)
Definition: qrsolv.f90:2