[8586] | 1 | MODULE bdyice |
---|
| 2 | !!====================================================================== |
---|
| 3 | !! *** MODULE bdyice *** |
---|
[9656] | 4 | !! Unstructured Open Boundary Cond. : Open boundary conditions for sea-ice (SI3) |
---|
[8586] | 5 | !!====================================================================== |
---|
| 6 | !! History : 3.3 ! 2010-09 (D. Storkey) Original code |
---|
[9656] | 7 | !! 3.4 ! 2012-01 (C. Rousset) add new sea ice model |
---|
| 8 | !! 4.0 ! 2018 (C. Rousset) SI3 compatibility |
---|
[8586] | 9 | !!---------------------------------------------------------------------- |
---|
[9570] | 10 | #if defined key_si3 |
---|
[8586] | 11 | !!---------------------------------------------------------------------- |
---|
[9656] | 12 | !! 'key_si3' SI3 sea ice model |
---|
[8586] | 13 | !!---------------------------------------------------------------------- |
---|
| 14 | !! bdy_ice : Application of open boundaries to ice |
---|
| 15 | !! bdy_ice_frs : Application of Flow Relaxation Scheme |
---|
| 16 | !!---------------------------------------------------------------------- |
---|
| 17 | USE oce ! ocean dynamics and tracers variables |
---|
[8637] | 18 | USE ice ! sea-ice: variables |
---|
| 19 | USE icevar ! sea-ice: operations |
---|
[9880] | 20 | USE icecor ! sea-ice: corrections |
---|
[8637] | 21 | USE icectl ! sea-ice: control prints |
---|
[8586] | 22 | USE phycst ! physical constant |
---|
| 23 | USE eosbn2 ! equation of state |
---|
| 24 | USE par_oce ! ocean parameters |
---|
| 25 | USE dom_oce ! ocean space and time domain variables |
---|
| 26 | USE sbc_oce ! Surface boundary condition: ocean fields |
---|
| 27 | USE bdy_oce ! ocean open boundary conditions |
---|
| 28 | ! |
---|
| 29 | USE lbclnk ! ocean lateral boundary conditions (or mpp link) |
---|
| 30 | USE in_out_manager ! write to numout file |
---|
| 31 | USE lib_mpp ! distributed memory computing |
---|
| 32 | USE lib_fortran ! to use key_nosignedzero |
---|
| 33 | USE timing ! Timing |
---|
| 34 | |
---|
| 35 | IMPLICIT NONE |
---|
| 36 | PRIVATE |
---|
| 37 | |
---|
| 38 | PUBLIC bdy_ice ! routine called in sbcmod |
---|
[9656] | 39 | PUBLIC bdy_ice_dyn ! routine called in icedyn_rhg_evp |
---|
[8586] | 40 | |
---|
| 41 | !!---------------------------------------------------------------------- |
---|
[9598] | 42 | !! NEMO/OCE 4.0 , NEMO Consortium (2018) |
---|
[10069] | 43 | !! $Id$ |
---|
[10068] | 44 | !! Software governed by the CeCILL license (see ./LICENSE) |
---|
[8586] | 45 | !!---------------------------------------------------------------------- |
---|
| 46 | CONTAINS |
---|
| 47 | |
---|
| 48 | SUBROUTINE bdy_ice( kt ) |
---|
| 49 | !!---------------------------------------------------------------------- |
---|
| 50 | !! *** SUBROUTINE bdy_ice *** |
---|
| 51 | !! |
---|
[9890] | 52 | !! ** Purpose : Apply open boundary conditions for sea ice |
---|
[8586] | 53 | !! |
---|
| 54 | !!---------------------------------------------------------------------- |
---|
| 55 | INTEGER, INTENT(in) :: kt ! Main time step counter |
---|
| 56 | ! |
---|
[9890] | 57 | INTEGER :: jbdy ! BDY set index |
---|
[8586] | 58 | !!---------------------------------------------------------------------- |
---|
[11042] | 59 | ! controls |
---|
| 60 | IF( ln_timing ) CALL timing_start('bdy_ice_thd') ! timing |
---|
| 61 | IF( ln_icediachk ) CALL ice_cons_hsm(0,'bdy_ice_thd', rdiag_v, rdiag_s, rdiag_t, rdiag_fv, rdiag_fs, rdiag_ft) ! conservation |
---|
[8586] | 62 | ! |
---|
| 63 | CALL ice_var_glo2eqv |
---|
| 64 | ! |
---|
[9890] | 65 | DO jbdy = 1, nb_bdy |
---|
[8586] | 66 | ! |
---|
[9890] | 67 | SELECT CASE( cn_ice(jbdy) ) |
---|
[9124] | 68 | CASE('none') ; CYCLE |
---|
[9890] | 69 | CASE('frs' ) ; CALL bdy_ice_frs( idx_bdy(jbdy), dta_bdy(jbdy), kt, jbdy ) |
---|
[8586] | 70 | CASE DEFAULT |
---|
| 71 | CALL ctl_stop( 'bdy_ice : unrecognised option for open boundaries for ice fields' ) |
---|
| 72 | END SELECT |
---|
| 73 | ! |
---|
| 74 | END DO |
---|
| 75 | ! |
---|
[9880] | 76 | CALL ice_cor( kt , 0 ) ! -- In case categories are out of bounds, do a remapping |
---|
[9885] | 77 | ! ! i.e. inputs have not the same ice thickness distribution (set by rn_himean) |
---|
| 78 | ! ! than the regional simulation |
---|
[9124] | 79 | CALL ice_var_agg(1) |
---|
| 80 | ! |
---|
[11042] | 81 | ! controls |
---|
| 82 | IF( ln_icediachk ) CALL ice_cons_hsm(1,'bdy_ice_thd', rdiag_v, rdiag_s, rdiag_t, rdiag_fv, rdiag_fs, rdiag_ft) ! conservation |
---|
| 83 | IF( ln_icectl ) CALL ice_prt ( kt, iiceprt, jiceprt, 1, ' - ice thermo bdy - ' ) ! prints |
---|
| 84 | IF( ln_timing ) CALL timing_stop ('bdy_ice_thd') ! timing |
---|
[8586] | 85 | ! |
---|
| 86 | END SUBROUTINE bdy_ice |
---|
| 87 | |
---|
| 88 | |
---|
[9890] | 89 | SUBROUTINE bdy_ice_frs( idx, dta, kt, jbdy ) |
---|
[8586] | 90 | !!------------------------------------------------------------------------------ |
---|
| 91 | !! *** SUBROUTINE bdy_ice_frs *** |
---|
| 92 | !! |
---|
[9890] | 93 | !! ** Purpose : Apply the Flow Relaxation Scheme for sea-ice fields |
---|
[8586] | 94 | !! |
---|
| 95 | !! Reference : Engedahl H., 1995: Use of the flow relaxation scheme in a three- |
---|
| 96 | !! dimensional baroclinic ocean model with realistic topography. Tellus, 365-382. |
---|
| 97 | !!------------------------------------------------------------------------------ |
---|
| 98 | TYPE(OBC_INDEX), INTENT(in) :: idx ! OBC indices |
---|
| 99 | TYPE(OBC_DATA), INTENT(in) :: dta ! OBC external data |
---|
| 100 | INTEGER, INTENT(in) :: kt ! main time-step counter |
---|
[9890] | 101 | INTEGER, INTENT(in) :: jbdy ! BDY set index |
---|
[8586] | 102 | ! |
---|
| 103 | INTEGER :: jpbound ! 0 = incoming ice |
---|
| 104 | ! ! 1 = outgoing ice |
---|
[9890] | 105 | INTEGER :: i_bdy, jgrd ! dummy loop indices |
---|
| 106 | INTEGER :: ji, jj, jk, jl, ib, jb |
---|
[8586] | 107 | REAL(wp) :: zwgt, zwgt1 ! local scalar |
---|
| 108 | REAL(wp) :: ztmelts, zdh |
---|
[11048] | 109 | REAL(wp), POINTER :: flagu, flagv ! short cuts |
---|
[8586] | 110 | !!------------------------------------------------------------------------------ |
---|
| 111 | ! |
---|
| 112 | jgrd = 1 ! Everything is at T-points here |
---|
| 113 | ! |
---|
| 114 | DO jl = 1, jpl |
---|
[9890] | 115 | DO i_bdy = 1, idx%nblenrim(jgrd) |
---|
| 116 | ji = idx%nbi(i_bdy,jgrd) |
---|
| 117 | jj = idx%nbj(i_bdy,jgrd) |
---|
[11048] | 118 | IF( ji == 1 .OR. ji == jpi .OR. jj == 1 .OR. jj == jpj ) CYCLE ! to remove |
---|
[9890] | 119 | zwgt = idx%nbw(i_bdy,jgrd) |
---|
| 120 | zwgt1 = 1.e0 - idx%nbw(i_bdy,jgrd) |
---|
| 121 | a_i(ji,jj,jl) = ( a_i(ji,jj,jl) * zwgt1 + dta%a_i(i_bdy,jl) * zwgt ) * tmask(ji,jj,1) ! Leads fraction |
---|
| 122 | h_i(ji,jj,jl) = ( h_i(ji,jj,jl) * zwgt1 + dta%h_i(i_bdy,jl) * zwgt ) * tmask(ji,jj,1) ! Ice depth |
---|
| 123 | h_s(ji,jj,jl) = ( h_s(ji,jj,jl) * zwgt1 + dta%h_s(i_bdy,jl) * zwgt ) * tmask(ji,jj,1) ! Snow depth |
---|
[8586] | 124 | |
---|
| 125 | ! ----------------- |
---|
| 126 | ! Pathological case |
---|
| 127 | ! ----------------- |
---|
| 128 | ! In case a) snow load would be in excess or b) ice is coming into a warmer environment that would lead to |
---|
| 129 | ! very large transformation from snow to ice (see icethd_dh.F90) |
---|
| 130 | |
---|
| 131 | ! Then, a) transfer the snow excess into the ice (different from icethd_dh) |
---|
[9935] | 132 | zdh = MAX( 0._wp, ( rhos * h_s(ji,jj,jl) + ( rhoi - rau0 ) * h_i(ji,jj,jl) ) * r1_rau0 ) |
---|
[8586] | 133 | ! Or, b) transfer all the snow into ice (if incoming ice is likely to melt as it comes into a warmer environment) |
---|
[9935] | 134 | !zdh = MAX( 0._wp, h_s(ji,jj,jl) * rhos / rhoi ) |
---|
[8586] | 135 | |
---|
| 136 | ! recompute h_i, h_s |
---|
| 137 | h_i(ji,jj,jl) = MIN( hi_max(jl), h_i(ji,jj,jl) + zdh ) |
---|
[9935] | 138 | h_s(ji,jj,jl) = MAX( 0._wp, h_s(ji,jj,jl) - zdh * rhoi / rhos ) |
---|
[8586] | 139 | |
---|
| 140 | ENDDO |
---|
| 141 | ENDDO |
---|
[10425] | 142 | CALL lbc_bdy_lnk( 'bdyice', a_i(:,:,:), 'T', 1., jbdy ) |
---|
| 143 | CALL lbc_bdy_lnk( 'bdyice', h_i(:,:,:), 'T', 1., jbdy ) |
---|
| 144 | CALL lbc_bdy_lnk( 'bdyice', h_s(:,:,:), 'T', 1., jbdy ) |
---|
[8586] | 145 | |
---|
| 146 | DO jl = 1, jpl |
---|
[9890] | 147 | DO i_bdy = 1, idx%nblenrim(jgrd) |
---|
| 148 | ji = idx%nbi(i_bdy,jgrd) |
---|
| 149 | jj = idx%nbj(i_bdy,jgrd) |
---|
[11048] | 150 | IF( ji == 1 .OR. ji == jpi .OR. jj == 1 .OR. jj == jpj ) CYCLE ! to remove |
---|
| 151 | flagu => idx%flagu(i_bdy,jgrd) |
---|
| 152 | flagv => idx%flagv(i_bdy,jgrd) |
---|
[8586] | 153 | ! condition on ice thickness depends on the ice velocity |
---|
| 154 | ! if velocity is outward (strictly), then ice thickness, volume... must be equal to adjacent values |
---|
[9890] | 155 | jpbound = 0 ; ib = ji ; jb = jj |
---|
[8586] | 156 | ! |
---|
[11048] | 157 | IF( u_ice(ji ,jj ) < 0. .AND. flagu == 1. ) jpbound = 1 ; ib = ji+1 |
---|
| 158 | IF( u_ice(ji-1,jj ) > 0. .AND. flagu == -1. ) jpbound = 1 ; ib = ji-1 |
---|
| 159 | IF( v_ice(ji ,jj ) < 0. .AND. flagv == 1. ) jpbound = 1 ; jb = jj+1 |
---|
| 160 | IF( v_ice(ji ,jj-1) > 0. .AND. flagv == -1. ) jpbound = 1 ; jb = jj-1 |
---|
[8586] | 161 | ! |
---|
[9890] | 162 | IF( nn_ice_dta(jbdy) == 0 ) jpbound = 0 ; ib = ji ; jb = jj ! case ice boundaries = initial conditions |
---|
| 163 | ! ! do not make state variables dependent on velocity |
---|
[8586] | 164 | ! |
---|
[9890] | 165 | IF( a_i(ib,jb,jl) > 0._wp ) THEN ! there is ice at the boundary |
---|
[8586] | 166 | ! |
---|
[9890] | 167 | a_i(ji,jj,jl) = a_i(ib,jb,jl) ! concentration |
---|
| 168 | h_i(ji,jj,jl) = h_i(ib,jb,jl) ! thickness ice |
---|
| 169 | h_s(ji,jj,jl) = h_s(ib,jb,jl) ! thickness snw |
---|
[8586] | 170 | ! |
---|
[9888] | 171 | SELECT CASE( jpbound ) |
---|
| 172 | ! |
---|
| 173 | CASE( 0 ) ! velocity is inward |
---|
| 174 | ! |
---|
[9890] | 175 | oa_i(ji,jj, jl) = rn_ice_age(jbdy) * a_i(ji,jj,jl) ! age |
---|
| 176 | a_ip(ji,jj, jl) = 0._wp ! pond concentration |
---|
| 177 | v_ip(ji,jj, jl) = 0._wp ! pond volume |
---|
| 178 | t_su(ji,jj, jl) = rn_ice_tem(jbdy) ! temperature surface |
---|
| 179 | t_s (ji,jj,:,jl) = rn_ice_tem(jbdy) ! temperature snw |
---|
| 180 | t_i (ji,jj,:,jl) = rn_ice_tem(jbdy) ! temperature ice |
---|
| 181 | s_i (ji,jj, jl) = rn_ice_sal(jbdy) ! salinity |
---|
| 182 | sz_i(ji,jj,:,jl) = rn_ice_sal(jbdy) ! salinity profile |
---|
[9888] | 183 | ! |
---|
| 184 | CASE( 1 ) ! velocity is outward |
---|
| 185 | ! |
---|
[9890] | 186 | oa_i(ji,jj, jl) = oa_i(ib,jb, jl) ! age |
---|
| 187 | a_ip(ji,jj, jl) = a_ip(ib,jb, jl) ! pond concentration |
---|
| 188 | v_ip(ji,jj, jl) = v_ip(ib,jb, jl) ! pond volume |
---|
| 189 | t_su(ji,jj, jl) = t_su(ib,jb, jl) ! temperature surface |
---|
| 190 | t_s (ji,jj,:,jl) = t_s (ib,jb,:,jl) ! temperature snw |
---|
| 191 | t_i (ji,jj,:,jl) = t_i (ib,jb,:,jl) ! temperature ice |
---|
| 192 | s_i (ji,jj, jl) = s_i (ib,jb, jl) ! salinity |
---|
| 193 | sz_i(ji,jj,:,jl) = sz_i(ib,jb,:,jl) ! salinity profile |
---|
[9888] | 194 | ! |
---|
| 195 | END SELECT |
---|
[9885] | 196 | ! |
---|
[9888] | 197 | IF( nn_icesal == 1 ) THEN ! if constant salinity |
---|
| 198 | s_i (ji,jj ,jl) = rn_icesal |
---|
| 199 | sz_i(ji,jj,:,jl) = rn_icesal |
---|
| 200 | ENDIF |
---|
[8586] | 201 | ! |
---|
[9888] | 202 | ! global fields |
---|
| 203 | v_i (ji,jj,jl) = h_i(ji,jj,jl) * a_i(ji,jj,jl) ! volume ice |
---|
| 204 | v_s (ji,jj,jl) = h_s(ji,jj,jl) * a_i(ji,jj,jl) ! volume snw |
---|
| 205 | sv_i(ji,jj,jl) = MIN( s_i(ji,jj,jl) , sss_m(ji,jj) ) * v_i(ji,jj,jl) ! salt content |
---|
[8586] | 206 | DO jk = 1, nlay_s |
---|
[9935] | 207 | e_s(ji,jj,jk,jl) = rhos * ( rcpi * ( rt0 - t_s(ji,jj,jk,jl) ) + rLfus ) ! enthalpy in J/m3 |
---|
[9888] | 208 | e_s(ji,jj,jk,jl) = e_s(ji,jj,jk,jl) * v_s(ji,jj,jl) * r1_nlay_s ! enthalpy in J/m2 |
---|
| 209 | END DO |
---|
[8586] | 210 | DO jk = 1, nlay_i |
---|
[9935] | 211 | ztmelts = - rTmlt * sz_i(ji,jj,jk,jl) ! Melting temperature in C |
---|
[9888] | 212 | t_i(ji,jj,jk,jl) = MIN( t_i(ji,jj,jk,jl), ztmelts + rt0 ) ! Force t_i to be lower than melting point => likely conservation issue |
---|
| 213 | ! |
---|
[9935] | 214 | e_i(ji,jj,jk,jl) = rhoi * ( rcpi * ( ztmelts - ( t_i(ji,jj,jk,jl) - rt0 ) ) & ! enthalpy in J/m3 |
---|
| 215 | & + rLfus * ( 1._wp - ztmelts / ( t_i(ji,jj,jk,jl) - rt0 ) ) & |
---|
| 216 | & - rcp * ztmelts ) |
---|
[9888] | 217 | e_i(ji,jj,jk,jl) = e_i(ji,jj,jk,jl) * v_i(ji,jj,jl) * r1_nlay_i ! enthalpy in J/m2 |
---|
[8586] | 218 | END DO |
---|
| 219 | ! |
---|
[9888] | 220 | ELSE ! no ice at the boundary |
---|
[9885] | 221 | ! |
---|
[9888] | 222 | a_i (ji,jj, jl) = 0._wp |
---|
| 223 | h_i (ji,jj, jl) = 0._wp |
---|
| 224 | h_s (ji,jj, jl) = 0._wp |
---|
| 225 | oa_i(ji,jj, jl) = 0._wp |
---|
| 226 | a_ip(ji,jj, jl) = 0._wp |
---|
| 227 | v_ip(ji,jj, jl) = 0._wp |
---|
| 228 | t_su(ji,jj, jl) = rt0 |
---|
| 229 | t_s (ji,jj,:,jl) = rt0 |
---|
| 230 | t_i (ji,jj,:,jl) = rt0 |
---|
| 231 | |
---|
| 232 | IF( nn_icesal == 1 ) THEN ! if constant salinity |
---|
| 233 | s_i (ji,jj ,jl) = rn_icesal |
---|
| 234 | sz_i(ji,jj,:,jl) = rn_icesal |
---|
| 235 | ELSE ! if variable salinity |
---|
| 236 | s_i (ji,jj,jl) = rn_simin |
---|
| 237 | sz_i(ji,jj,:,jl) = rn_simin |
---|
| 238 | ENDIF |
---|
| 239 | ! |
---|
| 240 | ! global fields |
---|
| 241 | v_i (ji,jj, jl) = 0._wp |
---|
| 242 | v_s (ji,jj, jl) = 0._wp |
---|
| 243 | sv_i(ji,jj, jl) = 0._wp |
---|
| 244 | e_s (ji,jj,:,jl) = 0._wp |
---|
| 245 | e_i (ji,jj,:,jl) = 0._wp |
---|
| 246 | |
---|
[8586] | 247 | ENDIF |
---|
[9888] | 248 | |
---|
[8586] | 249 | END DO |
---|
| 250 | ! |
---|
[9888] | 251 | END DO ! jl |
---|
[9890] | 252 | |
---|
[10425] | 253 | CALL lbc_bdy_lnk( 'bdyice', a_i (:,:,:) , 'T', 1., jbdy ) |
---|
| 254 | CALL lbc_bdy_lnk( 'bdyice', h_i (:,:,:) , 'T', 1., jbdy ) |
---|
| 255 | CALL lbc_bdy_lnk( 'bdyice', h_s (:,:,:) , 'T', 1., jbdy ) |
---|
| 256 | CALL lbc_bdy_lnk( 'bdyice', oa_i(:,:,:) , 'T', 1., jbdy ) |
---|
| 257 | CALL lbc_bdy_lnk( 'bdyice', a_ip(:,:,:) , 'T', 1., jbdy ) |
---|
| 258 | CALL lbc_bdy_lnk( 'bdyice', v_ip(:,:,:) , 'T', 1., jbdy ) |
---|
| 259 | CALL lbc_bdy_lnk( 'bdyice', s_i (:,:,:) , 'T', 1., jbdy ) |
---|
| 260 | CALL lbc_bdy_lnk( 'bdyice', t_su(:,:,:) , 'T', 1., jbdy ) |
---|
| 261 | CALL lbc_bdy_lnk( 'bdyice', v_i (:,:,:) , 'T', 1., jbdy ) |
---|
| 262 | CALL lbc_bdy_lnk( 'bdyice', v_s (:,:,:) , 'T', 1., jbdy ) |
---|
| 263 | CALL lbc_bdy_lnk( 'bdyice', sv_i(:,:,:) , 'T', 1., jbdy ) |
---|
| 264 | CALL lbc_bdy_lnk( 'bdyice', t_s (:,:,:,:), 'T', 1., jbdy ) |
---|
| 265 | CALL lbc_bdy_lnk( 'bdyice', e_s (:,:,:,:), 'T', 1., jbdy ) |
---|
| 266 | CALL lbc_bdy_lnk( 'bdyice', t_i (:,:,:,:), 'T', 1., jbdy ) |
---|
| 267 | CALL lbc_bdy_lnk( 'bdyice', e_i (:,:,:,:), 'T', 1., jbdy ) |
---|
[8586] | 268 | ! |
---|
| 269 | END SUBROUTINE bdy_ice_frs |
---|
| 270 | |
---|
| 271 | |
---|
| 272 | SUBROUTINE bdy_ice_dyn( cd_type ) |
---|
| 273 | !!------------------------------------------------------------------------------ |
---|
| 274 | !! *** SUBROUTINE bdy_ice_dyn *** |
---|
| 275 | !! |
---|
[9890] | 276 | !! ** Purpose : Apply dynamics boundary conditions for sea-ice. |
---|
[8586] | 277 | !! |
---|
[9890] | 278 | !! ** Method : if this adjacent grid point is not ice free, then u_ice and v_ice take its value |
---|
| 279 | !! if is ice free, then u_ice and v_ice are unchanged by BDY |
---|
| 280 | !! they keep values calculated in rheology |
---|
| 281 | !! |
---|
[8586] | 282 | !!------------------------------------------------------------------------------ |
---|
| 283 | CHARACTER(len=1), INTENT(in) :: cd_type ! nature of velocity grid-points |
---|
| 284 | ! |
---|
[9890] | 285 | INTEGER :: i_bdy, jgrd ! dummy loop indices |
---|
| 286 | INTEGER :: ji, jj ! local scalar |
---|
| 287 | INTEGER :: jbdy ! BDY set index |
---|
[8586] | 288 | REAL(wp) :: zmsk1, zmsk2, zflag |
---|
| 289 | !!------------------------------------------------------------------------------ |
---|
[9905] | 290 | IF( ln_timing ) CALL timing_start('bdy_ice_dyn') |
---|
[8586] | 291 | ! |
---|
[9890] | 292 | DO jbdy=1, nb_bdy |
---|
[8586] | 293 | ! |
---|
[9890] | 294 | SELECT CASE( cn_ice(jbdy) ) |
---|
[8586] | 295 | ! |
---|
| 296 | CASE('none') |
---|
| 297 | CYCLE |
---|
| 298 | ! |
---|
| 299 | CASE('frs') |
---|
| 300 | ! |
---|
[9890] | 301 | IF( nn_ice_dta(jbdy) == 0 ) CYCLE ! case ice boundaries = initial conditions |
---|
| 302 | ! ! do not change ice velocity (it is only computed by rheology) |
---|
[8586] | 303 | SELECT CASE ( cd_type ) |
---|
| 304 | ! |
---|
| 305 | CASE ( 'U' ) |
---|
| 306 | jgrd = 2 ! u velocity |
---|
[9890] | 307 | DO i_bdy = 1, idx_bdy(jbdy)%nblenrim(jgrd) |
---|
| 308 | ji = idx_bdy(jbdy)%nbi(i_bdy,jgrd) |
---|
| 309 | jj = idx_bdy(jbdy)%nbj(i_bdy,jgrd) |
---|
[11048] | 310 | IF( ji == 1 .OR. ji == jpi .OR. jj == 1 .OR. jj == jpj ) CYCLE ! to remove |
---|
[9890] | 311 | zflag = idx_bdy(jbdy)%flagu(i_bdy,jgrd) |
---|
[11048] | 312 | ! i-1 i i | ! i i i+1 | ! i i i+1 | |
---|
| 313 | ! > ice > | ! o > ice | ! o > o | |
---|
| 314 | ! => set at u_ice(i-1) ! => set to O ! => unchanged |
---|
| 315 | IF( zflag == -1. .AND. ji > 1 .AND. ji < jpi ) THEN |
---|
| 316 | IF ( vt_i(ji ,jj) > 0. ) THEN ; u_ice(ji,jj) = u_ice(ji-1,jj) |
---|
| 317 | ELSEIF( vt_i(ji+1,jj) > 0. ) THEN ; u_ice(ji,jj) = 0._wp |
---|
| 318 | END IF |
---|
| 319 | END IF |
---|
| 320 | ! | i i+1 i+1 ! | i i i+1 ! | i i i+1 |
---|
| 321 | ! | > ice > ! | ice > o ! | o > o |
---|
| 322 | ! => set at u_ice(i+1) ! => set to O ! => unchanged |
---|
| 323 | IF( zflag == 1. .AND. ji+1 < jpi+1 ) THEN |
---|
| 324 | IF ( vt_i(ji+1,jj) > 0. ) THEN ; u_ice(ji,jj) = u_ice(ji+1,jj) |
---|
| 325 | ELSEIF( vt_i(ji ,jj) > 0. ) THEN ; u_ice(ji,jj) = 0._wp |
---|
| 326 | END IF |
---|
| 327 | END IF |
---|
[8586] | 328 | ! |
---|
[11048] | 329 | IF( zflag == 0. ) u_ice(ji,jj) = 0._wp ! u_ice = 0 if north/south bdy |
---|
[8586] | 330 | ! |
---|
| 331 | END DO |
---|
[10425] | 332 | CALL lbc_bdy_lnk( 'bdyice', u_ice(:,:), 'U', -1., jbdy ) |
---|
[8586] | 333 | ! |
---|
| 334 | CASE ( 'V' ) |
---|
| 335 | jgrd = 3 ! v velocity |
---|
[9890] | 336 | DO i_bdy = 1, idx_bdy(jbdy)%nblenrim(jgrd) |
---|
| 337 | ji = idx_bdy(jbdy)%nbi(i_bdy,jgrd) |
---|
| 338 | jj = idx_bdy(jbdy)%nbj(i_bdy,jgrd) |
---|
[11048] | 339 | IF( ji == 1 .OR. ji == jpi .OR. jj == 1 .OR. jj == jpj ) CYCLE ! to remove |
---|
[9890] | 340 | zflag = idx_bdy(jbdy)%flagv(i_bdy,jgrd) |
---|
[11048] | 341 | ! ¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨ ! ¨¨¨¨ïce¨¨¨(jj+1)¨¨ ! ¨¨¨¨¨¨ö¨¨¨¨(jj+1) |
---|
| 342 | ! ^ (jj ) ! ^ (jj ) ! ^ (jj ) |
---|
| 343 | ! ice (jj ) ! o (jj ) ! o (jj ) |
---|
| 344 | ! ^ (jj-1) ! ! |
---|
| 345 | ! => set to u_ice(jj-1) ! => set to 0 ! => unchanged |
---|
| 346 | IF( zflag == -1. .AND. jj > 1 .AND. jj < jpj ) THEN |
---|
| 347 | IF ( vt_i(ji,jj ) > 0. ) THEN ; v_ice(ji,jj) = v_ice(ji,jj-1) |
---|
| 348 | ELSEIF( vt_i(ji,jj+1) > 0. ) THEN ; v_ice(ji,jj) = 0._wp |
---|
| 349 | END IF |
---|
| 350 | END IF |
---|
| 351 | ! ^ (jj+1) ! ! |
---|
| 352 | ! ice (jj+1) ! o (jj+1) ! o (jj+1) |
---|
| 353 | ! ^ (jj ) ! ^ (jj ) ! ^ (jj ) |
---|
| 354 | ! ________________ ! ____ice___(jj )_ ! _____o____(jj ) |
---|
| 355 | ! => set to u_ice(jj+1) ! => set to 0 ! => unchanged |
---|
| 356 | IF( zflag == 1. .AND. jj < jpj ) THEN |
---|
| 357 | IF ( vt_i(ji,jj+1) > 0. ) THEN ; v_ice(ji,jj) = v_ice(ji,jj+1) |
---|
| 358 | ELSEIF( vt_i(ji,jj ) > 0. ) THEN ; v_ice(ji,jj) = 0._wp |
---|
| 359 | END IF |
---|
| 360 | END IF |
---|
[8586] | 361 | ! |
---|
[11048] | 362 | IF( zflag == 0. ) v_ice(ji,jj) = 0._wp ! v_ice = 0 if west/east bdy |
---|
[8586] | 363 | ! |
---|
| 364 | END DO |
---|
[10425] | 365 | CALL lbc_bdy_lnk( 'bdyice', v_ice(:,:), 'V', -1., jbdy ) |
---|
[8586] | 366 | ! |
---|
| 367 | END SELECT |
---|
| 368 | ! |
---|
| 369 | CASE DEFAULT |
---|
| 370 | CALL ctl_stop( 'bdy_ice_dyn : unrecognised option for open boundaries for ice fields' ) |
---|
| 371 | END SELECT |
---|
| 372 | ! |
---|
| 373 | END DO |
---|
| 374 | ! |
---|
[9905] | 375 | IF( ln_timing ) CALL timing_stop('bdy_ice_dyn') |
---|
| 376 | ! |
---|
[8586] | 377 | END SUBROUTINE bdy_ice_dyn |
---|
| 378 | |
---|
| 379 | #else |
---|
| 380 | !!--------------------------------------------------------------------------------- |
---|
| 381 | !! Default option Empty module |
---|
| 382 | !!--------------------------------------------------------------------------------- |
---|
| 383 | CONTAINS |
---|
| 384 | SUBROUTINE bdy_ice( kt ) ! Empty routine |
---|
[9927] | 385 | IMPLICIT NONE |
---|
| 386 | INTEGER, INTENT( in ) :: kt |
---|
[8586] | 387 | WRITE(*,*) 'bdy_ice: You should not have seen this print! error?', kt |
---|
| 388 | END SUBROUTINE bdy_ice |
---|
| 389 | #endif |
---|
| 390 | |
---|
| 391 | !!================================================================================= |
---|
| 392 | END MODULE bdyice |
---|