1 | MODULE bdyice_lim |
---|
2 | !!====================================================================== |
---|
3 | !! *** MODULE bdyice_lim *** |
---|
4 | !! Unstructured Open Boundary Cond. : Open boundary conditions for sea-ice (LIM2 and LIM3) |
---|
5 | !!====================================================================== |
---|
6 | !! History : 3.3 ! 2010-09 (D. Storkey) Original code |
---|
7 | !! 3.4 ! 2011 (D. Storkey) rewrite in preparation for OBC-BDY merge |
---|
8 | !! - ! 2012-01 (C. Rousset) add lim3 and remove useless jk loop |
---|
9 | !!---------------------------------------------------------------------- |
---|
10 | #if defined key_bdy && ( defined key_lim2 || defined key_lim3 ) |
---|
11 | !!---------------------------------------------------------------------- |
---|
12 | !! 'key_bdy' and Unstructured Open Boundary Conditions |
---|
13 | !! 'key_lim2' LIM-2 sea ice model |
---|
14 | !! 'key_lim3' LIM-3 sea ice model |
---|
15 | !!---------------------------------------------------------------------- |
---|
16 | !! bdy_ice_lim : Application of open boundaries to ice |
---|
17 | !! bdy_ice_frs : Application of Flow Relaxation Scheme |
---|
18 | !!---------------------------------------------------------------------- |
---|
19 | USE timing ! Timing |
---|
20 | USE phycst ! physical constant |
---|
21 | USE eosbn2 ! equation of state |
---|
22 | USE oce ! ocean dynamics and tracers variables |
---|
23 | #if defined key_lim2 |
---|
24 | USE par_ice_2 |
---|
25 | USE ice_2 ! LIM_2 ice variables |
---|
26 | USE dom_ice_2 ! sea-ice domain |
---|
27 | #elif defined key_lim3 |
---|
28 | USE ice ! LIM_3 ice variables |
---|
29 | USE limvar |
---|
30 | USE limctl |
---|
31 | #endif |
---|
32 | USE par_oce ! ocean parameters |
---|
33 | USE dom_oce ! ocean space and time domain variables |
---|
34 | USE sbc_oce ! Surface boundary condition: ocean fields |
---|
35 | USE bdy_oce ! ocean open boundary conditions |
---|
36 | USE lbclnk ! ocean lateral boundary conditions (or mpp link) |
---|
37 | USE in_out_manager ! write to numout file |
---|
38 | USE lib_mpp ! distributed memory computing |
---|
39 | USE lib_fortran ! to use key_nosignedzero |
---|
40 | |
---|
41 | IMPLICIT NONE |
---|
42 | PRIVATE |
---|
43 | |
---|
44 | PUBLIC bdy_ice_lim ! routine called in sbcmod |
---|
45 | PUBLIC bdy_ice_lim_dyn ! routine called in limrhg |
---|
46 | |
---|
47 | !!---------------------------------------------------------------------- |
---|
48 | !! NEMO/OPA 3.3 , NEMO Consortium (2010) |
---|
49 | !! $Id$ |
---|
50 | !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) |
---|
51 | !!---------------------------------------------------------------------- |
---|
52 | CONTAINS |
---|
53 | |
---|
54 | SUBROUTINE bdy_ice_lim( kt ) |
---|
55 | !!---------------------------------------------------------------------- |
---|
56 | !! *** SUBROUTINE bdy_ice_lim *** |
---|
57 | !! |
---|
58 | !! ** Purpose : - Apply open boundary conditions for ice (LIM2 and LIM3) |
---|
59 | !! |
---|
60 | !!---------------------------------------------------------------------- |
---|
61 | INTEGER, INTENT( in ) :: kt ! Main time step counter |
---|
62 | INTEGER :: ib_bdy ! Loop index |
---|
63 | |
---|
64 | #if defined key_lim3 |
---|
65 | CALL lim_var_glo2eqv |
---|
66 | #endif |
---|
67 | |
---|
68 | DO ib_bdy=1, nb_bdy |
---|
69 | |
---|
70 | SELECT CASE( cn_ice_lim(ib_bdy) ) |
---|
71 | CASE('none') |
---|
72 | CYCLE |
---|
73 | CASE('frs') |
---|
74 | CALL bdy_ice_frs( idx_bdy(ib_bdy), dta_bdy(ib_bdy), kt, ib_bdy ) |
---|
75 | CASE DEFAULT |
---|
76 | CALL ctl_stop( 'bdy_ice_lim : unrecognised option for open boundaries for ice fields' ) |
---|
77 | END SELECT |
---|
78 | |
---|
79 | END DO |
---|
80 | |
---|
81 | #if defined key_lim3 |
---|
82 | CALL lim_var_zapsmall |
---|
83 | CALL lim_var_agg(1) |
---|
84 | IF( ln_limctl ) CALL lim_prt( kt, iiceprt, jiceprt, 1, ' - ice thermo bdy - ' ) |
---|
85 | #endif |
---|
86 | |
---|
87 | END SUBROUTINE bdy_ice_lim |
---|
88 | |
---|
89 | SUBROUTINE bdy_ice_frs( idx, dta, kt, ib_bdy ) |
---|
90 | !!------------------------------------------------------------------------------ |
---|
91 | !! *** SUBROUTINE bdy_ice_frs *** |
---|
92 | !! |
---|
93 | !! ** Purpose : Apply the Flow Relaxation Scheme for sea-ice fields in the case |
---|
94 | !! of unstructured open boundaries. |
---|
95 | !! |
---|
96 | !! Reference : Engedahl H., 1995: Use of the flow relaxation scheme in a three- |
---|
97 | !! dimensional baroclinic ocean model with realistic topography. Tellus, 365-382. |
---|
98 | !!------------------------------------------------------------------------------ |
---|
99 | TYPE(OBC_INDEX), INTENT(in) :: idx ! OBC indices |
---|
100 | TYPE(OBC_DATA), INTENT(in) :: dta ! OBC external data |
---|
101 | INTEGER, INTENT(in) :: kt ! main time-step counter |
---|
102 | INTEGER, INTENT(in) :: ib_bdy ! BDY set index |
---|
103 | |
---|
104 | INTEGER :: jpbound ! 0 = incoming ice |
---|
105 | ! 1 = outgoing ice |
---|
106 | INTEGER :: jb, jk, jgrd, jl ! dummy loop indices |
---|
107 | INTEGER :: ji, jj, ii, ij ! local scalar |
---|
108 | REAL(wp) :: zwgt, zwgt1 ! local scalar |
---|
109 | REAL(wp) :: ztmelts, zdh |
---|
110 | #if defined key_lim2 && ! defined key_lim2_vp && defined key_agrif |
---|
111 | USE ice_2, vt_s => hsnm |
---|
112 | USE ice_2, vt_i => hicm |
---|
113 | #endif |
---|
114 | |
---|
115 | !!------------------------------------------------------------------------------ |
---|
116 | ! |
---|
117 | IF( nn_timing == 1 ) CALL timing_start('bdy_ice_frs') |
---|
118 | ! |
---|
119 | jgrd = 1 ! Everything is at T-points here |
---|
120 | ! |
---|
121 | #if defined key_lim2 |
---|
122 | DO jb = 1, idx%nblenrim(jgrd) |
---|
123 | ji = idx%nbi(jb,jgrd) |
---|
124 | jj = idx%nbj(jb,jgrd) |
---|
125 | zwgt = idx%nbw(jb,jgrd) |
---|
126 | zwgt1 = 1.e0 - idx%nbw(jb,jgrd) |
---|
127 | frld (ji,jj) = ( frld (ji,jj) * zwgt1 + dta%frld (jb) * zwgt ) * tmask(ji,jj,1) ! Leads fraction |
---|
128 | hicif(ji,jj) = ( hicif(ji,jj) * zwgt1 + dta%hicif(jb) * zwgt ) * tmask(ji,jj,1) ! Ice depth |
---|
129 | hsnif(ji,jj) = ( hsnif(ji,jj) * zwgt1 + dta%hsnif(jb) * zwgt ) * tmask(ji,jj,1) ! Snow depth |
---|
130 | END DO |
---|
131 | |
---|
132 | CALL lbc_bdy_lnk( frld, 'T', 1., ib_bdy ) ! lateral boundary conditions |
---|
133 | CALL lbc_bdy_lnk( hicif, 'T', 1., ib_bdy ) |
---|
134 | CALL lbc_bdy_lnk( hsnif, 'T', 1., ib_bdy ) |
---|
135 | |
---|
136 | vt_i(:,:) = hicif(:,:) * frld(:,:) |
---|
137 | vt_s(:,:) = hsnif(:,:) * frld(:,:) |
---|
138 | ! |
---|
139 | #elif defined key_lim3 |
---|
140 | |
---|
141 | DO jl = 1, jpl |
---|
142 | DO jb = 1, idx%nblenrim(jgrd) |
---|
143 | ji = idx%nbi(jb,jgrd) |
---|
144 | jj = idx%nbj(jb,jgrd) |
---|
145 | zwgt = idx%nbw(jb,jgrd) |
---|
146 | zwgt1 = 1.e0 - idx%nbw(jb,jgrd) |
---|
147 | a_i (ji,jj,jl) = ( a_i (ji,jj,jl) * zwgt1 + dta%a_i (jb,jl) * zwgt ) * tmask(ji,jj,1) ! Leads fraction |
---|
148 | ht_i(ji,jj,jl) = ( ht_i(ji,jj,jl) * zwgt1 + dta%ht_i(jb,jl) * zwgt ) * tmask(ji,jj,1) ! Ice depth |
---|
149 | ht_s(ji,jj,jl) = ( ht_s(ji,jj,jl) * zwgt1 + dta%ht_s(jb,jl) * zwgt ) * tmask(ji,jj,1) ! Snow depth |
---|
150 | |
---|
151 | ! ----------------- |
---|
152 | ! Pathological case |
---|
153 | ! ----------------- |
---|
154 | ! In case a) snow load would be in excess or b) ice is coming into a warmer environment that would lead to |
---|
155 | ! very large transformation from snow to ice (see limthd_dh.F90) |
---|
156 | |
---|
157 | ! Then, a) transfer the snow excess into the ice (different from limthd_dh) |
---|
158 | zdh = MAX( 0._wp, ( rhosn * ht_s(ji,jj,jl) + ( rhoic - rau0 ) * ht_i(ji,jj,jl) ) * r1_rau0 ) |
---|
159 | ! Or, b) transfer all the snow into ice (if incoming ice is likely to melt as it comes into a warmer environment) |
---|
160 | !zdh = MAX( 0._wp, ht_s(ji,jj,jl) * rhosn / rhoic ) |
---|
161 | |
---|
162 | ! recompute ht_i, ht_s |
---|
163 | ht_i(ji,jj,jl) = MIN( hi_max(jl), ht_i(ji,jj,jl) + zdh ) |
---|
164 | ht_s(ji,jj,jl) = MAX( 0._wp, ht_s(ji,jj,jl) - zdh * rhoic / rhosn ) |
---|
165 | |
---|
166 | ENDDO |
---|
167 | CALL lbc_bdy_lnk( a_i(:,:,jl), 'T', 1., ib_bdy ) |
---|
168 | CALL lbc_bdy_lnk( ht_i(:,:,jl), 'T', 1., ib_bdy ) |
---|
169 | CALL lbc_bdy_lnk( ht_s(:,:,jl), 'T', 1., ib_bdy ) |
---|
170 | ENDDO |
---|
171 | ! retrieve at_i |
---|
172 | at_i(:,:) = 0._wp |
---|
173 | DO jl = 1, jpl |
---|
174 | at_i(:,:) = a_i(:,:,jl) + at_i(:,:) |
---|
175 | END DO |
---|
176 | |
---|
177 | DO jl = 1, jpl |
---|
178 | DO jb = 1, idx%nblenrim(jgrd) |
---|
179 | ji = idx%nbi(jb,jgrd) |
---|
180 | jj = idx%nbj(jb,jgrd) |
---|
181 | |
---|
182 | ! condition on ice thickness depends on the ice velocity |
---|
183 | ! if velocity is outward (strictly), then ice thickness, volume... must be equal to adjacent values |
---|
184 | jpbound = 0; ii = ji; ij = jj; |
---|
185 | |
---|
186 | IF( u_ice(ji+1,jj ) < 0. .AND. umask(ji-1,jj ,1) == 0. ) jpbound = 1; ii = ji+1; ij = jj |
---|
187 | IF( u_ice(ji-1,jj ) > 0. .AND. umask(ji+1,jj ,1) == 0. ) jpbound = 1; ii = ji-1; ij = jj |
---|
188 | IF( v_ice(ji ,jj+1) < 0. .AND. vmask(ji ,jj-1,1) == 0. ) jpbound = 1; ii = ji ; ij = jj+1 |
---|
189 | IF( v_ice(ji ,jj-1) > 0. .AND. vmask(ji ,jj+1,1) == 0. ) jpbound = 1; ii = ji ; ij = jj-1 |
---|
190 | |
---|
191 | IF( nn_ice_lim_dta(ib_bdy) == 0 ) jpbound = 0; ii = ji; ij = jj ! case ice boundaries = initial conditions |
---|
192 | ! do not make state variables dependent on velocity |
---|
193 | |
---|
194 | |
---|
195 | rswitch = MAX( 0.0_wp , SIGN ( 1.0_wp , at_i(ii,ij) - 0.01 ) ) ! 0 if no ice |
---|
196 | |
---|
197 | ! concentration and thickness |
---|
198 | a_i (ji,jj,jl) = a_i (ii,ij,jl) * rswitch |
---|
199 | ht_i(ji,jj,jl) = ht_i(ii,ij,jl) * rswitch |
---|
200 | ht_s(ji,jj,jl) = ht_s(ii,ij,jl) * rswitch |
---|
201 | |
---|
202 | ! Ice and snow volumes |
---|
203 | v_i(ji,jj,jl) = ht_i(ji,jj,jl) * a_i(ji,jj,jl) |
---|
204 | v_s(ji,jj,jl) = ht_s(ji,jj,jl) * a_i(ji,jj,jl) |
---|
205 | |
---|
206 | SELECT CASE( jpbound ) |
---|
207 | |
---|
208 | CASE( 0 ) ! velocity is inward |
---|
209 | |
---|
210 | ! Ice salinity, age, temperature |
---|
211 | sm_i(ji,jj,jl) = rswitch * rn_ice_sal(ib_bdy) + ( 1.0 - rswitch ) * rn_simin |
---|
212 | oa_i(ji,jj,jl) = rswitch * rn_ice_age(ib_bdy) * a_i(ji,jj,jl) |
---|
213 | t_su(ji,jj,jl) = rswitch * rn_ice_tem(ib_bdy) + ( 1.0 - rswitch ) * rn_ice_tem(ib_bdy) |
---|
214 | DO jk = 1, nlay_s |
---|
215 | t_s(ji,jj,jk,jl) = rswitch * rn_ice_tem(ib_bdy) + ( 1.0 - rswitch ) * rt0 |
---|
216 | END DO |
---|
217 | DO jk = 1, nlay_i |
---|
218 | t_i(ji,jj,jk,jl) = rswitch * rn_ice_tem(ib_bdy) + ( 1.0 - rswitch ) * rt0 |
---|
219 | s_i(ji,jj,jk,jl) = rswitch * rn_ice_sal(ib_bdy) + ( 1.0 - rswitch ) * rn_simin |
---|
220 | END DO |
---|
221 | |
---|
222 | CASE( 1 ) ! velocity is outward |
---|
223 | |
---|
224 | ! Ice salinity, age, temperature |
---|
225 | sm_i(ji,jj,jl) = rswitch * sm_i(ii,ij,jl) + ( 1.0 - rswitch ) * rn_simin |
---|
226 | oa_i(ji,jj,jl) = rswitch * oa_i(ii,ij,jl) |
---|
227 | t_su(ji,jj,jl) = rswitch * t_su(ii,ij,jl) + ( 1.0 - rswitch ) * rt0 |
---|
228 | DO jk = 1, nlay_s |
---|
229 | t_s(ji,jj,jk,jl) = rswitch * t_s(ii,ij,jk,jl) + ( 1.0 - rswitch ) * rt0 |
---|
230 | END DO |
---|
231 | DO jk = 1, nlay_i |
---|
232 | t_i(ji,jj,jk,jl) = rswitch * t_i(ii,ij,jk,jl) + ( 1.0 - rswitch ) * rt0 |
---|
233 | s_i(ji,jj,jk,jl) = rswitch * s_i(ii,ij,jk,jl) + ( 1.0 - rswitch ) * rn_simin |
---|
234 | END DO |
---|
235 | |
---|
236 | END SELECT |
---|
237 | |
---|
238 | ! if salinity is constant, then overwrite rn_ice_sal |
---|
239 | IF( nn_icesal == 1 ) THEN |
---|
240 | sm_i(ji,jj,jl) = rn_icesal |
---|
241 | s_i (ji,jj,:,jl) = rn_icesal |
---|
242 | ENDIF |
---|
243 | |
---|
244 | ! contents |
---|
245 | smv_i(ji,jj,jl) = MIN( sm_i(ji,jj,jl) , sss_m(ji,jj) ) * v_i(ji,jj,jl) |
---|
246 | DO jk = 1, nlay_s |
---|
247 | ! Snow energy of melting |
---|
248 | e_s(ji,jj,jk,jl) = rswitch * rhosn * ( cpic * ( rt0 - t_s(ji,jj,jk,jl) ) + lfus ) |
---|
249 | ! Multiply by volume, so that heat content in J/m2 |
---|
250 | e_s(ji,jj,jk,jl) = e_s(ji,jj,jk,jl) * v_s(ji,jj,jl) * r1_nlay_s |
---|
251 | END DO |
---|
252 | DO jk = 1, nlay_i |
---|
253 | ztmelts = - tmut * s_i(ji,jj,jk,jl) + rt0 !Melting temperature in K |
---|
254 | ! heat content per unit volume |
---|
255 | e_i(ji,jj,jk,jl) = rswitch * rhoic * & |
---|
256 | ( cpic * ( ztmelts - t_i(ji,jj,jk,jl) ) & |
---|
257 | + lfus * ( 1.0 - (ztmelts-rt0) / MIN((t_i(ji,jj,jk,jl)-rt0),-epsi20) ) & |
---|
258 | - rcp * ( ztmelts - rt0 ) ) |
---|
259 | ! Mutliply by ice volume, and divide by number of layers to get heat content in J/m2 |
---|
260 | e_i(ji,jj,jk,jl) = e_i(ji,jj,jk,jl) * a_i(ji,jj,jl) * ht_i(ji,jj,jl) * r1_nlay_i |
---|
261 | END DO |
---|
262 | |
---|
263 | END DO |
---|
264 | |
---|
265 | CALL lbc_bdy_lnk( a_i(:,:,jl), 'T', 1., ib_bdy ) |
---|
266 | CALL lbc_bdy_lnk( ht_i(:,:,jl), 'T', 1., ib_bdy ) |
---|
267 | CALL lbc_bdy_lnk( ht_s(:,:,jl), 'T', 1., ib_bdy ) |
---|
268 | CALL lbc_bdy_lnk( v_i(:,:,jl), 'T', 1., ib_bdy ) |
---|
269 | CALL lbc_bdy_lnk( v_s(:,:,jl), 'T', 1., ib_bdy ) |
---|
270 | |
---|
271 | CALL lbc_bdy_lnk( smv_i(:,:,jl), 'T', 1., ib_bdy ) |
---|
272 | CALL lbc_bdy_lnk( sm_i(:,:,jl), 'T', 1., ib_bdy ) |
---|
273 | CALL lbc_bdy_lnk( oa_i(:,:,jl), 'T', 1., ib_bdy ) |
---|
274 | CALL lbc_bdy_lnk( t_su(:,:,jl), 'T', 1., ib_bdy ) |
---|
275 | DO jk = 1, nlay_s |
---|
276 | CALL lbc_bdy_lnk(t_s(:,:,jk,jl), 'T', 1., ib_bdy ) |
---|
277 | CALL lbc_bdy_lnk(e_s(:,:,jk,jl), 'T', 1., ib_bdy ) |
---|
278 | END DO |
---|
279 | DO jk = 1, nlay_i |
---|
280 | CALL lbc_bdy_lnk(t_i(:,:,jk,jl), 'T', 1., ib_bdy ) |
---|
281 | CALL lbc_bdy_lnk(e_i(:,:,jk,jl), 'T', 1., ib_bdy ) |
---|
282 | END DO |
---|
283 | |
---|
284 | END DO !jl |
---|
285 | |
---|
286 | #endif |
---|
287 | ! |
---|
288 | IF( nn_timing == 1 ) CALL timing_stop('bdy_ice_frs') |
---|
289 | ! |
---|
290 | END SUBROUTINE bdy_ice_frs |
---|
291 | |
---|
292 | |
---|
293 | SUBROUTINE bdy_ice_lim_dyn( cd_type ) |
---|
294 | !!------------------------------------------------------------------------------ |
---|
295 | !! *** SUBROUTINE bdy_ice_lim_dyn *** |
---|
296 | !! |
---|
297 | !! ** Purpose : Apply dynamics boundary conditions for sea-ice in the cas of unstructured open boundaries. |
---|
298 | !! u_ice and v_ice are equal to the value of the adjacent grid point if this latter is not ice free |
---|
299 | !! if adjacent grid point is ice free, then u_ice and v_ice are equal to ocean velocities |
---|
300 | !! |
---|
301 | !! 2013-06 : C. Rousset |
---|
302 | !!------------------------------------------------------------------------------ |
---|
303 | !! |
---|
304 | CHARACTER(len=1), INTENT(in) :: cd_type ! nature of velocity grid-points |
---|
305 | INTEGER :: jb, jgrd ! dummy loop indices |
---|
306 | INTEGER :: ji, jj ! local scalar |
---|
307 | INTEGER :: ib_bdy ! Loop index |
---|
308 | REAL(wp) :: zmsk1, zmsk2, zflag |
---|
309 | !!------------------------------------------------------------------------------ |
---|
310 | ! |
---|
311 | IF( nn_timing == 1 ) CALL timing_start('bdy_ice_lim_dyn') |
---|
312 | ! |
---|
313 | DO ib_bdy=1, nb_bdy |
---|
314 | ! |
---|
315 | SELECT CASE( cn_ice_lim(ib_bdy) ) |
---|
316 | |
---|
317 | CASE('none') |
---|
318 | |
---|
319 | CYCLE |
---|
320 | |
---|
321 | CASE('frs') |
---|
322 | |
---|
323 | IF( nn_ice_lim_dta(ib_bdy) == 0 ) CYCLE ! case ice boundaries = initial conditions |
---|
324 | ! do not change ice velocity (it is only computed by rheology) |
---|
325 | |
---|
326 | SELECT CASE ( cd_type ) |
---|
327 | |
---|
328 | CASE ( 'U' ) |
---|
329 | |
---|
330 | jgrd = 2 ! u velocity |
---|
331 | DO jb = 1, idx_bdy(ib_bdy)%nblenrim(jgrd) |
---|
332 | ji = idx_bdy(ib_bdy)%nbi(jb,jgrd) |
---|
333 | jj = idx_bdy(ib_bdy)%nbj(jb,jgrd) |
---|
334 | zflag = idx_bdy(ib_bdy)%flagu(jb,jgrd) |
---|
335 | |
---|
336 | IF ( ABS( zflag ) == 1. ) THEN ! eastern and western boundaries |
---|
337 | ! one of the two zmsk is always 0 (because of zflag) |
---|
338 | zmsk1 = 1._wp - MAX( 0.0_wp, SIGN ( 1.0_wp , - vt_i(ji+1,jj) ) ) ! 0 if no ice |
---|
339 | zmsk2 = 1._wp - MAX( 0.0_wp, SIGN ( 1.0_wp , - vt_i(ji-1,jj) ) ) ! 0 if no ice |
---|
340 | |
---|
341 | ! u_ice = u_ice of the adjacent grid point except if this grid point is ice-free (then u_ice = u_oce) |
---|
342 | u_ice (ji,jj) = u_ice(ji+1,jj) * 0.5_wp * ABS( zflag + 1._wp ) * zmsk1 + & |
---|
343 | & u_ice(ji-1,jj) * 0.5_wp * ABS( zflag - 1._wp ) * zmsk2 + & |
---|
344 | & u_oce(ji ,jj) * ( 1._wp - MIN( 1._wp, zmsk1 + zmsk2 ) ) |
---|
345 | ELSE ! everywhere else |
---|
346 | !u_ice(ji,jj) = u_oce(ji,jj) |
---|
347 | u_ice(ji,jj) = 0._wp |
---|
348 | ENDIF |
---|
349 | ! mask ice velocities |
---|
350 | rswitch = MAX( 0.0_wp , SIGN ( 1.0_wp , at_i(ji,jj) - 0.01_wp ) ) ! 0 if no ice |
---|
351 | u_ice(ji,jj) = rswitch * u_ice(ji,jj) |
---|
352 | |
---|
353 | ENDDO |
---|
354 | |
---|
355 | CALL lbc_bdy_lnk( u_ice(:,:), 'U', -1., ib_bdy ) |
---|
356 | |
---|
357 | CASE ( 'V' ) |
---|
358 | |
---|
359 | jgrd = 3 ! v velocity |
---|
360 | DO jb = 1, idx_bdy(ib_bdy)%nblenrim(jgrd) |
---|
361 | ji = idx_bdy(ib_bdy)%nbi(jb,jgrd) |
---|
362 | jj = idx_bdy(ib_bdy)%nbj(jb,jgrd) |
---|
363 | zflag = idx_bdy(ib_bdy)%flagv(jb,jgrd) |
---|
364 | |
---|
365 | IF ( ABS( zflag ) == 1. ) THEN ! northern and southern boundaries |
---|
366 | ! one of the two zmsk is always 0 (because of zflag) |
---|
367 | zmsk1 = 1._wp - MAX( 0.0_wp, SIGN ( 1.0_wp , - vt_i(ji,jj+1) ) ) ! 0 if no ice |
---|
368 | zmsk2 = 1._wp - MAX( 0.0_wp, SIGN ( 1.0_wp , - vt_i(ji,jj-1) ) ) ! 0 if no ice |
---|
369 | |
---|
370 | ! u_ice = u_ice of the adjacent grid point except if this grid point is ice-free (then u_ice = u_oce) |
---|
371 | v_ice (ji,jj) = v_ice(ji,jj+1) * 0.5_wp * ABS( zflag + 1._wp ) * zmsk1 + & |
---|
372 | & v_ice(ji,jj-1) * 0.5_wp * ABS( zflag - 1._wp ) * zmsk2 + & |
---|
373 | & v_oce(ji,jj ) * ( 1._wp - MIN( 1._wp, zmsk1 + zmsk2 ) ) |
---|
374 | ELSE ! everywhere else |
---|
375 | !v_ice(ji,jj) = v_oce(ji,jj) |
---|
376 | v_ice(ji,jj) = 0._wp |
---|
377 | ENDIF |
---|
378 | ! mask ice velocities |
---|
379 | rswitch = MAX( 0.0_wp , SIGN ( 1.0_wp , at_i(ji,jj) - 0.01 ) ) ! 0 if no ice |
---|
380 | v_ice(ji,jj) = rswitch * v_ice(ji,jj) |
---|
381 | |
---|
382 | ENDDO |
---|
383 | |
---|
384 | CALL lbc_bdy_lnk( v_ice(:,:), 'V', -1., ib_bdy ) |
---|
385 | |
---|
386 | END SELECT |
---|
387 | |
---|
388 | CASE DEFAULT |
---|
389 | CALL ctl_stop( 'bdy_ice_lim_dyn : unrecognised option for open boundaries for ice fields' ) |
---|
390 | END SELECT |
---|
391 | |
---|
392 | ENDDO |
---|
393 | |
---|
394 | IF( nn_timing == 1 ) CALL timing_stop('bdy_ice_lim_dyn') |
---|
395 | |
---|
396 | END SUBROUTINE bdy_ice_lim_dyn |
---|
397 | |
---|
398 | #else |
---|
399 | !!--------------------------------------------------------------------------------- |
---|
400 | !! Default option Empty module |
---|
401 | !!--------------------------------------------------------------------------------- |
---|
402 | CONTAINS |
---|
403 | SUBROUTINE bdy_ice_lim( kt ) ! Empty routine |
---|
404 | WRITE(*,*) 'bdy_ice_lim: You should not have seen this print! error?', kt |
---|
405 | END SUBROUTINE bdy_ice_lim |
---|
406 | #endif |
---|
407 | |
---|
408 | !!================================================================================= |
---|
409 | END MODULE bdyice_lim |
---|