-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathbarotp.F+98ZA
executable file
·663 lines (662 loc) · 22.9 KB
/
barotp.F+98ZA
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
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
#if defined(ROW_LAND)
#define SEA_P .true.
#define SEA_U .true.
#define SEA_V .true.
#elif defined(ROW_ALLSEA)
#define SEA_P allip(j).or.ip(i,j).ne.0
#define SEA_U alliu(j).or.iu(i,j).ne.0
#define SEA_V alliv(j).or.iv(i,j).ne.0
#else
#define SEA_P ip(i,j).ne.0
#define SEA_U iu(i,j).ne.0
#define SEA_V iv(i,j).ne.0
#endif
subroutine barotp(m,n)
use mod_xc ! HYCOM communication interface
use mod_cb_arrays ! HYCOM saved arrays
use mod_pipe ! HYCOM debugging interface
use mod_tides ! HYCOM tides
#if defined(STOKES)
use mod_stokes ! HYCOM Stokes Drift
#endif
implicit none
c
integer m,n
c
c --- ------------------------------------------------------------------------
c --- advance barotropic equations.
c --- on entry: -n- is time t-dt, -m- is time t
c --- on exit: -m- is time t, -n- is time t+dt
c --- time level 3 is only used internally (n and m are always 1 or 2).
c
c --- LeapFrog version based on:
c --- Y. Morel, Baraille, R., Pichon A. (2008) "Time splitting and
c --- linear stability of the slow part of the barotropic component",
c --- Ocean Modeling, 23, pp 73-81.
c --- ------------------------------------------------------------------------
c
logical lpipe_barotp
parameter (lpipe_barotp=.false.)
logical ldebug_barotp
parameter (ldebug_barotp=.false.)
c
real q,pbudel,pbvdel,utndcy,vtndcy,wblpf
real d11,d12,d21,d22,ubp,vbp
real*8 sump
integer i,j,l,lll,ml,nl,mn,lstep1,margin,mbdy
c & ,iffstep
logical ldrag
c data iffstep/0/
c save iffstep
c
#if defined(RELO)
real, save, allocatable, dimension(:,:) ::
& pbavo,ubavo,vbavo,displd
c
if (.not.allocated(pbavo)) then
allocate(
& pbavo(1-nbdy:idm+nbdy,1-nbdy:jdm+nbdy),
& ubavo(1-nbdy:idm+nbdy,1-nbdy:jdm+nbdy),
& vbavo(1-nbdy:idm+nbdy,1-nbdy:jdm+nbdy),
& displd(1-nbdy:idm+nbdy,1-nbdy:jdm+nbdy) )
call mem_stat_add( 4*(idm+2*nbdy)*(jdm+2*nbdy) )
pbavo = r_init
ubavo = r_init
vbavo = r_init
displd = r_init
endif
#else
real, save, dimension(1-nbdy:idm+nbdy,1-nbdy:jdm+nbdy) ::
& pbavo,ubavo,vbavo,displd
#endif
c
mbdy = 6
c
c --- utotn,vtotn from momtum is time step t-1 to t+1 barotropic tendency
call xctilr(utotn( 1-nbdy,1-nbdy ),1, 1, 6,6, halo_uv)
call xctilr(vtotn( 1-nbdy,1-nbdy ),1, 1, 6,6, halo_vv)
c
if (lpipe .and. lpipe_barotp) then
c --- compare two model runs.
call pipe_compare_sym2(utotn, iu,'barotp:utotn',
& vtotn, iv,'barotp:vtotn')
call pipe_compare_sym1(pvtrop,iq,'barotp:pvtrp')
endif
c
c --- explicit time integration of barotropic flow (forward-backward scheme)
c --- in order to combine forward-backward scheme with leapfrog treatment of
c --- coriolis term, v-eqn must be solved before u-eqn every other time step
c
if (btrlfr) then
if (delt1.ne.baclin) then !not on very 1st time step
c --- start at time level t-dt and go to t+dt.
lstep1 = lstep + lstep !more stable, but also more expensive
!$OMP PARALLEL DO PRIVATE(j,i)
!$OMP& SCHEDULE(STATIC,jblk)
do j=1,jj
do i=1,ii
pbavo(i,j) = pbavg(i,j,n) !save t-1 for RA filter
ubavo(i,j) = ubavg(i,j,n) !save t-1 for RA filter
vbavo(i,j) = vbavg(i,j,n) !save t-1 for RA filter
c
pbavg(i,j,3) = pbavg(i,j,n)
ubavg(i,j,3) = ubavg(i,j,n)
vbavg(i,j,3) = vbavg(i,j,n)
enddo !i
enddo !j
else !1st time step
c --- start at time level t and go to t+dt.
lstep1 = lstep
!$OMP PARALLEL DO PRIVATE(j,i)
!$OMP& SCHEDULE(STATIC,jblk)
do j=1,jj
do i=1,ii
pbavo(i,j) = 0.0 !makes correct mean height safe
pbavg(i,j,n) = pbavg(i,j,m)
ubavg(i,j,n) = ubavg(i,j,m)
vbavg(i,j,n) = vbavg(i,j,m)
pbavg(i,j,3) = pbavg(i,j,m)
ubavg(i,j,3) = ubavg(i,j,m)
vbavg(i,j,3) = vbavg(i,j,m)
enddo !i
enddo !j
endif !usual:1st time step
else
c --- start at time level t and go to t+dt.
lstep1 = lstep !original, less stable, method
!$OMP PARALLEL DO PRIVATE(j,i)
!$OMP& SCHEDULE(STATIC,jblk)
do j=1,jj
do i=1,ii
pbavo(i,j) = 0.0 !makes correct mean height safe
pbavg(i,j,n) = pbavg(i,j,m)
ubavg(i,j,n) = ubavg(i,j,m)
vbavg(i,j,n) = vbavg(i,j,m)
enddo !i
enddo !j
endif !btrlfr
c
ldrag = tidflg.gt.0 .and. drgscl.ne.0.0 .and. thkdrg.eq.0.0
c
if (ldrag) then
displd(:,:) = 0.0
endif
c
c --- time step loop
c
if (btrlfr) then
wblpf = 0.0 !1st minor time step, lll=1, only
else
wblpf = wbaro
endif
c
do 840 lll=1,lstep1,2
c
call xctilr(pbavg( 1-nbdy,1-nbdy,1 ),1, 3, 6,6, halo_ps)
call xctilr(ubavg( 1-nbdy,1-nbdy,1 ),1, 3, 6,6, halo_uv)
call xctilr(vbavg( 1-nbdy,1-nbdy,1 ),1, 3, 6,6, halo_vv)
c
if (lpipe .and. lpipe_barotp) then
call pipe_compare_sym1(
& pbavg(1-nbdy,1-nbdy,nl),ip,'barot+:pbavn')
call pipe_compare_sym2(
& ubavg(1-nbdy,1-nbdy,nl),iu,'barot+:ubavn',
& vbavg(1-nbdy,1-nbdy,nl),iv,'barot+:vbavn')
call pipe_compare_sym1(
& pbavg(1-nbdy,1-nbdy,ml),ip,'barot+:pbavm')
call pipe_compare_sym2(
& ubavg(1-nbdy,1-nbdy,ml),iu,'barot+:ubavm',
& vbavg(1-nbdy,1-nbdy,ml),iv,'barot+:vbavm')
endif
c
c --- odd minor time step.
c
ml=n
nl=3
c
c --- continuity equation, and tidal drag on p-grid
c
c --- rhs: pbavg, ubavg+, vbavg+
c --- lhs: pbavg
c
margin = mbdy - 1
c
!$OMP PARALLEL DO PRIVATE(j,i,pbudel,pbvdel,
!$OMP& ubp,vbp,d11,d12,d21,d22,q)
!$OMP& SCHEDULE(STATIC,jblk)
do j=1-margin,jj+margin
do i=1-margin,ii+margin
if (SEA_P) then
#if defined(STOKES)
c
c Barotropic Stokes flow included here
c
pbudel = (ubavg(i+1,j,ml)+usdbavg(i+1,j))*
& (depthu(i+1,j)*scuy(i+1,j))
& -(ubavg(i, j,ml)+usdbavg(i, j))*
& (depthu(i, j)*scuy(i, j))
pbvdel = (vbavg(i,j+1,ml)+vsdbavg(i,j+1))*
& (depthv(i,j+1)*scvx(i,j+1))
& -(vbavg(i,j, ml)+vsdbavg(i,j ))*
& (depthv(i,j )*scvx(i,j ))
#else
pbudel = ubavg(i+1,j,ml)*(depthu(i+1,j)*scuy(i+1,j))
& -ubavg(i ,j,ml)*(depthu(i ,j)*scuy(i ,j))
pbvdel = vbavg(i,j+1,ml)*(depthv(i,j+1)*scvx(i,j+1))
& -vbavg(i,j ,ml)*(depthv(i,j )*scvx(i,j ))
#endif
pbavg(i,j,nl)=
& ((1.-wblpf)*pbavg(i,j,ml)+
& wblpf *pbavg(i,j,nl) )-
& (1.+wblpf)*dlt*(pbudel + pbvdel)*scp2i(i,j)
c
if (ldrag) then
c
c --- tidal drag tensor on p-grid:
c --- ub = ub - (dlt/H)*(t.11*ub + t.12*vb)
c --- vb = vb - (dlt/H)*(t.21*ub + t.22*vb)
c --- solve implicitly by inverting the matrix:
c --- 1+(dlt/H)*t.11 (dlt/H)*t.12
c --- (dlt/H)*t.21 1+(dlt/H)*t.22
c --- use depths (H) rather than onem*pbavg (h) for stability.
c
ubp = 0.5*(ubavg(i+1,j,nl)+ubavg(i,j,nl))
vbp = 0.5*(vbavg(i,j+1,nl)+vbavg(i,j,nl))
d11 = -dlt/depths(i,j) * drgten(1,1,i,j)
d12 = -dlt/depths(i,j) * drgten(1,2,i,j)
d21 = -dlt/depths(i,j) * drgten(2,1,i,j)
d22 = -dlt/depths(i,j) * drgten(2,2,i,j)
q = 1.0/((1.0-d11)*(1.0-d22)-d12*d21)
c --- set util5,util6 to the ubavg,vbavg drag increment
util5(i,j) = q*(ubp*(1.0-d22)+vbp*d12) - ubp
util6(i,j) = q*(ubp*d21+vbp*(1.0-d11)) - vbp
c --- add an explicit antidrag correction
* util5(i,j) = util5(i,j) - (d11*untide(i,j)+
* & d12*vntide(i,j) )
* util6(i,j) = util6(i,j) - (d21*untide(i,j)+
* & d22*vntide(i,j) )
c --- dissipation per m^2
displd(i,j) = displd(i,j) +
& (ubp*util5(i,j) + vbp*util6(i,j))*
& depths(i,j)*qthref/dlt
c
* if (ldebug_barotp .and. i.eq.itest.and.j.eq.jtest) then
* write (lp,'(i9,2i5,i3,3x,a,4g15.6)')
* & nstep,i+i0,j+j0,lll,
* & 'ubp,new,vbp,new =',
* & ubp,ubp+util5(i,j),
* & vbp,vbp+util6(i,j)
* endif !debug
else
util5(i,j) = 0.0
util6(i,j) = 0.0
endif !ldrag
endif !ip
enddo !i
enddo !j
c
mn=ml
c
c --- u momentum equation, 1st
c
c --- rhs: pbavg+, vbavg+, pvtrop+
c --- lhs: ubavg
c
margin = margin - 1
c
!$OMP PARALLEL DO PRIVATE(j,i,utndcy)
!$OMP& SCHEDULE(STATIC,jblk)
do j=1-margin,jj+margin
do i=1-margin,ii+margin
if (SEA_U) then
utndcy=-thref*(pbavg(i,j,nl)-pbavg(i-1,j,nl))*scuxi(i,j)+
& ((vbavg(i ,j, mn)*depthv(i ,j)
& +vbavg(i ,j+1,mn)*depthv(i ,j+1))+
& (vbavg(i-1,j, mn)*depthv(i-1,j)
& +vbavg(i-1,j+1,mn)*depthv(i-1,j+1)))*
& (0.125*(pvtrop(i,j)+pvtrop(i,j+1)))
c
ubavg(i,j,nl)=
& ((1.-wblpf)*ubavg(i,j,ml)+
& wblpf *ubavg(i,j,nl) )+
& (1.+wblpf)*dlt*(utndcy+utotn(i,j))+
& 0.5*(util5(i,j)+util5(i-1,j))
c
* if (ldebug_barotp .and. i.eq.itest.and.j.eq.jtest) then
* write (lp,'(i9,2i5,i3,3x,a,6g15.6)')
* & nstep,i+i0,j+j0,lll,
* & 'u_old,u_new,p_grad,corio,u_star,drag =',
* & ubavg(i,j,ml),ubavg(i,j,nl),
* & -thref*(pbavg(i,j,nl)-pbavg(i-1,j,nl))*scuxi(i,j)*dlt,
* & (vbavg(i ,j, mn)*depthv(i ,j)
* & +vbavg(i ,j+1,mn)*depthv(i ,j+1)
* & +vbavg(i-1,j, mn)*depthv(i-1,j)
* & +vbavg(i-1,j+1,mn)*depthv(i-1,j+1))
* & *(pvtrop(i,j)+pvtrop(i,j+1))
* & *.125 * dlt,utotn(i,j) * dlt,
* & 0.5*(util5(i,j)+util5(i-1,j))
* endif !debug
endif !iu
enddo !i
enddo !j
c
mn = nl
c
c --- v momentum equation, 2nd
c --- rhs: pbavg+, ubavg+, pvtrop+
c --- lhs: vbavg
c
margin = margin - 1
c
!$OMP PARALLEL DO PRIVATE(j,i,vtndcy)
!$OMP& SCHEDULE(STATIC,jblk)
do j=1-margin,jj+margin
do i=1-margin,ii+margin
if (SEA_V) then
vtndcy=-thref*(pbavg(i,j,nl)-pbavg(i,j-1,nl))*scvyi(i,j)-
& ((ubavg(i, j ,mn)*depthu(i, j )
& +ubavg(i+1,j ,mn)*depthu(i+1,j ))+
& (ubavg(i, j-1,mn)*depthu(i, j-1)
& +ubavg(i+1,j-1,mn)*depthu(i+1,j-1)))*
& (0.125*(pvtrop(i,j)+pvtrop(i+1,j)))
c
vbavg(i,j,nl)=
& ((1.-wblpf)*vbavg(i,j,ml)+
& wblpf *vbavg(i,j,nl) )+
& (1.+wblpf)*dlt*(vtndcy+vtotn(i,j))+
& 0.5*(util6(i,j)+util6(i,j-1))
c
* if (ldebug_barotp .and. i.eq.itest.and.j.eq.jtest) then
* write (lp,'(i9,2i5,i3,3x,a,6g15.6)')
* & nstep,i+i0,j+j0,lll,
* & 'v_old,v_new,p_grad,corio,v_star,drag =',
* & vbavg(i,j,ml),vbavg(i,j,nl),
* & -thref*(pbavg(i,j,nl)-pbavg(i,j-1,nl))*scvyi(i,j)*dlt,
* & -(ubavg(i, j ,mn)*depthu(i,j )
* & +ubavg(i+1,j ,mn)*depthu(i+1,j )
* & +ubavg(i, j-1,mn)*depthu(i,j-1)
* & +ubavg(i+1,j-1,mn)*depthu(i+1,j-1))
* & *(pvtrop(i,j)+pvtrop(i+1,j))
* & *.125 * dlt, vtotn(i,j) * dlt,
* & 0.5*(util6(i,j)+util6(i,j-1))
* endif !debug
endif !iv
enddo !i
enddo !j
c
if (ldebug_barotp) then
call xcsync(flush_lp)
endif
c
#if ! defined(RELO)
if (lbflag.eq.1) then
call latbdp( nl)
elseif (lbflag.eq.3) then
call latbdf( nl,lll)
endif
#endif
if (lbflag.eq.2) then
call latbdt( nl,lll)
elseif (lbflag.eq.4) then
call latbdtf(nl,lll)
endif
c
if (lpipe .and. lpipe_barotp) then
call pipe_compare_sym1(
& pbavg(1-nbdy,1-nbdy,nl),ip,'barot+:pbavn')
call pipe_compare_sym2(
& ubavg(1-nbdy,1-nbdy,nl),iu,'barot+:ubavn',
& vbavg(1-nbdy,1-nbdy,nl),iv,'barot+:vbavn')
call pipe_compare_sym1(
& pbavg(1-nbdy,1-nbdy,ml),ip,'barot+:pbavm')
call pipe_compare_sym2(
& ubavg(1-nbdy,1-nbdy,ml),iu,'barot+:ubavm',
& vbavg(1-nbdy,1-nbdy,ml),iv,'barot+:vbavm')
endif
c
c --- even minor time step.
c
ml=3
nl=n
wblpf = wbaro !used for all subsequent time steps: lll=2,lstep1
c
c --- continuity equation
c
c --- rhs: pbavg, ubavg+, vbavg+
c --- lhs: pbavg
c
margin = mbdy - 1
c
!$OMP PARALLEL DO PRIVATE(j,i,pbudel,pbvdel,
!$OMP& ubp,vbp,d11,d12,d21,d22,q)
!$OMP& SCHEDULE(STATIC,jblk)
do j=1-margin,jj+margin
do i=1-margin,ii+margin
if (SEA_P) then
#if defined(STOKES)
c
c Barotropic Stokes flow included here
c
pbudel = (ubavg(i+1,j,ml)+usdbavg(i+1,j))*
& (depthu(i+1,j)*scuy(i+1,j))
& -(ubavg(i, j,ml)+usdbavg(i, j))*
& (depthu(i ,j)*scuy(i, j))
pbvdel = (vbavg(i,j+1,ml)+vsdbavg(i,j+1))*
& (depthv(i,j+1)*scvx(i,j+1))
& -(vbavg(i,j, ml)+vsdbavg(i,j ))*
& (depthv(i,j )*scvx(i,j ))
#else
pbudel = ubavg(i+1,j,ml)*(depthu(i+1,j)*scuy(i+1,j))
& -ubavg(i ,j,ml)*(depthu(i ,j)*scuy(i ,j))
pbvdel = vbavg(i,j+1,ml)*(depthv(i,j+1)*scvx(i,j+1))
& -vbavg(i,j ,ml)*(depthv(i,j )*scvx(i,j ))
#endif
pbavg(i,j,nl)=
& ((1.-wblpf)*pbavg(i,j,ml)+
& wblpf *pbavg(i,j,nl) )-
& (1.+wblpf)*dlt*(pbudel + pbvdel)*scp2i(i,j)
c
if (ldrag) then
c --- tidal drag tensor on p-grid:
c --- ub = ub - (dlt/H)*(t.11*ub + t.12*vb)
c --- vb = vb - (dlt/H)*(t.21*ub + t.22*vb)
c --- solve implicitly by inverting the matrix:
c --- 1+(dlt/H)*t.11 (dlt/H)*t.12
c --- (dlt/H)*t.21 1+(dlt/H)*t.22
c --- use depths (H) rather than onem*pbavg (h) for stability.
c
ubp = 0.5*(ubavg(i+1,j,nl)+ubavg(i,j,nl))
vbp = 0.5*(vbavg(i,j+1,nl)+vbavg(i,j,nl))
d11 = -dlt/depths(i,j) * drgten(1,1,i,j)
d12 = -dlt/depths(i,j) * drgten(1,2,i,j)
d21 = -dlt/depths(i,j) * drgten(2,1,i,j)
d22 = -dlt/depths(i,j) * drgten(2,2,i,j)
q = 1.0/((1.0-d11)*(1.0-d22)-d12*d21)
c --- set util5,util6 to the ubavg,vbavg drag increment
util5(i,j) = q*(ubp*(1.0-d22)+vbp*d12) - ubp
util6(i,j) = q*(ubp*d21+vbp*(1.0-d11)) - vbp
c --- add an explicit antidrag correction
* util5(i,j) = util5(i,j) - (d11*untide(i,j)+
* & d12*vntide(i,j) )
* util6(i,j) = util6(i,j) - (d21*untide(i,j)+
* & d22*vntide(i,j) )
c --- dissipation per m^2
displd(i,j) = displd(i,j) +
& (ubp*util5(i,j) + vbp*util6(i,j))*
& depths(i,j)*qthref/dlt
c
* if (ldebug_barotp .and. i.eq.itest.and.j.eq.jtest) then
* write (lp,'(i9,2i5,i3,3x,a,4g15.6)')
* & nstep,i+i0,j+j0,lll+1,
* & 'ubp,new,vbp,new =',
* & ubp,ubp+util5(i,j),
* & vbp,vbp+util6(i,j)
* endif !debug
else
util5(i,j) = 0.0
util6(i,j) = 0.0
endif !ldrag
endif !ip
enddo !i
enddo !j
c
mn=ml
c
c --- v momentum equation, 1st
c
c --- rhs: pbavg+, ubavg+, pvtrop+
c --- lhs: vbavg
c
margin = margin - 1
c
!$OMP PARALLEL DO PRIVATE(j,i,vtndcy)
!$OMP& SCHEDULE(STATIC,jblk)
do j=1-margin,jj+margin
do i=1-margin,ii+margin
if (SEA_V) then
vtndcy=-thref*(pbavg(i,j,nl)-pbavg(i,j-1,nl))*scvyi(i,j)-
& ((ubavg(i, j ,mn)*depthu(i, j )
& +ubavg(i+1,j ,mn)*depthu(i+1,j ))+
& (ubavg(i, j-1,mn)*depthu(i, j-1)
& +ubavg(i+1,j-1,mn)*depthu(i+1,j-1)))*
& (0.125*(pvtrop(i,j)+pvtrop(i+1,j)))
c
vbavg(i,j,nl)=
& ((1.-wblpf)*vbavg(i,j,ml)+
& wblpf *vbavg(i,j,nl))+
& (1.+wblpf)*dlt*(vtndcy+vtotn(i,j))+
& 0.5*(util6(i,j)+util6(i,j-1))
c
* if (ldebug_barotp .and. i.eq.itest.and.j.eq.jtest) then
* write (lp,'(i9,2i5,i3,3x,a,6g15.6)')
* & nstep,i+i0,j+j0,lll+1,
* & 'v_old,v_new,p_grad,corio,v_star,drag =',
* & vbavg(i,j,ml),vbavg(i,j,nl),
* & -thref*(pbavg(i,j,nl)-pbavg(i,j-1,nl))*scvyi(i,j)*dlt,
* & -(ubavg(i, j ,mn)*depthu(i,j )
* & +ubavg(i+1,j ,mn)*depthu(i+1,j )
* & +ubavg(i, j-1,mn)*depthu(i,j-1)
* & +ubavg(i+1,j-1,mn)*depthu(i+1,j-1))
* & *(pvtrop(i,j)+pvtrop(i+1,j))
* & *.125 * dlt, vtotn(i,j) * dlt,
* & 0.5*(util6(i,j)+util6(i,j-1))
* endif !debug
endif !iv
enddo !i
enddo !j
c
mn=nl
c
c --- u momentum equation, 2nd
c
c --- rhs: pbavg+, vbavg+, pvtrop+
c --- lhs: ubavg
c
margin = margin - 1
c
!$OMP PARALLEL DO PRIVATE(j,i,utndcy)
!$OMP& SCHEDULE(STATIC,jblk)
do j=1-margin,jj+margin
do i=1-margin,ii+margin
if (SEA_U) then
utndcy=-thref*(pbavg(i,j,nl)-pbavg(i-1,j,nl))*scuxi(i,j)+
& ((vbavg(i ,j, mn)*depthv(i ,j)
& +vbavg(i ,j+1,mn)*depthv(i ,j+1))+
& (vbavg(i-1,j, mn)*depthv(i-1,j)
& +vbavg(i-1,j+1,mn)*depthv(i-1,j+1)))*
& (0.125*(pvtrop(i,j)+pvtrop(i,j+1)))
c
ubavg(i,j,nl)=
& ((1.-wblpf)*ubavg(i,j,ml)+
& wblpf *ubavg(i,j,nl) )+
& (1.+wblpf)*dlt*(utndcy+utotn(i,j))+
& 0.5*(util5(i,j)+util5(i-1,j))
c
* if (ldebug_barotp .and. i.eq.itest.and.j.eq.jtest) then
* write (lp,'(i9,2i5,i3,3x,a,6g15.6)')
* & nstep,i+i0,j+j0,lll+1,
* & 'u_old,u_new,p_grad,corio,u_star,drag =',
* & ubavg(i,j,ml),ubavg(i,j,nl),
* & -thref*(pbavg(i,j,nl)-pbavg(i-1,j,nl))*scuxi(i,j)*dlt,
* & (vbavg(i ,j, mn)*depthv(i ,j)
* & +vbavg(i ,j+1,mn)*depthv(i ,j+1)
* & +vbavg(i-1,j, mn)*depthv(i-1,j)
* & +vbavg(i-1,j+1,mn)*depthv(i-1,j+1))
* & *(pvtrop(i,j)+pvtrop(i,j+1))
* & *.125 * dlt,utotn(i,j) * dlt,
* & 0.5*(util5(i,j)+util5(i-1,j))
* endif !debug
endif !iu
enddo !i
enddo !j
c
if (ldebug_barotp) then
call xcsync(flush_lp)
endif
c
#if ! defined(RELO)
if (lbflag.eq.1) then
call latbdp( nl)
elseif (lbflag.eq.3) then
call latbdf( nl,lll+1)
endif
#endif
if (lbflag.eq.2) then
call latbdt( nl,lll+1)
elseif (lbflag.eq.4) then
call latbdtf(nl,lll+1)
endif
c
840 continue ! lll=1,lstep1,2
c
if (ldrag) then !disp_count updated in momtum
displd_mn(:,:) = displd_mn(:,:) + displd(:,:)/real(lstep1)
endif
c
if (lbflag.eq.1) then
c
c --- correct mean height.
c --- this should not be required - so there may be a bug in the bc.
c
!$OMP PARALLEL DO PRIVATE(j,i)
!$OMP& SCHEDULE(STATIC,jblk)
do j=1,jj
do i=1-margin,ii+margin
if (SEA_P) then
util1(i,j) = pbavg(i,j,n)*scp2(i,j)
endif !ip
enddo !i
enddo !j
call xcsum(sump, util1,ipa)
q = sump/area
c
c --- rhs: pbavg
c --- lhs: pbavg
c
margin = 0
c
!$OMP PARALLEL DO PRIVATE(j,i)
!$OMP& SCHEDULE(STATIC,jblk)
do j=1-margin,jj+margin
do i=1-margin,ii+margin
if (SEA_P) then
pbavo(i,j) = pbavo(i,j) - q
pbavg(i,j,1) = pbavg(i,j,1) - q
pbavg(i,j,2) = pbavg(i,j,2) - q
pbavg(i,j,3) = pbavg(i,j,3) - q
endif !ip
enddo !i
enddo !j
endif
if (lpipe .and. lpipe_barotp) then
call pipe_compare(pbavg(1-nbdy,1-nbdy,1), ip,'barotp:pbav1')
call pipe_compare(pbavg(1-nbdy,1-nbdy,2), ip,'barotp:pbav2')
call pipe_compare(pbavg(1-nbdy,1-nbdy,3), ip,'barotp:pbav3')
endif
c
if (btrlfr .and. delt1.ne.baclin) then !not on very 1st time step
c --- Robert-Asselin time filter
!$OMP PARALLEL DO PRIVATE(j,i,q)
!$OMP& SCHEDULE(STATIC,jblk)
do j=1,jj
do i=1,ii
q = 0.5*ra2fac*( pbavo(i,j) + !t-1
& pbavg(i,j,n) - !t+1
& 2.0*pbavg(i,j,m) ) !t
pbavg(i,j,m) = pbavg(i,j,m) + q
q = 0.5*ra2fac*( ubavo(i,j) + !t-1
& ubavg(i,j,n) - !t+1
& 2.0*ubavg(i,j,m) ) !t
ubavg(i,j,m) = ubavg(i,j,m) + q
q = 0.5*ra2fac*( vbavo(i,j) + !t-1
& vbavg(i,j,n) - !t+1
& 2.0*vbavg(i,j,m) ) !t
vbavg(i,j,m) = vbavg(i,j,m) + q
enddo !i
enddo !j
endif !btrlfr & not on very 1st time step
c
return
end subroutine barotp
c
c> Revision history:
c>
c> Mar. 1995 - changed vertical velocity averaging interval from 10 cm to 1 m
c> (loops 33,35)
c> Mar. 1995 - changed order of loop nesting in loop 842
c> July 1997 - eliminated 3-D arrays -uold,vold- (used in time smoothing)
c> Aug. 1997 - transferred loops preceding loop 840 to momeq2.f
c> Jan. 2000 - added latbdp for lateral boundary ports
c> Aug. 2001 - two barotropic time steps per loop, for halo efficiency
c> Nov. 2006 - added lbflag==3 (latbdf) and thref_bt (mod_tides)
c> Nov. 2006 - removed thref_bt (and mod_tides)
c> Apr. 2007 - added btrlfr: leapfrog time step; see also momtum
c> Apr. 2010 - bugfixes for 1st time step and 1st miner time step
c> Apr 2011 - added Robert-Asselin filtering for btrlfr
c> Aug 2011 - reworked Robert-Asselin filtering for btrlfr
c> Mar. 2012 - added latbdtf for nesting with Flather b.c.'s.
c> Jan. 2013 - added tidal drag tensor
c> June 2013 - added lbflag==6 for latbdtc
c> Apr. 2014 - replace ip with ipa for mass sums
c> May 2014 - use land/sea masks (e.g. ip) to skip land
c> May 2014 - removed lbflag==6 for latbdtc