1 | MODULE bdydyn2d |
---|
2 | !!====================================================================== |
---|
3 | !! *** MODULE bdydyn *** |
---|
4 | !! Unstructured Open Boundary Cond. : Apply boundary conditions to barotropic solution |
---|
5 | !!====================================================================== |
---|
6 | !! History : 3.4 ! 2011 (D. Storkey) new module as part of BDY rewrite |
---|
7 | !!---------------------------------------------------------------------- |
---|
8 | #if defined key_bdy |
---|
9 | !!---------------------------------------------------------------------- |
---|
10 | !! 'key_bdy' : Unstructured Open Boundary Condition |
---|
11 | !!---------------------------------------------------------------------- |
---|
12 | !! bdy_dyn2d : Apply open boundary conditions to barotropic variables. |
---|
13 | !! bdy_dyn2d_fla : Apply Flather condition |
---|
14 | !!---------------------------------------------------------------------- |
---|
15 | USE timing ! Timing |
---|
16 | USE oce ! ocean dynamics and tracers |
---|
17 | USE dom_oce ! ocean space and time domain |
---|
18 | USE bdy_oce ! ocean open boundary conditions |
---|
19 | USE dynspg_oce ! for barotropic variables |
---|
20 | USE phycst ! physical constants |
---|
21 | USE lbclnk ! ocean lateral boundary conditions (or mpp link) |
---|
22 | USE in_out_manager ! |
---|
23 | |
---|
24 | IMPLICIT NONE |
---|
25 | PRIVATE |
---|
26 | |
---|
27 | PUBLIC bdy_dyn2d ! routine called in dynspg_ts and bdy_dyn |
---|
28 | |
---|
29 | !!---------------------------------------------------------------------- |
---|
30 | !! NEMO/OPA 3.3 , NEMO Consortium (2010) |
---|
31 | !! $Id: bdydyn.F90 2528 2010-12-27 17:33:53Z rblod $ |
---|
32 | !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) |
---|
33 | !!---------------------------------------------------------------------- |
---|
34 | CONTAINS |
---|
35 | |
---|
36 | SUBROUTINE bdy_dyn2d( kt ) |
---|
37 | !!---------------------------------------------------------------------- |
---|
38 | !! *** SUBROUTINE bdy_dyn2d *** |
---|
39 | !! |
---|
40 | !! ** Purpose : - Apply open boundary conditions for barotropic variables |
---|
41 | !! |
---|
42 | !!---------------------------------------------------------------------- |
---|
43 | INTEGER, INTENT(in) :: kt ! Main time step counter |
---|
44 | !! |
---|
45 | INTEGER :: ib_bdy ! Loop counter |
---|
46 | |
---|
47 | DO ib_bdy=1, nb_bdy |
---|
48 | |
---|
49 | SELECT CASE( nn_dyn2d(ib_bdy) ) |
---|
50 | CASE(jp_none) |
---|
51 | CYCLE |
---|
52 | CASE(jp_frs) |
---|
53 | CALL bdy_dyn2d_frs( idx_bdy(ib_bdy), dta_bdy(ib_bdy) ) |
---|
54 | CASE(jp_flather) |
---|
55 | CALL bdy_dyn2d_fla( idx_bdy(ib_bdy), dta_bdy(ib_bdy) ) |
---|
56 | CASE DEFAULT |
---|
57 | CALL ctl_stop( 'bdy_dyn2d : unrecognised option for open boundaries for barotropic variables' ) |
---|
58 | END SELECT |
---|
59 | ENDDO |
---|
60 | |
---|
61 | END SUBROUTINE bdy_dyn2d |
---|
62 | |
---|
63 | SUBROUTINE bdy_dyn2d_frs( idx, dta ) |
---|
64 | !!---------------------------------------------------------------------- |
---|
65 | !! *** SUBROUTINE bdy_dyn2d_frs *** |
---|
66 | !! |
---|
67 | !! ** Purpose : - Apply the Flow Relaxation Scheme for barotropic velocities |
---|
68 | !! at open boundaries. |
---|
69 | !! |
---|
70 | !! References :- Engedahl H., 1995: Use of the flow relaxation scheme in |
---|
71 | !! a three-dimensional baroclinic ocean model with realistic |
---|
72 | !! topography. Tellus, 365-382. |
---|
73 | !!---------------------------------------------------------------------- |
---|
74 | TYPE(OBC_INDEX), INTENT(in) :: idx ! OBC indices |
---|
75 | TYPE(OBC_DATA), INTENT(in) :: dta ! OBC external data |
---|
76 | !! |
---|
77 | INTEGER :: jb, jk ! dummy loop indices |
---|
78 | INTEGER :: ii, ij, igrd ! local integers |
---|
79 | REAL(wp) :: zwgt ! boundary weight |
---|
80 | !!---------------------------------------------------------------------- |
---|
81 | ! |
---|
82 | IF( nn_timing == 1 ) CALL timing_start('bdy_dyn2d_frs') |
---|
83 | ! |
---|
84 | igrd = 2 ! Relaxation of zonal velocity |
---|
85 | DO jb = 1, idx%nblen(igrd) |
---|
86 | ii = idx%nbi(jb,igrd) |
---|
87 | ij = idx%nbj(jb,igrd) |
---|
88 | zwgt = idx%nbw(jb,igrd) |
---|
89 | pu2d(ii,ij) = ( pu2d(ii,ij) + zwgt * ( dta%u2d(jb) - pu2d(ii,ij) ) ) * umask(ii,ij,1) |
---|
90 | END DO |
---|
91 | ! |
---|
92 | igrd = 3 ! Relaxation of meridional velocity |
---|
93 | DO jb = 1, idx%nblen(igrd) |
---|
94 | ii = idx%nbi(jb,igrd) |
---|
95 | ij = idx%nbj(jb,igrd) |
---|
96 | zwgt = idx%nbw(jb,igrd) |
---|
97 | pv2d(ii,ij) = ( pv2d(ii,ij) + zwgt * ( dta%v2d(jb) - pv2d(ii,ij) ) ) * vmask(ii,ij,1) |
---|
98 | END DO |
---|
99 | CALL lbc_lnk( pu2d, 'U', -1. ) |
---|
100 | CALL lbc_lnk( pv2d, 'V', -1. ) ! Boundary points should be updated |
---|
101 | ! |
---|
102 | IF( nn_timing == 1 ) CALL timing_stop('bdy_dyn2d_frs') |
---|
103 | ! |
---|
104 | |
---|
105 | END SUBROUTINE bdy_dyn2d_frs |
---|
106 | |
---|
107 | |
---|
108 | SUBROUTINE bdy_dyn2d_fla( idx, dta ) |
---|
109 | !!---------------------------------------------------------------------- |
---|
110 | !! *** SUBROUTINE bdy_dyn2d_fla *** |
---|
111 | !! |
---|
112 | !! - Apply Flather boundary conditions on normal barotropic velocities |
---|
113 | !! |
---|
114 | !! ** WARNINGS about FLATHER implementation: |
---|
115 | !!1. According to Palma and Matano, 1998 "after ssh" is used. |
---|
116 | !! In ROMS and POM implementations, it is "now ssh". In the current |
---|
117 | !! implementation (tested only in the EEL-R5 conf.), both cases were unstable. |
---|
118 | !! So I use "before ssh" in the following. |
---|
119 | !! |
---|
120 | !!2. We assume that the normal ssh gradient at the bdy is zero. As a matter of |
---|
121 | !! fact, the model ssh just inside the dynamical boundary is used (the outside |
---|
122 | !! ssh in the code is not updated). |
---|
123 | !! |
---|
124 | !! References: Flather, R. A., 1976: A tidal model of the northwest European |
---|
125 | !! continental shelf. Mem. Soc. R. Sci. Liege, Ser. 6,10, 141-164. |
---|
126 | !!---------------------------------------------------------------------- |
---|
127 | TYPE(OBC_INDEX), INTENT(in) :: idx ! OBC indices |
---|
128 | TYPE(OBC_DATA), INTENT(in) :: dta ! OBC external data |
---|
129 | |
---|
130 | INTEGER :: jb, igrd ! dummy loop indices |
---|
131 | INTEGER :: ii, ij, iim1, iip1, ijm1, ijp1 ! 2D addresses |
---|
132 | REAL(wp) :: zcorr ! Flather correction |
---|
133 | REAL(wp) :: zforc ! temporary scalar |
---|
134 | !!---------------------------------------------------------------------- |
---|
135 | |
---|
136 | IF( nn_timing == 1 ) CALL timing_start('bdy_dyn2d_fla') |
---|
137 | |
---|
138 | ! ---------------------------------! |
---|
139 | ! Flather boundary conditions :! |
---|
140 | ! ---------------------------------! |
---|
141 | |
---|
142 | !!! REPLACE spgu with nemo_wrk work space |
---|
143 | |
---|
144 | ! Fill temporary array with ssh data (here spgu): |
---|
145 | igrd = 1 |
---|
146 | spgu(:,:) = 0.0 |
---|
147 | DO jb = 1, idx%nblenrim(igrd) |
---|
148 | ii = idx%nbi(jb,igrd) |
---|
149 | ij = idx%nbj(jb,igrd) |
---|
150 | spgu(ii, ij) = dta%ssh(jb) |
---|
151 | END DO |
---|
152 | ! |
---|
153 | igrd = 2 ! Flather bc on u-velocity; |
---|
154 | ! ! remember that flagu=-1 if normal velocity direction is outward |
---|
155 | ! ! I think we should rather use after ssh ? |
---|
156 | DO jb = 1, idx%nblenrim(igrd) |
---|
157 | ii = idx%nbi(jb,igrd) |
---|
158 | ij = idx%nbj(jb,igrd) |
---|
159 | iim1 = ii + MAX( 0, INT( idx%flagu(jb) ) ) ! T pts i-indice inside the boundary |
---|
160 | iip1 = ii - MIN( 0, INT( idx%flagu(jb) ) ) ! T pts i-indice outside the boundary |
---|
161 | ! |
---|
162 | zcorr = - idx%flagu(jb) * SQRT( grav * phur(ii, ij) ) * ( pssh(iim1, ij) - spgu(iip1,ij) ) |
---|
163 | zforc = dta%u2d(jb) |
---|
164 | pu2d(ii,ij) = zforc + zcorr * umask(ii,ij,1) |
---|
165 | END DO |
---|
166 | ! |
---|
167 | igrd = 3 ! Flather bc on v-velocity |
---|
168 | ! ! remember that flagv=-1 if normal velocity direction is outward |
---|
169 | DO jb = 1, idx%nblenrim(igrd) |
---|
170 | ii = idx%nbi(jb,igrd) |
---|
171 | ij = idx%nbj(jb,igrd) |
---|
172 | ijm1 = ij + MAX( 0, INT( idx%flagv(jb) ) ) ! T pts j-indice inside the boundary |
---|
173 | ijp1 = ij - MIN( 0, INT( idx%flagv(jb) ) ) ! T pts j-indice outside the boundary |
---|
174 | ! |
---|
175 | zcorr = - idx%flagv(jb) * SQRT( grav * phvr(ii, ij) ) * ( pssh(ii, ijm1) - spgu(ii,ijp1) ) |
---|
176 | zforc = dta%v2d(jb) |
---|
177 | pv2d(ii,ij) = zforc + zcorr * vmask(ii,ij,1) |
---|
178 | END DO |
---|
179 | CALL lbc_lnk( pu2d, 'U', -1. ) ! Boundary points should be updated |
---|
180 | CALL lbc_lnk( pv2d, 'V', -1. ) ! |
---|
181 | ! |
---|
182 | IF( nn_timing == 1 ) CALL timing_stop('bdy_dyn2d_fla') |
---|
183 | ! |
---|
184 | END SUBROUTINE bdy_dyn2d_fla |
---|
185 | #else |
---|
186 | !!---------------------------------------------------------------------- |
---|
187 | !! Dummy module NO Unstruct Open Boundary Conditions |
---|
188 | !!---------------------------------------------------------------------- |
---|
189 | CONTAINS |
---|
190 | SUBROUTINE bdy_dyn2d( kt ) ! Empty routine |
---|
191 | WRITE(*,*) 'bdy_dyn_frs: You should not have seen this print! error?', kt |
---|
192 | END SUBROUTINE bdy_dyn2d |
---|
193 | #endif |
---|
194 | |
---|
195 | !!====================================================================== |
---|
196 | END MODULE bdydyn2d |
---|