ROMS
Loading...
Searching...
No Matches
ad_set_vbc.F
Go to the documentation of this file.
1#include "cppdefs.h"
3#ifdef ADJOINT
4!
5!git $Id$
6!================================================== Hernan G. Arango ===
7! Copyright (c) 2002-2025 The ROMS Group Andrew M. Moore !
8! Licensed under a MIT/X style license !
9! See License_ROMS.md !
10!=======================================================================
11! !
12! This module sets vertical boundary conditons for adjoint momentum !
13! and tracers. !
14! !
15! BASIC STATE variables needed: btflx, stflx, dqdt, t, sss, !
16! z_r, v, u !
17! (btflx and stflx are over written) !
18! !
19!=======================================================================
20!
21 implicit none
22!
23 PRIVATE
24 PUBLIC :: ad_set_vbc
25!
26 CONTAINS
27
28# ifdef SOLVE3D
29!
30!***********************************************************************
31 SUBROUTINE ad_set_vbc (ng, tile)
32!***********************************************************************
33!
34 USE mod_param
35 USE mod_grid
36 USE mod_forces
37 USE mod_ocean
38 USE mod_stepping
39!
40! Imported variable declarations.
41!
42 integer, intent(in) :: ng, tile
43!
44! Local variable declarations.
45!
46 character (len=*), parameter :: myfile = &
47 & __FILE__
48!
49# include "tile.h"
50!
51# ifdef PROFILE
52 CALL wclock_on (ng, iadm, 6, __line__, myfile)
53# endif
54 CALL ad_set_vbc_tile (ng, tile, &
55 & lbi, ubi, lbj, ubj, &
56 & imins, imaxs, jmins, jmaxs, &
57# if defined SCORRECTION || defined SRELAXATION
58 & grid(ng) % Hz, &
59 & grid(ng) % ad_Hz, &
60# endif
61# if defined UV_LOGDRAG
62 & grid(ng) % ZoBot, &
63# elif defined UV_LDRAG
64 & grid(ng) % rdrag, &
65# elif defined UV_QDRAG
66 & grid(ng) % rdrag2, &
67# endif
68# if defined UV_LOGDRAG && (!defined BBL_MODEL_NOT_YET || defined ICESHELF)
69 & grid(ng) % z_r, &
70 & grid(ng) % ad_z_r, &
71 & grid(ng) % z_w, &
72 & grid(ng) % ad_z_w, &
73# endif
74# if defined ICESHELF
75 & grid(ng) % zice, &
76# endif
77# if defined QCORRECTION || defined SALINITY
78 & ocean(ng) % t, &
79 & ocean(ng) % ad_t, &
80# endif
81# if !defined BBL_MODEL_NOT_YET || defined ICESHELF
82 & ocean(ng) % u, &
83 & ocean(ng) % ad_u, &
84 & ocean(ng) % v, &
85 & ocean(ng) % ad_v, &
86# endif
87# ifdef QCORRECTION
88 & forces(ng) % dqdt, &
89 & forces(ng) % sst, &
90# endif
91# if defined SCORRECTION || defined SRELAXATION
92 & forces(ng) % sss, &
93# endif
94# if defined ICESHELF
95# ifdef SHORTWAVE
96 & forces(ng) % srflx, &
97# endif
98 & forces(ng) % ad_sustr, &
99 & forces(ng) % ad_svstr, &
100# endif
101# ifndef BBL_MODEL_NOT_YET
102 & forces(ng) % ad_bustr, &
103 & forces(ng) % ad_bustr_sol, &
104 & forces(ng) % ad_bvstr, &
105 & forces(ng) % ad_bvstr_sol, &
106# endif
107 & forces(ng) % stflux, &
108 & forces(ng) % btflux, &
109# if defined QCORRECTION || defined SALINITY
110 & forces(ng) % stflx, &
111 & forces(ng) % ad_stflx, &
112# endif
113# ifdef SALINITY
114 & forces(ng) % btflx, &
115 & forces(ng) % ad_btflx, &
116# endif
117 & nrhs(ng))
118# ifdef PROFILE
119 CALL wclock_off (ng, iadm, 6, __line__, myfile)
120# endif
121!
122 RETURN
123 END SUBROUTINE ad_set_vbc
124!
125!***********************************************************************
126 SUBROUTINE ad_set_vbc_tile (ng, tile, &
127 & LBi, UBi, LBj, UBj, &
128 & IminS, ImaxS, JminS, JmaxS, &
129# if defined SCORRECTION || defined SRELAXATION
130 & Hz, ad_Hz, &
131# endif
132# if defined UV_LOGDRAG
133 & ZoBot, &
134# elif defined UV_LDRAG
135 & rdrag, &
136# elif defined UV_QDRAG
137 & rdrag2, &
138# endif
139# if defined UV_LOGDRAG && (!defined BBL_MODEL_NOT_YET || defined ICESHELF)
140 & z_r, ad_z_r, &
141 & z_w, ad_z_w, &
142# endif
143# if defined ICESHELF
144 & zice, &
145# endif
146# if defined QCORRECTION || defined SALINITY
147 & t, ad_t, &
148# endif
149# if !defined BBL_MODEL_NOT_YET || defined ICESHELF
150 & u, ad_u, &
151 & v, ad_v, &
152# endif
153# ifdef QCORRECTION
154 & dqdt, sst, &
155# endif
156# if defined SCORRECTION || defined SRELAXATION
157 & sss, &
158# endif
159# if defined ICESHELF
160# ifdef SHORTWAVE
161 & srflx, &
162# endif
163 & ad_sustr, ad_svstr, &
164# endif
165# ifndef BBL_MODEL_NOT_YET
166 & ad_bustr, ad_bustr_sol, &
167 & ad_bvstr, ad_bvstr_sol, &
168# endif
169 & stflux, btflux, &
170# if defined QCORRECTION || defined SALINITY
171 & stflx, ad_stflx, &
172# endif
173# ifdef SALINITY
174 & btflx, ad_btflx, &
175# endif
176 & nrhs)
177!***********************************************************************
178!
179 USE mod_param
180 USE mod_scalars
181!
182 USE ad_bc_2d_mod
183# ifdef DISTRIBUTE
185# endif
186!
187! Imported variable declarations.
188!
189 integer, intent(in) :: ng, tile
190 integer, intent(in) :: LBi, UBi, LBj, UBj
191 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
192 integer, intent(in) :: nrhs
193!
194# ifdef ASSUMED_SHAPE
195# if defined SCORRECTION || defined SRELAXATION
196 real(r8), intent(in) :: Hz(LBi:,LBj:,:)
197# endif
198# if defined UV_LOGDRAG
199 real(r8), intent(in) :: ZoBot(LBi:,LBj:)
200# elif defined UV_LDRAG
201 real(r8), intent(in) :: rdrag(LBi:,LBj:)
202# elif defined UV_QDRAG
203 real(r8), intent(in) :: rdrag2(LBi:,LBj:)
204# endif
205# if defined UV_LOGDRAG && (!defined BBL_MODEL_NOT_YET || defined ICESHELF)
206 real(r8), intent(in) :: z_r(LBi:,LBj:,:)
207 real(r8), intent(in) :: z_w(LBi:,LBj:,0:)
208# endif
209# if defined ICESHELF
210 real(r8), intent(in) :: zice(LBi:,LBj:)
211# endif
212# if defined QCORRECTION || defined SALINITY
213 real(r8), intent(in) :: t(LBi:,LBj:,:,:,:)
214# endif
215# if !defined BBL_MODEL_NOT_YET || defined ICESHELF
216 real(r8), intent(in) :: u(LBi:,LBj:,:,:)
217 real(r8), intent(in) :: v(LBi:,LBj:,:,:)
218# endif
219# if defined QCORRECTION || defined SALINITY
220 real(r8), intent(in) :: stflx(LBi:,LBj:,:)
221# endif
222# ifdef SALINITY
223 real(r8), intent(in) :: btflx(LBi:,LBj:,:)
224# endif
225# ifdef QCORRECTION
226 real(r8), intent(in) :: dqdt(LBi:,LBj:)
227 real(r8), intent(in) :: sst(LBi:,LBj:)
228# endif
229# if defined SCORRECTION || defined SRELAXATION
230 real(r8), intent(in) :: sss(LBi:,LBj:)
231# endif
232 real(r8), intent(in) :: stflux(LBi:,LBj:,:)
233 real(r8), intent(in) :: btflux(LBi:,LBj:,:)
234# if defined SCORRECTION || defined SRELAXATION
235 real(r8), intent(inout) :: ad_Hz(LBi:,LBj:,:)
236# endif
237# if defined QCORRECTION || defined SALINITY
238 real(r8), intent(inout) :: ad_t(LBi:,LBj:,:,:,:)
239# endif
240# if !defined BBL_MODEL_NOT_YET || defined ICESHELF
241 real(r8), intent(inout) :: ad_u(LBi:,LBj:,:,:)
242 real(r8), intent(inout) :: ad_v(LBi:,LBj:,:,:)
243# endif
244# if defined UV_LOGDRAG && (!defined BBL_MODEL_NOT_YET || defined ICESHELF)
245 real(r8), intent(inout) :: ad_z_r(LBi:,LBj:,:)
246 real(r8), intent(inout) :: ad_z_w(LBi:,LBj:,0:)
247# endif
248# if defined ICESHELF
249# ifdef SHORTWAVE
250 real(r8), intent(inout) :: srflx(LBi:,LBj:)
251# endif
252 real(r8), intent(inout) :: ad_sustr(LBi:,LBj:)
253 real(r8), intent(inout) :: ad_svstr(LBi:,LBj:)
254# endif
255# ifndef BBL_MODEL_NOT_YET
256 real(r8), intent(inout) :: ad_bustr(LBi:,LBj:)
257 real(r8), intent(inout) :: ad_bvstr(LBi:,LBj:)
258# endif
259# if defined QCORRECTION || defined SALINITY
260 real(r8), intent(inout) :: ad_stflx(LBi:,LBj:,:)
261# endif
262# ifdef SALINITY
263 real(r8), intent(inout) :: ad_btflx(LBi:,LBj:,:)
264# endif
265# ifndef BBL_MODEL_NOT_YET
266 real(r8), intent(out) :: ad_bustr_sol(LBi:,LBj:)
267 real(r8), intent(out) :: ad_bvstr_sol(LBi:,LBj:)
268# endif
269# else
270# if defined SCORRECTION || defined SRELAXATION
271 real(r8), intent(in) :: Hz(LBi:UBi,LBj:UBj,N(ng))
272# endif
273# if defined UV_LOGDRAG
274 real(r8), intent(in) :: ZoBot(LBi:UBi,LBj:UBj)
275# elif defined UV_LDRAG
276 real(r8), intent(in) :: rdrag(LBi:UBi,LBj:UBj)
277# elif defined UV_QDRAG
278 real(r8), intent(in) :: rdrag2(LBi:UBi,LBj:UBj)
279# endif
280# if defined UV_LOGDRAG && (!defined BBL_MODEL_NOT_YET || defined ICESHELF)
281 real(r8), intent(in) :: z_r(LBi:UBi,LBj:UBj,N(ng))
282 real(r8), intent(in) :: z_w(LBi:UBi,LBj:UBj,0:N(ng))
283# endif
284# if defined ICESHELF
285 real(r8), intent(in) :: zice(LBi:UBi,LBj:UBj)
286# endif
287# if defined QCORRECTION || defined SALINITY
288 real(r8), intent(in) :: t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
289# endif
290# if !defined BBL_MODEL_NOT_YET || defined ICESHELF
291 real(r8), intent(in) :: u(LBi:UBi,LBj:UBj,N(ng),2)
292 real(r8), intent(in) :: v(LBi:UBi,LBj:UBj,N(ng),2)
293# endif
294# if defined QCORRECTION || defined SALINITY
295 real(r8), intent(in) :: stflx(LBi:UBi,LBj:UBj,NT(ng))
296# endif
297# ifdef SALINITY
298 real(r8), intent(in) :: btflx(LBi:UBi,LBj:UBj,NT(ng))
299# endif
300# ifdef QCORRECTION
301 real(r8), intent(in) :: dqdt(LBi:UBi,LBj:UBj)
302 real(r8), intent(in) :: sst(LBi:UBi,LBj:UBj)
303# endif
304# if defined SCORRECTION || defined SRELAXATION
305 real(r8), intent(in) :: sss(LBi:UBi,LBj:UBj)
306# endif
307 real(r8), intent(in) :: stflx(LBi:UBi,LBj:UBj,NT(ng))
308 real(r8), intent(in) :: btflx(LBi:UBi,LBj:UBj,NT(ng))
309# if defined SCORRECTION || defined SRELAXATION
310 real(r8), intent(inout) :: ad_Hz(LBi:UBi,LBj:UBj,N(ng))
311# endif
312# if defined QCORRECTION || defined SALINITY
313 real(r8), intent(inout) :: ad_t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
314# endif
315# if !defined BBL_MODEL_NOT_YET || defined ICESHELF
316 real(r8), intent(inout) :: ad_u(LBi:UBi,LBj:UBj,N(ng),2)
317 real(r8), intent(inout) :: ad_v(LBi:UBi,LBj:UBj,N(ng),2)
318# endif
319# if defined UV_LOGDRAG && (!defined BBL_MODEL_NOT_YET || defined ICESHELF)
320 real(r8), intent(inout) :: ad_z_r(LBi:UBi,LBj:UBj,N(ng))
321 real(r8), intent(inout) :: ad_z_w(LBi:UBi,LBj:UBj,0:N(ng))
322# endif
323# if defined ICESHELF
324# ifdef SHORTWAVE
325 real(r8), intent(inout) :: srflx(LBi:UBi,LBj:UBj)
326# endif
327 real(r8), intent(inout) :: ad_sustr(LBi:UBi,LBj:UBj)
328 real(r8), intent(inout) :: ad_svstr(LBi:UBi,LBj:UBj)
329# endif
330# ifndef BBL_MODEL_NOT_YET
331 real(r8), intent(inout) :: ad_bustr(LBi:UBi,LBj:UBj)
332 real(r8), intent(inout) :: ad_bvstr(LBi:UBi,LBj:UBj)
333# endif
334# if defined QCORRECTION || defined SALINITY
335 real(r8), intent(inout) :: ad_stflx(LBi:UBi,LBj:UBj,NT(ng))
336# endif
337# ifdef SALINITY
338 real(r8), intent(inout) :: ad_btflx(LBi:UBi,LBj:UBj,NT(ng))
339# endif
340# ifndef BBL_MODEL_NOT_YET
341 real(r8), intent(inout) :: ad_bustr_sol(LBi:UBi,LBj:UBj)
342 real(r8), intent(inout) :: ad_bvstr_sol(LBi:UBi,LBj:UBj)
343# endif
344# endif
345!
346! Local variable declarations.
347!
348 integer :: i, j, itrc
349!
350 real(r8) :: EmP, ad_EmP
351 real(r8) :: cff1, cff2, cff3
352 real(r8) :: ad_cff1, ad_cff2, ad_cff3, adfac, adfac1, adfac2
353!
354# if (!defined BBL_MODEL_NOT_YET || defined ICESHELF) && \
355 defined uv_logdrag
356 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: wrk
357 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_wrk
358# endif
359
360# include "set_bounds.h"
361!
362!-----------------------------------------------------------------------
363! Initialize adjoint private variables.
364!-----------------------------------------------------------------------
365!
366 ad_emp=0.0_r8
367 ad_cff1=0.0_r8
368 ad_cff2=0.0_r8
369 ad_cff3=0.0_r8
370
371# if (!defined BBL_MODEL_NOT_YET || defined ICESHELF) && \
372 defined uv_logdrag
373 DO j=jmins,jmaxs
374 DO i=imins,imaxs
375 ad_wrk(i,j)=0.0_r8
376 END DO
377 END DO
378# endif
379
380# ifndef BBL_MODEL_NOT_YET
381!
382!-----------------------------------------------------------------------
383! Set kinematic bottom momentum flux (m2/s2).
384!-----------------------------------------------------------------------
385!
386! Apply boundary conditions.
387!
388# ifdef DISTRIBUTE
389!^ CALL mp_exchange2d (ng, tile, iTLM, 2, &
390!^ & LBi, UBi, LBj, UBj, &
391!^ & NghostPoints, &
392!^ & EWperiodic(ng), NSperiodic(ng), &
393!^ & tl_bustr, tl_bvstr)
394!^
395 CALL ad_mp_exchange2d (ng, tile, iadm, 2, &
396 & lbi, ubi, lbj, ubj, &
397 & nghostpoints, &
398 & ewperiodic(ng), nsperiodic(ng), &
399 & ad_bustr, ad_bvstr)
400# endif
401!^ CALL bc_u2d_tile (ng, tile, &
402!^ & LBi, UBi, LBj, UBj, &
403!^ & tl_bustr)
404!^
405 CALL ad_bc_u2d_tile (ng, tile, &
406 & lbi, ubi, lbj, ubj, &
407 & ad_bustr)
408!^ CALL bc_v2d_tile (ng, tile, &
409!^ & LBi, UBi, LBj, UBj, &
410!^ & tl_bvstr)
411!^
412 CALL ad_bc_v2d_tile (ng, tile, &
413 & lbi, ubi, lbj, ubj, &
414 & ad_bvstr)
415!
416! Save adjoint bottom momentum flux for output purposes.
417!
418 DO j=jstrr,jendr
419 DO i=istr,iendr
420 ad_bustr_sol(i,j)=ad_bustr(i,j)
421 END DO
422 IF (j.ge.jstr) THEN
423 DO i=istr,iendr
424 ad_bvstr_sol(i,j)=ad_bvstr(i,j)
425 END DO
426 END IF
427 END DO
428!
429!-----------------------------------------------------------------------
430! Set kinematic bottom momentum flux (m2/s2).
431!-----------------------------------------------------------------------
432# if defined UV_LOGDRAG
433!
434! Set logarithmic bottom stress.
435!
436 DO j=jstrv-1,jend
437 DO i=istru-1,iend
438 cff1=1.0_r8/log((z_r(i,j,1)-z_w(i,j,0))/zobot(i,j))
439 cff2=vonkar*vonkar*cff1*cff1
440 wrk(i,j)=min(cdb_max,max(cdb_min,cff2))
441 END DO
442 END DO
443 DO j=jstrv,jend
444 DO i=istr,iend
445 cff1=0.25_r8*(u(i ,j ,1,nrhs)+ &
446 & u(i+1,j ,1,nrhs)+ &
447 & u(i ,j-1,1,nrhs)+ &
448 & u(i+1,j-1,1,nrhs))
449 cff2=sqrt(cff1*cff1+v(i,j,1,nrhs)*v(i,j,1,nrhs))
450!^ tl_bvstr(i,j)=0.5_r8* &
451!^ & ((tl_wrk(i,j-1)+tl_wrk(i,j))* &
452!^ & v(i,j,1,nrhs)*cff2+ &
453!^ & (wrk(i,j-1)+wrk(i,j))* &
454!^ & (tl_v(i,j,1,nrhs)*cff2+ &
455!^ & v(i,j,1,nrhs)*tl_cff2))
456!^
457 adfac=0.5_r8*ad_bvstr(i,j)
458 adfac1=adfac*v(i,j,1,nrhs)*cff2
459 adfac2=adfac*(wrk(i,j-1)+wrk(i,j))
460 ad_wrk(i,j-1)=ad_wrk(i,j-1)+adfac1
461 ad_wrk(i,j )=ad_wrk(i,j )+adfac1
462 ad_v(i,j,1,nrhs)=ad_v(i,j,1,nrhs)+adfac2*cff2
463 ad_cff2=ad_cff2+adfac2*v(i,j,1,nrhs)
464 ad_bvstr(i,j)=0.0_r8
465 IF (cff2.ne.0.0_r8) THEN
466!< tl_cff2=(cff1*tl_cff1+v(i,j,1,nrhs)*tl_v(i,j,1,nrhs))/cff2
467!<
468 adfac=ad_cff2/cff2
469 ad_v(i,j,1,nrhs)=ad_v(i,j,1,nrhs)+adfac*v(i,j,1,nrhs)
470 ad_cff1=ad_cff1+adfac*cff1
471 ad_cff2=0.0_r8
472 ELSE
473!^ tl_cff2=0.0_r8
474!^
475 ad_cff2=0.0_r8
476 END IF
477!^ tl_cff1=0.25_r8*(tl_u(i ,j ,1,nrhs)+ &
478!^ & tl_u(i+1,j ,1,nrhs)+ &
479!^ & tl_u(i ,j-1,1,nrhs)+ &
480!^ & tl_u(i+1,j-1,1,nrhs))
481!^
482 adfac=0.25_r8*ad_cff1
483 ad_u(i ,j-1,1,nrhs)=ad_u(i ,j-1,1,nrhs)+adfac
484 ad_u(i+1,j-1,1,nrhs)=ad_u(i+1,j-1,1,nrhs)+adfac
485 ad_u(i ,j ,1,nrhs)=ad_u(i ,j ,1,nrhs)+adfac
486 ad_u(i+1,j ,1,nrhs)=ad_u(i+1,j ,1,nrhs)+adfac
487 ad_cff1=0.0_r8
488 END DO
489 END DO
490 DO j=jstr,jend
491 DO i=istru,iend
492 cff1=0.25_r8*(v(i ,j ,1,nrhs)+ &
493 & v(i ,j+1,1,nrhs)+ &
494 & v(i-1,j ,1,nrhs)+ &
495 & v(i-1,j+1,1,nrhs))
496 cff2=sqrt(u(i,j,1,nrhs)*u(i,j,1,nrhs)+cff1*cff1)
497!^ tl_bustr(i,j)=0.5_r8* &
498!^ & ((tl_wrk(i-1,j)+tl_wrk(i,j))* &
499!^ & u(i,j,1,nrhs)*cff2+ &
500!^ & (wrk(i-1,j)+wrk(i,j))* &
501!^ & (tl_u(i,j,1,nrhs)*cff2+ &
502!^ & u(i,j,1,nrhs)*tl_cff2))
503!^
504 adfac=0.5_r8*ad_bustr(i,j)
505 adfac1=adfac*u(i,j,1,nrhs)*cff2
506 adfac2=adfac*(wrk(i-1,j)+wrk(i,j))
507 ad_wrk(i-1,j)=ad_wrk(i-1,j)+adfac1
508 ad_wrk(i ,j)=ad_wrk(i ,j)+adfac1
509 ad_u(i,j,1,nrhs)=ad_u(i,j,1,nrhs)+adfac2*cff2
510 ad_cff2=ad_cff2+adfac1*u(i,j,1,nrhs)
511 ad_bustr(i,j)=0.0_r8
512 IF (cff2.ne.0.0_r8) THEN
513!^ tl_cff2=(u(i,j,1,nrhs)*tl_u(i,j,1,nrhs)+cff1*tl_cff1)/cff2
514!^
515 adfac=ad_cff2/cff2
516 ad_u(i,j,1,nrhs)=ad_u(i,j,1,nrhs)+adfac*u(i,j,1,nrhs)
517 ad_cff1=ad_cff1+adfac*cff1
518 ad_cff2=0.0_r8
519 ELSE
520!^ tl_cff2=0.0_r8
521!^
522 ad_cff2=0.0_r8
523 END IF
524!^ tl_cff1=0.25_r8*(tl_v(i ,j ,1,nrhs)+ &
525!^ & tl_v(i ,j+1,1,nrhs)+ &
526!^ & tl_v(i-1,j ,1,nrhs)+ &
527!^ & tl_v(i-1,j+1,1,nrhs))
528!^
529 adfac=0.25_r8*ad_cff1
530 ad_v(i-1,j ,1,nrhs)=ad_v(i-1,j ,1,nrhs)+adfac
531 ad_v(i ,j ,1,nrhs)=ad_v(i ,j ,1,nrhs)+adfac
532 ad_v(i-1,j+1,1,nrhs)=ad_v(i-1,j+1,1,nrhs)+adfac
533 ad_v(i ,j+1,1,nrhs)=ad_v(i ,j+1,1,nrhs)+adfac
534 ad_cff1=0.0_r8
535 END DO
536 END DO
537 DO j=jstrv-1,jend
538 DO i=istru-1,iend
539 cff1=1.0_r8/log((z_r(i,j,1)-z_w(i,j,0))/zobot(i,j))
540 cff2=vonkar*vonkar*cff1*cff1
541 cff3=max(cdb_min,cff2)
542 wrk(i,j)=min(cdb_max,cff3)
543!^ tl_wrk(i,j)=(0.5_r8-SIGN(0.5_r8,cff3-Cdb_max))*tl_cff3
544!^
545 ad_cff3=ad_cff3+ &
546 & (0.5_r8-sign(0.5_r8,cff3-cdb_max))*ad_wrk(i,j)
547 ad_wrk(i,j)=0.0_r8
548!^ tl_cff3=(0.5_r8-SIGN(0.5_r8,Cdb_min-cff2))*tl_cff2
549!^
550 ad_cff2=ad_cff2+ &
551 & (0.5_r8-sign(0.5_r8,cdb_min-cff2))*ad_cff3
552 ad_cff3=0.0_r8
553!^ tl_cff2=vonKar*vonKar*2.0_r8*cff1*tl_cff1
554!^
555 ad_cff1=ad_cff1+vonkar*vonkar*2.0_r8*cff1*ad_cff2
556 ad_cff2=0.0_r8
557!^ tl_cff1=-cff1*cff1*(tl_z_r(i,j,1)-tl_z_w(i,j,0))/ &
558!^ & (z_r(i,j,1)-z_w(i,j,0))
559!^
560 adfac=-cff1*cff1*ad_cff1/(z_r(i,j,1)-z_w(i,j,0))
561 ad_z_r(i,j,1)=ad_z_r(i,j,1)+adfac
562 ad_z_w(i,j,0)=ad_z_w(i,j,0)-adfac
563 ad_cff1=0.0_r8
564 END DO
565 END DO
566# elif defined UV_QDRAG
567!
568! Set quadratic bottom stress.
569!
570 DO j=jstrv,jend
571 DO i=istr,iend
572 cff1=0.25_r8*(u(i ,j ,1,nrhs)+ &
573 & u(i+1,j ,1,nrhs)+ &
574 & u(i ,j-1,1,nrhs)+ &
575 & u(i+1,j-1,1,nrhs))
576 cff2=sqrt(cff1*cff1+v(i,j,1,nrhs)*v(i,j,1,nrhs))
577!^ tl_bvstr(i,j)=0.5_r8*(rdrag2(i,j-1)+rdrag2(i,j))* &
578!^ & (tl_v(i,j,1,nrhs)*cff2+ &
579!^ & v(i,j,1,nrhs)*tl_cff2)
580!^
581 adfac=0.5_r8*(rdrag2(i,j-1)+rdrag2(i,j))* &
582 & ad_bvstr(i,j)
583 ad_v(i,j,1,nrhs)=ad_v(i,j,1,nrhs)+adfac*cff2
584 ad_cff2=ad_cff2+adfac*v(i,j,1,nrhs)
585 ad_bvstr(i,j)=0.0_r8
586 IF (cff2.ne.0.0_r8) THEN
587!^ tl_cff2=(cff1*tl_cff1+v(i,j,1,nrhs)*tl_v(i,j,1,nrhs))/cff2
588!^
589 adfac=ad_cff2/cff2
590 ad_v(i,j,1,nrhs)=ad_v(i,j,1,nrhs)+adfac*v(i,j,1,nrhs)
591 ad_cff1=ad_cff1+adfac*cff1
592 ad_cff2=0.0_r8
593 ELSE
594!^ tl_cff2=0.0_r8
595!^
596 ad_cff2=0.0_r8
597 END IF
598!^ tl_cff1=0.25_r8*(tl_u(i ,j ,1,nrhs)+ &
599!^ & tl_u(i+1,j ,1,nrhs)+ &
600!^ & tl_u(i ,j-1,1,nrhs)+ &
601!^ & tl_u(i+1,j-1,1,nrhs))
602!^
603 adfac=0.25_r8*ad_cff1
604 ad_u(i ,j-1,1,nrhs)=ad_u(i ,j-1,1,nrhs)+adfac
605 ad_u(i+1,j-1,1,nrhs)=ad_u(i+1,j-1,1,nrhs)+adfac
606 ad_u(i ,j ,1,nrhs)=ad_u(i ,j ,1,nrhs)+adfac
607 ad_u(i+1,j ,1,nrhs)=ad_u(i+1,j ,1,nrhs)+adfac
608 ad_cff1=0.0_r8
609 END DO
610 END DO
611 DO j=jstr,jend
612 DO i=istru,iend
613 cff1=0.25_r8*(v(i ,j ,1,nrhs)+ &
614 & v(i ,j+1,1,nrhs)+ &
615 & v(i-1,j ,1,nrhs)+ &
616 & v(i-1,j+1,1,nrhs))
617 cff2=sqrt(u(i,j,1,nrhs)*u(i,j,1,nrhs)+cff1*cff1)
618!^ tl_bustr(i,j)=0.5_r8*(rdrag2(i-1,j)+rdrag2(i,j))* &
619!^ & (tl_u(i,j,1,nrhs)*cff2+ &
620!^ & u(i,j,1,nrhs)*tl_cff2)
621!^
622 adfac=0.5_r8*(rdrag2(i-1,j)+rdrag2(i,j))* &
623 & ad_bustr(i,j)
624 ad_u(i,j,1,nrhs)=ad_u(i,j,1,nrhs)+adfac*cff2
625 ad_cff2=ad_cff2+adfac*u(i,j,1,nrhs)
626 ad_bustr(i,j)=0.0_r8
627 IF (cff2.ne.0.0_r8) THEN
628!^ tl_cff2=(u(i,j,1,nrhs)*tl_u(i,j,1,nrhs)+cff1*tl_cff1)/cff2
629!^
630 adfac=ad_cff2/cff2
631 ad_u(i,j,1,nrhs)=ad_u(i,j,1,nrhs)+adfac*u(i,j,1,nrhs)
632 ad_cff1=ad_cff1+adfac*cff1
633 ad_cff2=0.0_r8
634 ELSE
635!^ tl_cff2=0.0_r8
636!^
637 ad_cff2=0.0_r8
638 END IF
639!^ tl_cff1=0.25_r8*(tl_v(i ,j ,1,nrhs)+ &
640!^ & tl_v(i ,j+1,1,nrhs)+ &
641!^ & tl_v(i-1,j ,1,nrhs)+ &
642!^ & tl_v(i-1,j+1,1,nrhs))
643!^
644 adfac=0.25_r8*ad_cff1
645 ad_v(i-1,j ,1,nrhs)=ad_v(i-1,j ,1,nrhs)+adfac
646 ad_v(i ,j ,1,nrhs)=ad_v(i ,j ,1,nrhs)+adfac
647 ad_v(i-1,j+1,1,nrhs)=ad_v(i-1,j+1,1,nrhs)+adfac
648 ad_v(i ,j+1,1,nrhs)=ad_v(i ,j+1,1,nrhs)+adfac
649 ad_cff1=0.0_r8
650 END DO
651 END DO
652# elif defined UV_LDRAG
653!
654! Set linear bottom stress.
655!
656 DO j=jstrv,jend
657 DO i=istr,iend
658!^ tl_bvstr(i,j)=0.5_r8*(rdrag(i,j-1)+rdrag(i,j))* &
659!^ & tl_v(i,j,1,nrhs)
660!^
661 ad_v(i,j,1,nrhs)=ad_v(i,j,1,nrhs)+ &
662 & 0.5_r8*(rdrag(i,j-1)+rdrag(i,j))* &
663 & ad_bvstr(i,j)
664 ad_bvstr(i,j)=0.0_r8
665 END DO
666 END DO
667 DO j=jstr,jend
668 DO i=istru,iend
669!^ tl_bustr(i,j)=0.5_r8*(rdrag(i-1,j)+rdrag(i,j))* &
670!^ & tl_u(i,j,1,nrhs)
671!^
672 ad_u(i,j,1,nrhs)=ad_u(i,j,1,nrhs)+ &
673 & 0.5_r8*(rdrag(i-1,j)+rdrag(i,j))* &
674 & ad_bustr(i,j)
675 ad_bustr(i,j)=0.0_r8
676 END DO
677 END DO
678# endif
679# endif
680
681# ifdef ICESHELF
682!
683!-----------------------------------------------------------------------
684! If ice shelf cavities, replace surface wind stress with ice shelf
685! cavity stress (m2/s2).
686!-----------------------------------------------------------------------
687!
688! Apply periodic or gradient boundary conditions for output
689! purposes only.
690!
691# ifdef DISTRIBUTE
692!^ CALL mp_exchange2d (ng, tile, iTLM, 2, &
693!^ & LBi, UBi, LBj, UBj, &
694!^ & NghostPoints, &
695!^ & EWperiodic(ng), NSperiodic(ng), &
696!^ & tl_sustr, tl_svstr)
697!^
698 CALL ad_mp_exchange2d (ng, tile, iadm, 2, &
699 & lbi, ubi, lbj, ubj, &
700 & nghostpoints, &
701 & ewperiodic(ng), nsperiodic(ng), &
702 & ad_sustr, ad_svstr)
703# endif
704!^ CALL bc_u2d_tile (ng, tile, &
705!^ & LBi, UBi, LBj, UBj, &
706!^ & tl_sustr)
707!^
708 CALL ad_bc_u2d_tile (ng, tile, &
709 & lbi, ubi, lbj, ubj, &
710 & ad_sustr)
711!^ CALL bc_v2d_tile (ng, tile, &
712!^ & LBi, UBi, LBj, UBj, &
713!^ & tl_svstr)
714!^
715 CALL ad_bc_v2d_tile (ng, tile, &
716 & lbi, ubi, lbj, ubj, &
717 & ad_svstr)
718
719# if defined UV_LOGDRAG
720!
721! Set logarithmic ice shelf cavity stress.
722!
723 DO j=jstrv-1,jend
724 DO i=istru-1,iend
725 cff1=1.0_r8/log((z_w(i,j,n(ng))-z_r(i,j,n(ng)))/zobot(i,j))
726 cff2=vonkar*vonkar*cff1*cff1
727 cff3=max(cdb_min,cff2)
728 wrk(i,j)=min(cdb_max,cff3)
729 END DO
730 END DO
731 DO j=jstrv,jend
732 DO i=istr,iend
733 IF (zice(i,j)*zice(i,j-1).ne.0.0_r8) THEN
734 cff1=0.25_r8*(u(i ,j ,n(ng),nrhs)+ &
735 & u(i+1,j ,n(ng),nrhs)+ &
736 & u(i ,j-1,n(ng),nrhs)+ &
737 & u(i+1,j-1,n(ng),nrhs))
738 & cff2=sqrt(cff1*cff1+v(i,j,n(ng),nrhs)*v(i,j,n(ng),nrhs))
739!^ tl_svstr(i,j)=-0.5_r8* &
740!^ & ((tl_wrk(i,j-1)+tl_wrk(i,j))* &
741!^ & v(i,j,N(ng),nrhs)*cff2+ &
742!^ & (wrk(i,j-1)+wrk(i,j))* &
743!^ & (tl_v(i,j,N(ng),nrhs)*cff2+ &
744!^ & v(i,j,N(ng),nrhs)*tl_cff2))
745!^
746 adfac=-0.5_r8*ad_svstr(i,j)
747 adfac1=adfac*v(i,j,n(ng),nrhs)*cff2
748 adfac2=adfac*(wrk(i,j-1)+wrk(i,j))
749 ad_wrk(i,j-1)=ad_wrk(i,j-1)+adfac1
750 ad_wrk(i,j )=ad_wrk(i,j )+adfac1
751 ad_v(i,j,n(ng),nrhs)=ad_v(i,j,n(ng),nrhs)+adfac2*cff2
752 ad_cff2=ad_cff2+adfac2*v(i,j,n(ng),nrhs)
753 ad_svstr(i,j)=0.0_r8
754 IF (cff2.ne.0.0_r8) THEN
755!^ tl_cff2=(cff1*tl_cff1+ &
756!^ & v(i,j,N(ng),nrhs)*tl_v(i,j,N(ng),nrhs))/cff2
757!^
758 adfac=ad_cff2/cff2
759 ad_v(i,j,n(ng),nrhs)=ad_v(i,j,n(ng),nrhs)+ &
760 & adfac*v(i,j,n(ng),nrhs)
761 ad_cff1=ad_cff1+adfac*cff1
762 ad_cff2=0.0_r8
763 ELSE
764!^ tl_cff2=0.0_r8
765!^
766 ad_cff2=0.0_r8
767 END IF
768!^ tl_cff1=0.25_r8*(tl_u(i ,j ,N(ng),nrhs)+ &
769!^ & tl_u(i+1,j ,N(ng),nrhs)+ &
770!^ & tl_u(i ,j-1,N(ng),nrhs)+ &
771!^ & tl_u(i+1,j-1,N(ng),nrhs))
772!^
773 adfac=0.25_r8*ad_cff1
774 ad_u(i ,j-1,n(ng),nrhs)=ad_u(i ,j-1,n(ng),nrhs)+adfac
775 ad_u(i+1,j-1,n(ng),nrhs)=ad_u(i+1,j-1,n(ng),nrhs)+adfac
776 ad_u(i ,j ,n(ng),nrhs)=ad_u(i ,j ,n(ng),nrhs)+adfac
777 ad_u(i+1,j ,n(ng),nrhs)=ad_u(i+1,j ,n(ng),nrhs)+adfac
778 ad_cff1=0.0_r8
779 END IF
780 END DO
781 END DO
782 DO j=jstr,jend
783 DO i=istru,iend
784 IF (zice(i,j)*zice(i-1,j).ne.0.0_r8) THEN
785 cff1=0.25_r8*(v(i ,j ,n(ng),nrhs)+ &
786 & v(i ,j+1,n(ng),nrhs)+ &
787 & v(i-1,j ,n(ng),nrhs)+ &
788 & v(i-1,j+1,n(ng),nrhs))
789 cff2=sqrt(u(i,j,n(ng),nrhs)*u(i,j,n(ng),nrhs)+cff1*cff1)
790!^ tl_sustr(i,j)=-0.5_r8* &
791!^ & ((tl_wrk(i-1,j)+tl_wrk(i,j))* &
792!^ & u(i,j,N(ng),nrhs)*cff2+ &
793!^ & (wrk(i-1,j)+wrk(i,j))* &
794!^ & (tl_u(i,j,N(ng),nrhs)*cff2+ &
795!^ & u(i,j,N(ng),nrhs)*tl_cff2))
796!^
797 adfac=-0.5_r8*ad_sustr(i,j)
798 adfac1=adfac*u(i,j,n(ng),nrhs)*cff2
799 adfac2=adfac*(wrk(i-1,j)+wrk(i,j))
800 ad_wrk(i-1,j)=ad_wrk(i-1,j)+adfac1
801 ad_wrk(i ,j)=ad_wrk(i ,j)+adfac1
802 ad_u(i,j,n(ng),nrhs)=ad_u(i,j,n(ng),nrhs)+adfac2*cff2
803 ad_cff2=ad_cff2+adfac2*u(i,j,n(ng),nrhs)
804 ad_sustr(i,j)=0.0_r8
805 IF (cff2.ne.0.0_r8) THEN
806!^ tl_cff2=(u(i,j,N(ng),nrhs)*tl_u(i,j,N(ng),nrhs)+ &
807!^ & cff1*tl_cff1)/cff2
808!^
809 adfac=ad_cff2/cff2
810 ad_u(i,j,n(ng),nrhs)=ad_u(i,j,n(ng),nrhs)+ &
811 & adfac*u(i,j,n(ng),nrhs)
812 ad_cff1=ad_cff1+adfac*cff1
813 ad_cff2=0.0_r8
814 ELSE
815!^ tl_cff2=0.0_r8
816!^
817 ad_cff2=0.0_r8
818 END IF
819!^ tl_cff1=0.25_r8*(tl_v(i ,j ,N(ng),nrhs)+ &
820!^ & tl_v(i ,j+1,N(ng),nrhs)+ &
821!^ & tl_v(i-1,j ,N(ng),nrhs)+ &
822!^ & tl_v(i-1,j+1,N(ng),nrhs))
823!^
824 adfac=0.25_r8*ad_cff1
825 ad_v(i-1,j ,n(ng),nrhs)=ad_v(i-1,j ,n(ng),nrhs)+adfac
826 ad_v(i ,j ,n(ng),nrhs)=ad_v(i ,j ,n(ng),nrhs)+adfac
827 ad_v(i-1,j+1,n(ng),nrhs)=ad_v(i-1,j+1,n(ng),nrhs)+adfac
828 ad_v(i ,j+1,n(ng),nrhs)=ad_v(i ,j+1,n(ng),nrhs)+adfac
829 ad_cff1=0.0_r8
830 END IF
831 END DO
832 END DO
833 DO j=jstrv-1,jend
834 DO i=istru-1,iend
835 cff1=1.0_r8/log((z_w(i,j,n(ng))-z_r(i,j,n(ng)))/zobot(i,j))
836 cff2=vonkar*vonkar*cff1*cff1
837 cff3=max(cdb_min,cff2)
838 wrk(i,j)=min(cdb_max,cff3)
839!^ tl_wrk(i,j)=(0.5_r8-SIGN(0.5_r8,cff3-Cdb_max))*tl_cff3
840!^
841 ad_cff3=ad_cff3+ &
842 & (0.5_r8-sign(0.5_r8,cff3-cdb_max))*ad_wrk(i,j)
843 ad_wrk(i,j)=0.0_r8
844!^ tl_cff3=(0.5_r8-SIGN(0.5_r8,Cdb_min-cff2))*tl_cff2
845!^
846 ad_cff2=ad_cff2+ &
847 & (0.5_r8-sign(0.5_r8,cdb_min-cff2))*ad_cff3
848 ad_cff3=0.0_r8
849!^ tl_cff2=vonKar*vonKar*2.0_r8*cff1*tl_cff1
850!^
851 ad_cff1=ad_cff1+vonkar*vonkar*2.0_r8*cff1*ad_cff2
852 ad_cff2=0.0_r8
853!^ tl_cff1=-cff1*cff1*(tl_z_w(i,j,N(ng))-tl_z_r(i,j,N(ng)))/ &
854!^ & (z_w(i,j,N(ng))-z_r(i,j,N(ng)))
855!^
856 adfac=-cff1*cff1*ad_cff1/(z_w(i,j,n(ng))-z_r(i,j,n(ng)))
857 ad_z_r(i,j,n(ng))=ad_z_r(i,j,n(ng))-adfac
858 ad_z_w(i,j,n(ng))=ad_z_w(i,j,n(ng))+adfac
859 ad_cff1=0.0_r8
860 END DO
861 END DO
862# elif defined UV_QDRAG
863!
864! Set quadratic ice shelf cavity stress.
865!
866 DO j=jstrv,jend
867 DO i=istr,iend
868 IF (zice(i,j)*zice(i,j-1).ne.0.0_r8) THEN
869 cff1=0.25_r8*(u(i ,j ,n(ng),nrhs)+ &
870 & u(i+1,j ,n(ng),nrhs)+ &
871 & u(i ,j-1,n(ng),nrhs)+ &
872 & u(i+1,j-1,n(ng),nrhs))
873 cff2=sqrt(cff1*cff1+v(i,j,n(ng),nrhs)*v(i,j,n(ng),nrhs))
874!^ tl_svstr(i,j)=-0.5_r8*(rdrag2(i,j-1)+rdrag2(i,j))* &
875!^ & (tl_v(i,j,N(ng),nrhs)*cff2+ &
876!^ & v(i,j,N(ng),nrhs)*tl_cff2)
877!^
878 adfac=-0.5_r8*(rdrag2(i,j-1)+rdrag2(i,j))*ad_svstr(i,j)
879 ad_v(i,j,n(ng),nrhs)=ad_v(i,j,n(ng),nrhs)+adfac*cff2
880 ad_cff2=ad_cff2+adfac*v(i,j,n(ng),nrhs)
881 ad_svstr(i,j)=0.0_r8
882 IF (cff2.ne.0.0_r8) THEN
883!^ tl_cff2=(cff1*tl_cff1+ &
884!^ & v(i,j,N(ng),nrhs)*tl_v(i,j,N(ng),nrhs))/cff2
885!^
886 adfac=ad_cff2/cff2
887 ad_v(i,j,n(ng),nrhs)=ad_v(i,j,n(ng),nrhs)+ &
888 & adfac*v(i,j,n(ng),nrhs)
889 ad_cff1=ad_cff1+adfac*cff1
890 ad_cff2=0.0_r8
891 ELSE
892!^ tl_cff2=0.0_r8
893!^
894 ad_cff2=0.0_r8
895 END IF
896!^ tl_cff1=0.25_r8*(tl_u(i ,j ,N(ng),nrhs)+ &
897!^ & tl_u(i+1,j ,N(ng),nrhs)+ &
898!^ & tl_u(i ,j-1,N(ng),nrhs)+ &
899!^ & tl_u(i+1,j-1,N(ng),nrhs))
900!^
901 adfac=0.25_r8*ad_cff1
902 ad_u(i ,j-1,n(ng),nrhs)=ad_u(i ,j-1,n(ng),nrhs)+adfac
903 ad_u(i+1,j-1,n(ng),nrhs)=ad_u(i+1,j-1,n(ng),nrhs)+adfac
904 ad_u(i ,j ,n(ng),nrhs)=ad_u(i ,j ,n(ng),nrhs)+adfac
905 ad_u(i+1,j ,n(ng),nrhs)=ad_u(i+1,j ,n(ng),nrhs)+adfac
906 ad_cff1=0.0_r8
907 END IF
908 END DO
909 END DO
910 DO j=jstr,jend
911 DO i=istru,iend
912 IF (zice(i,j)*zice(i-1,j).ne.0.0_r8) THEN
913 cff1=0.25_r8*(v(i ,j ,n(ng),nrhs)+ &
914 & v(i ,j+1,n(ng),nrhs)+ &
915 & v(i-1,j ,n(ng),nrhs)+ &
916 & v(i-1,j+1,n(ng),nrhs))
917 cff2=sqrt(u(i,j,n(ng),nrhs)*u(i,j,n(ng),nrhs)+cff1*cff1)
918!^ tl_sustr(i,j)=-0.5_r8*(rdrag2(i-1,j)+rdrag2(i,j))* &
919!^ & (tl_u(i,j,N(ng),nrhs)*cff2+ &
920!^ & u(i,j,N(ng),nrhs)*tl_cff2)
921!^
922 adfac=-0.5_r8*(rdrag2(i-1,j)+rdrag2(i,j))*ad_sustr(i,j)
923 ad_u(i,j,n(ng),nrhs)=ad_u(i,j,n(ng),nrhs)+adfac*cff2
924 ad_cff2=ad_cff2+adfac*u(i,j,n(ng),nrhs)
925 ad_sustr(i,j)=0.0_r8
926 IF (cff2.ne.0.0_r8) THEN
927!^ tl_cff2=(u(i,j,N(ng),nrhs)*tl_u(i,j,N(ng),nrhs)+ &
928!^ & cff1*tl_cff1)/cff2
929!^
930 adfac=ad_cff2/cff2
931 ad_u(i,j,n(ng),nrhs)=ad_u(i,j,n(ng),nrhs)+ &
932 & adfac*u(i,j,n(ng),nrhs)
933 ad_cff1=ad_cff1+adfac*cff1
934 ad_cff2=0.0_r8
935 ELSE
936!^ tl_cff2=0.0_r8
937!^
938 ad_cff2=0.0_r8
939 END IF
940!^ tl_cff1=0.25_r8*(tl_v(i ,j ,N(ng),nrhs)+ &
941!^ & tl_v(i ,j+1,N(ng),nrhs)+ &
942!^ & tl_v(i-1,j ,N(ng),nrhs)+ &
943!^ & tl_v(i-1,j+1,N(ng),nrhs))
944!^
945 adfac=0.25_r8*ad_cff1
946 ad_v(i-1,j ,n(ng),nrhs)=ad_v(i-1,j ,n(ng),nrhs)+adfac
947 ad_v(i ,j ,n(ng),nrhs)=ad_v(i ,j ,n(ng),nrhs)+adfac
948 ad_v(i-1,j+1,n(ng),nrhs)=ad_v(i-1,j+1,n(ng),nrhs)+adfac
949 ad_v(i ,j+1,n(ng),nrhs)=ad_v(i ,j+1,n(ng),nrhs)+adfac
950 ad_cff1=0.0_r8
951 END IF
952 END DO
953 END DO
954# elif defined UV_LDRAG
955!
956! Set linear ice shelf cavity stress.
957!
958 DO j=jstrv,jend
959 DO i=istr,iend
960 IF (zice(i,j)*zice(i,j-1).ne.0.0_r8) THEN
961!^ tl_svstr(i,j)=-0.5_r8*(rdrag(i,j-1)+rdrag(i,j))* &
962!^ & tl_v(i,j,N(ng),nrhs)
963!^
964 ad_v(i,j,n(ng),nrhs)=ad_v(i,j,n(ng),nrhs)- &
965 & 0.5_r8*(rdrag(i,j-1)+rdrag(i,j))* &
966 & ad_svstr(i,j)
967 ad_svstr(i,j)=0.0_r8
968 END IF
969 END DO
970 END DO
971 DO j=jstr,jend
972 DO i=istru,iend
973 IF (zice(i,j)*zice(i-1,j).ne.0.0_r8) THEN
974!^ tl_sustr(i,j)=-0.5_r8*(rdrag(i-1,j)+rdrag(i,j))* &
975!^ & tl_u(i,j,N(ng),nrhs)
976!^
977 ad_u(i,j,n(ng),nrhs)=ad_u(i,j,n(ng),nrhs)- &
978 & 0.5_r8*(rdrag(i-1,j)+rdrag(i,j))* &
979 & ad_sustr(i,j)
980 ad_sustr(i,j)=0.0_r8
981 END IF
982 END DO
983 END DO
984# else
985 DO j=jstrv,jend
986 DO i=istr,iend
987 IF (zice(i,j)*zice(i,j-1).ne.0.0_r8) THEN
988!^ tl_svstr(i,j)=0.0_r8
989!^
990 ad_svstr(i,j)=0.0_r8
991 END IF
992 END DO
993 END DO
994 DO j=jstr,jend
995 DO i=istru,iend
996 IF (zice(i,j)*zice(i-1,j).ne.0.0_r8) THEN
997!^ tl_sustr(i,j)=0.0_r8
998!^
999 ad_sustr(i,j)=0.0_r8
1000 END IF
1001 END DO
1002 END DO
1003# endif
1004!
1005!-----------------------------------------------------------------------
1006! If ice shelf cavities, zero out for now the surface tracer flux
1007! over the ice.
1008!-----------------------------------------------------------------------
1009!
1010# ifdef SHORTWAVE
1011 DO j=jstrr,jendr
1012 DO i=istrr,iendr
1013 IF (zice(i,j).ne.0.0_r8) THEN
1014!! srflx(i,j)=0.0_r8
1015 END IF
1016 END DO
1017 END DO
1018# endif
1019 DO itrc=1,nt(ng)
1020 DO j=jstrr,jendr
1021 DO i=istrr,iendr
1022 IF (zice(i,j).ne.0.0_r8) THEN
1023!^ tl_stflx(i,j,itrc)=0.0_r8
1024!^
1025 ad_stflx(i,j,itrc)=0.0_r8
1026 END IF
1027 END DO
1028 END DO
1029 END DO
1030# endif
1031
1032# if defined BIOLOGY || defined SEDIMENT || defined T_PASSIVE
1033!
1034!-----------------------------------------------------------------------
1035! Load surface and bottom passive tracer fluxes (T m/s).
1036!-----------------------------------------------------------------------
1037!
1038 DO itrc=nat+1,nt(ng)
1039 DO j=jstrr,jendr
1040 DO i=istrr,iendr
1041!^ tl_btflx(i,j,itrc)=0.0_r8
1042!^
1043 ad_btflx(i,j,itrc)=0.0_r8
1044!^ tl_stflx(i,j,itrc)=0.0_r8
1045!^
1046 ad_stflx(i,j,itrc)=0.0_r8
1047 END DO
1048 END DO
1049 END DO
1050# endif
1051
1052# ifdef SALINITY
1053!
1054!-----------------------------------------------------------------------
1055! Add salt flux correction or multiply flux by salinity.
1056!-----------------------------------------------------------------------
1057!
1058 DO j=jstrr,jendr
1059 DO i=istrr,iendr
1060 emp=stflux(i,j,isalt)
1061!^ tl_btflx(i,j,isalt)=btflx(i,j,isalt)* &
1062!^ & tl_t(i,j,1,nrhs,isalt)
1063!^
1064 ad_t(i,j,1,nrhs,isalt)=ad_t(i,j,1,nrhs,isalt)+ &
1065 & btflx(i,j,isalt)*ad_btflx(i,j,isalt)
1066# if !(defined AD_SENSITIVITY || defined I4DVAR_ANA_SENSITIVITY || \
1067 defined opt_observations || defined so_semi)
1068 ad_btflx(i,j,isalt)=0.0_r8
1069# endif
1070# if defined SCORRECTION
1071!^ tl_stflx(i,j,isalt)=(EmP*tl_t(i,j,N(ng),nrhs,isalt)+ &
1072!^ & tl_EmP*t(i,j,N(ng),nrhs,isalt))- &
1073!^ & Tnudg(isalt,ng)* &
1074!^ & (tl_Hz(i,j,N(ng))* &
1075!^ & (t(i,j,N(ng),nrhs,isalt)-sss(i,j))+ &
1076!^ & Hz(i,j,N(ng))* &
1077!^ & tl_t(i,j,N(ng),nrhs,isalt))
1078!^
1079 adfac=tnudg(isalt,ng)*ad_stflx(i,j,isalt)
1080 ad_hz(i,j,n(ng))=ad_hz(i,j,n(ng))- &
1081 & (t(i,j,n(ng),nrhs,isalt)-sss(i,j))*adfac
1082 ad_t(i,j,n(ng),nrhs,isalt)=ad_t(i,j,n(ng),nrhs,isalt)- &
1083 & hz(i,j,n(ng))*adfac
1084 ad_t(i,j,n(ng),nrhs,isalt)=ad_t(i,j,n(ng),nrhs,isalt)+ &
1085 & emp*ad_stflx(i,j,isalt)
1086 ad_emp=ad_emp+ &
1087 t(i,j,n(ng),nrhs,isalt)*ad_stflx(i,j,isalt)
1088# if !(defined ADJUST_STFLUX || defined AD_SENSITIVITY || \
1089 defined i4dvar_ana_sensitivity || defined opt_observations || \
1090 defined so_semi)
1091 ad_stflx(i,j,isalt)=0.0_r8
1092# endif
1093# elif defined SRELAXATION
1094!^ tl_stflx(i,j,isalt)=-Tnudg(isalt,ng)* &
1095!^ & (tl_Hz(i,j,N(ng))* &
1096!^ & (t(i,j,N(ng),nrhs,isalt)-sss(i,j))+ &
1097!^ & Hz(i,j,N(ng))* &
1098!^ & tl_t(i,j,N(ng),nrhs,isalt))
1099!^
1100 adfac=-tnudg(isalt,ng)*ad_stflx(i,j,isalt)
1101 ad_hz(i,j,n(ng))=ad_hz(i,j,n(ng))+ &
1102 & (t(i,j,n(ng),nrhs,isalt)-sss(i,j))*adfac
1103 ad_t(i,j,n(ng),nrhs,isalt)=ad_t(i,j,n(ng),nrhs,isalt)+ &
1104 & hz(i,j,n(ng))*adfac
1105# if !(defined ADJUST_STFLUX || defined AD_SENSITIVITY || \
1106 defined i4dvar_ana_sensitivity || defined opt_observations || \
1107 defined so_semi)
1108 ad_stflx(i,j,isalt)=0.0_r8
1109# endif
1110# else
1111!^ tl_stflx(i,j,isalt)=EmP*tl_t(i,j,N(ng),nrhs,isalt)+ &
1112!^ & tl_EmP*t(i,j,N(ng),nrhs,isalt)
1113!^
1114 ad_t(i,j,n(ng),nrhs,isalt)=ad_t(i,j,n(ng),nrhs,isalt)+ &
1115 & emp*ad_stflx(i,j,isalt)
1116 ad_emp=ad_emp+ &
1117 t(i,j,n(ng),nrhs,isalt)*ad_stflx(i,j,isalt)
1118# if !(defined ADJUST_STFLUX || defined AD_SENSITIVITY || \
1119 defined i4dvar_ana_sensitivity || defined opt_observations || \
1120 defined so_semi)
1121 ad_stflx(i,j,isalt)=0.0_r8
1122# endif
1123# endif
1124!^ tl_EmP=0.0_r8
1125!^
1126 ad_emp=0.0_r8
1127 END DO
1128 END DO
1129# endif
1130
1131# ifdef LIMIT_STFLX_COOLING
1132!
1133!-----------------------------------------------------------------------
1134! If net heat flux is cooling and SST is at freezing point or below
1135! then suppress further cooling. Note: stflx sign convention is that
1136! positive means heating the ocean (J Wilkin).
1137!-----------------------------------------------------------------------
1138!
1139! Below the surface heat flux stflx(:,:,itemp) is ZERO if cooling AND
1140! the SST is cooler that the threshold. The value is retained if
1141! warming.
1142!
1143! cff3 = 0 if SST warmer than threshold (cff1) - change nothing
1144! cff3 = 1 if SST colder than threshold (cff1)
1145!
1146! 0.5*(cff2-ABS(cff2)) = 0 if flux is warming
1147! = stflx(:,:,itemp) if flux is cooling
1148!
1149 cff1=-2.0_r8 ! nominal SST threshold to cease cooling
1150 DO j=jstrr,jendr
1151 DO i=istrr,iendr
1152 cff2=stflx(i,j,itemp)
1153 cff3=0.5_r8*(1.0_r8+sign(1.0_r8,cff1-t(i,j,n(ng),nrhs,itemp)))
1154!^ tl_stflx(i,j,itemp)=(1.0_r8- &
1155!^ & cff3*0.5_r8*(1.0_r8-SIGN(1.0_r8,cff2)))* &
1156!^ & tl_cff2
1157!^
1158 ad_cff2=ad_cff2+ &
1159 & (1.0_r8-cff3*0.5_r8*(1.0_r8-sign(1.0_r8,cff2)))* &
1160 & ad_stflx(i,j,itemp)
1161 ad_stflx(i,j,itemp)=0.0_r8
1162!^ tl_cff3=0.5_r8*SIGN(1.0_r8, cff1-t(i,j,N(ng),nrhs,itemp))*0.0
1163!^ tl_cff3=0.0_r8
1164!^
1165!^ tl_cff2=tl_stflx(i,j,itemp)
1166!^
1167 ad_stflx(i,j,itemp)=ad_stflx(i,j,itemp)+ad_cff2
1168 ad_cff2=0.0_r8
1169 END DO
1170 END DO
1171# endif
1172
1173# ifdef QCORRECTION
1174!
1175!-----------------------------------------------------------------------
1176! Add in flux correction to surface net heat flux (degC m/s).
1177!-----------------------------------------------------------------------
1178!
1179! Add in net heat flux correction.
1180!
1181 DO j=jstrr,jendr
1182 DO i=istrr,iendr
1183!^ tl_stflx(i,j,itemp)=tl_stflx(i,j,itemp)+ &
1184!^ & dqdt(i,j)*tl_t(i,j,N(ng),nrhs,itemp)
1185!^
1186 ad_t(i,j,n(ng),nrhs,itemp)=ad_t(i,j,n(ng),nrhs,itemp)+ &
1187 & dqdt(i,j)*ad_stflx(i,j,itemp)
1188# if !(defined ADJUST_STFLUX || defined AD_SENSITIVITY || \
1189 defined i4dvar_ana_sensitivity || defined opt_observations || \
1190 defined so_semi)
1191 ad_stflx(i,j,itemp)=0.0_r8
1192# endif
1193 END DO
1194 END DO
1195# endif
1196!
1197!-----------------------------------------------------------------------
1198! Adjoint of loading surface and bottom net heat flux (degC m/s) into
1199! state variables "stflx" and "btflx".
1200!
1201! Notice that the forcing net heat flux stflux(:,:,itemp) is processed
1202! elsewhere from data, coupling, bulk flux parameterization,
1203! or analytical formulas.
1204!
1205! During model coupling, we need to make sure that this forcing is
1206! unaltered during the coupling interval when ROMS timestep size is
1207! smaller. Notice that further manipulations to state variable
1208! stflx(:,:,itemp) are allowed below.
1209!-----------------------------------------------------------------------
1210!
1211 DO j=jstrr,jendr
1212 DO i=istrr,iendr
1213!^ tl_btflx(i,j,itemp)=0.0_r8 ! not needed in TLM
1214!^
1215!! ad_btflx(i,j,itemp)=0.0_r8 ! not needed in ADM
1216!^ tl_stflx(i,j,itemp)=0.0_r8 ! not needed in TLM
1217!^
1218!! ad_stflx(i,j,itemp)=0.0_r8 ! not needed in ADM
1219 END DO
1220 END DO
1221!
1222 RETURN
1223 END SUBROUTINE ad_set_vbc_tile
1224
1225# else
1226
1227!
1228!***********************************************************************
1229 SUBROUTINE ad_set_vbc (ng, tile)
1230!***********************************************************************
1231!
1232 USE mod_param
1233 USE mod_forces
1234 USE mod_grid
1235 USE mod_ocean
1236 USE mod_stepping
1237!
1238! Imported variable declarations.
1239!
1240 integer, intent(in) :: ng, tile
1241!
1242! Local variable declarations.
1243!
1244 character (len=*), parameter :: MyFile = &
1245 & __FILE__
1246!
1247# include "tile.h"
1248!
1249# ifdef PROFILE
1250 CALL wclock_on (ng, iadm, 6, __line__, myfile)
1251# endif
1252 CALL ad_set_vbc_tile (ng, tile, &
1253 & lbi, ubi, lbj, ubj, &
1254 & imins, imaxs, jmins, jmaxs, &
1255 & krhs(ng), kstp(ng), knew(ng), &
1256# if defined UV_LDRAG
1257 & grid(ng) % rdrag, &
1258# elif defined UV_QDRAG
1259 & grid(ng) % rdrag2, &
1260# endif
1261 & ocean(ng) % ubar, &
1262 & ocean(ng) % vbar, &
1263 & ocean(ng) % ad_ubar, &
1264 & ocean(ng) % ad_vbar, &
1265 & forces(ng) % ad_bustr, &
1266 & forces(ng) % ad_bustr_sol, &
1267 & forces(ng) % ad_bvstr, &
1268 & forces(ng) % ad_bvstr_sol)
1269# ifdef PROFILE
1270 CALL wclock_off (ng, iadm, 6, __line__, myfile)
1271# endif
1272!
1273 RETURN
1274 END SUBROUTINE ad_set_vbc
1275!
1276!***********************************************************************
1277 SUBROUTINE ad_set_vbc_tile (ng, tile, &
1278 & LBi, UBi, LBj, UBj, &
1279 & IminS, ImaxS, JminS, JmaxS, &
1280 & krhs, kstp, knew, &
1281# if defined UV_LDRAG
1282 & rdrag, &
1283# elif defined UV_QDRAG
1284 & rdrag2, &
1285# endif
1286 & ubar, vbar, &
1287 & ad_ubar, ad_vbar, &
1288 & ad_bustr, ad_bustr_sol, &
1289 & ad_bvstr, ad_bvstr_sol)
1290!***********************************************************************
1291!
1292 USE mod_param
1293 USE mod_scalars
1294!
1295 USE ad_bc_2d_mod
1296# ifdef DISTRIBUTE
1298# endif
1299!
1300! Imported variable declarations.
1301!
1302 integer, intent(in) :: ng, tile
1303 integer, intent(in) :: LBi, UBi, LBj, UBj
1304 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
1305 integer, intent(in) :: krhs, kstp, knew
1306!
1307# ifdef ASSUMED_SHAPE
1308# if defined UV_LDRAG
1309 real(r8), intent(in) :: rdrag(LBi:,LBj:)
1310# elif defined UV_QDRAG
1311 real(r8), intent(in) :: rdrag2(LBi:,LBj:)
1312# endif
1313 real(r8), intent(in) :: ubar(LBi:,LBj:,:)
1314 real(r8), intent(in) :: vbar(LBi:,LBj:,:)
1315 real(r8), intent(inout) :: ad_ubar(LBi:,LBj:,:)
1316 real(r8), intent(inout) :: ad_vbar(LBi:,LBj:,:)
1317 real(r8), intent(inout) :: ad_bustr(LBi:,LBj:)
1318 real(r8), intent(inout) :: ad_bvstr(LBi:,LBj:)
1319 real(r8), intent(out) :: ad_bustr_sol(LBi:,LBj:)
1320 real(r8), intent(out) :: ad_bvstr_sol(LBi:,LBj:)
1321# else
1322# if defined UV_LDRAG
1323 real(r8), intent(in) :: rdrag(LBi:UBi,LBj:UBj)
1324# elif defined UV_QDRAG
1325 real(r8), intent(in) :: rdrag2(LBi:UBi,LBj:UBj)
1326# endif
1327 real(r8), intent(in) :: ubar(LBi:UBi,LBj:UBj,:)
1328 real(r8), intent(in) :: vbar(LBi:UBi,LBj:UBj,:)
1329 real(r8), intent(inout) :: ad_ubar(LBi:UBi,LBj:UBj,:)
1330 real(r8), intent(inout) :: ad_vbar(LBi:UBi,LBj:UBj,:)
1331 real(r8), intent(inout) :: ad_bustr(LBi:UBi,LBj:UBj)
1332 real(r8), intent(inout) :: ad_bvstr(LBi:UBi,LBj:UBj)
1333 real(r8), intent(out) :: ad_bustr_sol(LBi:UBi,LBj:UBj)
1334 real(r8), intent(out) :: ad_bvstr_sol(LBi:UBi,LBj:UBj)
1335# endif
1336!
1337! Local variable declarations.
1338!
1339 integer :: i, j
1340!
1341 real(r8) :: cff1, cff2
1342 real(r8) :: ad_cff1, ad_cff2
1343 real(r8) :: adfac
1344
1345# include "set_bounds.h"
1346!
1347!-----------------------------------------------------------------------
1348! Initialize adjoint private variables.
1349!-----------------------------------------------------------------------
1350!
1351 ad_cff1=0.0_r8
1352 ad_cff2=0.0_r8
1353!
1354!-----------------------------------------------------------------------
1355! Set kinematic barotropic bottom momentum stress (m2/s2).
1356!-----------------------------------------------------------------------
1357!
1358! Apply boundary conditions.
1359!
1360# ifdef DISTRIBUTE
1361!^ CALL mp_exchange2d (ng, tile, iTLM, 2, &
1362!^ & LBi, UBi, LBj, UBj, &
1363!^ & NghostPoints, &
1364!^ & EWperiodic(ng), NSperiodic(ng), &
1365!^ & tl_bustr, tl_bvstr)
1366!^
1367 CALL ad_mp_exchange2d (ng, tile, iadm, 2, &
1368 & lbi, ubi, lbj, ubj, &
1369 & nghostpoints, &
1370 & ewperiodic(ng), nsperiodic(ng), &
1371 & ad_bustr, ad_bvstr)
1372# endif
1373!^ CALL bc_v2d_tile (ng, tile, &
1374!^ & LBi, UBi, LBj, UBj, &
1375!^ & tl_bvstr)
1376!^
1377 CALL ad_bc_v2d_tile (ng, tile, &
1378 & lbi, ubi, lbj, ubj, &
1379 & ad_bvstr)
1380!^ CALL bc_u2d_tile (ng, tile, &
1381!^ & LBi, UBi, LBj, UBj, &
1382!^ & tl_bustr)
1383!^
1384 CALL ad_bc_u2d_tile (ng, tile, &
1385 & lbi, ubi, lbj, ubj, &
1386 & ad_bustr)
1387!
1388! Save adjoint bottom momentum flux for output purposes.
1389!
1390 DO j=jstrr,jendr
1391 DO i=istr,iendr
1392 ad_bustr_sol(i,j)=ad_bustr(i,j)
1393 END DO
1394 IF (j.ge.jstr) THEN
1395 DO i=istr,iendr
1396 ad_bvstr_sol(i,j)=ad_bvstr(i,j)
1397 END DO
1398 END IF
1399 END DO
1400
1401# if defined UV_LDRAG
1402!
1403! Set linear bottom stress.
1404!
1405 DO j=jstrv,jend
1406 DO i=istr,iend
1407!^ tl_bvstr(i,j)=0.5_r8*(rdrag(i,j-1)+rdrag(i,j))* &
1408!^ & tl_vbar(i,j,krhs)
1409!^
1410 ad_vbar(i,j,krhs)=ad_vbar(i,j,krhs)+ &
1411 & 0.5_r8*(rdrag(i,j-1)+rdrag(i,j))* &
1412 & ad_bvstr(i,j)
1413 ad_bvstr(i,j)=0.0_r8
1414 END DO
1415 END DO
1416 DO j=jstr,jend
1417 DO i=istru,iend
1418!^ tl_bustr(i,j)=0.5_r8*(rdrag(i,j-1)+rdrag(i,j))* &
1419!^ & tl_ubar(i,j,krhs)
1420!^
1421 ad_ubar(i,j,krhs)=ad_ubar(i,j,krhs)+ &
1422 & 0.5_r8*(rdrag(i,j-1)+rdrag(i,j))* &
1423 & ad_bustr(i,j)
1424 ad_bustr(i,j)=0.0_r8
1425 END DO
1426 END DO
1427# elif defined UV_QDRAG
1428!
1429! Set quadratic bottom stress.
1430!
1431 DO j=jstrv,jend
1432 DO i=istr,iend
1433 cff1=0.25_r8*(ubar(i ,j ,krhs)+ &
1434 & ubar(i+1,j ,krhs)+ &
1435 & ubar(i ,j-1,krhs)+ &
1436 & ubar(i+1,j-1,krhs))
1437 cff2=sqrt(cff1*cff1+vbar(i,j,krhs)*vbar(i,j,krhs))
1438!^ tl_bvstr(i,j)=0.5_r8*(rdrag2(i,j-1)+rdrag2(i,j))* &
1439!^ & (tl_vbar(i,j,krhs)*cff2+ &
1440!^ & vbar(i,j,krhs)*tl_cff2)
1441!^
1442 adfac=0.5_r8*(rdrag2(i,j-1)+rdrag2(i,j))*ad_bvstr(i,j)
1443 ad_vbar(i,j,krhs)=ad_vbar(i,j,krhs)+adfac*cff2
1444 ad_cff2=ad_cff2+vbar(i,j,krhs)*adfac
1445 ad_bvstr(i,j)=0.0_r8
1446 IF (cff2.ne.0.0_r8) THEN
1447!^ tl_cff2=(cff1*tl_cff1+vbar(i,j,krhs)*tl_vbar(i,j,krhs))/cff2
1448!^
1449 adfac=ad_cff2/cff2
1450 ad_vbar(i,j,krhs)=ad_vbar(i,j,krhs)+ &
1451 & adfac*vbar(i,j,krhs)
1452 ad_cff1=cff1*adfac
1453 ad_cff2=0.0_r8
1454 ELSE
1455!^ tl_cff2=0.0_r8
1456!^
1457 ad_cff2=0.0_r8
1458 END IF
1459!^ tl_cff1=0.25_r8*(tl_ubar(i ,j ,krhs)+ &
1460!^ & tl_ubar(i+1,j ,krhs)+ &
1461!^ & tl_ubar(i ,j-1,krhs)+ &
1462!^ & tl_ubar(i+1,j-1,krhs))
1463!^
1464 adfac=0.25_r8*ad_cff1
1465 ad_ubar(i ,j-1,krhs)=ad_ubar(i ,j-1,krhs)+adfac
1466 ad_ubar(i+1,j-1,krhs)=ad_ubar(i+1,j-1,krhs)+adfac
1467 ad_ubar(i ,j ,krhs)=ad_ubar(i ,j ,krhs)+adfac
1468 ad_ubar(i+1,j ,krhs)=ad_ubar(i+1,j ,krhs)+adfac
1469 ad_cff1=0.0_r8
1470 END DO
1471 END DO
1472 DO j=jstr,jend
1473 DO i=istru,iend
1474 cff1=0.25_r8*(vbar(i ,j ,krhs)+ &
1475 & vbar(i ,j+1,krhs)+ &
1476 & vbar(i-1,j ,krhs)+ &
1477 & vbar(i-1,j+1,krhs))
1478 cff2=sqrt(ubar(i,j,krhs)*ubar(i,j,krhs)+cff1*cff1)
1479!^ tl_bustr(i,j)=0.5_r8*(rdrag2(i-1,j)+rdrag2(i,j))* &
1480!^ & (tl_ubar(i,j,krhs)*cff2+ &
1481!^ & ubar(i,j,krhs)*tl_cff2)
1482!^
1483 adfac=0.5_r8*(rdrag2(i-1,j)+rdrag2(i,j))*ad_bustr(i,j)
1484 ad_ubar(i,j,krhs)=ad_ubar(i,j,krhs)+adfac*cff2
1485 ad_cff2=ad_cff2+adfac*ubar(i,j,krhs)
1486 ad_bustr(i,j)=0.0_r8
1487 IF (cff2.ne.0.0_r8) THEN
1488!^ tl_cff2=(ubar(i,j,krhs)*tl_ubar(i,j,krhs)+cff1*tl_cff1)/cff2
1489!^
1490 adfac=ad_cff2/cff2
1491 ad_ubar(i,j,krhs)=ad_ubar(i,j,krhs)+ &
1492 & adfac*ubar(i,j,krhs)
1493 ad_cff1=ad_cff1+adfac*cff1
1494 ad_cff2=0.0_r8
1495 ELSE
1496!^ tl_cff2=0.0_r8
1497!^
1498 ad_cff2=0.0_r8
1499 END IF
1500!^ tl_cff1=0.25_r8*(tl_vbar(i ,j ,krhs)+ &
1501!^ & tl_vbar(i ,j+1,krhs)+ &
1502!^ & tl_vbar(i-1,j ,krhs)+ &
1503!^ & tl_vbar(i-1,j+1,krhs))
1504!^
1505 adfac=0.25_r8*ad_cff1
1506 ad_vbar(i-1,j ,krhs)=ad_vbar(i-1,j ,krhs)+adfac
1507 ad_vbar(i ,j ,krhs)=ad_vbar(i ,j ,krhs)+adfac
1508 ad_vbar(i-1,j+1,krhs)=ad_vbar(i-1,j+1,krhs)+adfac
1509 ad_vbar(i ,j+1,krhs)=ad_vbar(i ,j+1,krhs)+adfac
1510 ad_cff1=0.0_r8
1511 END DO
1512 END DO
1513# endif
1514!
1515 RETURN
1516 END SUBROUTINE ad_set_vbc_tile
1517
1518# endif
1519#endif
1520 END MODULE ad_set_vbc_mod
subroutine ad_bc_v2d_tile(ng, tile, lbi, ubi, lbj, ubj, ad_a)
Definition ad_bc_2d.F:428
subroutine ad_bc_u2d_tile(ng, tile, lbi, ubi, lbj, ubj, ad_a)
Definition ad_bc_2d.F:202
subroutine ad_set_vbc_tile(ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, hz, ad_hz, zobot, z_r, ad_z_r, z_w, ad_z_w, zice, t, ad_t, u, ad_u, v, ad_v, dqdt, sst, sss, ad_sustr, ad_svstr, stflux, btflux, stflx, ad_stflx, btflx, ad_btflx, nrhs)
Definition ad_set_vbc.F:177
subroutine, public ad_set_vbc(ng, tile)
Definition ad_set_vbc.F:32
type(t_forces), dimension(:), allocatable forces
Definition mod_forces.F:554
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
integer nat
Definition mod_param.F:499
integer nghostpoints
Definition mod_param.F:710
integer, parameter iadm
Definition mod_param.F:665
real(dp) cdb_min
real(dp) vonkar
logical, dimension(:), allocatable ewperiodic
logical, dimension(:), allocatable nsperiodic
real(dp), dimension(:,:), allocatable tnudg
integer isalt
integer itemp
real(dp) cdb_max
integer, dimension(:), allocatable kstp
integer, dimension(:), allocatable knew
integer, dimension(:), allocatable nrhs
integer, dimension(:), allocatable krhs
subroutine ad_mp_exchange2d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, nghost, ew_periodic, ns_periodic, ad_a, ad_b, ad_c, ad_d)
recursive subroutine wclock_off(ng, model, region, line, routine)
Definition timers.F:148
recursive subroutine wclock_on(ng, model, region, line, routine)
Definition timers.F:3