1 |
C $Header: /u/gcmpack/MITgcm_contrib/mlosch/optim_m1qn3/optim_sub.F,v 1.2 2012/04/27 09:50:46 mlosch Exp $ |
2 |
C $Name: $ |
3 |
|
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 |
#ifdef DYNAMIC |
52 |
_RL, dimension(:), allocatable :: xx, adxx |
53 |
#else |
54 |
integer nmax |
55 |
parameter( nmax = MAX_INDEPEND ) |
56 |
_RL xx(nmax) |
57 |
_RL adxx(nmax) |
58 |
#endif |
59 |
|
60 |
logical coldStart |
61 |
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 |
integer ndz |
68 |
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 |
allocate( xx(nn) ) |
89 |
allocate( adxx(nn) ) |
90 |
#endif |
91 |
|
92 |
#ifndef DYNAMIC |
93 |
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 |
coldStart=.false. |
142 |
if ( optimcycle .eq. 0 ) coldStart=.true. |
143 |
if ( coldStart ) then |
144 |
c-- cold start |
145 |
print *, ' OPTIM_SUB: cold start, optimcycle =', optimcycle |
146 |
imode(2) = 0 |
147 |
c this variable is the only one of the m1qn3-related common blocks |
148 |
c that needs to be initialized here to make sure that we have a |
149 |
c clean start |
150 |
reentry = 0 |
151 |
c ff has be read in optim_readparms, so we do not read it here again |
152 |
objf = ff |
153 |
df1 = objf-fmin |
154 |
c open output file for m1qn3 |
155 |
open(io,file=fname_m1qn3,status='unknown') |
156 |
else |
157 |
c-- warm restart |
158 |
c requires restoring the state of m1qn3 with pickup file |
159 |
print *, ' OPTIM_SUB: warm start, optimcycle =', optimcycle |
160 |
imode(2) = 1 |
161 |
call optim_store_m1qn3(ndz,iz,dz,niter,nsim,epsg,df1, |
162 |
I optimcycle, |
163 |
I .false.) |
164 |
c re-open output file for m1qn3 |
165 |
open(io,file=fname_m1qn3,status='old',position='append') |
166 |
endif |
167 |
c-- read the model output into xx,adxx |
168 |
if ( indic .eq. 4 ) then |
169 |
do i = 1,nn |
170 |
xx(i) = 0. |
171 |
adxx(i) = 0. |
172 |
enddo |
173 |
c |
174 |
print *, ' OPTIM_SUB: read model state' |
175 |
call optim_readdata( nn, ctrlname, .false., objf, xx ) |
176 |
call optim_readdata( nn, costname, .false., objf, adxx ) |
177 |
print *, ' OPTIM_SUB after reading nn, objf = ', nn, objf, |
178 |
& xx(1), adxx(1) |
179 |
else |
180 |
print *, ' OPTIM_SUB: indic = ', indic, ' is not possible' |
181 |
stop 'ABNORMAL in S/R OPTIM_SUB' |
182 |
endif |
183 |
|
184 |
c-- call the minimizer, a slightly modified version of m1qn3 v3.3 |
185 |
c (Gilbert & Lemarechal, 1989), downloaded in April 2012. |
186 |
c simul_rc,euclid,ctonbe,ctcabe are external subroutines that |
187 |
c are provided within m1qn3. euclid, ctonbe, ctcabe can be replaced |
188 |
c by something more efficient, simul_rc is a dummy routine for |
189 |
c the reverse communication mode and should not be changed. |
190 |
print *, ' OPTIM_SUB: call m1qn3_offline' |
191 |
call m1qn3_offline (simul_rc,euclid,ctonbe,ctcabe, |
192 |
& nn,xx,objf,adxx,dxmin,df1, |
193 |
& epsg,normtype,impres,io,imode,omode,niter,nsim, |
194 |
& iz,dz,ndz,reverse,indic,izs,rzs,dzs) |
195 |
close(io) |
196 |
print *, ' OPTIM_SUB: returned from m1qn3_offline' |
197 |
|
198 |
c write state of m1qn3 into pickup file for warm restart |
199 |
call optim_store_m1qn3(ndz,iz,dz,niter,nsim,epsg,df1, |
200 |
I optimcycle, |
201 |
I .true.) |
202 |
c write model control vector |
203 |
print *,' OPTIMS_SUB: writing ', nn,' sized control to file ', |
204 |
& ctrlname |
205 |
c give the cost function a funny value to make sure that nobody |
206 |
c mistakes it for the real one |
207 |
call optim_writedata( nn, ctrlname, .false., -9999., xx ) |
208 |
|
209 |
c clean up |
210 |
#ifdef DYNAMIC |
211 |
deallocate(xx, adxx) |
212 |
#endif /* DYNAMIC */ |
213 |
deallocate(dz) |
214 |
|
215 |
c stopping criterion |
216 |
if ( reverse .lt. 0 ) then |
217 |
print *, ' OPTIM_SUB: reverse = ', reverse |
218 |
print *, ' OPTIM_SUB: optimization terminated with omode = ', |
219 |
& omode |
220 |
stop 'ABNORMAL in S/R OPTIM_SUB' |
221 |
endif |
222 |
|
223 |
return |
224 |
end |