/[MITgcm]/MITgcm/lsopt/lsupdxx.F
ViewVC logotype

Diff of /MITgcm/lsopt/lsupdxx.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph | View Patch Patch

revision 1.1 by heimbach, Tue Feb 5 20:34:34 2002 UTC revision 1.2 by heimbach, Fri Nov 15 04:03:24 2002 UTC
# Line 0  Line 1 
1    
2          subroutine lsupdxx(
3         &     nn, ifail, lphprint
4         &     , jmin, jmax, nupdate
5         &     , ff, fmin, fold, gnorm0, dotdg
6         &     , gg, dd, xx, xdiff
7         &     , tmin, tmax, tact, epsx
8         &     )
9    
10    c     ==================================================================
11    c     SUBROUTINE lsupdxx
12    c     ==================================================================
13    c
14    c     o conceived for variable online/offline version
15    c       computes - new descent direction dd based on latest
16    c                  available gradient
17    c                - new tact based on new dd
18    c                - new control vector xx needed for offline run
19    c
20    c     o started: Patrick Heimbach, MIT/EAPS
21    c                29-Feb-2000:
22    c
23    c     o Version 2.1.0, 02-Mar-2000: Patrick Heimbach, MIT/EAPS
24    c
25    c     ==================================================================
26    c     SUBROUTINE lsupdxx
27    c     ==================================================================
28    c
29    
30    #include <blas1.h>
31    
32          implicit none
33    
34    c-----------------------------------------
35    c declare arguments
36    c-----------------------------------------
37          integer nn, jmin, jmax, nupdate, ifail
38          double precision    ff, fmin, fold, gnorm0, dotdg
39          double precision    gg(nn), dd(nn), xx(nn), xdiff(nn)
40          double precision    tmin, tmax, tact, epsx
41          logical lphprint
42    
43    c-----------------------------------------
44    C declare local variables
45    c-----------------------------------------
46          integer i
47          double precision    fdiff, preco
48    
49          double precision     SDOT
50          external SDOT
51    
52    c     ==================================================================
53    
54    c-----------------------------------------
55    c use Fletchers scaling
56    c and initialize diagional to 1.
57    c-----------------------------------------
58    c
59          if ( ( jmax .eq. 0 ) .or. (nupdate .eq. 0 ) ) then
60    
61             if (jmax .eq. 0) then
62                fold = fmin
63                if (lphprint)
64         &           print *, 'pathei-lsopt: using fold = fmin = ', fmin
65             end if
66             fdiff = fold - ff
67             if (jmax .eq. 0) fdiff = ABS(fdiff)
68            
69             preco = 2. * fdiff / (gnorm0*gnorm0)
70             do i = 1, nn
71                dd(i)    = -gg(i)*preco
72             end do
73    
74             if (lphprint)
75         &        print *, 'pathei-lsopt: first estimate of dd via ',
76         &        'fold - ff'
77    
78    c-----------------------------------------
79    c use the matrix stored in [diag]
80    c and the (y,s) pairs
81    c-----------------------------------------
82    
83             else
84    
85                do i = 1, nn
86                   dd(i) = -gg(i)
87                end do
88    
89                if (jmax .gt. 0) then
90                   call hessupd( nn, nupdate, dd, jmin, jmax, xdiff,
91         &              lphprint )
92                else
93                   if (lphprint)
94         &              print *, 'pathei-lsopt: no hessupd for first optim.'
95                end if
96    
97             endif
98    
99    c-----------------------------------------
100    c check whether new direction is a descent one
101    c-----------------------------------------
102             dotdg = SDOT( nn, dd, 1, gg, 1 )
103             if (dotdg .ge. 0.0) then
104                ifail = 4
105                goto 999
106             end if
107    
108    c----------------------------------
109    c declare arguments
110    c----------------------------------
111    
112          tmin = 0.
113          do i = 1, nn
114             tmin = max( tmin, abs(dd(i)) )
115          end do
116          tmin = epsx/tmin
117    
118    c----------------------------------
119    c make sure that t is between
120    c tmin and tmax
121    c----------------------------------
122    
123          tact  = 1.0
124          tmax = 1.0e+10
125          if (tact.le.tmin) then
126             tact = tmin
127             if (tact.gt.tmax) then
128                tmin = tmax
129             endif
130          endif
131    
132          if (tact .gt. tmax) then
133              tact = tmax
134              ifail = 7
135          endif
136    
137    c----------------------------------
138    c compute new x
139    c----------------------------------
140          do i = 1, nn
141             xdiff(i) = xx(i) + tact*dd(i)
142          end do
143    
144    c----------------------------------
145    c save new x to file for offline version
146    c----------------------------------
147    
148     999  continue
149    
150          return
151    
152          end

Legend:
Removed from v.1.1  
changed lines
  Added in v.1.2

  ViewVC Help
Powered by ViewVC 1.1.22