/[MITgcm]/MITgcm_contrib/mlosch/optim_m1qn3/optim_sub.F
ViewVC logotype

Annotation of /MITgcm_contrib/mlosch/optim_m1qn3/optim_sub.F

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


Revision 1.6 - (hide annotations) (download)
Mon Jun 4 12:23:34 2012 UTC (11 years, 10 months ago) by mlosch
Branch: MAIN
Changes since 1.5: +10 -2 lines
add extra output for extra debuggung

1 mlosch 1.6 C $Header: /u/gcmpack/MITgcm_contrib/mlosch/optim_m1qn3/optim_sub.F,v 1.5 2012/05/10 10:25:14 mlosch Exp $
2 mlosch 1.2 C $Name: $
3 mlosch 1.1
4     c Include ECCO_CPPOPTIONS because the ecco_ctrl,cost files
5     c have headers with options for OBCS masks.
6     #include "ECCO_CPPOPTIONS.h"
7    
8     subroutine optim_sub(
9     I nn, ff
10     & )
11    
12     c ==================================================================
13     c subroutine optim_sub
14     c ==================================================================
15     c
16     c o This is the main driver routine for the offline version of
17     c m1qn3 (m1qn3_offline). It
18     c - sets all m1qn3 relevant parameters
19     c - reads the model state (control vector, cost function,
20     c and gradient)
21     c - reads the state of m1qn3_offline
22     c - calls m1qn3_offline
23     c - saves model state (control vector) and state of m1qn3_offline
24     c The routine is somewhat lengthy and could be split into separate
25     c subroutines, but I felt that it was easier to write and test it
26     c in this form
27     c Martin Losch (Martin.Losch@awi.de), Apr, 2012
28     c
29     c ==================================================================
30    
31     implicit none
32    
33     c == global variables ==
34    
35     #include "EEPARAMS.h"
36     #include "SIZE.h"
37    
38     #include "ctrl.h"
39     #include "optim.h"
40     #include "m1qn3_common.h"
41    
42     c == routine arguments ==
43    
44     integer nn
45     _RL ff
46    
47     c == local variables ==
48    
49     _RL objf
50    
51 mlosch 1.2 #ifdef DYNAMIC
52     _RL, dimension(:), allocatable :: xx, adxx
53 mlosch 1.1 #else
54     integer nmax
55     parameter( nmax = MAX_INDEPEND )
56     _RL xx(nmax)
57     _RL adxx(nmax)
58     #endif
59    
60 mlosch 1.4 CML logical coldStart
61 mlosch 1.1 c formal parameters of m1qn3
62     integer reverse
63     integer impres,imode(3),omode,niter,nsim,iz(5),indic
64     _RL dxmin,df1
65     character*3 normtype
66     c work arrays
67 mlosch 1.2 integer ndz
68 mlosch 1.1 CML _RL dz(ndz)
69     double precision, dimension(:), allocatable :: dz
70     c extra dummy variables
71     integer izs(1)
72     _RS rzs(1)
73     _RL dzs(1)
74     integer, parameter :: io = 60
75     character*(*), parameter :: fname_m1qn3 = 'm1qn3_output.txt'
76     c end of m1qn3 parameters
77    
78     integer i
79    
80     c == external ==
81    
82     external simul_rc,euclid,ctonbe,ctcabe
83    
84     c == end of interface ==
85    
86     c-- Allocate memory for the control variables and the gradient vector.
87     #if defined(DYNAMIC)
88 mlosch 1.2 allocate( xx(nn) )
89     allocate( adxx(nn) )
90 mlosch 1.1 #endif
91    
92 mlosch 1.2 #ifndef DYNAMIC
93 mlosch 1.1 if (nn .gt. nmax) then
94     print*,' OPTIMUM: Not enough space.'
95     print*,' nmax = ',nmax
96     print*,' nn = ',nn
97     print*
98     print*,' Set MAX_INDEPEND in Makefile .ge. ',nn
99     print*
100     stop ' ... stopped in OPTIMUM.'
101     endif
102     #endif
103    
104     print*, ' OPTIM_SUB: Calling m1qn3_optim for iteration: ',
105     & optimcycle
106     print*, ' OPTIM_SUB: with nn, REAL_BYTE = ', nn, REAL_BYTE
107    
108     c can be 'two','sup','dfn', see m1qn3 documentation for details
109     normtype='two'
110     c after reading data.optim some of these parameter values can be guessed
111     c impres=6, impres determines the amount of m1qn3-output see documentation
112     impres=iprint
113     c these should strictly be different (nsim>niter), but in practice
114     c it does not matter
115     niter = numiter
116     nsim = nfunc*niter
117     c epsg=1.d-8
118     dxmin=epsx
119     c will be set later
120     df1=-UNSET_RL
121     c
122     imode=(/0,1,0/)
123     omode=-1
124     c initialise work array
125     ndz = 3*nn+nupdate*(3*nn+1)
126     do i=1,5
127     iz(i)=0
128     enddo
129     allocate(dz(ndz))
130     do i=1,ndz
131     dz(i) = 0.
132     enddo
133     c these alway have to be set like this
134     reverse=1
135     indic=4
136     c initialise the dummy arguments that are not used
137     izs(1)=UNSET_I
138     rzs(1)=UNSET_RS
139     dzs(1)=UNSET_RL
140    
141 mlosch 1.4 c-- first read the model output into xx, adxx, and cost function
142     c value into objf
143     do i = 1,nn
144     xx(i) = 0.
145     adxx(i) = 0.
146     enddo
147     c
148     print *, ' OPTIM_SUB: read model state'
149     call optim_readdata( nn, ctrlname, .false., objf, xx )
150     call optim_readdata( nn, costname, .false., objf, adxx )
151     print *, ' OPTIM_SUB after reading ',
152     & ctrlname, ' and ', costname, ':'
153     print *, ' OPTIM_SUB nn = ', nn
154     print *, ' OPTIM_SUB objf = ', objf
155     print *, ' OPTIM_SUB xx(1) = ', xx(1)
156     print *, ' OPTIM_SUB adxx(1) = ', adxx(1)
157    
158 mlosch 1.3 if ( coldStart ) then
159 mlosch 1.1 c-- cold start
160     print *, ' OPTIM_SUB: cold start, optimcycle =', optimcycle
161     imode(2) = 0
162     c this variable is the only one of the m1qn3-related common blocks
163     c that needs to be initialized here to make sure that we have a
164     c clean start
165     reentry = 0
166 mlosch 1.5 c compute expected decrease of cost function from objf and fmin;
167     c this value is only used for a cold start of m1qn3_offline, for a
168     c warm start df1 is overwritten with data from a restart file
169     df1=objf-fmin
170     if ( df1 .le. 0. ) then
171     print *, ' OPTIM_SUB: df1 = objf-fmin = ', df1,
172     & ' but should be > 0.'
173     stop 'ABNORMAL in S/R OPTIM_SUB'
174     endif
175    
176 mlosch 1.1 c open output file for m1qn3
177     open(io,file=fname_m1qn3,status='unknown')
178     else
179     c-- warm restart
180     c requires restoring the state of m1qn3 with pickup file
181     print *, ' OPTIM_SUB: warm start, optimcycle =', optimcycle
182     imode(2) = 1
183     call optim_store_m1qn3(ndz,iz,dz,niter,nsim,epsg,df1,
184     I optimcycle,
185     I .false.)
186     c re-open output file for m1qn3
187     open(io,file=fname_m1qn3,status='old',position='append')
188     endif
189    
190     c-- call the minimizer, a slightly modified version of m1qn3 v3.3
191     c (Gilbert & Lemarechal, 1989), downloaded in April 2012.
192     c simul_rc,euclid,ctonbe,ctcabe are external subroutines that
193     c are provided within m1qn3. euclid, ctonbe, ctcabe can be replaced
194     c by something more efficient, simul_rc is a dummy routine for
195     c the reverse communication mode and should not be changed.
196 mlosch 1.6 print *, ' OPTIM_SUB: call m1qn3_offline ........'
197 mlosch 1.1 call m1qn3_offline (simul_rc,euclid,ctonbe,ctcabe,
198     & nn,xx,objf,adxx,dxmin,df1,
199     & epsg,normtype,impres,io,imode,omode,niter,nsim,
200     & iz,dz,ndz,reverse,indic,izs,rzs,dzs)
201     close(io)
202 mlosch 1.6 print *, ' OPTIM_SUB: ...........................'
203 mlosch 1.1 print *, ' OPTIM_SUB: returned from m1qn3_offline'
204 mlosch 1.6 print *, ' OPTIM_SUB: nn = ', nn
205     print *, ' OPTIM_SUB: xx(1) = ', xx(1), xx(2)
206     print *, ' OPTIM_SUB: adxx(1) = ', adxx(1), adxx(2)
207     print *, ' OPTIM_SUB: omode = ', omode
208     print *, ' OPTIM_SUB: niter = ', niter
209     print *, ' OPTIM_SUB: nsim = ', nsim
210     print *, ' OPTIM_SUB: reverse = ', reverse
211 mlosch 1.1
212     c write state of m1qn3 into pickup file for warm restart
213     call optim_store_m1qn3(ndz,iz,dz,niter,nsim,epsg,df1,
214     I optimcycle,
215     I .true.)
216     c write model control vector
217     print *,' OPTIMS_SUB: writing ', nn,' sized control to file ',
218     & ctrlname
219     c give the cost function a funny value to make sure that nobody
220     c mistakes it for the real one
221     call optim_writedata( nn, ctrlname, .false., -9999., xx )
222    
223     c clean up
224 mlosch 1.2 #ifdef DYNAMIC
225     deallocate(xx, adxx)
226     #endif /* DYNAMIC */
227 mlosch 1.1 deallocate(dz)
228    
229     c stopping criterion
230     if ( reverse .lt. 0 ) then
231     print *, ' OPTIM_SUB: reverse = ', reverse
232     print *, ' OPTIM_SUB: optimization terminated with omode = ',
233     & omode
234     stop 'ABNORMAL in S/R OPTIM_SUB'
235     endif
236    
237     return
238     end

  ViewVC Help
Powered by ViewVC 1.1.22