OOFiles: bvls.f

File bvls.f, 19.6 kB (added by dmitrey.kroshko, 9 months ago)
Line 
1 c=======================================================================
2       subroutine bvls(key, m, n, a, b, bl, bu, x, w, act, zz, istate,
3      +  loopA, istop, msg)
4 c=======================================================================
5       implicit double precision (a-h, o-z)
6       character*64 msg
7
8 c
9 c$$$$ calls qr
10 c--------------------Bounded Variable Least Squares---------------------
11 c
12 c        Robert L. Parker and Philip B. Stark    Version 3/19/90
13 c
14 c  Robert L. Parker                           Philip B. Stark
15 c  Scripps Institution of Oceanography        Department of Statistics
16 c  University of California, San Diego        University of California
17 c  La Jolla CA 92093                          Berkeley CA 94720-3860
18 c  rlparker@ucsd.edu                          stark@stat.berkeley.edu
19 c
20 c  Copyright of this software is reserved by the authors; however, this
21 c  algorithm and subroutine may be used freely for non-commercial
22 c  purposes, and may be distributed freely for non-commercial purposes. 
23 c
24 c  The authors do not warrant this software in any way: use it at your
25 c  own risk.
26 c
27 c
28 c  See the article ``Bounded Variable Least Squares:  An Algorithm and
29 c  Applications'' by P.B. Stark and R.L. Parker, in the journal
30 c  Computational Statistics, in press (1995) for further description
31 c  and applications to minimum l-1, l-2 and l-infinity fitting problems,
32 c  as well as finding bounds on linear functionals subject to bounds on
33 c  variables and fitting linear data within l-1, l-2 or l-infinity
34 c  measures of misfit.
35 c
36 c  BVLS solves the problem:
37 c
38 c          min  || a.x - b ||     such that   bl <= x <= bu
39 c                            2
40 c    where 
41 c               x  is an unknown n-vector
42 c               a  is a given m by n matrix
43 c               b  is a given  m-vector
44 c               bl is a given n-vector of lower bounds on the
45 c                                components of x.
46 c               bu is a given n-vector of upper bounds on the
47 c                                components of x.
48 c
49 c               
50 c-----------------------------------------------------------------------
51 c    Input parameters:
52 c
53 c  m, n, a, b, bl, bu   see above.   Let mm=min(m,n).
54 c
55 c  If key = 0, the subroutine solves the problem from scratch.
56 c
57 c  If key > 0 the routine initializes using the user's guess about
58 c   which components of  x  are `active', i.e. are stricly within their
59 c   bounds, which are at their lower bounds, and which are at their
60 c   upper bounds.  This information is supplied through the array 
61 c   istate.  istate(n+1) should contain the total number of components
62 c   at their bounds (the `bound variables').  The absolute values of the
63 c   first nbound=istate(n+1) entries of  istate  are the indices
64 c   of these `bound' components of  x.  The sign of istate(j), j=1,...,
65 c   nbound, indicates whether  x(|istate(j)|) is at its upper or lower
66 c   bound.  istate(j) is positive if the component is at its upper
67 c   bound, negative if the component is at its lower bound.
68 c   istate(j), j=nbound+1,...,n  contain the indices of the components
69 c   of  x  that are active (i.e. are expected to lie strictly within
70 c   their bounds).  When key > 0, the routine initially sets the active
71 c   components to the averages of their upper and lower bounds:
72 c   x(j)=(bl(j)+bu(j))/2, for j in the active set. 
73 c
74 c-----------------------------------------------------------------------
75 c    Output parameters:
76 c
77 c  x       the solution vector.
78 c
79 c  w(1)    the minimum 2-norm || a.x-b ||.
80 c
81 c  istate  vector indicating which components of  x  are active and
82 c          which are at their bounds (see the previous paragraph). 
83 c          istate can be supplied to the routine to give it a good
84 c          starting guess for the solution.
85 c
86 c  loopA   number of iterations taken in the main loop, Loop A.
87 c
88 c-----------------------------------------------------------------------
89 c    Working  arrays:
90 c
91 c  w      dimension n.               act      dimension m*(mm+2).
92 c  zz     dimension m.               istate   dimension n+1.
93 c
94 c-----------------------------------------------------------------------
95 c  Method: active variable method along the general plan of NNLS by
96 c  Lawson & Hanson, "Solving Least Squares Problems," 1974.  See
97 c  Algorithm 23.10.  Step numbers in comment statements refer to their
98 c  scheme.
99 c  For more details and further uses, see the article
100 c  "Bounded Variable Least Squares:  An Algorithm and Applications"
101 c  by Stark and Parker in 1995 Computational Statistics.
102 c
103 c-----------------------------------------------------------------------
104 c  A number of measures are taken to enhance numerical reliability:
105 c
106 c 1. As noted by Lawson and Hanson, roundoff errors in the computation
107 c   of the gradient of the misfit may cause a component on the bounds
108 c   to appear to want to become active, yet when the component is added
109 c   to the active set, it moves away from the feasible region.  In this
110 c   case the component is not made active, the gradient of the misfit
111 c   with respect to a change in that component is set to zero, and the
112 c   program returns to the Kuhn-Tucker test.  Flag  ifrom5  is used in
113 c   this test, which occurs at the end of Step 6.
114 c
115 c
116 c 2. When the least-squares minimizer after Step 6 is infeasible, it
117 c   is used in a convex interpolation with the previous solution to
118 c   obtain a feasible vector.  The constant in this interpolation is
119 c   supposed to put at least one component of  x   on a bound. There can
120 c   be difficulties:
121 c
122 c 2a. Sometimes, due to roundoff, no interpolated component ends up on
123 c   a bound.  The code in Step 11 uses the flag  jj, computed in Step 8,
124 c   to ensure that at least the component that determined the
125 c   interpolation constant  alpha  is moved to the appropriate bound. 
126 c   This guarantees that what Lawson and Hanson call `Loop B' is finite.
127 c
128 c 2b. The code in Step 11 also incorporates Lawson and Hanson's feature
129 c   that any components remaining infeasible at this stage (which must
130 c   be due to roundoff) are moved to their nearer bound.
131 c
132 c
133 c 3. If the columns of  a  passed to qr are linearly dependent, the new
134 c   potentially active component is not introduced: the gradient of the
135 c   misfit with respect to that component is set to zero, and control
136 c   returns to the Kuhn-Tucker test.
137 c
138 c
139 c 4. When some of the columns of  a  are approximately linearly
140 c   dependent, we have observed cycling of active components: a
141 c   component just moved to a bound desires immediately to become
142 c   active again; qr allows it to become active and a different
143 c   component is moved to its bound.   This component immediately wants
144 c   to become active, which qr allows, and the original component is
145 c   moved back to its bound.  We have taken two steps to avoid this
146 c   problem:
147 c
148 c 4a. First, the column of the matrix  a  corresponding to the new
149 c   potentially active component is passed to qr as the last column of
150 c   its matrix.  This ordering tends to make a component recently moved
151 c   to a bound fail the test mentioned in (1), above.
152 c
153 c 4b. Second, we have incorporated a test that prohibits short cycles.
154 c   If the most recent successful change to the active set was to move
155 c   the component x(jj) to a bound, x(jj) is not permitted to reenter
156 c   the solution at this stage.  This test occurs just after checking
157 c   the Kuhn-Tucker conditions, and uses the flag  jj, set in Step 8.
158 c   The flag  jj  is reset after Step 6 if Step 6 was entered from
159 c   Step 5 indicating that a new component has successfully entered the
160 c   active set. The test for resetting  jj  uses the flag  ifrom5,
161 c   which will not equal zero in case Step 6 was entered from Step 5.
162 c
163 c
164 c      dimension istop(1)
165       dimension a(m,n), b(m), x(n), bl(n), bu(n)
166 c     dimension w(n), act(m,min(m,n)+2), zz(m), istate(n+1)
167       dimension w(n), act(m,m+2), zz(m), istate(n+1)
168 c
169       data eps/1.0d-11/
170 c
171 c----------------------First Executable Statement-----------------------
172 c
173 c  Step 1.  Initialize everything--active and bound sets, initial
174 c   values, etc.
175 c
176 c  Initialize flags, etc.
177       istop = 1000
178       msg = 'Converged'
179
180
181       mm=min(m,n)
182       mm1 = mm + 1
183       jj = 0
184       ifrom5 = 0
185 c  Check consistency of given bounds  bl, bu.
186       bdiff = 0.0
187       do 1005 j=1, n
188         bdiff=max(bdiff, bu(j)-bl(j))
189         if (bl(j) .gt. bu(j)) then
190           istop = -1000
191           msg = 'Inconsistent bounds in BVLS.'
192           return
193 c          stop
194         endif
195  1005 continue
196       if (bdiff .eq. 0.0) then
197         istop = -1000
198         msg = 'No free variables in BVLS--check input bounds'
199         return
200 c        stop
201       endif
202 c
203 c  In a fresh initialization (key = 0) bind all variables at their lower
204 c   bounds.  If (key != 0), use the supplied  istate  vector to
205 c   initialize the variables.  istate(n+1) contains the number of
206 c   bound variables.  The absolute values of the first
207 c   nbound=istate(n+1) entries of  istate  are the indices of the bound
208 c   variables.  The sign of each entry determines whether the indicated
209 c   variable is at its upper (positive) or lower (negative) bound.
210       if (key .eq. 0) then
211         nbound=n
212         nact=0
213         do 1010 j=1, nbound
214           istate(j)=-j
215  1010   continue
216       else
217         nbound=istate(n+1)
218       endif
219       nact=n - nbound
220       if ( nact .gt. mm ) then
221         istop = -1000
222         msg = 'Too many active variables in BVLS starting solution!'
223         return
224 c        stop
225       endif
226       do 1100 k=1, nbound
227         j=abs(istate(k))
228         if (istate(k) .lt. 0) x(j)=bl(j)
229         if (istate(k) .gt. 0) x(j)=bu(j)
230  1100 continue
231 c
232 c  In a warm start (key != 0) initialize the active variables to
233 c   (bl+bu)/2.  This is needed in case the initial qr results in
234 c   active variables out-of-bounds and Steps 8-11 get executed the
235 c   first time through.
236       do 1150 k=nbound+1,n
237         kk=istate(k)
238         x(kk)=(bu(kk)+bl(kk))/2
239  1150 continue
240 c
241 c  Compute bnorm, the norm of the data vector b, for reference.
242       bsq=0.0
243       do 1200 i=1, m
244         bsq=bsq + b(i)**2
245  1200 continue
246       bnorm=sqrt(bsq)
247 c
248 c-----------------------------Main Loop---------------------------------
249 c
250 c  Initialization complete.  Begin major loop (Loop A).
251       do 15000 loopA=1, 3*n
252 c
253 c  Step 2.
254 c  Initialize the negative gradient vector w(*).
255  2000 obj=0.0
256       do 2050 j=1, n
257         w(j)=0.0
258  2050 continue
259 c
260 c  Compute the residual vector b-a.x , the negative gradient vector
261 c   w(*), and the current objective value obj = || a.x - b ||.
262 c   The residual vector is stored in the mm+1'st column of act(*,*).
263       do 2300 i=1, m
264         ri=b(i)
265         do 2100 j=1, n
266           ri=ri - a(i,j)*x(j)
267  2100   continue
268         obj=obj + ri**2
269         do 2200 j=1, n
270           w(j)=w(j) + a(i,j)*ri
271  2200   continue
272         act(i,mm1)=ri
273  2300 continue
274 c
275 c  Converged?  Stop if the misfit << || b ||, or if all components are
276 c   active (unless this is the first iteration from a warm start).
277       if (sqrt(obj) .le. bnorm*eps .or.
278      +  (loopA .gt. 1 .and. nbound .eq. 0)) then
279          istate(n+1)=nbound
280          w(1)=sqrt(obj)
281          return
282       endif
283 c
284 c  Add the contribution of the active components back into the residual.
285       do 2500 k=nbound+1, n
286         j=istate(k)
287         do 2400 i=1, m
288           act(i,mm1)=act(i,mm1) + a(i,j)*x(j)
289  2400   continue
290  2500 continue
291 c
292 c  The first iteration in a warm start requires immediate qr.
293       if (loopA .eq. 1 .and. key .ne. 0) goto 6000
294 c
295 c  Steps 3, 4.
296 c  Find the bound element that most wants to be active.
297  3000 worst=0.0
298       it=1
299       do 3100 j=1, nbound
300          ks=abs(istate(j))
301          bad=w(ks)*sign(1, istate(j))
302          if (bad .lt. worst) then
303             it=j
304             worst=bad
305             iact=ks
306          endif
307  3100 continue
308 c
309 c  Test whether the Kuhn-Tucker condition is met.
310       if (worst .ge. 0.0 ) then
311          istate(n+1)=nbound
312          w(1)=sqrt(obj)
313          return
314       endif
315 c
316 c  The component  x(iact)  is the one that most wants to become active.
317 c   If the last successful change in the active set was to move x(iact)
318 c   to a bound, don't let x(iact) in now: set the derivative of the
319 c   misfit with respect to x(iact) to zero and return to the Kuhn-Tucker
320 c   test.
321       if ( iact .eq. jj ) then
322         w(jj)=0.0
323         goto 3000
324       endif
325 c
326 c  Step 5.
327 c  Undo the effect of the new (potentially) active variable on the
328 c   residual vector.
329       if (istate(it) .gt. 0) bound=bu(iact)
330       if (istate(it) .lt. 0) bound=bl(iact)
331       do 5100 i=1, m
332         act(i,mm1)=act(i,mm1) + bound*a(i,iact)
333  5100 continue
334 c
335 c  Set flag ifrom5, indicating that Step 6 was entered from Step 5.
336 c   This forms the basis of a test for instability: the gradient
337 c   calculation shows that x(iact) wants to join the active set; if
338 c   qr puts x(iact) beyond the bound from which it came, the gradient
339 c   calculation was in error and the variable should not have been
340 c   introduced.
341       ifrom5=istate(it)
342 c
343 c  Swap the indices (in istate) of the new active variable and the
344 c   rightmost bound variable; `unbind' that location by decrementing
345 c   nbound.
346       istate(it)=istate(nbound)
347       nbound=nbound - 1
348       nact=nact + 1
349       istate(nbound+1)=iact
350 c
351       if (mm .lt. nact) then
352         istop = -1000
353         msg = 'Too many free variables in BVLS!'
354         return
355 c        stop
356       endif
357 c
358 c  Step 6.
359 c  Load array  act  with the appropriate columns of  a  for qr.  For
360 c   added stability, reverse the column ordering so that the most
361 c   recent addition to the active set is in the last column.  Also
362 c   copy the residual vector from act(., mm1) into act(., mm1+1).
363  6000 do 6200 i=1, m
364         act(i,mm1+1)=act(i,mm1)
365         do 6100 k=nbound+1, n
366           j=istate(k)
367           act(i,nact+1-k+nbound)=a(i,j)
368  6100   continue
369  6200 continue
370 c
371       call qr(m, nact, act, act(1,mm1+1), zz, resq)
372 c
373 c  Test for linear dependence in qr, and for an instability that moves
374 c   the variable just introduced away from the feasible region
375 c   (rather than into the region or all the way through it).
376 c   In either case, remove the latest vector introduced from the
377 c   active set and adjust the residual vector accordingly. 
378 c   Set the gradient component (w(iact)) to zero and return to
379 c   the Kuhn-Tucker test.
380       if (resq .lt. 0.0
381      +   .or. (ifrom5 .gt. 0 .and. zz(nact) .gt. bu(iact))
382      +   .or. (ifrom5 .lt. 0 .and. zz(nact) .lt. bl(iact))) then
383          nbound=nbound + 1
384          istate(nbound)=istate(nbound)*sign(1.0d0, x(iact)-bu(iact))
385          nact=nact - 1
386          do 6500 i=1, m
387            act(i,mm1)=act(i,mm1) - x(iact)*a(i,iact)
388  6500    continue
389          ifrom5=0
390          w(iact)=0.0
391          goto 3000
392       endif
393 c
394 c  If Step 6 was entered from Step 5 and we are here, a new variable
395 c   has been successfully introduced into the active set; the last
396 c   variable that was fixed at a bound is again permitted to become
397 c   active.
398       if ( ifrom5 .ne. 0 ) jj=0
399       ifrom5=0
400 c
401 c   Step 7.  Check for strict feasibility of the new qr solution.
402       do 7100 k=1, nact
403         k1=k
404         j=istate(k+nbound)
405         if (zz(nact+1-k).lt.bl(j) .or. zz(nact+1-k).gt.bu(j)) goto 8000
406  7100 continue
407       do 7200 k=1, nact
408         j=istate(k+nbound)
409         x(j)=zz(nact+1-k)
410  7200 continue
411 c  New iterate is feasible; back to the top.
412       goto 15000
413 c
414 c  Steps 8, 9.
415  8000 alpha=2.0
416       alf=alpha
417       do 8200 k=k1, nact
418         j=istate(k+nbound)
419         if (zz(nact+1-k) .gt. bu(j))
420      +    alf=(bu(j)-x(j))/(zz(nact+1-k)-x(j))
421         if (zz(nact+1-k) .lt. bl(j))
422      +    alf=(bl(j)-x(j))/(zz(nact+1-k)-x(j))
423         if (alf .lt. alpha) then
424           alpha=alf
425           jj=j
426           sj=sign(1.0, zz(nact+1-k)-bl(j))
427         endif
428  8200 continue
429 c
430 c  Step 10
431       do 10000 k=1, nact
432         j=istate(k+nbound)
433         x(j)=x(j) + alpha*(zz(nact+1-k)-x(j))
434 10000 continue
435 c
436 c  Step 11. 
437 c  Move the variable that determined alpha to the appropriate bound. 
438 c   (jj is its index; sj is + if zz(jj)> bu(jj), - if zz(jj)<bl(jj) ).
439 c   If any other component of  x  is infeasible at this stage, it must
440 c   be due to roundoff.  Bind every infeasible component and every
441 c   component at a bound to the appropriate bound.  Correct the
442 c   residual vector for any variables moved to bounds.  Since at least
443 c   one variable is removed from the active set in this step, Loop B
444 c   (Steps 6-11) terminates after at most  nact  steps.
445       noldb=nbound
446       do 11200 k=1, nact
447         j=istate(k+noldb)
448         if (((bu(j)-x(j)) .le. 0.0) .or.
449      +    (j .eq. jj .and. sj .gt. 0.0)) then
450 c  Move x(j) to its upper bound.
451           x(j)=bu(j)
452           istate(k+noldb)=istate(nbound+1)
453           istate(nbound+1)=j
454           nbound=nbound+1
455           do 11100 i=1, m
456              act(i,mm1)=act(i,mm1) - bu(j)*a(i,j)
457 11100     continue
458         else if (((x(j)-bl(j)) .le. 0.0) .or.
459      +    (j .eq. jj .and. sj .lt. 0.0)) then
460 c  Move x(j) to its lower bound.
461           x(j)=bl(j)
462           istate(k+noldb)=istate(nbound+1)
463           istate(nbound+1)=-j
464           nbound=nbound+1
465           do 11150 i=1, m
466              act(i,mm1)=act(i,mm1) - bl(j)*a(i,j)
467 11150     continue
468         endif
469 11200 continue
470       nact=n - nbound
471 c
472 c  If there are still active variables left repeat the qr; if not,
473 c    go back to the top.
474       if (nact .gt. 0 ) goto 6000
475      
476 c
477 15000 continue
478 c
479       msg = 'BVLS fails to converge'
480       istop = -1000
481       return
482 c      stop
483       end
484 c======================================================================
485       subroutine qr(m, n, a, b, x, resq)
486 c======================================================================
487       implicit double precision (a-h, o-z)                              *doub
488 c$$$$ calls no other routines
489 c  Relies on FORTRAN77 do-loop conventions!
490 c  Solves over-determined least-squares problem  ax ~ b
491 c  where  a  is an  m by n  matrix,  b  is an m-vector .
492 c  resq  is the sum of squared residuals of optimal solution.  Also used
493 c  to signal error conditions - if -2 , system is underdetermined,  if
494 c  -1,  system is singular.
495 c  Method - successive Householder rotations.  See Lawson & Hanson -
496 c  Solving Least Squares Problems (1974).
497 c  Routine will also work when m=n.
498 c*****   CAUTION -  a and b  are overwritten by this routine.
499       dimension a(m,n),b(m),x(n)
500       double precision sum, dot
501 c
502       resq=-2.0
503       if (m .lt. n) return
504       resq=-1.0
505 c   Loop ending on 1800 rotates  a  into upper triangular form.
506       do 1800 j=1, n
507 c   Find constants for rotation and diagonal entry.
508         sq=0.0
509         do 1100 i=j, m
510           sq=a(i,j)**2 + sq
511  1100   continue
512         if (sq .eq. 0.0) return
513         qv1=-sign(sqrt(sq), a(j,j))
514         u1=a(j,j) - qv1
515         a(j,j)=qv1
516         j1=j + 1
517 c  Rotate remaining columns of sub-matrix.
518         do 1400 jj=j1, n
519           dot=u1*a(j,jj)
520           do 1200 i=j1, m
521             dot=a(i,jj)*a(i,j) + dot
522  1200     continue
523           const=dot/abs(qv1*u1)
524           do 1300 i=j1, m
525             a(i,jj)=a(i,jj) - const*a(i,j)
526  1300     continue
527           a(j,jj)=a(j,jj) - const*u1
528  1400   continue
529 c  Rotate  b  vector.
530         dot=u1*b(j)
531         do 1600 i=j1, m
532           dot=b(i)*a(i,j) + dot
533  1600   continue
534         const=dot/abs(qv1*u1)
535         b(j)=b(j) - const*u1
536         do 1700 i=j1, m
537           b(i)=b(i) - const*a(i,j)
538  1700   continue
539  1800 continue
540 c  Solve triangular system by back-substitution.
541       do 2200 ii=1, n
542         i=n-ii+1
543         sum=b(i)
544         do 2100 j=i+1, n
545           sum=sum - a(i,j)*x(j)
546  2100   continue
547         if (a(i,i).eq. 0.0) return
548          x(i)=sum/a(i,i)
549  2200 continue
550 c  Find residual in overdetermined case.
551       resq=0.0
552       do 2300 i=n+1, m
553         resq=b(i)**2 + resq
554  2300 continue
555       return
556       end                                                               qr
557 c______________________________________________________________________
558