ROMS
Loading...
Searching...
No Matches
ad_obc_volcons_mod Module Reference

Functions/Subroutines

subroutine ad_obc_flux (ng, tile, kinp)
 
subroutine, public ad_obc_flux_tile (ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, kinp, umask, vmask, h, ad_h, om_v, on_u, ubar, vbar, zeta, ad_ubar, ad_vbar, ad_zeta)
 
subroutine, public ad_set_duv_bc_tile (ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, kinp, umask, vmask, om_v, on_u, ubar, vbar, ad_ubar, ad_vbar, drhs, duon, dvom, ad_drhs, ad_duon, ad_dvom)
 
subroutine ad_conserve_mass_tile (ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, kinp, umask, vmask, ad_ubar, ad_vbar)
 

Function/Subroutine Documentation

◆ ad_conserve_mass_tile()

subroutine ad_obc_volcons_mod::ad_conserve_mass_tile ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) imins,
integer, intent(in) imaxs,
integer, intent(in) jmins,
integer, intent(in) jmaxs,
integer, intent(in) kinp,
real(r8), dimension(lbi:,lbj:), intent(in) umask,
real(r8), dimension(lbi:,lbj:), intent(in) vmask,
real(r8), dimension(lbi:,lbj:,:), intent(inout) ad_ubar,
real(r8), dimension(lbi:,lbj:,:), intent(inout) ad_vbar )
private

Definition at line 617 of file ad_obc_volcons.F.

625!***********************************************************************
626!
627 USE mod_param
628 USE mod_scalars
629!
630! Imported variable declarations.
631!
632 integer, intent(in) :: ng, tile
633 integer, intent(in) :: LBi, UBi, LBj, UBj
634 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
635 integer, intent(in) :: kinp
636!
637# ifdef ASSUMED_SHAPE
638# ifdef MASKING
639 real(r8), intent(in) :: umask(LBi:,LBj:)
640 real(r8), intent(in) :: vmask(LBi:,LBj:)
641# endif
642 real(r8), intent(inout) :: ad_ubar(LBi:,LBj:,:)
643 real(r8), intent(inout) :: ad_vbar(LBi:,LBj:,:)
644# else
645# ifdef MASKING
646 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
647 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
648# endif
649 real(r8), intent(inout) :: ad_ubar(LBi:UBi,LBj:UBj,:)
650 real(r8), intent(inout) :: ad_vbar(LBi:UBi,LBj:UBj,:)
651# endif
652!
653! Local variable declarations.
654!
655 integer :: i, j
656
657# include "set_bounds.h"
658!
659!-----------------------------------------------------------------------
660! Corrects velocities across the open boundaries to enforce global
661! mass conservation constraint.
662!-----------------------------------------------------------------------
663!
664 IF (ad_volcons(inorth,ng)) THEN
665 IF (domain(ng)%Northern_Edge(tile)) THEN
666 DO i=istr,iend
667# ifdef MASKING
668!^ vbar(i,Jend+1,kinp)=vbar(i,Jend+1,kinp)* &
669!^ & vmask(i,Jend+1)
670!^
671 ad_vbar(i,jend+1,kinp)=ad_vbar(i,jend+1,kinp)* &
672 & vmask(i,jend+1)
673# endif
674!^ tl_vbar(i,Jend+1,kinp)=(tl_vbar(i,Jend+1,kinp)+tl_ubar_xs)
675!^
676 ad_ubar_xs=ad_ubar_xs+ad_vbar(i,jend+1,kinp)
677 END DO
678 END IF
679 END IF
680
681 IF (ad_volcons(isouth,ng)) THEN
682 IF (domain(ng)%Southern_Edge(tile)) THEN
683 DO i=istr,iend
684# ifdef MASKING
685!^ tl_vbar(i,Jstr,kinp)=tl_vbar(i,Jstr,kinp)* &
686!^ & vmask(i,Jstr)
687!^
688 ad_vbar(i,jstr,kinp)=ad_vbar(i,jstr,kinp)* &
689 & vmask(i,jstr)
690# endif
691!^ tl_vbar(i,Jstr,kinp)=(tl_vbar(i,Jstr,kinp)-tl_ubar_xs)
692!^
693 ad_ubar_xs=ad_ubar_xs-ad_vbar(i,jstr,kinp)
694 END DO
695 END IF
696 END IF
697
698 IF (ad_volcons(ieast,ng)) THEN
699 IF (domain(ng)%Eastern_Edge(tile)) THEN
700 DO j=jstr,jend
701# ifdef MASKING
702!^ tl_ubar(Iend+1,j,kinp)=tl_ubar(Iend+1,j,kinp)* &
703!^ & umask(Iend+1,j)
704!^
705 ad_ubar(iend+1,j,kinp)=ad_ubar(iend+1,j,kinp)* &
706 & umask(iend+1,j)
707# endif
708!^ tl_ubar(Iend+1,j,kinp)=tl_ubar(Iend+1,j,kinp)+tl_ubar_xs
709!^
710 ad_ubar_xs=ad_ubar_xs+ad_ubar(iend+1,j,kinp)
711 END DO
712 END IF
713 END IF
714
715 IF (ad_volcons(iwest,ng)) THEN
716 IF (domain(ng)%Western_Edge(tile)) THEN
717 DO j=jstr,jend
718# ifdef MASKING
719!^ tl_ubar(Istr,j,kinp)=tl_ubar(Istr,j,kinp)* &
720!^ & umask(Istr,j)
721!^
722 ad_ubar(istr,j,kinp)=ad_ubar(istr,j,kinp)* &
723 & umask(istr,j)
724# endif
725!^ tl_ubar(Istr,j,kinp)=tl_ubar(Istr,j,kinp)-tl_ubar_xs
726!^
727 ad_ubar_xs=ad_ubar_xs-ad_ubar(istr,j,kinp)
728 END DO
729 END IF
730 END IF
731
732 RETURN
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329
logical, dimension(:,:), allocatable ad_volcons
integer, parameter iwest
real(dp) ad_ubar_xs
integer, parameter isouth
integer, parameter ieast
integer, parameter inorth

References mod_scalars::ad_ubar_xs, mod_scalars::ad_volcons, mod_param::domain, mod_scalars::ieast, mod_scalars::inorth, mod_scalars::isouth, and mod_scalars::iwest.

◆ ad_obc_flux()

subroutine ad_obc_volcons_mod::ad_obc_flux ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) kinp )
private

Definition at line 26 of file ad_obc_volcons.F.

27!***********************************************************************
28!
29 USE mod_param
30 USE mod_grid
31 USE mod_ocean
32 USE mod_stepping
33!
34! Imported variable declarations.
35!
36 integer, intent(in) :: ng, tile, kinp
37!
38! Local variable declarations.
39!
40# include "tile.h"
41!
42 CALL ad_obc_flux_tile (ng, tile, &
43 & lbi, ubi, lbj, ubj, &
44 & imins, imaxs, jmins, jmaxs, &
45 & kinp, &
46# ifdef MASKING
47 & grid(ng) % umask, &
48 & grid(ng) % vmask, &
49# endif
50 & grid(ng) % h, &
51 & grid(ng) % ad_h, &
52 & grid(ng) % om_v, &
53 & grid(ng) % on_u, &
54 & ocean(ng) % ubar, &
55 & ocean(ng) % vbar, &
56 & ocean(ng) % zeta, &
57 & ocean(ng) % ad_ubar, &
58 & ocean(ng) % ad_vbar, &
59 & ocean(ng) % ad_zeta)
60 RETURN
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351

References ad_obc_flux_tile(), mod_grid::grid, and mod_ocean::ocean.

Here is the call graph for this function:

◆ ad_obc_flux_tile()

subroutine, public ad_obc_volcons_mod::ad_obc_flux_tile ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) imins,
integer, intent(in) imaxs,
integer, intent(in) jmins,
integer, intent(in) jmaxs,
integer, intent(in) kinp,
real(r8), dimension(lbi:,lbj:), intent(in) umask,
real(r8), dimension(lbi:,lbj:), intent(in) vmask,
real(r8), dimension(lbi:,lbj:), intent(in) h,
real(r8), dimension(lbi:,lbj:), intent(inout) ad_h,
real(r8), dimension(lbi:,lbj:), intent(in) om_v,
real(r8), dimension(lbi:,lbj:), intent(in) on_u,
real(r8), dimension(lbi:,lbj:,:), intent(in) ubar,
real(r8), dimension(lbi:,lbj:,:), intent(in) vbar,
real(r8), dimension(lbi:,lbj:,:), intent(in) zeta,
real(r8), dimension(lbi:,lbj:,:), intent(inout) ad_ubar,
real(r8), dimension(lbi:,lbj:,:), intent(inout) ad_vbar,
real(r8), dimension(lbi:,lbj:,:), intent(inout) ad_zeta )

Definition at line 64 of file ad_obc_volcons.F.

75!***********************************************************************
76!
77 USE mod_param
78 USE mod_parallel
79 USE mod_scalars
80!
81! Imported variable declarations.
82!
83 integer, intent(in) :: ng, tile
84 integer, intent(in) :: LBi, UBi, LBj, UBj
85 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
86 integer, intent(in) :: kinp
87!
88# ifdef ASSUMED_SHAPE
89# ifdef MASKING
90 real(r8), intent(in) :: umask(LBi:,LBj:)
91 real(r8), intent(in) :: vmask(LBi:,LBj:)
92# endif
93 real(r8), intent(in) :: h(LBi:,LBj:)
94 real(r8), intent(in) :: om_v(LBi:,LBj:)
95 real(r8), intent(in) :: on_u(LBi:,LBj:)
96
97 real(r8), intent(in) :: ubar(LBi:,LBj:,:)
98 real(r8), intent(in) :: vbar(LBi:,LBj:,:)
99 real(r8), intent(in) :: zeta(LBi:,LBj:,:)
100
101 real(r8), intent(inout) :: ad_h(LBi:,LBj:)
102 real(r8), intent(inout) :: ad_ubar(LBi:,LBj:,:)
103 real(r8), intent(inout) :: ad_vbar(LBi:,LBj:,:)
104 real(r8), intent(inout) :: ad_zeta(LBi:,LBj:,:)
105# else
106# ifdef MASKING
107 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
108 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
109# endif
110 real(r8), intent(in) :: h(LBi:UBi,LBj:UBj)
111 real(r8), intent(in) :: om_v(LBi:UBi,LBj:UBj)
112 real(r8), intent(in) :: on_u(LBi:UBi,LBj:UBj)
113
114 real(r8), intent(in) :: ubar(LBi:UBi,LBj:UBj,:)
115 real(r8), intent(in) :: vbar(LBi:UBi,LBj:UBj,:)
116 real(r8), intent(in) :: zeta(LBi:UBi,LBj:UBj,:)
117
118 real(r8), intent(inout) :: ad_h(LBi:UBi,LBj:UBj)
119 real(r8), intent(inout) :: ad_ubar(LBi:UBi,LBj:UBj,:)
120 real(r8), intent(inout) :: ad_vbar(LBi:UBi,LBj:UBj,:)
121 real(r8), intent(inout) :: ad_zeta(LBi:UBi,LBj:UBj,:)
122# endif
123!
124! Local variable declarations.
125!
126 integer :: i, j
127
128 real(r8) :: cff, my_area, my_flux
129 real(r8) :: adfac, ad_cff, ad_my_area, ad_my_flux
130
131# ifdef DISTRIBUTE
132 real(r8), dimension(2) :: rbuffer
133 character (len=3), dimension(2) :: op_handle
134# endif
135!
136# include "set_bounds.h"
137!
138!-----------------------------------------------------------------------
139! Initialize adjoint private variables.
140!-----------------------------------------------------------------------
141!
142 ad_cff=0.0_r8
143!
144!-----------------------------------------------------------------------
145! Perform adjoint global summation and compute correction velocity.
146!-----------------------------------------------------------------------
147!
148 IF (any(ad_volcons(:,ng))) THEN
149!^ tl_ubar_xs=tl_bc_flux/bc_area- &
150!^ & tl_bc_area*ubar_xs/bc_area
151!^
154 ad_ubar_xs=0.0_r8
155!^ tl_bc_flux=tl_bc_flux+tl_my_flux
156!^
157 ad_my_flux=ad_bc_flux
158!^ tl_bc_area=tl_bc_area+tl_my_area
159!^
160 ad_my_area=ad_bc_area
161 END IF
162!
163!-----------------------------------------------------------------------
164! Compute open segments cross-section area and mass flux.
165!-----------------------------------------------------------------------
166!
167 IF (ad_volcons(inorth,ng)) THEN
168 IF (domain(ng)%Northern_Edge(tile)) THEN
169 DO i=istr,iend
170 cff=0.5_r8*om_v(i,jend+1)* &
171 & (zeta(i,jend ,kinp)+h(i,jend )+ &
172 & zeta(i,jend+1,kinp)+h(i,jend+1))
173# ifdef MASKING
174 cff=cff*vmask(i,jend+1)
175# endif
176!^ tl_my_flux=tl_my_flux- &
177!^ & tl_cff*vbar(i,Jend+1,kinp)- &
178!^ & cff*tl_vbar(i,Jend+1,kinp)
179!^
180 ad_vbar(i,jend+1,kinp)=ad_vbar(i,jend+1,kinp)- &
181 & cff*ad_my_flux
182 ad_cff=ad_cff-ad_my_flux*vbar(i,jend+1,kinp)
183!^ tl_my_area=tl_my_area+tl_cff
184!^
185 ad_cff=ad_cff+ad_my_area
186# ifdef MASKING
187!^ tl_cff=tl_cff*vmask(i,Jend+1)
188!^
189 ad_cff=ad_cff*vmask(i,jend+1)
190# endif
191!^ tl_cff=0.5_r8*om_v(i,Jend+1)* &
192!^ & (tl_zeta(i,Jend ,kinp)+tl_h(i,Jend )+ &
193!^ & tl_zeta(i,Jend+1,kinp)+tl_h(i,Jend+1))
194!^
195 adfac=0.5_r8*om_v(i,jend+1)*ad_cff
196 ad_zeta(i,jend ,kinp)=ad_zeta(i,jend ,kinp)+adfac
197 ad_zeta(i,jend+1,kinp)=ad_zeta(i,jend+1,kinp)+adfac
198 ad_h(i,jend )=ad_h(i,jend )+adfac
199 ad_h(i,jend+1)=ad_h(i,jend+1)+adfac
200 ad_cff=0.0_r8
201 END DO
202 END IF
203 END IF
204
205 IF (ad_volcons(isouth,ng)) THEN
206 IF (domain(ng)%Southern_Edge(tile)) THEN
207 DO i=istr,iend
208 cff=0.5_r8*om_v(i,jstr)* &
209 & (zeta(i,jstr-1,kinp)+h(i,jstr-1)+ &
210 & zeta(i,jstr ,kinp)+h(i,jstr ))
211# ifdef MASKING
212 cff=cff*vmask(i,jstr)
213# endif
214!^ tl_my_flux=tl_my_flux+ &
215!^ & tl_cff*vbar(i,JstrV-1,kinp)+ &
216!^ & cff*tl_vbar(i,JstrV-1,kinp)
217!^
218 ad_vbar(i,jstrv-1,kinp)=ad_vbar(i,jstrv-1,kinp)+ &
219 & cff*ad_my_flux
220 ad_cff=ad_cff+ad_my_flux*vbar(i,jstrv-1,kinp)
221!^ tl_my_area=tl_my_area+tl_cff
222!^
223 ad_cff=ad_cff+ad_my_area
224# ifdef MASKING
225!^ tl_cff=tl_cff*vmask(i,Jstr)
226!^
227 ad_cff=ad_cff*vmask(i,jstr)
228# endif
229!^ tl_cff=0.5_r8*om_v(i,Jstr)* &
230!^ & (tl_zeta(i,Jstr-1,kinp)+tl_h(i,Jstr-1)+ &
231!^ & tl_zeta(i,Jstr ,kinp)+tl_h(i,Jstr ))
232!^
233 adfac=0.5_r8*om_v(i,jstr)*ad_cff
234 ad_zeta(i,jstr-1,kinp)=ad_zeta(i,jstr-1,kinp)+adfac
235 ad_zeta(i,jstr ,kinp)=ad_zeta(i,jstr ,kinp)+adfac
236 ad_h(i,jstr-1)=ad_h(i,jstr-1)+adfac
237 ad_h(i,jstr )=ad_h(i,jstr )+adfac
238 ad_cff=0.0_r8
239 END DO
240 END IF
241 END IF
242
243 IF (ad_volcons(ieast,ng)) THEN
244 IF (domain(ng)%Eastern_Edge(tile)) THEN
245 DO j=jstr,jend
246 cff=0.5_r8*on_u(iend+1,j)* &
247 & (zeta(iend ,j,kinp)+h(iend ,j)+ &
248 & zeta(iend+1,j,kinp)+h(iend+1,j))
249# ifdef MASKING
250 cff=cff*umask(iend+1,j)
251# endif
252!^ tl_my_flux=tl_my_flux- &
253!^ & tl_cff*ubar(Iend+1,j,kinp)- &
254!^ & cff*tl_ubar(Iend+1,j,kinp)
255!^
256 ad_ubar(iend+1,j,kinp)=ad_ubar(iend+1,j,kinp)- &
257 & cff*ad_my_flux
258 ad_cff=ad_cff-ad_my_flux*ubar(iend+1,j,kinp)
259!^ tl_my_area=tl_my_area+tl_cff
260!^
261 ad_cff=ad_cff+ad_my_area
262# ifdef MASKING
263!^ tl_cff=tl_cff*umask(Iend+1,j)
264!^
265 ad_cff=ad_cff*umask(iend+1,j)
266# endif
267!^ tl_cff=0.5_r8*on_u(Iend+1,j)* &
268!^ & (tl_zeta(Iend ,j,kinp)+tl_h(Iend ,j)+ &
269!^ & tl_zeta(Iend+1,j,kinp)+tl_h(Iend+1,j))
270!^
271 adfac=0.5_r8*on_u(iend+1,j)*ad_cff
272 ad_zeta(iend ,j,kinp)=ad_zeta(iend ,j,kinp)+adfac
273 ad_zeta(iend+1,j,kinp)=ad_zeta(iend+1,j,kinp)+adfac
274 ad_h(iend ,j)=ad_h(iend ,j)+adfac
275 ad_h(iend+1,j)=ad_h(iend+1,j)+adfac
276 ad_cff=0.0_r8
277 END DO
278 END IF
279 END IF
280
281 IF (ad_volcons(iwest,ng)) THEN
282 IF (domain(ng)%Western_Edge(tile)) THEN
283 DO j=jstr,jend
284 cff=0.5_r8*on_u(istr,j)* &
285 & (zeta(istr-1,j,kinp)+h(istr-1,j)+ &
286 & zeta(istr ,j,kinp)+h(istr ,j))
287# ifdef MASKING
288 cff=cff*umask(istr,j)
289# endif
290!^ tl_my_flux=tl_my_flux+ &
291!^ & tl_cff*ubar(Istr,j,kinp)+ &
292!^ & cff*tl_ubar(Istr,j,kinp)
293!^
294 ad_ubar(istr,j,kinp)=ad_ubar(istr,j,kinp)+ &
295 & cff*ad_my_flux
296 ad_cff=ad_cff+ad_my_flux*ubar(istr,j,kinp)
297!^ tl_my_area=tl_my_area+tl_cff
298!^
299 ad_cff=ad_cff+ad_my_area
300# ifdef MASKING
301!^ tl_cff=tl_cff*umask(Istr,j)
302!^
303 ad_cff=ad_cff*umask(istr,j)
304# endif
305!^ tl_cff=0.5_r8*on_u(Istr,j)* &
306!^ & (tl_zeta(Istr-1,j,kinp)+tl_h(Istr-1,j)+ &
307!^ & tl_zeta(Istr ,j,kinp)+tl_h(Istr ,j))
308!^
309 adfac=0.5_r8*on_u(istr,j)*ad_cff
310 ad_zeta(istr-1,j,kinp)=ad_zeta(istr-1,j,kinp)+adfac
311 ad_zeta(istr ,j,kinp)=ad_zeta(istr ,j,kinp)+adfac
312 ad_h(istr-1,j)=ad_h(istr-1,j)+adfac
313 ad_h(istr-1,j)=ad_h(istr-1,j)+adfac
314 ad_cff=0.0_r8
315 END DO
316 END IF
317 END IF
318
319!^ tl_my_area=0.0_r8
320!^ tl_my_flux=0.0_r8
321!^
322 ad_my_area=0.0_r8
323 ad_my_flux=0.0_r8
324
325 RETURN
real(dp) ad_bc_flux
real(dp) ubar_xs
real(dp) bc_area
real(dp) ad_bc_area

References mod_scalars::ad_bc_area, mod_scalars::ad_bc_flux, mod_scalars::ad_ubar_xs, mod_scalars::ad_volcons, mod_scalars::bc_area, mod_param::domain, mod_scalars::ieast, mod_scalars::inorth, mod_scalars::isouth, mod_scalars::iwest, and mod_scalars::ubar_xs.

Referenced by ad_obc_flux(), ad_step2d_mod::ad_step2d_tile(), ad_step2d_mod::ad_step2d_tile(), and ad_step2d_mod::ad_step2d_tile().

Here is the caller graph for this function:

◆ ad_set_duv_bc_tile()

subroutine, public ad_obc_volcons_mod::ad_set_duv_bc_tile ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) imins,
integer, intent(in) imaxs,
integer, intent(in) jmins,
integer, intent(in) jmaxs,
integer, intent(in) kinp,
real(r8), dimension(lbi:,lbj:), intent(in) umask,
real(r8), dimension(lbi:,lbj:), intent(in) vmask,
real(r8), dimension(lbi:,lbj:), intent(in) om_v,
real(r8), dimension(lbi:,lbj:), intent(in) on_u,
real(r8), dimension(lbi:,lbj:,:), intent(in) ubar,
real(r8), dimension(lbi:,lbj:,:), intent(in) vbar,
real(r8), dimension(lbi:,lbj:,:), intent(inout) ad_ubar,
real(r8), dimension(lbi:,lbj:,:), intent(inout) ad_vbar,
real(r8), dimension(imins:,jmins:), intent(in) drhs,
real(r8), dimension(imins:,jmins:), intent(in) duon,
real(r8), dimension(imins:,jmins:), intent(in) dvom,
real(r8), dimension(imins:,jmins:), intent(inout) ad_drhs,
real(r8), dimension(imins:,jmins:), intent(inout) ad_duon,
real(r8), dimension(imins:,jmins:), intent(inout) ad_dvom )

Definition at line 329 of file ad_obc_volcons.F.

341!***********************************************************************
342!
343 USE mod_param
344 USE mod_scalars
345 USE mod_parallel
346
347# ifdef DISTRIBUTE
348!
350 USE distribute_mod, ONLY : mp_reduce
351# endif
352!
353! Imported variable declarations.
354!
355 integer, intent(in) :: ng, tile
356 integer, intent(in) :: LBi, UBi, LBj, UBj
357 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
358 integer, intent(in) :: kinp
359!
360# ifdef ASSUMED_SHAPE
361# ifdef MASKING
362 real(r8), intent(in) :: umask(LBi:,LBj:)
363 real(r8), intent(in) :: vmask(LBi:,LBj:)
364# endif
365 real(r8), intent(in) :: om_v(LBi:,LBj:)
366 real(r8), intent(in) :: on_u(LBi:,LBj:)
367
368 real(r8), intent(in) :: ubar(LBi:,LBj:,:)
369 real(r8), intent(in) :: vbar(LBi:,LBj:,:)
370 real(r8), intent(in) :: Drhs(IminS:,JminS:)
371 real(r8), intent(in) :: Duon(IminS:,JminS:)
372 real(r8), intent(in) :: Dvom(IminS:,JminS:)
373
374 real(r8), intent(inout) :: ad_ubar(LBi:,LBj:,:)
375 real(r8), intent(inout) :: ad_vbar(LBi:,LBj:,:)
376 real(r8), intent(inout) :: ad_Drhs(IminS:,JminS:)
377 real(r8), intent(inout) :: ad_Duon(IminS:,JminS:)
378 real(r8), intent(inout) :: ad_Dvom(IminS:,JminS:)
379# else
380# ifdef MASKING
381 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
382 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
383# endif
384 real(r8), intent(in) :: om_v(LBi:UBi,LBj:UBj)
385 real(r8), intent(in) :: on_u(LBi:UBi,LBj:UBj)
386
387 real(r8), intent(in) :: ubar(LBi:UBi,LBj:UBj,:)
388 real(r8), intent(in) :: vbar(LBi:UBi,LBj:UBj,:)
389 real(r8), intent(in) :: Drhs(IminS:ImaxS,JminS:JmaxS)
390 real(r8), intent(in) :: Duon(IminS:ImaxS,JminS:JmaxS)
391 real(r8), intent(in) :: Dvom(IminS:ImaxS,JminS:JmaxS)
392
393 real(r8), intent(inout) :: ad_ubar(LBi:UBi,LBj:UBj,:)
394 real(r8), intent(inout) :: ad_vbar(LBi:UBi,LBj:UBj,:)
395 real(r8), intent(inout) :: ad_Drhs(IminS:ImaxS,JminS:JmaxS)
396 real(r8), intent(inout) :: ad_Duon(IminS:ImaxS,JminS:JmaxS)
397 real(r8), intent(inout) :: ad_Dvom(IminS:ImaxS,JminS:JmaxS)
398# endif
399!
400! Local variable declarations.
401!
402 integer :: NSUB, i, j
403
404 real(r8) :: adfac, adfac1, adfac2, adfac3
405 real(r8) :: ad_my_ubar_xs
406# ifdef DISTRIBUTE
407 real(r8) :: rbuffer
408 character (len=3) :: op_handle
409# endif
410
411# include "set_bounds.h"
412
413# ifdef DISTRIBUTE
414# define I_RANGE IstrU,MIN(Iend+1,Lm(ng))
415# define J_RANGE JstrV,MIN(Jend+1,Mm(ng))
416# else
417# define I_RANGE MAX(2,IstrU-1),MIN(Iend+1,Lm(ng))
418# define J_RANGE MAX(2,JstrV-1),MIN(Jend+1,Mm(ng))
419# endif
420!
421!-----------------------------------------------------------------------
422! Set vertically integrated mass fluxes "Duon" and "Dvom" along
423! the open boundaries in such a way that the integral volume is
424! conserved. This is done by applying "ubar_xs" correction to
425! the velocities.
426!-----------------------------------------------------------------------
427!
428 ad_my_ubar_xs=0.0_r8
429!
430# ifdef DISTRIBUTE
431
432! Do a special exchange to avoid having three ghost points for high
433! order numerical stencil.
434!
435 IF (ad_volcons(isouth,ng).or.ad_volcons(inorth,ng)) THEN
436!^ CALL mp_exchange2d (ng, tile, iTLM, 1, &
437!^ & IminS, ImaxS, JminS, JmaxS, &
438!^ & NghostPoints, &
439!^ & EWperiodic(ng), NSperiodic(ng), &
440!^ & tl_Dvom)
441!^
442 CALL ad_mp_exchange2d (ng, tile, iadm, 1, &
443 & imins, imaxs, jmins, jmaxs, &
444 & nghostpoints, &
445 & ewperiodic(ng), nsperiodic(ng), &
446 & ad_dvom)
447 END IF
448
449 IF (ad_volcons(iwest,ng).or.ad_volcons(ieast,ng)) THEN
450!^ CALL mp_exchange2d (ng, tile, iTLM, 1, &
451!^ & IminS, ImaxS, JminS, JmaxS, &
452!^ & NghostPoints, &
453!^ & EWperiodic(ng), NSperiodic(ng), &
454!^ & tl_Duon)
455!^
456 CALL ad_mp_exchange2d (ng, tile, iadm, 1, &
457 & imins, imaxs, jmins, jmaxs, &
458 & nghostpoints, &
459 & ewperiodic(ng), nsperiodic(ng), &
460 & ad_duon)
461 END IF
462# endif
463
464 IF (ad_volcons(inorth,ng)) THEN
465 IF (domain(ng)%Northern_Edge(tile)) THEN
466 DO i=-2+i_range+1
467# ifdef MASKING
468!^ tl_Dvom(i,Jend+1)=tl_Dvom(i,Jend+1)*vmask(i,Jend+1)
469!^
470 ad_dvom(i,jend+1)=ad_dvom(i,jend+1)*vmask(i,jend+1)
471# endif
472!^ tl_Dvom(i,Jend+1)=0.5_r8* &
473!^ & ((tl_Drhs(i,Jend+1)+tl_Drhs(i,Jend))* &
474!^ & (vbar(i,Jend+1,kinp)+ubar_xs)+ &
475!^ & (Drhs(i,Jend+1)+Drhs(i,Jend))* &
476!^ & (tl_vbar(i,Jend+1,kinp)+tl_ubar_xs))* &
477!^ & om_v(i,Jend+1)
478!^
479 adfac=0.5_r8*om_v(i,jend+1)*ad_dvom(i,jend+1)
480 adfac1=adfac*(vbar(i,jend+1,kinp)+ubar_xs)
481 adfac2=adfac*(drhs(i,jend+1)+drhs(i,jend))
482 ad_drhs(i,jend+1)=ad_drhs(i,jend+1)+adfac1
483 ad_drhs(i,jend )=ad_drhs(i,jend )+adfac1
484 ad_vbar(i,jend+1,kinp)=ad_vbar(i,jend+1,kinp)+adfac2
485 ad_my_ubar_xs=ad_my_ubar_xs+adfac2
486 ad_dvom(i,jend+1)=0.0_r8
487 END DO
488 END IF
489 END IF
490
491 IF (ad_volcons(isouth,ng)) THEN
492 IF (domain(ng)%Southern_Edge(tile)) THEN
493 DO i=-2+i_range+1
494# ifdef MASKING
495!^ tl_Dvom(i,Jstr)=tl_Dvom(i,Jstr)*vmask(i,Jstr)
496!^
497 ad_dvom(i,jstr)=ad_dvom(i,jstr)*vmask(i,jstr)
498# endif
499!^ tl_Dvom(i,Jstr)=0.5_r8* &
500!^ & ((tl_Drhs(i,Jstr)+tl_Drhs(i,Jstr-1))* &
501!^ & (vbar(i,Jstr,kinp)-ubar_xs)+ &
502!^ & (Drhs(i,Jstr)+Drhs(i,Jstr-1))* &
503!^ & (tl_vbar(i,Jstr,kinp)-tl_ubar_xs))* &
504!^ & om_v(i,Jstr)
505!^
506 adfac=0.5_r8*om_v(i,jstr)*ad_dvom(i,jstr)
507 adfac1=adfac*(vbar(i,jstr,kinp)-ubar_xs)
508 adfac2=adfac*(drhs(i,jstr)+drhs(i,jstr-1))
509 ad_drhs(i,jstr-1)=ad_drhs(i,jstr-1)+adfac1
510 ad_drhs(i,jstr )=ad_drhs(i,jstr )+adfac1
511 ad_vbar(i,jstr,kinp)=ad_vbar(i,jstr,kinp)+adfac2
512 ad_my_ubar_xs=ad_my_ubar_xs-adfac2
513 ad_dvom(i,jstr)=0.0_r8
514 END DO
515 END IF
516 END IF
517
518 IF (ad_volcons(ieast,ng)) THEN
519 IF (domain(ng)%Eastern_Edge(tile)) THEN
520 DO j=-2+j_range+1
521# ifdef MASKING
522!^ tl_Duon(Iend+1,j)=tl_Duon(Iend+1,j)*umask(Iend+1,j)
523!^
524 ad_duon(iend+1,j)=ad_duon(iend+1,j)*umask(iend+1,j)
525# endif
526!^ tl_Duon(Iend+1,j)=0.5_r8* &
527!^ & ((tl_Drhs(Iend+1,j)+tl_Drhs(Iend,j))* &
528!^ & (ubar(Iend+1,j,kinp)+ubar_xs)+ &
529!^ & (Drhs(Iend+1,j)+Drhs(Iend,j))* &
530!^ & (tl_ubar(Iend+1,j,kinp)+tl_ubar_xs))* &
531!^ & on_u(Iend+1,j)
532!^
533 adfac=0.5_r8*on_u(iend+1,j)*ad_duon(iend+1,j)
534 adfac1=adfac*(ubar(iend+1,j,kinp)+ubar_xs)
535 adfac2=adfac*(drhs(iend+1,j)+drhs(iend,j))
536 ad_drhs(iend ,j)=ad_drhs(iend ,j)+adfac1
537 ad_drhs(iend+1,j)=ad_drhs(iend+1,j)+adfac1
538 ad_ubar(iend+1,j,kinp)=ad_ubar(iend+1,j,kinp)+adfac2
539 ad_my_ubar_xs=ad_my_ubar_xs+adfac2
540 ad_duon(iend+1,j)=0.0_r8
541 END DO
542 END IF
543 END IF
544
545 IF (ad_volcons(iwest,ng)) THEN
546 IF (domain(ng)%Western_Edge(tile)) THEN
547 DO j=-2+j_range+1
548# ifdef MASKING
549!^ tl_Duon(Istr,j)=tl_Duon(Istr,j)*umask(Istr,j)
550!^
551 ad_duon(istr,j)=ad_duon(istr,j)*umask(istr,j)
552# endif
553!^ tl_Duon(Istr,j)=0.5_r8* &
554!^ & ((tl_Drhs(Istr,j)+tl_Drhs(Istr-1,j))* &
555!^ & (ubar(Istr,j,kinp)-ubar_xs)+ &
556!^ & (Drhs(Istr,j)+Drhs(Istr-1,j))* &
557!^ & (tl_ubar(Istr,j,kinp)-tl_ubar_xs))* &
558!^ & on_u(Istr,j)
559!^
560 adfac=0.5_r8*on_u(istr,j)*ad_duon(istr,j)
561 adfac1=adfac*(ubar(istr,j,kinp)-ubar_xs)
562 adfac2=adfac*(drhs(istr,j)+drhs(istr-1,j))
563 ad_drhs(istr-1,j)=ad_drhs(istr-1,j)+adfac1
564 ad_drhs(istr ,j)=ad_drhs(istr ,j)+adfac1
565 ad_ubar(istr,j,kinp)=ad_ubar(istr,j,kinp)+adfac2
566 ad_my_ubar_xs=ad_my_ubar_xs-adfac2
567 ad_duon(istr,j)=0.0_r8
568 END DO
569 END IF
570 END IF
571
572# undef I_RANGE
573# undef J_RANGE
574!
575!-----------------------------------------------------------------------
576! Perform global summation and compute correction velocity.
577!-----------------------------------------------------------------------
578!
579 IF (any(ad_volcons(:,ng))) THEN
580# ifdef DISTRIBUTE
581 nsub=1 ! distributed-memory
582# else
583 IF (domain(ng)%SouthWest_Corner(tile).and. &
584 & domain(ng)%NorthEast_Corner(tile)) THEN
585 nsub=1 ! non-tiled application
586 ELSE
587 nsub=ntilex(ng)*ntilee(ng) ! tiled application
588 END IF
589# endif
590!$OMP CRITICAL (AD_OBC_VOLUME)
591 IF (tile_count.eq.0) THEN
592 adfac3=0.0_r8
593 END IF
594 adfac3=adfac3+ad_my_ubar_xs
596 IF (tile_count.eq.nsub) THEN
597 tile_count=0
598# ifdef DISTRIBUTE
599 rbuffer=adfac3
600 op_handle='SUM'
601 CALL mp_reduce (ng, iadm, 1, rbuffer, op_handle)
602 adfac3=rbuffer
603# endif
604 IF (iif(ng).eq.nfast(ng)+1) THEN
606 ELSE
607 ad_ubar_xs=adfac3
608 ENDIF
609 END IF
610!$OMP END CRITICAL (AD_OBC_VOLUME)
611 END IF
612
613 RETURN
integer tile_count
integer, dimension(:), allocatable ntilex
Definition mod_param.F:685
integer nghostpoints
Definition mod_param.F:710
integer, parameter iadm
Definition mod_param.F:665
integer, dimension(:), allocatable ntilee
Definition mod_param.F:686
logical, dimension(:), allocatable ewperiodic
logical, dimension(:), allocatable nsperiodic
integer, dimension(:), allocatable nfast
integer, dimension(:), allocatable iif
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)

References mp_exchange_mod::ad_mp_exchange2d(), mod_scalars::ad_ubar_xs, mod_scalars::ad_volcons, mod_param::domain, mod_scalars::ewperiodic, mod_param::iadm, mod_scalars::ieast, mod_scalars::iif, mod_scalars::inorth, mod_scalars::isouth, mod_scalars::iwest, mod_scalars::nfast, mod_param::nghostpoints, mod_scalars::nsperiodic, mod_param::ntilee, mod_param::ntilex, mod_parallel::tile_count, and mod_scalars::ubar_xs.

Referenced by ad_step2d_mod::ad_step2d_tile(), ad_step2d_mod::ad_step2d_tile(), and ad_step2d_mod::ad_step2d_tile().

Here is the call graph for this function:
Here is the caller graph for this function: