ROMS
Loading...
Searching...
No Matches
ana_m3obc.h
Go to the documentation of this file.
1!!
2 SUBROUTINE ana_m3obc (ng, tile, model)
3!
4!! git $Id$
5!!======================================================================
6!! Copyright (c) 2002-2025 The ROMS Group !
7!! Licensed under a MIT/X style license !
8!! See License_ROMS.md !
9!=======================================================================
10! !
11! This routine sets 3D momentum open boundary conditions using !
12! analytical expressions. !
13! !
14!=======================================================================
15!
16 USE mod_param
17 USE mod_boundary
18 USE mod_ncparam
19!
20! Imported variable declarations.
21!
22 integer, intent(in) :: ng, tile, model
23!
24! Local variable declarations.
25!
26 character (len=*), parameter :: MyFile = &
27 & __FILE__
28!
29#include "tile.h"
30!
31 CALL ana_m3obc_tile (ng, tile, model, &
32 & lbi, ubi, lbj, ubj, &
33 & imins, imaxs, jmins, jmaxs)
34!
35! Set analytical header file name used.
36!
37#ifdef DISTRIBUTE
38 IF (lanafile) THEN
39#else
40 IF (lanafile.and.(tile.eq.0)) THEN
41#endif
42 ananame(14)=myfile
43 END IF
44!
45 RETURN
46 END SUBROUTINE ana_m3obc
47!
48!***********************************************************************
49 SUBROUTINE ana_m3obc_tile (ng, tile, model, &
50 & LBi, UBi, LBj, UBj, &
51 & IminS, ImaxS, JminS, JmaxS)
52!***********************************************************************
53!
54 USE mod_param
55 USE mod_boundary
56 USE mod_grid
57 USE mod_ncparam
58 USE mod_ocean
59 USE mod_scalars
60!
61! Imported variable declarations.
62!
63 integer, intent(in) :: ng, tile, model
64 integer, intent(in) :: LBi, UBi, LBj, UBj
65 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
66!
67! Local variable declarations.
68!
69 integer :: i, j, k
70!
71 real(r8) :: fac, val
72
73#include "set_bounds.h"
74!
75!-----------------------------------------------------------------------
76! 3D momentum open boundary conditions.
77!-----------------------------------------------------------------------
78!
79#if defined SED_TEST1
80 IF (lbc(iwest,isuvel,ng)%acquire.and. &
81 & lbc(iwest,isvvel,ng)%acquire.and. &
82 & domain(ng)%Western_Edge(tile)) THEN
83 fac=5.0e-06_r8
84 DO k=1,n(ng)
85 DO j=jstrt,jendt
86 val=0.5_r8*(ocean(ng)%zeta(0 ,j,knew)+ &
87 & grid(ng)%h(0 ,j)+ &
88 & ocean(ng)%zeta(1 ,j,knew)+ &
89 & grid(ng)%h(1 ,j))
90 boundary(ng)%u_west(j,k)=-log((val+0.5*
91 % (grid(ng)%z_r(istr-1,j,k)+ &
92 & grid(ng)%z_r(istr ,j,k)))/ &
93 & fac)/ &
94 & (log(val/fac)-1.0_r8+fac/val)
95 END DO
96 DO j=jstrp,jendt
97 boundary(ng)%v_west(j,k)=0.0_r8
98 END DO
99 END DO
100 END IF
101
102 IF (lbc(ieast,isuvel,ng)%acquire.and. &
103 & lbc(ieast,isvvel,ng)%acquire.and. &
104 & domain(ng)%Eastern_Edge(tile)) THEN
105 fac=5.0e-06_r8
106 DO k=1,n(ng)
107 DO j=jstrt,jendt
108 val=0.5_r8*(ocean(ng)%zeta(iend ,j,knew)+ &
109 & grid(ng)%h(iend ,j)+ &
110 & ocean(ng)%zeta(iend+1,j,knew)+ &
111 % GRID(ng)%h(iend+1,j))
112 boundary(ng)%u_east(j,k)=-log((val+0.5*
113 & (grid(ng)%z_r(iend ,j,k)+ &
114 & grid(ng)%z_r(iend+1,j,k)))/ &
115 & fac)/ &
116 & (log(val/fac)-1.0_r8+fac/val)
117 END DO
118 DO j=jstrp,jendt
119 boundary(ng)%v_east(j,k)=0.0_r8
120 END DO
121 END DO
122 END IF
123#else
124 IF (lbc(ieast,isuvel,ng)%acquire.and. &
125 & lbc(ieast,isvvel,ng)%acquire.and. &
126 & domain(ng)%Eastern_Edge(tile)) THEN
127 DO k=1,n(ng)
128 DO j=jstrt,jendt
129 boundary(ng)%u_east(j,k)=0.0_r8
130 END DO
131 DO j=jstrp,jendt
132 boundary(ng)%v_east(j,k)=0.0_r8
133 END DO
134 END DO
135 END IF
136
137 IF (lbc(iwest,isuvel,ng)%acquire.and. &
138 & lbc(iwest,isvvel,ng)%acquire.and. &
139 & domain(ng)%Western_Edge(tile)) THEN
140 DO k=1,n(ng)
141 DO j=jstrt,jendt
142 boundary(ng)%u_west(j,k)=0.0_r8
143 END DO
144 DO j=jstrp,jendt
145 boundary(ng)%v_west(j,k)=0.0_r8
146 END DO
147 END DO
148 END IF
149
150 IF (lbc(isouth,isuvel,ng)%acquire.and. &
151 & lbc(isouth,isvvel,ng)%acquire.and. &
152 & domain(ng)%Southern_Edge(tile)) THEN
153 DO k=1,n(ng)
154 DO i=istrp,iendt
155 boundary(ng)%u_south(i,k)=0.0_r8
156 END DO
157 DO i=istrt,iendt
158 boundary(ng)%v_south(i,k)=0.0_r8
159 END DO
160 END DO
161 END IF
162
163 IF (lbc(inorth,isuvel,ng)%acquire.and. &
164 & lbc(inorth,isvvel,ng)%acquire.and. &
165 & domain(ng)%Northern_Edge(tile)) THEN
166 DO k=1,n(ng)
167 DO i=istrp,iendt
168 boundary(ng)%u_north(i,k)=0.0_r8
169 END DO
170 DO i=istrt,iendt
171 boundary(ng)%v_north(i,k)=0.0_r8
172 END DO
173 END DO
174 END IF
175#endif
176!
177 RETURN
178 END SUBROUTINE ana_m3obc_tile
subroutine ana_m3obc(ng, tile, model)
Definition ana_m3obc.h:3
subroutine ana_m3obc_tile(ng, tile, model, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs)
Definition ana_m3obc.h:52
type(t_boundary), dimension(:), allocatable boundary
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
integer isvvel
integer isuvel
logical lanafile
character(len=256), dimension(39) ananame
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
integer, dimension(:), allocatable n
Definition mod_param.F:479
type(t_lbc), dimension(:,:,:), allocatable lbc
Definition mod_param.F:375
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329
integer, parameter iwest
integer, parameter isouth
integer, parameter ieast
integer, parameter inorth