[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 | ! |
---|
[11536] | 57 | INTEGER :: jbdy, ir ! BDY set index, rim index |
---|
| 58 | INTEGER :: ibeg, iend ! length of rim to be treated (rim 0 or rim 1) |
---|
| 59 | LOGICAL :: llrim0 ! indicate if rim 0 is treated |
---|
| 60 | LOGICAL, DIMENSION(4) :: llsend1, llrecv1 ! indicate how communications are to be carried out |
---|
[8586] | 61 | !!---------------------------------------------------------------------- |
---|
[11041] | 62 | ! controls |
---|
[13601] | 63 | IF( ln_timing ) CALL timing_start('bdy_ice_thd') ! timing |
---|
[8586] | 64 | ! |
---|
| 65 | CALL ice_var_glo2eqv |
---|
| 66 | ! |
---|
[11536] | 67 | llsend1(:) = .false. ; llrecv1(:) = .false. |
---|
| 68 | DO ir = 1, 0, -1 ! treat rim 1 before rim 0 |
---|
| 69 | IF( ir == 0 ) THEN ; llrim0 = .TRUE. |
---|
| 70 | ELSE ; llrim0 = .FALSE. |
---|
| 71 | END IF |
---|
| 72 | DO jbdy = 1, nb_bdy |
---|
| 73 | ! |
---|
| 74 | SELECT CASE( cn_ice(jbdy) ) |
---|
| 75 | CASE('none') ; CYCLE |
---|
| 76 | CASE('frs' ) ; CALL bdy_ice_frs( idx_bdy(jbdy), dta_bdy(jbdy), kt, jbdy, llrim0 ) |
---|
| 77 | CASE DEFAULT |
---|
| 78 | CALL ctl_stop( 'bdy_ice : unrecognised option for open boundaries for ice fields' ) |
---|
| 79 | END SELECT |
---|
| 80 | ! |
---|
| 81 | END DO |
---|
[8586] | 82 | ! |
---|
[11536] | 83 | ! Update bdy points |
---|
| 84 | IF( nn_hls > 1 .AND. ir == 1 ) CYCLE ! at least 2 halos will be corrected -> no need to correct rim 1 before rim 0 |
---|
| 85 | IF( nn_hls == 1 ) THEN ; llsend1(:) = .false. ; llrecv1(:) = .false. ; END IF |
---|
| 86 | DO jbdy = 1, nb_bdy |
---|
| 87 | IF( cn_ice(jbdy) == 'frs' ) THEN |
---|
| 88 | llsend1(:) = llsend1(:) .OR. lsend_bdyint(jbdy,1,:,ir) ! possibly every direction, T points |
---|
| 89 | llrecv1(:) = llrecv1(:) .OR. lrecv_bdyint(jbdy,1,:,ir) ! possibly every direction, T points |
---|
| 90 | END IF |
---|
| 91 | END DO ! jbdy |
---|
| 92 | IF( ANY(llsend1) .OR. ANY(llrecv1) ) THEN ! if need to send/recv in at least one direction |
---|
| 93 | ! exchange 3d arrays |
---|
[13472] | 94 | CALL lbc_lnk_multi('bdyice', a_i , 'T', 1._wp, h_i , 'T', 1._wp, h_s , 'T', 1._wp, oa_i, 'T', 1._wp & |
---|
| 95 | & , s_i , 'T', 1._wp, t_su, 'T', 1._wp, v_i , 'T', 1._wp, v_s , 'T', 1._wp, sv_i, 'T', 1._wp & |
---|
| 96 | & , a_ip, 'T', 1._wp, v_ip, 'T', 1._wp, v_il, 'T', 1._wp & |
---|
| 97 | & , kfillmode=jpfillnothing ,lsend=llsend1, lrecv=llrecv1 ) |
---|
[11536] | 98 | ! exchange 4d arrays : third dimension = 1 and then third dimension = jpk |
---|
[13472] | 99 | CALL lbc_lnk_multi('bdyice', t_s , 'T', 1._wp, e_s , 'T', 1._wp, kfillmode=jpfillnothing ,lsend=llsend1, lrecv=llrecv1 ) |
---|
| 100 | CALL lbc_lnk_multi('bdyice', t_i , 'T', 1._wp, e_i , 'T', 1._wp, kfillmode=jpfillnothing ,lsend=llsend1, lrecv=llrecv1 ) |
---|
[11536] | 101 | END IF |
---|
| 102 | END DO ! ir |
---|
[8586] | 103 | ! |
---|
[9880] | 104 | CALL ice_cor( kt , 0 ) ! -- In case categories are out of bounds, do a remapping |
---|
[9885] | 105 | ! ! i.e. inputs have not the same ice thickness distribution (set by rn_himean) |
---|
| 106 | ! ! than the regional simulation |
---|
[9124] | 107 | CALL ice_var_agg(1) |
---|
| 108 | ! |
---|
[11041] | 109 | ! controls |
---|
[13601] | 110 | IF( ln_icectl ) CALL ice_prt ( kt, iiceprt, jiceprt, 1, ' - ice thermo bdy - ' ) ! prints |
---|
| 111 | IF( ln_timing ) CALL timing_stop ('bdy_ice_thd') ! timing |
---|
[8586] | 112 | ! |
---|
| 113 | END SUBROUTINE bdy_ice |
---|
| 114 | |
---|
| 115 | |
---|
[11536] | 116 | SUBROUTINE bdy_ice_frs( idx, dta, kt, jbdy, llrim0 ) |
---|
[8586] | 117 | !!------------------------------------------------------------------------------ |
---|
| 118 | !! *** SUBROUTINE bdy_ice_frs *** |
---|
| 119 | !! |
---|
[9890] | 120 | !! ** Purpose : Apply the Flow Relaxation Scheme for sea-ice fields |
---|
[8586] | 121 | !! |
---|
| 122 | !! Reference : Engedahl H., 1995: Use of the flow relaxation scheme in a three- |
---|
| 123 | !! dimensional baroclinic ocean model with realistic topography. Tellus, 365-382. |
---|
| 124 | !!------------------------------------------------------------------------------ |
---|
[11536] | 125 | TYPE(OBC_INDEX), INTENT(in) :: idx ! OBC indices |
---|
| 126 | TYPE(OBC_DATA), INTENT(in) :: dta ! OBC external data |
---|
| 127 | INTEGER, INTENT(in) :: kt ! main time-step counter |
---|
| 128 | INTEGER, INTENT(in) :: jbdy ! BDY set index |
---|
| 129 | LOGICAL, INTENT(in) :: llrim0 ! indicate if rim 0 is treated |
---|
[8586] | 130 | ! |
---|
| 131 | INTEGER :: jpbound ! 0 = incoming ice |
---|
| 132 | ! ! 1 = outgoing ice |
---|
[11536] | 133 | INTEGER :: ibeg, iend ! length of rim to be treated (rim 0 or rim 1) |
---|
[9890] | 134 | INTEGER :: i_bdy, jgrd ! dummy loop indices |
---|
| 135 | INTEGER :: ji, jj, jk, jl, ib, jb |
---|
[8586] | 136 | REAL(wp) :: zwgt, zwgt1 ! local scalar |
---|
| 137 | REAL(wp) :: ztmelts, zdh |
---|
[11536] | 138 | REAL(wp), POINTER :: flagu, flagv ! short cuts |
---|
[8586] | 139 | !!------------------------------------------------------------------------------ |
---|
| 140 | ! |
---|
| 141 | jgrd = 1 ! Everything is at T-points here |
---|
[11536] | 142 | IF( llrim0 ) THEN ; ibeg = 1 ; iend = idx%nblenrim0(jgrd) |
---|
| 143 | ELSE ; ibeg = idx%nblenrim0(jgrd)+1 ; iend = idx%nblenrim(jgrd) |
---|
| 144 | END IF |
---|
[8586] | 145 | ! |
---|
| 146 | DO jl = 1, jpl |
---|
[11536] | 147 | DO i_bdy = ibeg, iend |
---|
[9890] | 148 | ji = idx%nbi(i_bdy,jgrd) |
---|
| 149 | jj = idx%nbj(i_bdy,jgrd) |
---|
| 150 | zwgt = idx%nbw(i_bdy,jgrd) |
---|
| 151 | zwgt1 = 1.e0 - idx%nbw(i_bdy,jgrd) |
---|
[11536] | 152 | a_i (ji,jj, jl) = ( a_i (ji,jj, jl) * zwgt1 + dta%a_i(i_bdy,jl) * zwgt ) * tmask(ji,jj,1) ! Ice concentration |
---|
| 153 | h_i (ji,jj, jl) = ( h_i (ji,jj, jl) * zwgt1 + dta%h_i(i_bdy,jl) * zwgt ) * tmask(ji,jj,1) ! Ice depth |
---|
| 154 | h_s (ji,jj, jl) = ( h_s (ji,jj, jl) * zwgt1 + dta%h_s(i_bdy,jl) * zwgt ) * tmask(ji,jj,1) ! Snow depth |
---|
| 155 | t_i (ji,jj,:,jl) = ( t_i (ji,jj,:,jl) * zwgt1 + dta%t_i(i_bdy,jl) * zwgt ) * tmask(ji,jj,1) ! Ice temperature |
---|
| 156 | t_s (ji,jj,:,jl) = ( t_s (ji,jj,:,jl) * zwgt1 + dta%t_s(i_bdy,jl) * zwgt ) * tmask(ji,jj,1) ! Snow temperature |
---|
| 157 | t_su(ji,jj, jl) = ( t_su(ji,jj, jl) * zwgt1 + dta%tsu(i_bdy,jl) * zwgt ) * tmask(ji,jj,1) ! Surf temperature |
---|
| 158 | s_i (ji,jj, jl) = ( s_i (ji,jj, jl) * zwgt1 + dta%s_i(i_bdy,jl) * zwgt ) * tmask(ji,jj,1) ! Ice salinity |
---|
| 159 | a_ip(ji,jj, jl) = ( a_ip(ji,jj, jl) * zwgt1 + dta%aip(i_bdy,jl) * zwgt ) * tmask(ji,jj,1) ! Ice pond concentration |
---|
| 160 | h_ip(ji,jj, jl) = ( h_ip(ji,jj, jl) * zwgt1 + dta%hip(i_bdy,jl) * zwgt ) * tmask(ji,jj,1) ! Ice pond depth |
---|
[13472] | 161 | h_il(ji,jj, jl) = ( h_il(ji,jj, jl) * zwgt1 + dta%hil(i_bdy,jl) * zwgt ) * tmask(ji,jj,1) ! Ice pond lid depth |
---|
[11536] | 162 | ! |
---|
| 163 | sz_i(ji,jj,:,jl) = s_i(ji,jj,jl) |
---|
| 164 | ! |
---|
| 165 | ! make sure ponds = 0 if no ponds scheme |
---|
| 166 | IF( .NOT.ln_pnd ) THEN |
---|
| 167 | a_ip(ji,jj,jl) = 0._wp |
---|
| 168 | h_ip(ji,jj,jl) = 0._wp |
---|
[13472] | 169 | h_il(ji,jj,jl) = 0._wp |
---|
[11536] | 170 | ENDIF |
---|
[13472] | 171 | |
---|
| 172 | IF( .NOT.ln_pnd_lids ) THEN |
---|
| 173 | h_il(ji,jj,jl) = 0._wp |
---|
| 174 | ENDIF |
---|
[11536] | 175 | ! |
---|
[8586] | 176 | ! ----------------- |
---|
| 177 | ! Pathological case |
---|
| 178 | ! ----------------- |
---|
| 179 | ! In case a) snow load would be in excess or b) ice is coming into a warmer environment that would lead to |
---|
| 180 | ! very large transformation from snow to ice (see icethd_dh.F90) |
---|
| 181 | |
---|
| 182 | ! Then, a) transfer the snow excess into the ice (different from icethd_dh) |
---|
[12489] | 183 | zdh = MAX( 0._wp, ( rhos * h_s(ji,jj,jl) + ( rhoi - rho0 ) * h_i(ji,jj,jl) ) * r1_rho0 ) |
---|
[8586] | 184 | ! Or, b) transfer all the snow into ice (if incoming ice is likely to melt as it comes into a warmer environment) |
---|
[9935] | 185 | !zdh = MAX( 0._wp, h_s(ji,jj,jl) * rhos / rhoi ) |
---|
[8586] | 186 | |
---|
| 187 | ! recompute h_i, h_s |
---|
| 188 | h_i(ji,jj,jl) = MIN( hi_max(jl), h_i(ji,jj,jl) + zdh ) |
---|
[9935] | 189 | h_s(ji,jj,jl) = MAX( 0._wp, h_s(ji,jj,jl) - zdh * rhoi / rhos ) |
---|
[11536] | 190 | ! |
---|
[8586] | 191 | ENDDO |
---|
| 192 | ENDDO |
---|
| 193 | |
---|
| 194 | DO jl = 1, jpl |
---|
[11536] | 195 | DO i_bdy = ibeg, iend |
---|
[9890] | 196 | ji = idx%nbi(i_bdy,jgrd) |
---|
| 197 | jj = idx%nbj(i_bdy,jgrd) |
---|
[11536] | 198 | flagu => idx%flagu(i_bdy,jgrd) |
---|
| 199 | flagv => idx%flagv(i_bdy,jgrd) |
---|
[8586] | 200 | ! condition on ice thickness depends on the ice velocity |
---|
| 201 | ! if velocity is outward (strictly), then ice thickness, volume... must be equal to adjacent values |
---|
[9890] | 202 | jpbound = 0 ; ib = ji ; jb = jj |
---|
[8586] | 203 | ! |
---|
[11536] | 204 | IF( flagu == 1. ) THEN |
---|
| 205 | IF( ji+1 > jpi ) CYCLE |
---|
| 206 | IF( u_ice(ji ,jj ) < 0. ) jpbound = 1 ; ib = ji+1 |
---|
| 207 | END IF |
---|
| 208 | IF( flagu == -1. ) THEN |
---|
| 209 | IF( ji-1 < 1 ) CYCLE |
---|
| 210 | IF( u_ice(ji-1,jj ) < 0. ) jpbound = 1 ; ib = ji-1 |
---|
| 211 | END IF |
---|
| 212 | IF( flagv == 1. ) THEN |
---|
| 213 | IF( jj+1 > jpj ) CYCLE |
---|
| 214 | IF( v_ice(ji ,jj ) < 0. ) jpbound = 1 ; jb = jj+1 |
---|
| 215 | END IF |
---|
| 216 | IF( flagv == -1. ) THEN |
---|
| 217 | IF( jj-1 < 1 ) CYCLE |
---|
| 218 | IF( v_ice(ji ,jj-1) < 0. ) jpbound = 1 ; jb = jj-1 |
---|
| 219 | END IF |
---|
[8586] | 220 | ! |
---|
[9890] | 221 | IF( nn_ice_dta(jbdy) == 0 ) jpbound = 0 ; ib = ji ; jb = jj ! case ice boundaries = initial conditions |
---|
| 222 | ! ! do not make state variables dependent on velocity |
---|
[8586] | 223 | ! |
---|
[9890] | 224 | IF( a_i(ib,jb,jl) > 0._wp ) THEN ! there is ice at the boundary |
---|
[8586] | 225 | ! |
---|
[11536] | 226 | a_i (ji,jj, jl) = a_i (ib,jb, jl) |
---|
| 227 | h_i (ji,jj, jl) = h_i (ib,jb, jl) |
---|
| 228 | h_s (ji,jj, jl) = h_s (ib,jb, jl) |
---|
| 229 | t_i (ji,jj,:,jl) = t_i (ib,jb,:,jl) |
---|
| 230 | t_s (ji,jj,:,jl) = t_s (ib,jb,:,jl) |
---|
| 231 | t_su(ji,jj, jl) = t_su(ib,jb, jl) |
---|
| 232 | s_i (ji,jj, jl) = s_i (ib,jb, jl) |
---|
| 233 | a_ip(ji,jj, jl) = a_ip(ib,jb, jl) |
---|
| 234 | h_ip(ji,jj, jl) = h_ip(ib,jb, jl) |
---|
[13472] | 235 | h_il(ji,jj, jl) = h_il(ib,jb, jl) |
---|
[8586] | 236 | ! |
---|
[11536] | 237 | sz_i(ji,jj,:,jl) = sz_i(ib,jb,:,jl) |
---|
[9885] | 238 | ! |
---|
[11536] | 239 | ! ice age |
---|
| 240 | IF ( jpbound == 0 ) THEN ! velocity is inward |
---|
| 241 | oa_i(ji,jj,jl) = rice_age(jbdy) * a_i(ji,jj,jl) |
---|
| 242 | ELSEIF( jpbound == 1 ) THEN ! velocity is outward |
---|
| 243 | oa_i(ji,jj,jl) = oa_i(ib,jb,jl) |
---|
| 244 | ENDIF |
---|
| 245 | ! |
---|
[9888] | 246 | IF( nn_icesal == 1 ) THEN ! if constant salinity |
---|
| 247 | s_i (ji,jj ,jl) = rn_icesal |
---|
| 248 | sz_i(ji,jj,:,jl) = rn_icesal |
---|
| 249 | ENDIF |
---|
[8586] | 250 | ! |
---|
[9888] | 251 | ! global fields |
---|
| 252 | v_i (ji,jj,jl) = h_i(ji,jj,jl) * a_i(ji,jj,jl) ! volume ice |
---|
| 253 | v_s (ji,jj,jl) = h_s(ji,jj,jl) * a_i(ji,jj,jl) ! volume snw |
---|
| 254 | sv_i(ji,jj,jl) = MIN( s_i(ji,jj,jl) , sss_m(ji,jj) ) * v_i(ji,jj,jl) ! salt content |
---|
[8586] | 255 | DO jk = 1, nlay_s |
---|
[9935] | 256 | e_s(ji,jj,jk,jl) = rhos * ( rcpi * ( rt0 - t_s(ji,jj,jk,jl) ) + rLfus ) ! enthalpy in J/m3 |
---|
[9888] | 257 | e_s(ji,jj,jk,jl) = e_s(ji,jj,jk,jl) * v_s(ji,jj,jl) * r1_nlay_s ! enthalpy in J/m2 |
---|
| 258 | END DO |
---|
[8586] | 259 | DO jk = 1, nlay_i |
---|
[9935] | 260 | ztmelts = - rTmlt * sz_i(ji,jj,jk,jl) ! Melting temperature in C |
---|
[9888] | 261 | 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 |
---|
| 262 | ! |
---|
[9935] | 263 | e_i(ji,jj,jk,jl) = rhoi * ( rcpi * ( ztmelts - ( t_i(ji,jj,jk,jl) - rt0 ) ) & ! enthalpy in J/m3 |
---|
| 264 | & + rLfus * ( 1._wp - ztmelts / ( t_i(ji,jj,jk,jl) - rt0 ) ) & |
---|
| 265 | & - rcp * ztmelts ) |
---|
[9888] | 266 | 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] | 267 | END DO |
---|
| 268 | ! |
---|
[11536] | 269 | ! melt ponds |
---|
| 270 | v_ip(ji,jj,jl) = h_ip(ji,jj,jl) * a_ip(ji,jj,jl) |
---|
[13472] | 271 | v_il(ji,jj,jl) = h_il(ji,jj,jl) * a_ip(ji,jj,jl) |
---|
[11536] | 272 | ! |
---|
[9888] | 273 | ELSE ! no ice at the boundary |
---|
[9885] | 274 | ! |
---|
[9888] | 275 | a_i (ji,jj, jl) = 0._wp |
---|
| 276 | h_i (ji,jj, jl) = 0._wp |
---|
| 277 | h_s (ji,jj, jl) = 0._wp |
---|
| 278 | oa_i(ji,jj, jl) = 0._wp |
---|
| 279 | t_su(ji,jj, jl) = rt0 |
---|
| 280 | t_s (ji,jj,:,jl) = rt0 |
---|
| 281 | t_i (ji,jj,:,jl) = rt0 |
---|
[11536] | 282 | |
---|
[13472] | 283 | a_ip(ji,jj,jl) = 0._wp |
---|
| 284 | h_ip(ji,jj,jl) = 0._wp |
---|
| 285 | h_il(ji,jj,jl) = 0._wp |
---|
[9888] | 286 | |
---|
| 287 | IF( nn_icesal == 1 ) THEN ! if constant salinity |
---|
| 288 | s_i (ji,jj ,jl) = rn_icesal |
---|
| 289 | sz_i(ji,jj,:,jl) = rn_icesal |
---|
| 290 | ELSE ! if variable salinity |
---|
| 291 | s_i (ji,jj,jl) = rn_simin |
---|
| 292 | sz_i(ji,jj,:,jl) = rn_simin |
---|
| 293 | ENDIF |
---|
| 294 | ! |
---|
| 295 | ! global fields |
---|
| 296 | v_i (ji,jj, jl) = 0._wp |
---|
| 297 | v_s (ji,jj, jl) = 0._wp |
---|
| 298 | sv_i(ji,jj, jl) = 0._wp |
---|
| 299 | e_s (ji,jj,:,jl) = 0._wp |
---|
| 300 | e_i (ji,jj,:,jl) = 0._wp |
---|
[13472] | 301 | v_ip(ji,jj, jl) = 0._wp |
---|
| 302 | v_il(ji,jj, jl) = 0._wp |
---|
[9888] | 303 | |
---|
[8586] | 304 | ENDIF |
---|
[9888] | 305 | |
---|
[8586] | 306 | END DO |
---|
| 307 | ! |
---|
[9888] | 308 | END DO ! jl |
---|
[8586] | 309 | ! |
---|
| 310 | END SUBROUTINE bdy_ice_frs |
---|
| 311 | |
---|
| 312 | |
---|
| 313 | SUBROUTINE bdy_ice_dyn( cd_type ) |
---|
| 314 | !!------------------------------------------------------------------------------ |
---|
| 315 | !! *** SUBROUTINE bdy_ice_dyn *** |
---|
| 316 | !! |
---|
[9890] | 317 | !! ** Purpose : Apply dynamics boundary conditions for sea-ice. |
---|
[8586] | 318 | !! |
---|
[9890] | 319 | !! ** Method : if this adjacent grid point is not ice free, then u_ice and v_ice take its value |
---|
| 320 | !! if is ice free, then u_ice and v_ice are unchanged by BDY |
---|
| 321 | !! they keep values calculated in rheology |
---|
| 322 | !! |
---|
[8586] | 323 | !!------------------------------------------------------------------------------ |
---|
| 324 | CHARACTER(len=1), INTENT(in) :: cd_type ! nature of velocity grid-points |
---|
| 325 | ! |
---|
[11536] | 326 | INTEGER :: i_bdy, jgrd ! dummy loop indices |
---|
| 327 | INTEGER :: ji, jj ! local scalar |
---|
| 328 | INTEGER :: jbdy, ir ! BDY set index, rim index |
---|
| 329 | INTEGER :: ibeg, iend ! length of rim to be treated (rim 0 or rim 1) |
---|
[8586] | 330 | REAL(wp) :: zmsk1, zmsk2, zflag |
---|
[11536] | 331 | LOGICAL, DIMENSION(4) :: llsend2, llrecv2, llsend3, llrecv3 ! indicate how communications are to be carried out |
---|
[8586] | 332 | !!------------------------------------------------------------------------------ |
---|
[9905] | 333 | IF( ln_timing ) CALL timing_start('bdy_ice_dyn') |
---|
[8586] | 334 | ! |
---|
[11536] | 335 | llsend2(:) = .false. ; llrecv2(:) = .false. |
---|
| 336 | llsend3(:) = .false. ; llrecv3(:) = .false. |
---|
| 337 | DO ir = 1, 0, -1 |
---|
| 338 | DO jbdy = 1, nb_bdy |
---|
[8586] | 339 | ! |
---|
[11536] | 340 | SELECT CASE( cn_ice(jbdy) ) |
---|
[8586] | 341 | ! |
---|
[11536] | 342 | CASE('none') |
---|
| 343 | CYCLE |
---|
| 344 | ! |
---|
| 345 | CASE('frs') |
---|
| 346 | ! |
---|
| 347 | IF( nn_ice_dta(jbdy) == 0 ) CYCLE ! case ice boundaries = initial conditions |
---|
| 348 | ! ! do not change ice velocity (it is only computed by rheology) |
---|
| 349 | SELECT CASE ( cd_type ) |
---|
| 350 | ! |
---|
| 351 | CASE ( 'U' ) |
---|
| 352 | jgrd = 2 ! u velocity |
---|
| 353 | IF( ir == 0 ) THEN ; ibeg = 1 ; iend = idx_bdy(jbdy)%nblenrim0(jgrd) |
---|
| 354 | ELSE ; ibeg = idx_bdy(jbdy)%nblenrim0(jgrd)+1 ; iend = idx_bdy(jbdy)%nblenrim(jgrd) |
---|
| 355 | END IF |
---|
| 356 | DO i_bdy = ibeg, iend |
---|
| 357 | ji = idx_bdy(jbdy)%nbi(i_bdy,jgrd) |
---|
| 358 | jj = idx_bdy(jbdy)%nbj(i_bdy,jgrd) |
---|
| 359 | zflag = idx_bdy(jbdy)%flagu(i_bdy,jgrd) |
---|
| 360 | ! i-1 i i | ! i i i+1 | ! i i i+1 | |
---|
| 361 | ! > ice > | ! o > ice | ! o > o | |
---|
| 362 | ! => set at u_ice(i-1) ! => set to O ! => unchanged |
---|
| 363 | IF( zflag == -1. .AND. ji > 1 .AND. ji < jpi ) THEN |
---|
| 364 | IF ( vt_i(ji ,jj) > 0. ) THEN ; u_ice(ji,jj) = u_ice(ji-1,jj) |
---|
| 365 | ELSEIF( vt_i(ji+1,jj) > 0. ) THEN ; u_ice(ji,jj) = 0._wp |
---|
| 366 | END IF |
---|
| 367 | END IF |
---|
| 368 | ! | i i+1 i+1 ! | i i i+1 ! | i i i+1 |
---|
| 369 | ! | > ice > ! | ice > o ! | o > o |
---|
| 370 | ! => set at u_ice(i+1) ! => set to O ! => unchanged |
---|
| 371 | IF( zflag == 1. .AND. ji+1 < jpi+1 ) THEN |
---|
| 372 | IF ( vt_i(ji+1,jj) > 0. ) THEN ; u_ice(ji,jj) = u_ice(ji+1,jj) |
---|
| 373 | ELSEIF( vt_i(ji ,jj) > 0. ) THEN ; u_ice(ji,jj) = 0._wp |
---|
| 374 | END IF |
---|
| 375 | END IF |
---|
| 376 | ! |
---|
| 377 | IF( zflag == 0. ) u_ice(ji,jj) = 0._wp ! u_ice = 0 if north/south bdy |
---|
| 378 | ! |
---|
| 379 | END DO |
---|
[8586] | 380 | ! |
---|
[11536] | 381 | CASE ( 'V' ) |
---|
| 382 | jgrd = 3 ! v velocity |
---|
| 383 | IF( ir == 0 ) THEN ; ibeg = 1 ; iend = idx_bdy(jbdy)%nblenrim0(jgrd) |
---|
| 384 | ELSE ; ibeg = idx_bdy(jbdy)%nblenrim0(jgrd)+1 ; iend = idx_bdy(jbdy)%nblenrim(jgrd) |
---|
| 385 | END IF |
---|
| 386 | DO i_bdy = ibeg, iend |
---|
| 387 | ji = idx_bdy(jbdy)%nbi(i_bdy,jgrd) |
---|
| 388 | jj = idx_bdy(jbdy)%nbj(i_bdy,jgrd) |
---|
| 389 | zflag = idx_bdy(jbdy)%flagv(i_bdy,jgrd) |
---|
| 390 | ! ! ice (jj+1) ! o (jj+1) |
---|
| 391 | ! ^ (jj ) ! ^ (jj ) ! ^ (jj ) |
---|
| 392 | ! ice (jj ) ! o (jj ) ! o (jj ) |
---|
| 393 | ! ^ (jj-1) ! ! |
---|
| 394 | ! => set to u_ice(jj-1) ! => set to 0 ! => unchanged |
---|
| 395 | IF( zflag == -1. .AND. jj > 1 .AND. jj < jpj ) THEN |
---|
| 396 | IF ( vt_i(ji,jj ) > 0. ) THEN ; v_ice(ji,jj) = v_ice(ji,jj-1) |
---|
| 397 | ELSEIF( vt_i(ji,jj+1) > 0. ) THEN ; v_ice(ji,jj) = 0._wp |
---|
| 398 | END IF |
---|
| 399 | END IF |
---|
| 400 | ! ^ (jj+1) ! ! |
---|
| 401 | ! ice (jj+1) ! o (jj+1) ! o (jj+1) |
---|
| 402 | ! ^ (jj ) ! ^ (jj ) ! ^ (jj ) |
---|
| 403 | ! ________________ ! ____ice___(jj )_ ! _____o____(jj ) |
---|
| 404 | ! => set to u_ice(jj+1) ! => set to 0 ! => unchanged |
---|
| 405 | IF( zflag == 1. .AND. jj < jpj ) THEN |
---|
| 406 | IF ( vt_i(ji,jj+1) > 0. ) THEN ; v_ice(ji,jj) = v_ice(ji,jj+1) |
---|
| 407 | ELSEIF( vt_i(ji,jj ) > 0. ) THEN ; v_ice(ji,jj) = 0._wp |
---|
| 408 | END IF |
---|
| 409 | END IF |
---|
| 410 | ! |
---|
| 411 | IF( zflag == 0. ) v_ice(ji,jj) = 0._wp ! v_ice = 0 if west/east bdy |
---|
| 412 | ! |
---|
| 413 | END DO |
---|
[8586] | 414 | ! |
---|
[11536] | 415 | END SELECT |
---|
[8586] | 416 | ! |
---|
[11536] | 417 | CASE DEFAULT |
---|
| 418 | CALL ctl_stop( 'bdy_ice_dyn : unrecognised option for open boundaries for ice fields' ) |
---|
[8586] | 419 | END SELECT |
---|
| 420 | ! |
---|
[11536] | 421 | END DO ! jbdy |
---|
| 422 | ! |
---|
| 423 | SELECT CASE ( cd_type ) |
---|
| 424 | CASE ( 'U' ) |
---|
| 425 | IF( nn_hls > 1 .AND. ir == 1 ) CYCLE ! at least 2 halos will be corrected -> no need to correct rim 1 before rim 0 |
---|
| 426 | IF( nn_hls == 1 ) THEN ; llsend2(:) = .false. ; llrecv2(:) = .false. ; END IF |
---|
| 427 | DO jbdy = 1, nb_bdy |
---|
| 428 | IF( cn_ice(jbdy) == 'frs' .AND. nn_ice_dta(jbdy) /= 0 ) THEN |
---|
| 429 | llsend2(:) = llsend2(:) .OR. lsend_bdyint(jbdy,2,:,ir) ! possibly every direction, U points |
---|
| 430 | llsend2(1) = llsend2(1) .OR. lsend_bdyext(jbdy,2,1,ir) ! neighbour might search point towards its west bdy |
---|
| 431 | llrecv2(:) = llrecv2(:) .OR. lrecv_bdyint(jbdy,2,:,ir) ! possibly every direction, U points |
---|
| 432 | llrecv2(2) = llrecv2(2) .OR. lrecv_bdyext(jbdy,2,2,ir) ! might search point towards east bdy |
---|
| 433 | END IF |
---|
| 434 | END DO |
---|
| 435 | IF( ANY(llsend2) .OR. ANY(llrecv2) ) THEN ! if need to send/recv in at least one direction |
---|
[13226] | 436 | CALL lbc_lnk( 'bdyice', u_ice, 'U', -1.0_wp, kfillmode=jpfillnothing ,lsend=llsend2, lrecv=llrecv2 ) |
---|
[11536] | 437 | END IF |
---|
| 438 | CASE ( 'V' ) |
---|
| 439 | IF( nn_hls > 1 .AND. ir == 1 ) CYCLE ! at least 2 halos will be corrected -> no need to correct rim 1 before rim 0 |
---|
| 440 | IF( nn_hls == 1 ) THEN ; llsend3(:) = .false. ; llrecv3(:) = .false. ; END IF |
---|
| 441 | DO jbdy = 1, nb_bdy |
---|
| 442 | IF( cn_ice(jbdy) == 'frs' .AND. nn_ice_dta(jbdy) /= 0 ) THEN |
---|
| 443 | llsend3(:) = llsend3(:) .OR. lsend_bdyint(jbdy,3,:,ir) ! possibly every direction, V points |
---|
| 444 | llsend3(3) = llsend3(3) .OR. lsend_bdyext(jbdy,3,3,ir) ! neighbour might search point towards its south bdy |
---|
| 445 | llrecv3(:) = llrecv3(:) .OR. lrecv_bdyint(jbdy,3,:,ir) ! possibly every direction, V points |
---|
| 446 | llrecv3(4) = llrecv3(4) .OR. lrecv_bdyext(jbdy,3,4,ir) ! might search point towards north bdy |
---|
| 447 | END IF |
---|
| 448 | END DO |
---|
| 449 | IF( ANY(llsend3) .OR. ANY(llrecv3) ) THEN ! if need to send/recv in at least one direction |
---|
[13226] | 450 | CALL lbc_lnk( 'bdyice', v_ice, 'V', -1.0_wp, kfillmode=jpfillnothing ,lsend=llsend3, lrecv=llrecv3 ) |
---|
[11536] | 451 | END IF |
---|
[8586] | 452 | END SELECT |
---|
[11536] | 453 | END DO ! ir |
---|
[8586] | 454 | ! |
---|
[9905] | 455 | IF( ln_timing ) CALL timing_stop('bdy_ice_dyn') |
---|
| 456 | ! |
---|
[8586] | 457 | END SUBROUTINE bdy_ice_dyn |
---|
| 458 | |
---|
| 459 | #else |
---|
| 460 | !!--------------------------------------------------------------------------------- |
---|
| 461 | !! Default option Empty module |
---|
| 462 | !!--------------------------------------------------------------------------------- |
---|
| 463 | CONTAINS |
---|
| 464 | SUBROUTINE bdy_ice( kt ) ! Empty routine |
---|
[9927] | 465 | IMPLICIT NONE |
---|
| 466 | INTEGER, INTENT( in ) :: kt |
---|
[8586] | 467 | WRITE(*,*) 'bdy_ice: You should not have seen this print! error?', kt |
---|
| 468 | END SUBROUTINE bdy_ice |
---|
| 469 | #endif |
---|
| 470 | |
---|
| 471 | !!================================================================================= |
---|
| 472 | END MODULE bdyice |
---|