1 |
C $Header: /u/gcmpack/MITgcm/pkg/ecco/cost_obcsw.F,v 1.3 2007/10/09 00:02:50 jmc Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "COST_CPPOPTIONS.h" |
5 |
#ifdef ALLOW_OBCS |
6 |
# include "OBCS_OPTIONS.h" |
7 |
#endif |
8 |
|
9 |
subroutine cost_obcsw( |
10 |
I myiter, |
11 |
I mytime, |
12 |
I mythid |
13 |
& ) |
14 |
|
15 |
c ================================================================== |
16 |
c SUBROUTINE cost_obcsw |
17 |
c ================================================================== |
18 |
c |
19 |
c o cost function contribution obc |
20 |
c |
21 |
c ================================================================== |
22 |
c SUBROUTINE cost_obcsw |
23 |
c ================================================================== |
24 |
|
25 |
implicit none |
26 |
|
27 |
c == global variables == |
28 |
|
29 |
#include "EEPARAMS.h" |
30 |
#include "SIZE.h" |
31 |
#include "PARAMS.h" |
32 |
#include "GRID.h" |
33 |
#include "DYNVARS.h" |
34 |
#ifdef ALLOW_OBCS |
35 |
# include "OBCS.h" |
36 |
#endif |
37 |
|
38 |
#include "cal.h" |
39 |
#include "ecco_cost.h" |
40 |
#include "ctrl.h" |
41 |
#include "ctrl_dummy.h" |
42 |
#include "optim.h" |
43 |
|
44 |
c == routine arguments == |
45 |
|
46 |
integer myiter |
47 |
_RL mytime |
48 |
integer mythid |
49 |
|
50 |
c == local variables == |
51 |
|
52 |
integer bi,bj |
53 |
integer i,j,k |
54 |
integer itlo,ithi |
55 |
integer jtlo,jthi |
56 |
integer jmin,jmax |
57 |
integer imin,imax |
58 |
integer irec |
59 |
integer il |
60 |
integer iobcs |
61 |
integer ip1 |
62 |
|
63 |
_RL fctile |
64 |
_RL fcthread |
65 |
_RL dummy |
66 |
|
67 |
character*(80) fnametheta |
68 |
character*(80) fnamesalt |
69 |
character*(80) fnameuvel |
70 |
character*(80) fnamevvel |
71 |
|
72 |
logical doglobalread |
73 |
logical ladinit |
74 |
|
75 |
#ifdef ECCO_VERBOSE |
76 |
character*(MAX_LEN_MBUF) msgbuf |
77 |
#endif |
78 |
|
79 |
c == external functions == |
80 |
|
81 |
integer ilnblnk |
82 |
external ilnblnk |
83 |
|
84 |
c == end of interface == |
85 |
|
86 |
jtlo = mybylo(mythid) |
87 |
jthi = mybyhi(mythid) |
88 |
itlo = mybxlo(mythid) |
89 |
ithi = mybxhi(mythid) |
90 |
jmin = 1 |
91 |
jmax = sny |
92 |
imin = 1 |
93 |
imax = snx |
94 |
|
95 |
c-- Read tiled data. |
96 |
doglobalread = .false. |
97 |
ladinit = .false. |
98 |
|
99 |
#ifdef ALLOW_OBCSW_COST_CONTRIBUTION |
100 |
|
101 |
ip1 = 1 |
102 |
fcthread = 0. _d 0 |
103 |
|
104 |
c-- Loop over records. |
105 |
do irec = 1,nmonsrec |
106 |
|
107 |
c-- temperature |
108 |
iobcs = 1 |
109 |
c-- Read time averages and the monthly mean data. |
110 |
il = ilnblnk( tbarfile ) |
111 |
write(fnametheta(1:80),'(2a,i10.10)') |
112 |
& tbarfile(1:il),'.',optimcycle |
113 |
call active_read_xyz( fnametheta, tbar, irec, |
114 |
& doglobalread, ladinit, |
115 |
& optimcycle, mythid, |
116 |
& xx_tbar_mean_dummy ) |
117 |
|
118 |
call mdsreadfieldyz( OBWtFile, readBinaryPrec, 'RS', |
119 |
& nr, OBWt, irec, mythid) |
120 |
|
121 |
do bj = jtlo,jthi |
122 |
do bi = itlo,ithi |
123 |
c-- Loop over the model layers |
124 |
fctile = 0. _d 0 |
125 |
do k = 1,nr |
126 |
c-- Compute model data misfit and cost function term for |
127 |
c the temperature field. |
128 |
do j = jmin,jmax |
129 |
i = OB_Iw(J,bi,bj) |
130 |
if (maskS(i+ip1,j,k,bi,bj) .ne. 0.) then |
131 |
fctile = fctile + |
132 |
& wobcsw(k,iobcs)*cosphi(i,j,bi,bj)* |
133 |
& (tbar(i,j,k,bi,bj) - OBWt(j,k,bi,bj))* |
134 |
& (tbar(i,j,k,bi,bj) - OBWt(j,k,bi,bj)) |
135 |
endif |
136 |
enddo |
137 |
enddo |
138 |
c-- End of loop over layers. |
139 |
fcthread = fcthread + fctile |
140 |
objf_obcsw(bi,bj) = objf_obcsw(bi,bj) + fctile |
141 |
enddo |
142 |
enddo |
143 |
|
144 |
c-- salt |
145 |
iobcs = 2 |
146 |
c-- Read time averages and the monthly mean data. |
147 |
il = ilnblnk( sbarfile ) |
148 |
write(fnamesalt(1:80),'(2a,i10.10)') |
149 |
& sbarfile(1:il),'.',optimcycle |
150 |
call active_read_xyz( fnamesalt, sbar, irec, |
151 |
& doglobalread, ladinit, |
152 |
& optimcycle, mythid, |
153 |
& xx_sbar_mean_dummy ) |
154 |
|
155 |
call mdsreadfieldyz( OBWsFile, readBinaryPrec, 'RS', |
156 |
& nr, OBWs, irec, mythid) |
157 |
|
158 |
do bj = jtlo,jthi |
159 |
do bi = itlo,ithi |
160 |
c-- Loop over the model layers |
161 |
fctile = 0. _d 0 |
162 |
do k = 1,nr |
163 |
c-- Compute model data misfit and cost function term for |
164 |
c the temperature field. |
165 |
do j = jmin,jmax |
166 |
i = OB_Iw(J,bi,bj) |
167 |
if (maskS(i+ip1,j,k,bi,bj) .ne. 0.) then |
168 |
fctile = fctile + |
169 |
& wobcsw(k,iobcs)*cosphi(i,j,bi,bj)* |
170 |
& (sbar(i,j,k,bi,bj) - OBWs(j,k,bi,bj))* |
171 |
& (sbar(i,j,k,bi,bj) - OBWs(j,k,bi,bj)) |
172 |
endif |
173 |
enddo |
174 |
enddo |
175 |
c-- End of loop over layers. |
176 |
fcthread = fcthread + fctile |
177 |
objf_obcsw(bi,bj) = objf_obcsw(bi,bj) + fctile |
178 |
enddo |
179 |
enddo |
180 |
|
181 |
c-- uvel |
182 |
iobcs = 3 |
183 |
c-- Read time averages and the monthly mean data. |
184 |
il = ilnblnk( ubarfile ) |
185 |
write(fnameuvel(1:80),'(2a,i10.10)') |
186 |
& ubarfile(1:il),'.',optimcycle |
187 |
call active_read_xyz( fnameuvel, ubar, irec, |
188 |
& doglobalread, ladinit, |
189 |
& optimcycle, mythid, |
190 |
& dummy ) |
191 |
|
192 |
call mdsreadfieldyz( OBWuFile, readBinaryPrec, 'RS', |
193 |
& nr, OBWu, irec, mythid) |
194 |
do bj = jtlo,jthi |
195 |
do bi = itlo,ithi |
196 |
c-- Loop over the model layers |
197 |
fctile = 0. _d 0 |
198 |
do k = 1,nr |
199 |
c-- Compute model data misfit and cost function term for |
200 |
c the temperature field. |
201 |
do j = jmin,jmax |
202 |
i = OB_Iw(J,bi,bj) |
203 |
if (maskW(i+ip1,j,k,bi,bj) .ne. 0.) then |
204 |
fctile = fctile + |
205 |
& wobcsw(k,iobcs)*cosphi(i,j,bi,bj)* |
206 |
& (ubar(i,j,k,bi,bj) - OBWu(j,k,bi,bj))* |
207 |
& (ubar(i,j,k,bi,bj) - OBWu(j,k,bi,bj)) |
208 |
endif |
209 |
enddo |
210 |
enddo |
211 |
c-- End of loop over layers. |
212 |
fcthread = fcthread + fctile |
213 |
objf_obcsw(bi,bj) = objf_obcsw(bi,bj) + fctile |
214 |
enddo |
215 |
enddo |
216 |
|
217 |
c-- vvel |
218 |
iobcs = 4 |
219 |
c-- Read time averages and the monthly mean data. |
220 |
il = ilnblnk( vbarfile ) |
221 |
write(fnamevvel(1:80),'(2a,i10.10)') |
222 |
& vbarfile(1:il),'.',optimcycle |
223 |
call active_read_xyz( fnamevvel, vbar, irec, |
224 |
& doglobalread, ladinit, |
225 |
& optimcycle, mythid, |
226 |
& dummy ) |
227 |
|
228 |
call mdsreadfieldyz( OBWvFile, readBinaryPrec, 'RS', |
229 |
& nr, OBWv, irec, mythid) |
230 |
|
231 |
do bj = jtlo,jthi |
232 |
do bi = itlo,ithi |
233 |
c-- Loop over the model layers |
234 |
fctile = 0. _d 0 |
235 |
do k = 1,nr |
236 |
c-- Compute model data misfit and cost function term for |
237 |
c the temperature field. |
238 |
do j = jmin,jmax |
239 |
i = OB_Iw(J,bi,bj) |
240 |
if (maskS(i,j,k,bi,bj) .ne. 0.) then |
241 |
fctile = fctile + |
242 |
& wobcsw(k,iobcs)*cosphi(i,j,bi,bj)* |
243 |
& (vbar(i,j,k,bi,bj) - OBWv(j,k,bi,bj))* |
244 |
& (vbar(i,j,k,bi,bj) - OBWv(j,k,bi,bj)) |
245 |
endif |
246 |
enddo |
247 |
enddo |
248 |
c-- End of loop over layers. |
249 |
fcthread = fcthread + fctile |
250 |
objf_obcsw(bi,bj) = objf_obcsw(bi,bj) + fctile |
251 |
enddo |
252 |
enddo |
253 |
|
254 |
#ifdef ECCO_VERBOSE |
255 |
c-- Print cost function for all tiles. |
256 |
_GLOBAL_SUM_RL( fcthread , myThid ) |
257 |
write(msgbuf,'(a)') ' ' |
258 |
call print_message( msgbuf, standardmessageunit, |
259 |
& SQUEEZE_RIGHT , mythid) |
260 |
write(msgbuf,'(a,i8.8)') |
261 |
& ' cost_obcsw: irec = ',irec |
262 |
call print_message( msgbuf, standardmessageunit, |
263 |
& SQUEEZE_RIGHT , mythid) |
264 |
write(msgbuf,'(a,a,d22.15)') |
265 |
& ' global cost function value', |
266 |
& ' (obcsw) = ',fcthread |
267 |
call print_message( msgbuf, standardmessageunit, |
268 |
& SQUEEZE_RIGHT , mythid) |
269 |
write(msgbuf,'(a)') ' ' |
270 |
call print_message( msgbuf, standardmessageunit, |
271 |
& SQUEEZE_RIGHT , mythid) |
272 |
#endif |
273 |
|
274 |
enddo |
275 |
c-- End of loop over records. |
276 |
|
277 |
#endif |
278 |
|
279 |
return |
280 |
end |
281 |
|
282 |
|
283 |
|
284 |
|
285 |
|
286 |
|
287 |
|