[3117] | 1 | MODULE bdydyn3d |
---|
| 2 | !!====================================================================== |
---|
[3182] | 3 | !! *** MODULE bdydyn3d *** |
---|
[3191] | 4 | !! Unstructured Open Boundary Cond. : Flow relaxation scheme on baroclinic velocities |
---|
[3117] | 5 | !!====================================================================== |
---|
[3191] | 6 | !! History : 3.4 ! 2011 (D. Storkey) new module as part of BDY rewrite |
---|
[3680] | 7 | !! 3.5 ! 2012 (S. Mocavero, I. Epicoco) Optimization of BDY communications |
---|
[3117] | 8 | !!---------------------------------------------------------------------- |
---|
| 9 | !! bdy_dyn3d : apply open boundary conditions to baroclinic velocities |
---|
| 10 | !! bdy_dyn3d_frs : apply Flow Relaxation Scheme |
---|
| 11 | !!---------------------------------------------------------------------- |
---|
[3182] | 12 | USE timing ! Timing |
---|
[3117] | 13 | USE oce ! ocean dynamics and tracers |
---|
| 14 | USE dom_oce ! ocean space and time domain |
---|
| 15 | USE bdy_oce ! ocean open boundary conditions |
---|
[4292] | 16 | USE bdylib ! for orlanski library routines |
---|
[3117] | 17 | USE lbclnk ! ocean lateral boundary conditions (or mpp link) |
---|
| 18 | USE in_out_manager ! |
---|
[3651] | 19 | Use phycst |
---|
[3117] | 20 | |
---|
| 21 | IMPLICIT NONE |
---|
| 22 | PRIVATE |
---|
| 23 | |
---|
[3191] | 24 | PUBLIC bdy_dyn3d ! routine called by bdy_dyn |
---|
[3651] | 25 | PUBLIC bdy_dyn3d_dmp ! routine called by step |
---|
[3117] | 26 | |
---|
| 27 | !!---------------------------------------------------------------------- |
---|
| 28 | !! NEMO/OPA 3.3 , NEMO Consortium (2010) |
---|
[5215] | 29 | !! $Id$ |
---|
[3117] | 30 | !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) |
---|
| 31 | !!---------------------------------------------------------------------- |
---|
| 32 | CONTAINS |
---|
| 33 | |
---|
| 34 | SUBROUTINE bdy_dyn3d( kt ) |
---|
| 35 | !!---------------------------------------------------------------------- |
---|
| 36 | !! *** SUBROUTINE bdy_dyn3d *** |
---|
| 37 | !! |
---|
| 38 | !! ** Purpose : - Apply open boundary conditions for baroclinic velocities |
---|
| 39 | !! |
---|
| 40 | !!---------------------------------------------------------------------- |
---|
[6140] | 41 | INTEGER, INTENT(in) :: kt ! Main time step counter |
---|
| 42 | ! |
---|
| 43 | INTEGER :: ib_bdy ! loop index |
---|
| 44 | !!---------------------------------------------------------------------- |
---|
| 45 | ! |
---|
[3117] | 46 | DO ib_bdy=1, nb_bdy |
---|
[6140] | 47 | ! |
---|
[4292] | 48 | SELECT CASE( cn_dyn3d(ib_bdy) ) |
---|
[6140] | 49 | CASE('none') ; CYCLE |
---|
| 50 | CASE('frs' ) ; CALL bdy_dyn3d_frs( idx_bdy(ib_bdy), dta_bdy(ib_bdy), kt, ib_bdy ) |
---|
| 51 | CASE('specified') ; CALL bdy_dyn3d_spe( idx_bdy(ib_bdy), dta_bdy(ib_bdy), kt, ib_bdy ) |
---|
| 52 | CASE('zero') ; CALL bdy_dyn3d_zro( idx_bdy(ib_bdy), dta_bdy(ib_bdy), kt, ib_bdy ) |
---|
| 53 | CASE('orlanski' ) ; CALL bdy_dyn3d_orlanski( idx_bdy(ib_bdy), dta_bdy(ib_bdy), ib_bdy, ll_npo=.false. ) |
---|
| 54 | CASE('orlanski_npo'); CALL bdy_dyn3d_orlanski( idx_bdy(ib_bdy), dta_bdy(ib_bdy), ib_bdy, ll_npo=.true. ) |
---|
| 55 | CASE DEFAULT ; CALL ctl_stop( 'bdy_dyn3d : unrecognised option for open boundaries for baroclinic velocities' ) |
---|
[3117] | 56 | END SELECT |
---|
[6140] | 57 | END DO |
---|
| 58 | ! |
---|
[3117] | 59 | END SUBROUTINE bdy_dyn3d |
---|
| 60 | |
---|
[6140] | 61 | |
---|
[3680] | 62 | SUBROUTINE bdy_dyn3d_spe( idx, dta, kt , ib_bdy ) |
---|
[3651] | 63 | !!---------------------------------------------------------------------- |
---|
| 64 | !! *** SUBROUTINE bdy_dyn3d_spe *** |
---|
| 65 | !! |
---|
| 66 | !! ** Purpose : - Apply a specified value for baroclinic velocities |
---|
| 67 | !! at open boundaries. |
---|
| 68 | !! |
---|
| 69 | !!---------------------------------------------------------------------- |
---|
[6140] | 70 | INTEGER , INTENT(in) :: kt ! time step index |
---|
| 71 | TYPE(OBC_INDEX), INTENT(in) :: idx ! OBC indices |
---|
| 72 | TYPE(OBC_DATA) , INTENT(in) :: dta ! OBC external data |
---|
| 73 | INTEGER , INTENT(in) :: ib_bdy ! BDY set index |
---|
| 74 | ! |
---|
[3651] | 75 | INTEGER :: jb, jk ! dummy loop indices |
---|
| 76 | INTEGER :: ii, ij, igrd ! local integers |
---|
| 77 | REAL(wp) :: zwgt ! boundary weight |
---|
| 78 | !!---------------------------------------------------------------------- |
---|
| 79 | ! |
---|
| 80 | IF( nn_timing == 1 ) CALL timing_start('bdy_dyn3d_spe') |
---|
| 81 | ! |
---|
| 82 | igrd = 2 ! Relaxation of zonal velocity |
---|
| 83 | DO jb = 1, idx%nblenrim(igrd) |
---|
| 84 | DO jk = 1, jpkm1 |
---|
| 85 | ii = idx%nbi(jb,igrd) |
---|
| 86 | ij = idx%nbj(jb,igrd) |
---|
| 87 | ua(ii,ij,jk) = dta%u3d(jb,jk) * umask(ii,ij,jk) |
---|
| 88 | END DO |
---|
| 89 | END DO |
---|
| 90 | ! |
---|
| 91 | igrd = 3 ! Relaxation of meridional velocity |
---|
| 92 | DO jb = 1, idx%nblenrim(igrd) |
---|
| 93 | DO jk = 1, jpkm1 |
---|
| 94 | ii = idx%nbi(jb,igrd) |
---|
| 95 | ij = idx%nbj(jb,igrd) |
---|
| 96 | va(ii,ij,jk) = dta%v3d(jb,jk) * vmask(ii,ij,jk) |
---|
| 97 | END DO |
---|
| 98 | END DO |
---|
[4292] | 99 | CALL lbc_bdy_lnk( ua, 'U', -1., ib_bdy ) ! Boundary points should be updated |
---|
| 100 | CALL lbc_bdy_lnk( va, 'V', -1., ib_bdy ) |
---|
[3651] | 101 | ! |
---|
[6140] | 102 | IF( kt == nit000 ) CLOSE( unit = 102 ) |
---|
| 103 | ! |
---|
[3651] | 104 | IF( nn_timing == 1 ) CALL timing_stop('bdy_dyn3d_spe') |
---|
[6140] | 105 | ! |
---|
[3651] | 106 | END SUBROUTINE bdy_dyn3d_spe |
---|
| 107 | |
---|
[6140] | 108 | |
---|
[3680] | 109 | SUBROUTINE bdy_dyn3d_zro( idx, dta, kt, ib_bdy ) |
---|
[3651] | 110 | !!---------------------------------------------------------------------- |
---|
| 111 | !! *** SUBROUTINE bdy_dyn3d_zro *** |
---|
| 112 | !! |
---|
| 113 | !! ** Purpose : - baroclinic velocities = 0. at open boundaries. |
---|
| 114 | !! |
---|
| 115 | !!---------------------------------------------------------------------- |
---|
[6140] | 116 | INTEGER , INTENT(in) :: kt ! time step index |
---|
| 117 | TYPE(OBC_INDEX), INTENT(in) :: idx ! OBC indices |
---|
| 118 | TYPE(OBC_DATA) , INTENT(in) :: dta ! OBC external data |
---|
[3680] | 119 | INTEGER, INTENT(in) :: ib_bdy ! BDY set index |
---|
[6140] | 120 | ! |
---|
[3651] | 121 | INTEGER :: ib, ik ! dummy loop indices |
---|
[6140] | 122 | INTEGER :: ii, ij, igrd ! local integers |
---|
[3651] | 123 | REAL(wp) :: zwgt ! boundary weight |
---|
| 124 | !!---------------------------------------------------------------------- |
---|
| 125 | ! |
---|
| 126 | IF( nn_timing == 1 ) CALL timing_start('bdy_dyn3d_zro') |
---|
| 127 | ! |
---|
| 128 | igrd = 2 ! Everything is at T-points here |
---|
| 129 | DO ib = 1, idx%nblenrim(igrd) |
---|
| 130 | ii = idx%nbi(ib,igrd) |
---|
| 131 | ij = idx%nbj(ib,igrd) |
---|
| 132 | DO ik = 1, jpkm1 |
---|
| 133 | ua(ii,ij,ik) = 0._wp |
---|
| 134 | END DO |
---|
| 135 | END DO |
---|
| 136 | |
---|
| 137 | igrd = 3 ! Everything is at T-points here |
---|
| 138 | DO ib = 1, idx%nblenrim(igrd) |
---|
| 139 | ii = idx%nbi(ib,igrd) |
---|
| 140 | ij = idx%nbj(ib,igrd) |
---|
| 141 | DO ik = 1, jpkm1 |
---|
| 142 | va(ii,ij,ik) = 0._wp |
---|
| 143 | END DO |
---|
| 144 | END DO |
---|
| 145 | ! |
---|
[3680] | 146 | CALL lbc_bdy_lnk( ua, 'U', -1., ib_bdy ) ; CALL lbc_bdy_lnk( va, 'V', -1.,ib_bdy ) ! Boundary points should be updated |
---|
[3651] | 147 | ! |
---|
[6140] | 148 | IF( kt == nit000 ) CLOSE( unit = 102 ) |
---|
| 149 | ! |
---|
| 150 | IF( nn_timing == 1 ) CALL timing_stop('bdy_dyn3d_zro') |
---|
| 151 | ! |
---|
| 152 | END SUBROUTINE bdy_dyn3d_zro |
---|
[3651] | 153 | |
---|
| 154 | |
---|
[3680] | 155 | SUBROUTINE bdy_dyn3d_frs( idx, dta, kt, ib_bdy ) |
---|
[3117] | 156 | !!---------------------------------------------------------------------- |
---|
| 157 | !! *** SUBROUTINE bdy_dyn3d_frs *** |
---|
| 158 | !! |
---|
| 159 | !! ** Purpose : - Apply the Flow Relaxation Scheme for baroclinic velocities |
---|
| 160 | !! at open boundaries. |
---|
| 161 | !! |
---|
| 162 | !! References :- Engedahl H., 1995: Use of the flow relaxation scheme in |
---|
| 163 | !! a three-dimensional baroclinic ocean model with realistic |
---|
| 164 | !! topography. Tellus, 365-382. |
---|
| 165 | !!---------------------------------------------------------------------- |
---|
[6140] | 166 | INTEGER , INTENT(in) :: kt ! time step index |
---|
| 167 | TYPE(OBC_INDEX), INTENT(in) :: idx ! OBC indices |
---|
| 168 | TYPE(OBC_DATA) , INTENT(in) :: dta ! OBC external data |
---|
[3680] | 169 | INTEGER, INTENT(in) :: ib_bdy ! BDY set index |
---|
[6140] | 170 | ! |
---|
[3117] | 171 | INTEGER :: jb, jk ! dummy loop indices |
---|
| 172 | INTEGER :: ii, ij, igrd ! local integers |
---|
| 173 | REAL(wp) :: zwgt ! boundary weight |
---|
| 174 | !!---------------------------------------------------------------------- |
---|
| 175 | ! |
---|
[3182] | 176 | IF( nn_timing == 1 ) CALL timing_start('bdy_dyn3d_frs') |
---|
[3117] | 177 | ! |
---|
| 178 | igrd = 2 ! Relaxation of zonal velocity |
---|
| 179 | DO jb = 1, idx%nblen(igrd) |
---|
| 180 | DO jk = 1, jpkm1 |
---|
| 181 | ii = idx%nbi(jb,igrd) |
---|
| 182 | ij = idx%nbj(jb,igrd) |
---|
| 183 | zwgt = idx%nbw(jb,igrd) |
---|
| 184 | ua(ii,ij,jk) = ( ua(ii,ij,jk) + zwgt * ( dta%u3d(jb,jk) - ua(ii,ij,jk) ) ) * umask(ii,ij,jk) |
---|
| 185 | END DO |
---|
| 186 | END DO |
---|
| 187 | ! |
---|
| 188 | igrd = 3 ! Relaxation of meridional velocity |
---|
| 189 | DO jb = 1, idx%nblen(igrd) |
---|
| 190 | DO jk = 1, jpkm1 |
---|
| 191 | ii = idx%nbi(jb,igrd) |
---|
| 192 | ij = idx%nbj(jb,igrd) |
---|
| 193 | zwgt = idx%nbw(jb,igrd) |
---|
| 194 | va(ii,ij,jk) = ( va(ii,ij,jk) + zwgt * ( dta%v3d(jb,jk) - va(ii,ij,jk) ) ) * vmask(ii,ij,jk) |
---|
| 195 | END DO |
---|
| 196 | END DO |
---|
[4292] | 197 | CALL lbc_bdy_lnk( ua, 'U', -1., ib_bdy ) ! Boundary points should be updated |
---|
| 198 | CALL lbc_bdy_lnk( va, 'V', -1., ib_bdy ) |
---|
[3117] | 199 | ! |
---|
[6140] | 200 | IF( kt == nit000 ) CLOSE( unit = 102 ) |
---|
| 201 | ! |
---|
| 202 | IF( nn_timing == 1 ) CALL timing_stop('bdy_dyn3d_frs') |
---|
| 203 | ! |
---|
| 204 | END SUBROUTINE bdy_dyn3d_frs |
---|
[3117] | 205 | |
---|
[3182] | 206 | |
---|
[4292] | 207 | SUBROUTINE bdy_dyn3d_orlanski( idx, dta, ib_bdy, ll_npo ) |
---|
| 208 | !!---------------------------------------------------------------------- |
---|
| 209 | !! *** SUBROUTINE bdy_dyn3d_orlanski *** |
---|
| 210 | !! |
---|
| 211 | !! - Apply Orlanski radiation to baroclinic velocities. |
---|
| 212 | !! - Wrapper routine for bdy_orlanski_3d |
---|
| 213 | !! |
---|
| 214 | !! |
---|
| 215 | !! References: Marchesiello, McWilliams and Shchepetkin, Ocean Modelling vol. 3 (2001) |
---|
| 216 | !!---------------------------------------------------------------------- |
---|
| 217 | TYPE(OBC_INDEX), INTENT(in) :: idx ! OBC indices |
---|
| 218 | TYPE(OBC_DATA), INTENT(in) :: dta ! OBC external data |
---|
| 219 | INTEGER, INTENT(in) :: ib_bdy ! BDY set index |
---|
| 220 | LOGICAL, INTENT(in) :: ll_npo ! switch for NPO version |
---|
| 221 | |
---|
| 222 | INTEGER :: jb, igrd ! dummy loop indices |
---|
| 223 | !!---------------------------------------------------------------------- |
---|
| 224 | |
---|
| 225 | IF( nn_timing == 1 ) CALL timing_start('bdy_dyn3d_orlanski') |
---|
| 226 | ! |
---|
| 227 | !! Note that at this stage the ub and ua arrays contain the baroclinic velocities. |
---|
| 228 | ! |
---|
| 229 | igrd = 2 ! Orlanski bc on u-velocity; |
---|
| 230 | ! |
---|
| 231 | CALL bdy_orlanski_3d( idx, igrd, ub, ua, dta%u3d, ll_npo ) |
---|
| 232 | |
---|
| 233 | igrd = 3 ! Orlanski bc on v-velocity |
---|
| 234 | ! |
---|
| 235 | CALL bdy_orlanski_3d( idx, igrd, vb, va, dta%v3d, ll_npo ) |
---|
| 236 | ! |
---|
| 237 | CALL lbc_bdy_lnk( ua, 'U', -1., ib_bdy ) ! Boundary points should be updated |
---|
| 238 | CALL lbc_bdy_lnk( va, 'V', -1., ib_bdy ) |
---|
| 239 | ! |
---|
| 240 | IF( nn_timing == 1 ) CALL timing_stop('bdy_dyn3d_orlanski') |
---|
| 241 | ! |
---|
| 242 | END SUBROUTINE bdy_dyn3d_orlanski |
---|
| 243 | |
---|
| 244 | |
---|
[3651] | 245 | SUBROUTINE bdy_dyn3d_dmp( kt ) |
---|
| 246 | !!---------------------------------------------------------------------- |
---|
| 247 | !! *** SUBROUTINE bdy_dyn3d_dmp *** |
---|
| 248 | !! |
---|
| 249 | !! ** Purpose : Apply damping for baroclinic velocities at open boundaries. |
---|
| 250 | !! |
---|
| 251 | !!---------------------------------------------------------------------- |
---|
[6140] | 252 | INTEGER, INTENT(in) :: kt ! time step index |
---|
| 253 | ! |
---|
[3651] | 254 | INTEGER :: jb, jk ! dummy loop indices |
---|
[6140] | 255 | INTEGER :: ib_bdy ! loop index |
---|
[3651] | 256 | INTEGER :: ii, ij, igrd ! local integers |
---|
| 257 | REAL(wp) :: zwgt ! boundary weight |
---|
| 258 | !!---------------------------------------------------------------------- |
---|
| 259 | ! |
---|
[6140] | 260 | IF( nn_timing == 1 ) CALL timing_start('bdy_dyn3d_dmp') |
---|
[3651] | 261 | ! |
---|
| 262 | DO ib_bdy=1, nb_bdy |
---|
[4292] | 263 | IF ( ln_dyn3d_dmp(ib_bdy) .and. cn_dyn3d(ib_bdy) /= 'none' ) THEN |
---|
[3651] | 264 | igrd = 2 ! Relaxation of zonal velocity |
---|
| 265 | DO jb = 1, idx_bdy(ib_bdy)%nblen(igrd) |
---|
| 266 | ii = idx_bdy(ib_bdy)%nbi(jb,igrd) |
---|
| 267 | ij = idx_bdy(ib_bdy)%nbj(jb,igrd) |
---|
| 268 | zwgt = idx_bdy(ib_bdy)%nbd(jb,igrd) |
---|
| 269 | DO jk = 1, jpkm1 |
---|
[3703] | 270 | ua(ii,ij,jk) = ( ua(ii,ij,jk) + zwgt * ( dta_bdy(ib_bdy)%u3d(jb,jk) - & |
---|
[4354] | 271 | ub(ii,ij,jk) + ub_b(ii,ij)) ) * umask(ii,ij,jk) |
---|
[3651] | 272 | END DO |
---|
| 273 | END DO |
---|
| 274 | ! |
---|
| 275 | igrd = 3 ! Relaxation of meridional velocity |
---|
| 276 | DO jb = 1, idx_bdy(ib_bdy)%nblen(igrd) |
---|
| 277 | ii = idx_bdy(ib_bdy)%nbi(jb,igrd) |
---|
| 278 | ij = idx_bdy(ib_bdy)%nbj(jb,igrd) |
---|
| 279 | zwgt = idx_bdy(ib_bdy)%nbd(jb,igrd) |
---|
| 280 | DO jk = 1, jpkm1 |
---|
[3703] | 281 | va(ii,ij,jk) = ( va(ii,ij,jk) + zwgt * ( dta_bdy(ib_bdy)%v3d(jb,jk) - & |
---|
[4354] | 282 | vb(ii,ij,jk) + vb_b(ii,ij)) ) * vmask(ii,ij,jk) |
---|
[3651] | 283 | END DO |
---|
| 284 | END DO |
---|
| 285 | ENDIF |
---|
[6140] | 286 | END DO |
---|
[3651] | 287 | ! |
---|
| 288 | CALL lbc_lnk( ua, 'U', -1. ) ; CALL lbc_lnk( va, 'V', -1. ) ! Boundary points should be updated |
---|
| 289 | ! |
---|
[6140] | 290 | IF( nn_timing == 1 ) CALL timing_stop('bdy_dyn3d_dmp') |
---|
| 291 | ! |
---|
[3651] | 292 | END SUBROUTINE bdy_dyn3d_dmp |
---|
| 293 | |
---|
[3117] | 294 | !!====================================================================== |
---|
| 295 | END MODULE bdydyn3d |
---|