[367] | 1 | MODULE obcdyn_bt |
---|
[1151] | 2 | #if ( defined key_dynspg_ts || defined key_dynspg_exp ) && defined key_obc |
---|
[367] | 3 | !!================================================================================= |
---|
| 4 | !! *** MODULE obcdyn_bt *** |
---|
| 5 | !! Ocean dynamics: Radiation/prescription of sea surface heights |
---|
| 6 | !! on each open boundary |
---|
| 7 | !!================================================================================= |
---|
| 8 | |
---|
| 9 | !!--------------------------------------------------------------------------------- |
---|
| 10 | !! obc_dyn_bt : call the subroutine for each open boundary |
---|
| 11 | !! obc_dyn_bt_east : Flather's algorithm at the east open boundary |
---|
| 12 | !! obc_dyn_bt_west : Flather's algorithm at the west open boundary |
---|
| 13 | !! obc_dyn_bt_north : Flather's algorithm at the north open boundary |
---|
| 14 | !! obc_dyn_bt_south : Flather's algorithm at the south open boundary |
---|
| 15 | !!---------------------------------------------------------------------------------- |
---|
| 16 | |
---|
| 17 | !!---------------------------------------------------------------------------------- |
---|
| 18 | !! * Modules used |
---|
| 19 | USE oce ! ocean dynamics and tracers |
---|
| 20 | USE dom_oce ! ocean space and time domain |
---|
| 21 | USE phycst ! physical constants |
---|
| 22 | USE obc_oce ! ocean open boundary conditions |
---|
| 23 | USE lbclnk ! ocean lateral boundary conditions (or mpp link) |
---|
| 24 | USE lib_mpp ! distributed memory computing |
---|
[2200] | 25 | USE obcdta ! ocean open boundary conditions |
---|
[367] | 26 | USE in_out_manager ! I/O manager |
---|
| 27 | USE dynspg_oce ! surface pressure gradient (free surface with time-splitting) |
---|
| 28 | |
---|
| 29 | IMPLICIT NONE |
---|
| 30 | PRIVATE |
---|
| 31 | |
---|
| 32 | !! * Accessibility |
---|
| 33 | PUBLIC obc_dyn_bt ! routine called in dynnxt (explicit free surface case) |
---|
| 34 | |
---|
| 35 | !!--------------------------------------------------------------------------------- |
---|
[2287] | 36 | !! NEMO/OPA 3.3 , NEMO Consortium (2010) |
---|
[1152] | 37 | !! $Id$ |
---|
[2287] | 38 | !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) |
---|
[367] | 39 | !!---------------------------------------------------------------------- |
---|
| 40 | |
---|
| 41 | CONTAINS |
---|
| 42 | |
---|
| 43 | SUBROUTINE obc_dyn_bt ( kt ) |
---|
| 44 | !!------------------------------------------------------------------------------ |
---|
| 45 | !! SUBROUTINE obc_dyn_bt |
---|
| 46 | !! *********************** |
---|
| 47 | !! ** Purpose : |
---|
| 48 | !! Apply Flather's algorithm at open boundaries for the explicit |
---|
| 49 | !! free surface case and free surface case with time-splitting |
---|
| 50 | !! |
---|
| 51 | !! This routine is called in dynnxt.F routine and updates ua, va and sshn. |
---|
| 52 | !! |
---|
| 53 | !! The logical variable lp_obc_east, and/or lp_obc_west, and/or lp_obc_north, |
---|
| 54 | !! and/or lp_obc_south allow the user to determine which boundary is an |
---|
| 55 | !! open one (must be done in the param_obc.h90 file). |
---|
| 56 | !! |
---|
| 57 | !! ** Reference : |
---|
| 58 | !! Flather, R. A., 1976, Mem. Soc. R. Sci. Liege, Ser. 6, 10, 141-164 |
---|
| 59 | !! |
---|
| 60 | !! History : |
---|
| 61 | !! 9.0 ! 05-12 (V. Garnier) original |
---|
| 62 | !!---------------------------------------------------------------------- |
---|
| 63 | !! * Arguments |
---|
| 64 | INTEGER, INTENT( in ) :: kt |
---|
| 65 | |
---|
| 66 | !!---------------------------------------------------------------------- |
---|
| 67 | |
---|
| 68 | IF( lp_obc_east ) CALL obc_dyn_bt_east |
---|
| 69 | IF( lp_obc_west ) CALL obc_dyn_bt_west |
---|
| 70 | IF( lp_obc_north ) CALL obc_dyn_bt_north |
---|
| 71 | IF( lp_obc_south ) CALL obc_dyn_bt_south |
---|
| 72 | |
---|
| 73 | IF( lk_mpp ) THEN |
---|
| 74 | IF( kt >= nit000+3 .AND. ln_rstart ) THEN |
---|
| 75 | CALL lbc_lnk( sshb, 'T', 1. ) |
---|
| 76 | CALL lbc_lnk( ub , 'U', -1. ) |
---|
| 77 | CALL lbc_lnk( vb , 'V', -1. ) |
---|
| 78 | END IF |
---|
| 79 | CALL lbc_lnk( sshn, 'T', 1. ) |
---|
| 80 | CALL lbc_lnk( ua , 'U', -1. ) |
---|
| 81 | CALL lbc_lnk( va , 'V', -1. ) |
---|
| 82 | ENDIF |
---|
| 83 | |
---|
| 84 | END SUBROUTINE obc_dyn_bt |
---|
| 85 | |
---|
| 86 | # if defined key_dynspg_exp |
---|
| 87 | SUBROUTINE obc_dyn_bt_east |
---|
| 88 | !!------------------------------------------------------------------------------ |
---|
| 89 | !! *** SUBROUTINE obc_dyn_bt_east *** |
---|
| 90 | !! |
---|
| 91 | !! ** Purpose : |
---|
| 92 | !! Apply Flather's algorithm on east OBC velocities ua, va |
---|
| 93 | !! Fix sea surface height (sshn) on east open boundary |
---|
| 94 | !! The logical lfbceast must be .TRUE. |
---|
| 95 | !! |
---|
| 96 | !! History : |
---|
| 97 | !! 9.0 ! 05-12 (V. Garnier) original |
---|
| 98 | !!------------------------------------------------------------------------------ |
---|
| 99 | !! * Local declaration |
---|
| 100 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
| 101 | !!------------------------------------------------------------------------------ |
---|
| 102 | |
---|
| 103 | DO ji = nie0, nie1 |
---|
| 104 | DO jk = 1, jpkm1 |
---|
| 105 | DO jj = 1, jpj |
---|
| 106 | ua(ji,jj,jk) = ua(ji,jj,jk) + sqrt( grav*hur (ji,jj) ) & |
---|
| 107 | & * ( ( sshn(ji,jj) + sshn(ji+1,jj) ) * 0.5 & |
---|
| 108 | & - sshfoe(jj) ) * uemsk(jj,jk) |
---|
| 109 | END DO |
---|
| 110 | END DO |
---|
| 111 | END DO |
---|
| 112 | DO ji = nie0p1, nie1p1 |
---|
| 113 | DO jj = 1, jpj |
---|
| 114 | sshn(ji,jj) = sshn(ji,jj) * (1.-temsk(jj,1)) + temsk(jj,1)*sshfoe(jj) |
---|
| 115 | END DO |
---|
| 116 | END DO |
---|
| 117 | |
---|
| 118 | END SUBROUTINE obc_dyn_bt_east |
---|
| 119 | |
---|
| 120 | |
---|
| 121 | SUBROUTINE obc_dyn_bt_west |
---|
| 122 | !!------------------------------------------------------------------------------ |
---|
| 123 | !! *** SUBROUTINE obc_dyn_bt_west *** |
---|
| 124 | !! |
---|
| 125 | !! ** Purpose : |
---|
| 126 | !! Apply Flather algorithm on west OBC velocities ua, va |
---|
| 127 | !! Fix sea surface height (sshn) on west open boundary |
---|
| 128 | !! The logical lfbcwest must be .TRUE. |
---|
| 129 | !! |
---|
| 130 | !! History : |
---|
| 131 | !! 9.0 ! 05-12 (V. Garnier) original |
---|
| 132 | !!------------------------------------------------------------------------------ |
---|
| 133 | !! * Local declaration |
---|
| 134 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
| 135 | !!------------------------------------------------------------------------------ |
---|
| 136 | |
---|
| 137 | DO ji = niw0, niw1 |
---|
| 138 | DO jk = 1, jpkm1 |
---|
| 139 | DO jj = 1, jpj |
---|
| 140 | ua(ji,jj,jk) = ua(ji,jj,jk) - sqrt( grav*hur (ji,jj) ) & |
---|
| 141 | & * ( ( sshn(ji,jj) + sshn(ji+1,jj) ) * 0.5 & |
---|
| 142 | & - sshfow(jj) ) * uwmsk(jj,jk) |
---|
| 143 | END DO |
---|
| 144 | END DO |
---|
| 145 | DO jj = 1, jpj |
---|
| 146 | sshn(ji,jj) = sshn(ji,jj) * (1.-twmsk(jj,1)) + twmsk(jj,1)*sshfow(jj) |
---|
| 147 | END DO |
---|
| 148 | END DO |
---|
| 149 | |
---|
| 150 | END SUBROUTINE obc_dyn_bt_west |
---|
| 151 | |
---|
| 152 | SUBROUTINE obc_dyn_bt_north |
---|
| 153 | !!------------------------------------------------------------------------------ |
---|
| 154 | !! *** SUBROUTINE obc_dyn_bt_north *** |
---|
| 155 | !! |
---|
| 156 | !! ** Purpose : |
---|
| 157 | !! Apply Flather algorithm on north OBC velocities ua, va |
---|
| 158 | !! Fix sea surface height (sshn) on north open boundary |
---|
| 159 | !! The logical lfbcnorth must be .TRUE. |
---|
| 160 | !! |
---|
| 161 | !! History : |
---|
| 162 | !! 9.0 ! 05-12 (V. Garnier) original |
---|
| 163 | !!------------------------------------------------------------------------------ |
---|
| 164 | !! * Local declaration |
---|
| 165 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
| 166 | !!------------------------------------------------------------------------------ |
---|
| 167 | |
---|
| 168 | DO jj = njn0, njn1 |
---|
| 169 | DO jk = 1, jpkm1 |
---|
| 170 | DO ji = 1, jpi |
---|
| 171 | va(ji,jj,jk) = va(ji,jj,jk) + sqrt( grav*hvr (ji,jj) ) & |
---|
| 172 | & * ( ( sshn(ji,jj) + sshn(ji,jj+1) ) * 0.5 & |
---|
| 173 | & - sshfon(ji) ) * vnmsk(ji,jk) |
---|
| 174 | END DO |
---|
| 175 | END DO |
---|
| 176 | END DO |
---|
| 177 | DO jj = njn0p1, njn1p1 |
---|
| 178 | DO ji = 1, jpi |
---|
| 179 | sshn(ji,jj)= sshn(ji,jj) * (1.-tnmsk(ji,1)) + sshfon(ji)*tnmsk(ji,1) |
---|
| 180 | END DO |
---|
| 181 | END DO |
---|
| 182 | |
---|
| 183 | END SUBROUTINE obc_dyn_bt_north |
---|
| 184 | |
---|
| 185 | SUBROUTINE obc_dyn_bt_south |
---|
| 186 | !!------------------------------------------------------------------------------ |
---|
| 187 | !! *** SUBROUTINE obc_dyn_bt_south *** |
---|
| 188 | !! |
---|
| 189 | !! ** Purpose : |
---|
| 190 | !! Apply Flather algorithm on south OBC velocities ua, va |
---|
| 191 | !! Fix sea surface height (sshn) on south open boundary |
---|
| 192 | !! The logical lfbcsouth must be .TRUE. |
---|
| 193 | !! |
---|
| 194 | !! History : |
---|
| 195 | !! 9.0 ! 05-12 (V. Garnier) original |
---|
| 196 | !!------------------------------------------------------------------------------ |
---|
| 197 | !! * Local declaration |
---|
| 198 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
| 199 | |
---|
| 200 | !!------------------------------------------------------------------------------ |
---|
| 201 | |
---|
| 202 | DO jj = njs0, njs1 |
---|
| 203 | DO jk = 1, jpkm1 |
---|
| 204 | DO ji = 1, jpi |
---|
| 205 | va(ji,jj,jk) = va(ji,jj,jk) - sqrt( grav*hvr (ji,jj) ) & |
---|
| 206 | & * ( ( sshn(ji,jj) + sshn(ji,jj+1) ) * 0.5 & |
---|
| 207 | & - sshfos(ji) ) * vsmsk(ji,jk) |
---|
| 208 | END DO |
---|
| 209 | END DO |
---|
| 210 | DO ji = 1, jpi |
---|
| 211 | sshn(ji,jj)= sshn(ji,jj) * (1.-tsmsk(ji,1)) + tsmsk(ji,1) * sshfos(ji) |
---|
| 212 | END DO |
---|
| 213 | END DO |
---|
| 214 | |
---|
| 215 | END SUBROUTINE obc_dyn_bt_south |
---|
| 216 | |
---|
| 217 | # elif defined key_dynspg_ts |
---|
| 218 | |
---|
| 219 | SUBROUTINE obc_dyn_bt_east |
---|
| 220 | !!------------------------------------------------------------------------------ |
---|
| 221 | !! *** SUBROUTINE obc_dyn_bt_east *** |
---|
| 222 | !! |
---|
| 223 | !! ** Purpose : |
---|
| 224 | !! Apply Flather's algorithm on east OBC velocities ua, va |
---|
| 225 | !! Fix sea surface height (sshn) on east open boundary |
---|
| 226 | !! The logical lfbceast must be .TRUE. |
---|
| 227 | !! |
---|
| 228 | !! History : |
---|
| 229 | !! 9.0 ! 05-12 (V. Garnier) original |
---|
| 230 | !!------------------------------------------------------------------------------ |
---|
| 231 | !! * Local declaration |
---|
| 232 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
| 233 | !!------------------------------------------------------------------------------ |
---|
| 234 | |
---|
| 235 | DO ji = nie0, nie1 |
---|
| 236 | DO jk = 1, jpkm1 |
---|
| 237 | DO jj = 1, jpj |
---|
| 238 | ua(ji,jj,jk) = ( ua(ji,jj,jk) + sshfoe_b(ji,jj) ) * uemsk(jj,jk) |
---|
| 239 | END DO |
---|
| 240 | END DO |
---|
| 241 | END DO |
---|
| 242 | DO ji = nie0p1, nie1p1 |
---|
| 243 | DO jj = 1, jpj |
---|
| 244 | sshn(ji,jj) = sshn(ji,jj) * (1.-temsk(jj,1)) + temsk(jj,1)*sshn_b(ji,jj) |
---|
| 245 | END DO |
---|
| 246 | END DO |
---|
| 247 | |
---|
| 248 | END SUBROUTINE obc_dyn_bt_east |
---|
| 249 | |
---|
| 250 | SUBROUTINE obc_dyn_bt_west |
---|
| 251 | !!------------------------------------------------------------------------------ |
---|
| 252 | !! *** SUBROUTINE obc_dyn_bt_west *** |
---|
| 253 | !! |
---|
| 254 | !! ** Purpose : |
---|
| 255 | !! ** Purpose : |
---|
| 256 | !! Apply Flather algorithm on west OBC velocities ua, va |
---|
| 257 | !! Fix sea surface height (sshn) on west open boundary |
---|
| 258 | !! The logical lfbcwest must be .TRUE. |
---|
| 259 | !! |
---|
| 260 | !! History : |
---|
| 261 | !! 9.0 ! 05-12 (V. Garnier) original |
---|
| 262 | !!------------------------------------------------------------------------------ |
---|
| 263 | !! * Local declaration |
---|
| 264 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
| 265 | !!------------------------------------------------------------------------------ |
---|
| 266 | |
---|
| 267 | DO ji = niw0, niw1 |
---|
| 268 | DO jk = 1, jpkm1 |
---|
| 269 | DO jj = 1, jpj |
---|
| 270 | ua(ji,jj,jk) = ( ua(ji,jj,jk) + sshfow_b(ji,jj) ) * uwmsk(jj,jk) |
---|
| 271 | END DO |
---|
| 272 | END DO |
---|
| 273 | DO jj = 1, jpj |
---|
| 274 | sshn(ji,jj) = sshn(ji,jj) * (1.-twmsk(jj,1)) + twmsk(jj,1)*sshn_b(ji,jj) |
---|
| 275 | END DO |
---|
| 276 | END DO |
---|
| 277 | |
---|
| 278 | END SUBROUTINE obc_dyn_bt_west |
---|
| 279 | |
---|
| 280 | SUBROUTINE obc_dyn_bt_north |
---|
| 281 | !!------------------------------------------------------------------------------ |
---|
| 282 | !! SUBROUTINE obc_dyn_bt_north |
---|
| 283 | !! ************************* |
---|
| 284 | !! ** Purpose : |
---|
| 285 | !! Apply Flather algorithm on north OBC velocities ua, va |
---|
| 286 | !! Fix sea surface height (sshn) on north open boundary |
---|
| 287 | !! The logical lfbcnorth must be .TRUE. |
---|
| 288 | !! |
---|
| 289 | !! History : |
---|
| 290 | !! 9.0 ! 05-12 (V. Garnier) original |
---|
| 291 | !!------------------------------------------------------------------------------ |
---|
| 292 | !! * Local declaration |
---|
| 293 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
| 294 | !!------------------------------------------------------------------------------ |
---|
| 295 | |
---|
| 296 | DO jj = njn0, njn1 |
---|
| 297 | DO jk = 1, jpkm1 |
---|
| 298 | DO ji = 1, jpi |
---|
| 299 | va(ji,jj,jk) = ( va(ji,jj,jk) + sshfon_b(ji,jj) ) * vnmsk(jj,jk) |
---|
| 300 | END DO |
---|
| 301 | END DO |
---|
| 302 | END DO |
---|
| 303 | DO jj = njn0p1, njn1p1 |
---|
| 304 | DO ji = 1, jpi |
---|
| 305 | sshn(ji,jj)= sshn(ji,jj) * (1.-tnmsk(ji,1)) + sshn_b(ji,jj)*tnmsk(ji,1) |
---|
| 306 | END DO |
---|
| 307 | END DO |
---|
| 308 | |
---|
| 309 | END SUBROUTINE obc_dyn_bt_north |
---|
| 310 | |
---|
| 311 | SUBROUTINE obc_dyn_bt_south |
---|
| 312 | !!------------------------------------------------------------------------------ |
---|
| 313 | !! SUBROUTINE obc_dyn_bt_south |
---|
| 314 | !! ************************* |
---|
| 315 | !! ** Purpose : |
---|
| 316 | !! Apply Flather algorithm on south OBC velocities ua, va |
---|
| 317 | !! Fix sea surface height (sshn) on south open boundary |
---|
| 318 | !! The logical lfbcsouth must be .TRUE. |
---|
| 319 | !! |
---|
| 320 | !! History : |
---|
| 321 | !! 9.0 ! 05-12 (V. Garnier) original |
---|
| 322 | !!------------------------------------------------------------------------------ |
---|
| 323 | !! * Local declaration |
---|
| 324 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
| 325 | |
---|
| 326 | !!------------------------------------------------------------------------------ |
---|
| 327 | |
---|
| 328 | DO jj = njs0, njs1 |
---|
| 329 | DO jk = 1, jpkm1 |
---|
| 330 | DO ji = 1, jpi |
---|
| 331 | va(ji,jj,jk) = ( va(ji,jj,jk) + sshfos_b(ji,jj) ) * vsmsk(jj,jk) |
---|
| 332 | END DO |
---|
| 333 | END DO |
---|
| 334 | DO ji = 1, jpi |
---|
| 335 | sshn(ji,jj)= sshn(ji,jj) * (1.-tsmsk(ji,1)) + tsmsk(ji,1) * sshn_b(ji,jj) |
---|
| 336 | END DO |
---|
| 337 | END DO |
---|
| 338 | |
---|
| 339 | END SUBROUTINE obc_dyn_bt_south |
---|
| 340 | |
---|
| 341 | # endif |
---|
| 342 | #else |
---|
| 343 | !!================================================================================= |
---|
| 344 | !! *** MODULE obcdyn_bt *** |
---|
| 345 | !! Ocean dynamics: Radiation of velocities on each open boundary |
---|
| 346 | !!================================================================================= |
---|
| 347 | CONTAINS |
---|
| 348 | |
---|
| 349 | SUBROUTINE obc_dyn_bt |
---|
| 350 | ! No open boundaries ==> empty routine |
---|
| 351 | END SUBROUTINE obc_dyn_bt |
---|
| 352 | #endif |
---|
| 353 | |
---|
| 354 | END MODULE obcdyn_bt |
---|