forked from altMITgcm/MITgcm
-
Notifications
You must be signed in to change notification settings - Fork 0
/
cost_obcs_ageos.F
551 lines (474 loc) · 18.2 KB
/
cost_obcs_ageos.F
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
C $Header: /u/gcmpack/MITgcm/pkg/ecco/cost_obcs_ageos.F,v 1.10 2014/10/20 03:16:12 gforget Exp $
C $Name: $
#include "ECCO_OPTIONS.h"
#ifdef ALLOW_CTRL
# include "CTRL_OPTIONS.h"
#endif
subroutine cost_obcs_ageos( myiter, mytime, mythid )
c ==================================================================
c SUBROUTINE cost_obcs_ageos
c ==================================================================
c
c o cost function contribution obc -- Ageostrophic boundary flow.
c
c started: G. Gebbie gebbie@mit.edu 4-Feb-2003
c
c warning: masks redundantly assume no gradient of topography at
c boundary.
c
c ==================================================================
c SUBROUTINE cost_obcs_ageos
c ==================================================================
implicit none
c == global variables ==
#include "EEPARAMS.h"
#include "SIZE.h"
#include "GRID.h"
#include "DYNVARS.h"
#include "PARAMS.h"
#ifdef ALLOW_OBCS
# include "OBCS_GRID.h"
#endif
#ifdef ALLOW_CAL
# include "cal.h"
#endif
#ifdef ALLOW_ECCO
# include "ecco_cost.h"
#endif
#ifdef ALLOW_CTRL
# include "CTRL_SIZE.h"
# include "ctrl.h"
# include "ctrl_dummy.h"
# include "optim.h"
# include "CTRL_OBCS.h"
#endif
c == routine arguments ==
integer myiter
_RL mytime
integer mythid
#if (defined (ALLOW_CTRL) && \
defined (ALLOW_OBCS) && \
defined (ALLOW_ECCO))
c == local variables ==
integer bi,bj
integer i,j,k
integer itlo,ithi
integer jtlo,jthi
integer jmin,jmax
integer imin,imax
integer irec
integer levmon
integer levoff
integer iltheta
integer ilsalt
integer iluvel
integer ilvvel
integer ip1, jp1
_RL fctile
_RL fcthread
_RL rholoc (1-olx:snx+olx,1-oly:sny+oly,nsx,nsy)
_RL xzgrdrho(1-olx:snx+olx,Nr,nsx,nsy)
_RL yzgrdrho(1-oly:sny+oly,Nr,nsx,nsy)
_RL xzdvel1 (1-olx:snx+olx,nr,nsx,nsy)
_RL xzdvel2 (1-olx:snx+olx,nr,nsx,nsy)
_RL yzdvel1 (1-oly:sny+oly,nr,nsx,nsy)
_RL yzdvel2 (1-oly:sny+oly,nr,nsx,nsy)
_RL maskxzageos (1-olx:snx+olx,nr,nsx,nsy)
_RL maskyzageos (1-oly:sny+oly,nr,nsx,nsy)
_RL dummy
character*(80) fnametheta
character*(80) fnamesalt
character*(80) fnameuvel
character*(80) fnamevvel
logical doglobalread
logical ladinit
character*(MAX_LEN_MBUF) msgbuf
c == external functions ==
integer ilnblnk
external ilnblnk
c == end of interface ==
jtlo = mybylo(mythid)
jthi = mybyhi(mythid)
itlo = mybxlo(mythid)
ithi = mybxhi(mythid)
jmin = 1
jmax = sny
imin = 1
imax = snx
c-- Read tiled data.
doglobalread = .false.
ladinit = .false.
#ifdef OBCS_AGEOS_COST_CONTRIBUTION
#ifdef ECCO_VERBOSE
_BEGIN_MASTER( mythid )
write(msgbuf,'(a)') ' '
call print_message( msgbuf, standardmessageunit,
& SQUEEZE_RIGHT , mythid)
write(msgbuf,'(a,i8.8)')
& ' cost_obcs_ageos: number of records to process = ',nmonsrec
call print_message( msgbuf, standardmessageunit,
& SQUEEZE_RIGHT , mythid)
write(msgbuf,'(a)') ' '
call print_message( msgbuf, standardmessageunit,
& SQUEEZE_RIGHT , mythid)
_END_MASTER( mythid )
#endif
if (optimcycle .ge. 0) then
iltheta = ilnblnk( tbarfile )
write(fnametheta(1:80),'(2a,i10.10)')
& tbarfile(1:iltheta),'.',optimcycle
ilsalt = ilnblnk( sbarfile )
write(fnamesalt(1:80),'(2a,i10.10)')
& sbarfile(1:ilsalt),'.',optimcycle
iluvel = ilnblnk( ubarfile )
write(fnameuvel(1:80),'(2a,i10.10)')
& ubarfile(1:iluvel),'.',optimcycle
ilvvel = ilnblnk( vbarfile )
write(fnamevvel(1:80),'(2a,i10.10)')
& vbarfile(1:ilvvel),'.',optimcycle
endif
fcthread = 0. _d 0
fctile = 0. _d 0
cgg Code safe: always initialize to zero.
do bj = jtlo,jthi
do bi = itlo,ithi
do k = 1,nr
do i = 1-olx,snx+olx
maskxzageos(i,k,bi,bj) = 0. _d 0
xzdvel1(i,k,bi,bj) = 0. _d 0
xzdvel2(i,k,bi,bj) = 0. _d 0
xzgrdrho(i,k,bi,bj) = 0. _d 0
enddo
do j = 1-oly,sny+oly
maskyzageos(j,k,bi,bj) = 0. _d 0
yzdvel1(j,k,bi,bj) = 0. _d 0
yzdvel2(j,k,bi,bj) = 0. _d 0
yzgrdrho(j,k,bi,bj) = 0. _d 0
enddo
enddo
enddo
enddo
do bj = jtlo,jthi
do bi = itlo,ithi
do j = 1-oly,sny+oly
do i = 1-olx,snx+olx
rholoc(i,j,bi,bj) = 0. _d 0
enddo
enddo
enddo
enddo
c-- Loop over records.
do irec = 1,nmonsrec
c-- Read time averages and the monthly mean data.
call active_read_xyz( fnametheta, tbar, irec,
& doglobalread, ladinit,
& optimcycle, mythid,
& xx_tbar_mean_dummy )
c-- Read time averages and the monthly mean data.
call active_read_xyz( fnamesalt, sbar, irec,
& doglobalread, ladinit,
& optimcycle, mythid,
& xx_sbar_mean_dummy )
c-- Read time averages and the monthly mean data.
call active_read_xyz( fnameuvel, ubar, irec,
& doglobalread, ladinit,
& optimcycle, mythid,
& xx_ubar_mean_dummy )
c-- Read time averages and the monthly mean data.
call active_read_xyz( fnamevvel, vbar, irec,
& doglobalread, ladinit,
& optimcycle, mythid,
& xx_vbar_mean_dummy )
cgg Minor problem : grad T,S needs overlap values.
_BARRIER
_EXCH_XYZ_RL(tbar,myThid)
_EXCH_XYZ_RL(sbar,myThid)
do bj = jtlo,jthi
do bi = itlo,ithi
#ifdef ALLOW_OBCSN_CONTROL
jp1 = 0
cgg Make a mask for the velocity shear comparison.
do k = 1,nr-1
do i = imin, imax
j = OB_Jn(i,bi,bj)
cgg All these points need to be wet.
if ( j.eq.OB_indexNone ) then
maskxzageos(i,k,bi,bj) = 0.
else
maskxzageos(i,k,bi,bj) =
& hfacC(i,j+jp1,k,bi,bj)*hfacC(i+1,j+jp1,k,bi,bj) *
& hfacC(i-1,j+jp1,k,bi,bj)*hfacC(i,j+jp1,k+1,bi,bj)*
& hfacC(i-1,j+jp1,k+1,bi,bj)*hfacC(i+1,j+jp1,k+1,bi,bj)*
& hfacS(i,j+jp1,k,bi,bj)*hfacS(i,j+jp1,k+1,bi,bj)
endif
enddo
enddo
do k = 1,nr
C-- jmc: both calls below are wrong if more than 1 tile => stop here
IF ( bi.NE.1 .OR. bj.NE.1 )
& STOP 'COST_OBCS_AGEOS wrong with more than 1 tile/proc'
call find_rho_2d(
I iMin, iMax, jMin, jMax, k,
I tbar(1-OLx,1-OLy,k,bi,bj),
I sbar(1-OLx,1-OLy,k,bi,bj),
O rholoc,
I k, bi, bj, myThid )
_EXCH_XY_RL(rholoc , myThid)
cgg Compute centered difference horizontal gradient on bdy.
do i = imin, imax
j = OB_Jn(i,bi,bj)
if ( j.eq.OB_indexNone ) j = 1
xzgrdrho(i,k,bi,bj) =
& (rholoc(i-1,j+jp1,bi,bj)-rholoc(i+1,j+jp1,bi,bj))
& /(2.*dxc(i,j+jp1,bi,bj))
enddo
enddo
cgg Compute vertical shear from geostrophy/thermal wind.
cgg Above level 4 needs not be geostrophic.
cgg Please get rid of the "4" ASAP. Ridiculous.
do k = 4,Nr-1
do i = imin,imax
j = OB_Jn(i,bi,bj)
if ( j.eq.OB_indexNone ) j = 1
xzdvel1(i,k,bi,bj) = vbar(i,j+jp1,k ,bi,bj)
& - vbar(i,j+jp1,k+1,bi,bj)
xzdvel2(i,k,bi,bj)=((xzgrdrho(i,k,bi,bj)*delz(k)/2.)+
& (xzgrdrho(i,k+1,bi,bj)*delz(k+1)/2.))
& * gravity / f0 / rhonil
fctile = fctile + 100*wcurrent(k,bi,bj) *
& maskxzageos(i,k,bi,bj)*
& (xzdvel2(i,k,bi,bj) - xzdvel1(i,k,bi,bj))*
& (xzdvel2(i,k,bi,bj) - xzdvel1(i,k,bi,bj))
if (maskxzageos(i,k,bi,bj) .ne. 0) then
cgg print*,'XX shear,geos shear',xzdvel1(i,k,bi,bj),xzdvel2(i,k,bi,bj)
cgg print*,'i,j,k,fctile N',i,j,k,fctile
endif
enddo
enddo
c-- End of loop over layers.
#endif /* ALLOW_OBCSN_CONTROL */
#ifdef ALLOW_OBCSS_CONTROL
jp1 = 1
cgg Make a mask for the velocity shear comparison.
do k = 1,nr-1
do i = imin, imax
j = OB_Js(i,bi,bj)
if ( j.eq.OB_indexNone ) then
maskxzageos(i,k,bi,bj) = 0.
else
cgg All these points need to be wet.
maskxzageos(i,k,bi,bj) =
& hfacC(i,j+jp1,k,bi,bj)*hfacC(i+1,j+jp1,k,bi,bj) *
& hfacC(i-1,j+jp1,k,bi,bj)*hfacC(i,j+jp1,k+1,bi,bj)*
& hfacC(i-1,j+jp1,k+1,bi,bj)*hfacC(i+1,j+jp1,k+1,bi,bj)*
& hfacS(i,j+jp1,k,bi,bj)*hfacS(i,j+jp1,k+1,bi,bj)
endif
enddo
enddo
do k = 1,nr
C-- jmc: both calls below are wrong if more than 1 tile => stop here
IF ( bi.NE.1 .OR. bj.NE.1 )
& STOP 'COST_OBCS_AGEOS wrong with more than 1 tile/proc'
call find_rho_2d(
I iMin, iMax, jMin, jMax, k,
I tbar(1-OLx,1-OLy,k,bi,bj),
I sbar(1-OLx,1-OLy,k,bi,bj),
O rholoc,
I k, bi, bj, myThid )
_EXCH_XY_RL(rholoc , myThid)
cgg Compute centered difference horizontal gradient on bdy.
do i = imin, imax
j = OB_Js(i,bi,bj)
if ( j.eq.OB_indexNone ) j = 1
xzgrdrho(i,k,bi,bj) =
& (rholoc(i-1,j+jp1,bi,bj)-rholoc(i+1,j+jp1,bi,bj))
& /(2.*dxc(i,j+jp1,bi,bj))
enddo
enddo
cgg Compute vertical shear from geostrophy/thermal wind.
do k = 4,Nr-1
do i = imin,imax
j = OB_Js(i,bi,bj)
if ( j.eq.OB_indexNone ) j = 1
cgg Retrieve the model vertical shear.
xzdvel1(i,k,bi,bj) = vbar(i,j+jp1,k ,bi,bj)
& - vbar(i,j+jp1,k+1,bi,bj)
cgg Compute vertical shear from geostrophy/thermal wind.
xzdvel2(i,k,bi,bj) =((xzgrdrho(i,k,bi,bj)*delz(k)/2.)+
& (xzgrdrho(i,k+1,bi,bj)*delz(k+1)/2.))
& * gravity / f0 /rhonil
cgg Make a comparison.
fctile = fctile + 100*wcurrent(k,bi,bj) *
& maskxzageos(i,k,bi,bj)*
& (xzdvel2(i,k,bi,bj) - xzdvel1(i,k,bi,bj))*
& (xzdvel2(i,k,bi,bj) - xzdvel1(i,k,bi,bj))
cgg print*,'fctile S',fctile
enddo
enddo
c-- End of loop over layers.
#endif
#ifdef ALLOW_OBCSW_CONTROL
ip1 = 1
cgg Make a mask for the velocity shear comparison.
do k = 1,nr-1
do j = jmin, jmax
i = OB_Iw(j,bi,bj)
cgg All these points need to be wet.
if ( i.eq.OB_indexNone ) then
maskyzageos(j,k,bi,bj) = 0.
else
maskyzageos(j,k,bi,bj) =
& hfacC(i+ip1,j,k,bi,bj)*hfacC(i+ip1,j+1,k,bi,bj) *
& hfacC(i+ip1,j-1,k,bi,bj)*hfacC(i+ip1,j,k+1,bi,bj)*
& hfacC(i+ip1,j-1,k+1,bi,bj)*hfacC(i+ip1,j+1,k+1,bi,bj)*
& hfacW(i+ip1,j,k,bi,bj)*hfacW(i+ip1,j,k+1,bi,bj)
endif
enddo
enddo
do k = 1,nr
IF ( bi.NE.1 .OR. bj.NE.1 )
& STOP 'COST_OBCS_AGEOS wrong with more than 1 tile/proc'
call find_rho_2d(
I iMin, iMax, jMin, jMax, k,
I tbar(1-OLx,1-OLy,k,bi,bj),
I sbar(1-OLx,1-OLy,k,bi,bj),
O rholoc,
I k, bi, bj, myThid )
_EXCH_XY_RL(rholoc , myThid)
cgg Compute centered difference horizontal gradient on bdy.
do j = jmin, jmax
i = OB_Iw(j,bi,bj)
if ( i.eq.OB_indexNone ) i = 1
cgg Negative sign due to geostrophy.
yzgrdrho(j,k,bi,bj) =
& (rholoc(i+ip1,j+1,bi,bj)-rholoc(i+ip1,j-1,bi,bj))
& /(2.*dyc(i+ip1,j,bi,bj))
enddo
enddo
cgg Compute vertical shear from geostrophy/thermal wind.
do k = 4,Nr-1
do j = jmin,jmax
i = OB_Iw(j,bi,bj)
if ( i.eq.OB_indexNone ) i = 1
cgg Retrieve the model vertical shear.
yzdvel1(j,k,bi,bj) = ubar(i+ip1,j,k ,bi,bj)
& - ubar(i+ip1,j,k+1,bi,bj)
cgg Compute vertical shear from geostrophy/thermal wind.
yzdvel2(j,k,bi,bj) =((yzgrdrho(j,k ,bi,bj)*delz(k)/2.)+
& (yzgrdrho(j,k+1,bi,bj)*delz(k+1)/2.))
& * gravity / f0 / rhonil
cgg Make a comparison.
fctile = fctile + 100*wcurrent(k,bi,bj) *
& maskyzageos(j,k,bi,bj) *
& (yzdvel2(j,k,bi,bj) - yzdvel1(j,k,bi,bj))*
& (yzdvel2(j,k,bi,bj) - yzdvel1(j,k,bi,bj))
enddo
enddo
c-- End of loop over layers.
#endif
#ifdef ALLOW_OBCSE_CONTROL
ip1 = 0
cgg Make a mask for the velocity shear comparison.
do k = 1,nr-1
do j = jmin, jmax
i = OB_Ie(j,bi,bj)
if ( i.eq.OB_indexNone ) then
maskyzageos(j,k,bi,bj) =0.
else
cgg All these points need to be wet.
maskyzageos(j,k,bi,bj) =
& hfacC(i+ip1,j,k,bi,bj)*hfacC(i+ip1,j+1,k,bi,bj) *
& hfacC(i+ip1,j-1,k,bi,bj)*hfacC(i+ip1,j,k+1,bi,bj)*
& hfacC(i+ip1,j-1,k+1,bi,bj)*hfacC(i+ip1,j+1,k+1,bi,bj)*
& hfacW(i+ip1,j,k,bi,bj)*hfacW(i+ip1,j,k+1,bi,bj)
endif
enddo
enddo
do k = 1,nr
IF ( bi.NE.1 .OR. bj.NE.1 )
& STOP 'COST_OBCS_AGEOS wrong with more than 1 tile/proc'
call find_rho_2d(
I iMin, iMax, jMin, jMax, k,
I tbar(1-OLx,1-OLy,k,bi,bj),
I sbar(1-OLx,1-OLy,k,bi,bj),
O rholoc,
I k, bi, bj, myThid )
_EXCH_XY_RL(rholoc , myThid)
cgg Compute centered difference horizontal gradient on bdy.
do j = jmin, jmax
i = OB_Ie(j,bi,bj)
if ( i.eq.OB_indexNone ) i = 1
cgg Negative sign due to geostrophy.
yzgrdrho(j,k,bi,bj) =
& (rholoc(i+ip1,,j+1,bi,bj)-rholoc(i+ip1,j-1,bi,bj))
& /(2.*dyc(i+ip1,j,bi,bj))
enddo
enddo
cgg Compute vertical shear from geostrophy/thermal wind.
do k = 4,Nr-1
do j = jmin,jmax
i = OB_Ie(j,bi,bj)
if ( i.eq.OB_indexNone ) i = 1
cgg Retrieve the model vertical shear.
yzdvel1(j,k,bi,bj) = ubar(i+ip1,j,k ,bi,bj)
& - ubar(i+ip1,j,k+1,bi,bj)
cgg Compute vertical shear from geostrophy/thermal wind.
yzdvel2(j,k,bi,bj) =((yzgrdrho(j,k ,bi,bj)*delz(k)/2.)+
& (yzgrdrho(j,k+1,bi,bj)*delz(k+1)/2.))
& * gravity / f0 /rhonil
cgg Make a comparison.
fctile = fctile + 100*wcurrent(k,bi,bj) *
& maskyzageos(j,k,bi,bj) *
& (yzdvel2(j,k,bi,bj) - yzdvel1(j,k,bi,bj))*
& (yzdvel2(j,k,bi,bj) - yzdvel1(j,k,bi,bj))
enddo
enddo
c-- End of loop over layers.
#endif
fcthread = fcthread + fctile
objf_ageos(bi,bj) = objf_ageos(bi,bj) + fctile
#ifdef ECCO_VERBOSE
c-- Print cost function for each tile in each thread.
write(msgbuf,'(a)') ' '
call print_message( msgbuf, standardmessageunit,
& SQUEEZE_RIGHT , mythid)
write(msgbuf,'(a,i8.8,1x,i3.3,1x,i3.3)')
& ' cost_Theta: irec,bi,bj = ',irec,bi,bj
call print_message( msgbuf, standardmessageunit,
& SQUEEZE_RIGHT , mythid)
write(msgbuf,'(a,d22.15)')
& ' cost function (temperature) = ',
& fctile
call print_message( msgbuf, standardmessageunit,
& SQUEEZE_RIGHT , mythid)
write(msgbuf,'(a)') ' '
call print_message( msgbuf, standardmessageunit,
& SQUEEZE_RIGHT , mythid)
#endif
enddo
enddo
#ifdef ECCO_VERBOSE
c-- Print cost function for all tiles.
_GLOBAL_SUM_RL( fcthread , myThid )
write(msgbuf,'(a)') ' '
call print_message( msgbuf, standardmessageunit,
& SQUEEZE_RIGHT , mythid)
write(msgbuf,'(a,i8.8)')
& ' cost_Theta: irec = ',irec
call print_message( msgbuf, standardmessageunit,
& SQUEEZE_RIGHT , mythid)
write(msgbuf,'(a,a,d22.15)')
& ' global cost function value',
& ' (temperature) = ',fcthread
call print_message( msgbuf, standardmessageunit,
& SQUEEZE_RIGHT , mythid)
write(msgbuf,'(a)') ' '
call print_message( msgbuf, standardmessageunit,
& SQUEEZE_RIGHT , mythid)
#endif
enddo
c-- End of loop over records.
#endif
#endif /* ALLOW_CTRL, ALLOW_OBCS and ALLOW_ECCO */
return
end