1 | MODULE bdyice |
---|
2 | !!====================================================================== |
---|
3 | !! *** MODULE bdyice *** |
---|
4 | !! Unstructured Open Boundary Cond. : Open boundary conditions for sea-ice (SI3) |
---|
5 | !!====================================================================== |
---|
6 | !! History : 3.3 ! 2010-09 (D. Storkey) Original code |
---|
7 | !! 3.4 ! 2012-01 (C. Rousset) add new sea ice model |
---|
8 | !! 4.0 ! 2018 (C. Rousset) SI3 compatibility |
---|
9 | !!---------------------------------------------------------------------- |
---|
10 | #if defined key_si3 |
---|
11 | !!---------------------------------------------------------------------- |
---|
12 | !! 'key_si3' SI3 sea ice model |
---|
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 |
---|
18 | USE ice ! sea-ice: variables |
---|
19 | USE icevar ! sea-ice: operations |
---|
20 | USE icecor ! sea-ice: corrections |
---|
21 | USE icectl ! sea-ice: control prints |
---|
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 |
---|
39 | PUBLIC bdy_ice_dyn ! routine called in icedyn_rhg_evp |
---|
40 | |
---|
41 | !!---------------------------------------------------------------------- |
---|
42 | !! NEMO/OCE 4.0 , NEMO Consortium (2018) |
---|
43 | !! $Id$ |
---|
44 | !! Software governed by the CeCILL license (see ./LICENSE) |
---|
45 | !!---------------------------------------------------------------------- |
---|
46 | CONTAINS |
---|
47 | |
---|
48 | SUBROUTINE bdy_ice( kt ) |
---|
49 | !!---------------------------------------------------------------------- |
---|
50 | !! *** SUBROUTINE bdy_ice *** |
---|
51 | !! |
---|
52 | !! ** Purpose : Apply open boundary conditions for sea ice |
---|
53 | !! |
---|
54 | !!---------------------------------------------------------------------- |
---|
55 | INTEGER, INTENT(in) :: kt ! Main time step counter |
---|
56 | ! |
---|
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(8) :: llsend1, llrecv1 ! indicate how communications are to be carried out |
---|
61 | !!---------------------------------------------------------------------- |
---|
62 | ! controls |
---|
63 | IF( ln_timing ) CALL timing_start('bdy_ice_thd') ! timing |
---|
64 | ! |
---|
65 | CALL ice_var_glo2eqv |
---|
66 | ! |
---|
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 |
---|
82 | ! |
---|
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 |
---|
94 | CALL lbc_lnk('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 ) |
---|
98 | ! exchange 4d arrays : third dimension = 1 and then third dimension = jpk |
---|
99 | CALL lbc_lnk('bdyice', t_s , 'T', 1._wp, e_s , 'T', 1._wp, kfillmode=jpfillnothing ,lsend=llsend1, lrecv=llrecv1 ) |
---|
100 | CALL lbc_lnk('bdyice', t_i , 'T', 1._wp, e_i , 'T', 1._wp, kfillmode=jpfillnothing ,lsend=llsend1, lrecv=llrecv1 ) |
---|
101 | END IF |
---|
102 | END DO ! ir |
---|
103 | ! |
---|
104 | CALL ice_cor( kt , 0 ) ! -- In case categories are out of bounds, do a remapping |
---|
105 | ! ! i.e. inputs have not the same ice thickness distribution (set by rn_himean) |
---|
106 | ! ! than the regional simulation |
---|
107 | CALL ice_var_agg(1) |
---|
108 | ! |
---|
109 | ! controls |
---|
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 |
---|
112 | ! |
---|
113 | END SUBROUTINE bdy_ice |
---|
114 | |
---|
115 | |
---|
116 | SUBROUTINE bdy_ice_frs( idx, dta, kt, jbdy, llrim0 ) |
---|
117 | !!------------------------------------------------------------------------------ |
---|
118 | !! *** SUBROUTINE bdy_ice_frs *** |
---|
119 | !! |
---|
120 | !! ** Purpose : Apply the Flow Relaxation Scheme for sea-ice fields |
---|
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 | !!------------------------------------------------------------------------------ |
---|
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 |
---|
130 | ! |
---|
131 | INTEGER :: jpbound ! 0 = incoming ice |
---|
132 | ! ! 1 = outgoing ice |
---|
133 | INTEGER :: ibeg, iend ! length of rim to be treated (rim 0 or rim 1) |
---|
134 | INTEGER :: i_bdy, jgrd ! dummy loop indices |
---|
135 | INTEGER :: ji, jj, jk, jl, ib, jb |
---|
136 | REAL(wp) :: zwgt, zwgt1 ! local scalar |
---|
137 | REAL(wp) :: ztmelts, zdh |
---|
138 | REAL(wp), POINTER :: flagu, flagv ! short cuts |
---|
139 | !!------------------------------------------------------------------------------ |
---|
140 | ! |
---|
141 | jgrd = 1 ! Everything is at T-points here |
---|
142 | IF( llrim0 ) THEN ; ibeg = 1 ; iend = idx%nblenrim0(jgrd) |
---|
143 | ELSE ; ibeg = idx%nblenrim0(jgrd)+1 ; iend = idx%nblenrim(jgrd) |
---|
144 | END IF |
---|
145 | ! |
---|
146 | DO jl = 1, jpl |
---|
147 | DO i_bdy = ibeg, iend |
---|
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) |
---|
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) = dta%t_i(i_bdy,jl) * tmask(ji,jj,1) ! Ice temperature |
---|
156 | t_s (ji,jj,:,jl) = dta%t_s(i_bdy,jl) * tmask(ji,jj,1) ! Snow temperature |
---|
157 | t_su(ji,jj, jl) = dta%tsu(i_bdy,jl) * tmask(ji,jj,1) ! Surf temperature |
---|
158 | s_i (ji,jj, jl) = dta%s_i(i_bdy,jl) * 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 |
---|
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 |
---|
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 |
---|
169 | h_il(ji,jj,jl) = 0._wp |
---|
170 | ENDIF |
---|
171 | |
---|
172 | IF( .NOT.ln_pnd_lids ) THEN |
---|
173 | h_il(ji,jj,jl) = 0._wp |
---|
174 | ENDIF |
---|
175 | ! |
---|
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) |
---|
183 | zdh = MAX( 0._wp, ( rhos * h_s(ji,jj,jl) + ( rhoi - rho0 ) * h_i(ji,jj,jl) ) * r1_rho0 ) |
---|
184 | ! Or, b) transfer all the snow into ice (if incoming ice is likely to melt as it comes into a warmer environment) |
---|
185 | !zdh = MAX( 0._wp, h_s(ji,jj,jl) * rhos / rhoi ) |
---|
186 | |
---|
187 | ! recompute h_i, h_s |
---|
188 | h_i(ji,jj,jl) = MIN( hi_max(jl), h_i(ji,jj,jl) + zdh ) |
---|
189 | h_s(ji,jj,jl) = MAX( 0._wp, h_s(ji,jj,jl) - zdh * rhoi / rhos ) |
---|
190 | ! |
---|
191 | ENDDO |
---|
192 | ENDDO |
---|
193 | |
---|
194 | DO jl = 1, jpl |
---|
195 | DO i_bdy = ibeg, iend |
---|
196 | ji = idx%nbi(i_bdy,jgrd) |
---|
197 | jj = idx%nbj(i_bdy,jgrd) |
---|
198 | flagu => idx%flagu(i_bdy,jgrd) |
---|
199 | flagv => idx%flagv(i_bdy,jgrd) |
---|
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 |
---|
202 | jpbound = 0 ; ib = ji ; jb = jj |
---|
203 | ! |
---|
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 |
---|
220 | ! |
---|
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 |
---|
223 | ! |
---|
224 | IF( a_i(ib,jb,jl) > 0._wp ) THEN ! there is ice at the boundary |
---|
225 | ! |
---|
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) |
---|
235 | h_il(ji,jj, jl) = h_il(ib,jb, jl) |
---|
236 | ! |
---|
237 | sz_i(ji,jj,:,jl) = sz_i(ib,jb,:,jl) |
---|
238 | ! |
---|
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 | ! |
---|
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 |
---|
250 | ! |
---|
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 |
---|
255 | DO jk = 1, nlay_s |
---|
256 | t_s(ji,jj,jk,jl) = MIN( t_s(ji,jj,jk,jl), -0.15_wp + rt0 ) ! Force t_s to be lower than -0.15deg (arbitrary) => likely conservation issue |
---|
257 | ! ! otherwise instant melting can occur |
---|
258 | e_s(ji,jj,jk,jl) = rhos * ( rcpi * ( rt0 - t_s(ji,jj,jk,jl) ) + rLfus ) ! enthalpy in J/m3 |
---|
259 | e_s(ji,jj,jk,jl) = e_s(ji,jj,jk,jl) * v_s(ji,jj,jl) * r1_nlay_s ! enthalpy in J/m2 |
---|
260 | END DO |
---|
261 | t_su(ji,jj,jl) = MIN( t_su(ji,jj,jl), -0.15_wp + rt0 ) ! Force t_su to be lower than -0.15deg (arbitrary) |
---|
262 | DO jk = 1, nlay_i |
---|
263 | ztmelts = - rTmlt * sz_i(ji,jj,jk,jl) ! Melting temperature in C |
---|
264 | t_i(ji,jj,jk,jl) = MIN( t_i(ji,jj,jk,jl), (ztmelts-0.15_wp) + rt0 ) ! Force t_i to be lower than melting point (-0.15) => likely conservation issue |
---|
265 | ! |
---|
266 | e_i(ji,jj,jk,jl) = rhoi * ( rcpi * ( ztmelts - ( t_i(ji,jj,jk,jl) - rt0 ) ) & ! enthalpy in J/m3 |
---|
267 | & + rLfus * ( 1._wp - ztmelts / ( t_i(ji,jj,jk,jl) - rt0 ) ) & |
---|
268 | & - rcp * ztmelts ) |
---|
269 | e_i(ji,jj,jk,jl) = e_i(ji,jj,jk,jl) * v_i(ji,jj,jl) * r1_nlay_i ! enthalpy in J/m2 |
---|
270 | END DO |
---|
271 | ! |
---|
272 | ! melt ponds |
---|
273 | v_ip(ji,jj,jl) = h_ip(ji,jj,jl) * a_ip(ji,jj,jl) |
---|
274 | v_il(ji,jj,jl) = h_il(ji,jj,jl) * a_ip(ji,jj,jl) |
---|
275 | ! |
---|
276 | ELSE ! no ice at the boundary |
---|
277 | ! |
---|
278 | a_i (ji,jj, jl) = 0._wp |
---|
279 | h_i (ji,jj, jl) = 0._wp |
---|
280 | h_s (ji,jj, jl) = 0._wp |
---|
281 | oa_i(ji,jj, jl) = 0._wp |
---|
282 | t_su(ji,jj, jl) = rt0 |
---|
283 | t_s (ji,jj,:,jl) = rt0 |
---|
284 | t_i (ji,jj,:,jl) = rt0 |
---|
285 | |
---|
286 | a_ip(ji,jj,jl) = 0._wp |
---|
287 | h_ip(ji,jj,jl) = 0._wp |
---|
288 | h_il(ji,jj,jl) = 0._wp |
---|
289 | |
---|
290 | IF( nn_icesal == 1 ) THEN ! if constant salinity |
---|
291 | s_i (ji,jj ,jl) = rn_icesal |
---|
292 | sz_i(ji,jj,:,jl) = rn_icesal |
---|
293 | ELSE ! if variable salinity |
---|
294 | s_i (ji,jj,jl) = rn_simin |
---|
295 | sz_i(ji,jj,:,jl) = rn_simin |
---|
296 | ENDIF |
---|
297 | ! |
---|
298 | ! global fields |
---|
299 | v_i (ji,jj, jl) = 0._wp |
---|
300 | v_s (ji,jj, jl) = 0._wp |
---|
301 | sv_i(ji,jj, jl) = 0._wp |
---|
302 | e_s (ji,jj,:,jl) = 0._wp |
---|
303 | e_i (ji,jj,:,jl) = 0._wp |
---|
304 | v_ip(ji,jj, jl) = 0._wp |
---|
305 | v_il(ji,jj, jl) = 0._wp |
---|
306 | |
---|
307 | ENDIF |
---|
308 | |
---|
309 | END DO |
---|
310 | ! |
---|
311 | END DO ! jl |
---|
312 | ! |
---|
313 | END SUBROUTINE bdy_ice_frs |
---|
314 | |
---|
315 | |
---|
316 | SUBROUTINE bdy_ice_dyn( cd_type ) |
---|
317 | !!------------------------------------------------------------------------------ |
---|
318 | !! *** SUBROUTINE bdy_ice_dyn *** |
---|
319 | !! |
---|
320 | !! ** Purpose : Apply dynamics boundary conditions for sea-ice. |
---|
321 | !! |
---|
322 | !! ** Method : if this adjacent grid point is not ice free, then u_ice and v_ice take its value |
---|
323 | !! if is ice free, then u_ice and v_ice are unchanged by BDY |
---|
324 | !! they keep values calculated in rheology |
---|
325 | !! |
---|
326 | !!------------------------------------------------------------------------------ |
---|
327 | CHARACTER(len=1), INTENT(in) :: cd_type ! nature of velocity grid-points |
---|
328 | ! |
---|
329 | INTEGER :: i_bdy, jgrd ! dummy loop indices |
---|
330 | INTEGER :: ji, jj ! local scalar |
---|
331 | INTEGER :: jbdy, ir ! BDY set index, rim index |
---|
332 | INTEGER :: ibeg, iend ! length of rim to be treated (rim 0 or rim 1) |
---|
333 | INTEGER, DIMENSION(3) :: idir3 |
---|
334 | REAL(wp) :: zmsk1, zmsk2, zflag |
---|
335 | LOGICAL, DIMENSION(8) :: llsend2, llrecv2, llsend3, llrecv3 ! indicate how communications are to be carried out |
---|
336 | !!------------------------------------------------------------------------------ |
---|
337 | IF( ln_timing ) CALL timing_start('bdy_ice_dyn') |
---|
338 | ! |
---|
339 | llsend2(:) = .false. ; llrecv2(:) = .false. |
---|
340 | llsend3(:) = .false. ; llrecv3(:) = .false. |
---|
341 | DO ir = 1, 0, -1 |
---|
342 | DO jbdy = 1, nb_bdy |
---|
343 | ! |
---|
344 | SELECT CASE( cn_ice(jbdy) ) |
---|
345 | ! |
---|
346 | CASE('none') |
---|
347 | CYCLE |
---|
348 | ! |
---|
349 | CASE('frs') |
---|
350 | ! |
---|
351 | IF( nn_ice_dta(jbdy) == 0 ) CYCLE ! case ice boundaries = initial conditions |
---|
352 | ! ! do not change ice velocity (it is only computed by rheology) |
---|
353 | SELECT CASE ( cd_type ) |
---|
354 | ! |
---|
355 | CASE ( 'U' ) |
---|
356 | jgrd = 2 ! u velocity |
---|
357 | IF( ir == 0 ) THEN ; ibeg = 1 ; iend = idx_bdy(jbdy)%nblenrim0(jgrd) |
---|
358 | ELSE ; ibeg = idx_bdy(jbdy)%nblenrim0(jgrd)+1 ; iend = idx_bdy(jbdy)%nblenrim(jgrd) |
---|
359 | END IF |
---|
360 | DO i_bdy = ibeg, iend |
---|
361 | ji = idx_bdy(jbdy)%nbi(i_bdy,jgrd) |
---|
362 | jj = idx_bdy(jbdy)%nbj(i_bdy,jgrd) |
---|
363 | zflag = idx_bdy(jbdy)%flagu(i_bdy,jgrd) |
---|
364 | ! i-1 i i | ! i i i+1 | ! i i i+1 | |
---|
365 | ! > ice > | ! o > ice | ! o > o | |
---|
366 | ! => set at u_ice(i-1) ! => set to O ! => unchanged |
---|
367 | IF( zflag == -1. .AND. ji > 1 .AND. ji < jpi ) THEN |
---|
368 | IF ( vt_i(ji ,jj) > 0. ) THEN ; u_ice(ji,jj) = u_ice(ji-1,jj) |
---|
369 | ELSEIF( vt_i(ji+1,jj) > 0. ) THEN ; u_ice(ji,jj) = u_oce(ji,jj) |
---|
370 | END IF |
---|
371 | END IF |
---|
372 | ! | i i+1 i+1 ! | i i i+1 ! | i i i+1 |
---|
373 | ! | > ice > ! | ice > o ! | o > o |
---|
374 | ! => set at u_ice(i+1) ! => set to O ! => unchanged |
---|
375 | IF( zflag == 1. .AND. ji+1 < jpi+1 ) THEN |
---|
376 | IF ( vt_i(ji+1,jj) > 0. ) THEN ; u_ice(ji,jj) = u_ice(ji+1,jj) |
---|
377 | ELSEIF( vt_i(ji ,jj) > 0. ) THEN ; u_ice(ji,jj) = u_oce(ji,jj) |
---|
378 | END IF |
---|
379 | END IF |
---|
380 | ! |
---|
381 | IF( zflag == 0. ) u_ice(ji,jj) = 0._wp ! u_ice = 0 if north/south bdy |
---|
382 | ! |
---|
383 | END DO |
---|
384 | ! |
---|
385 | CASE ( 'V' ) |
---|
386 | jgrd = 3 ! v velocity |
---|
387 | IF( ir == 0 ) THEN ; ibeg = 1 ; iend = idx_bdy(jbdy)%nblenrim0(jgrd) |
---|
388 | ELSE ; ibeg = idx_bdy(jbdy)%nblenrim0(jgrd)+1 ; iend = idx_bdy(jbdy)%nblenrim(jgrd) |
---|
389 | END IF |
---|
390 | DO i_bdy = ibeg, iend |
---|
391 | ji = idx_bdy(jbdy)%nbi(i_bdy,jgrd) |
---|
392 | jj = idx_bdy(jbdy)%nbj(i_bdy,jgrd) |
---|
393 | zflag = idx_bdy(jbdy)%flagv(i_bdy,jgrd) |
---|
394 | ! ! ice (jj+1) ! o (jj+1) |
---|
395 | ! ^ (jj ) ! ^ (jj ) ! ^ (jj ) |
---|
396 | ! ice (jj ) ! o (jj ) ! o (jj ) |
---|
397 | ! ^ (jj-1) ! ! |
---|
398 | ! => set to u_ice(jj-1) ! => set to 0 ! => unchanged |
---|
399 | IF( zflag == -1. .AND. jj > 1 .AND. jj < jpj ) THEN |
---|
400 | IF ( vt_i(ji,jj ) > 0. ) THEN ; v_ice(ji,jj) = v_ice(ji,jj-1) |
---|
401 | ELSEIF( vt_i(ji,jj+1) > 0. ) THEN ; v_ice(ji,jj) = v_oce(ji,jj) |
---|
402 | END IF |
---|
403 | END IF |
---|
404 | ! ^ (jj+1) ! ! |
---|
405 | ! ice (jj+1) ! o (jj+1) ! o (jj+1) |
---|
406 | ! ^ (jj ) ! ^ (jj ) ! ^ (jj ) |
---|
407 | ! ________________ ! ____ice___(jj )_ ! _____o____(jj ) |
---|
408 | ! => set to u_ice(jj+1) ! => set to 0 ! => unchanged |
---|
409 | IF( zflag == 1. .AND. jj < jpj ) THEN |
---|
410 | IF ( vt_i(ji,jj+1) > 0. ) THEN ; v_ice(ji,jj) = v_ice(ji,jj+1) |
---|
411 | ELSEIF( vt_i(ji,jj ) > 0. ) THEN ; v_ice(ji,jj) = v_oce(ji,jj) |
---|
412 | END IF |
---|
413 | END IF |
---|
414 | ! |
---|
415 | IF( zflag == 0. ) v_ice(ji,jj) = 0._wp ! v_ice = 0 if west/east bdy |
---|
416 | ! |
---|
417 | END DO |
---|
418 | ! |
---|
419 | END SELECT |
---|
420 | ! |
---|
421 | CASE DEFAULT |
---|
422 | CALL ctl_stop( 'bdy_ice_dyn : unrecognised option for open boundaries for ice fields' ) |
---|
423 | END SELECT |
---|
424 | ! |
---|
425 | END DO ! jbdy |
---|
426 | ! |
---|
427 | SELECT CASE ( cd_type ) |
---|
428 | CASE ( 'U' ) |
---|
429 | IF( nn_hls > 1 .AND. ir == 1 ) CYCLE ! at least 2 halos will be corrected -> no need to correct rim 1 before rim 0 |
---|
430 | IF( nn_hls == 1 ) THEN ; llsend2(:) = .false. ; llrecv2(:) = .false. ; END IF |
---|
431 | DO jbdy = 1, nb_bdy |
---|
432 | IF( cn_ice(jbdy) == 'frs' .AND. nn_ice_dta(jbdy) /= 0 ) THEN |
---|
433 | llsend2( : ) = llsend2( : ) .OR. lsend_bdyint(jbdy,2, : ,ir) ! possibly every direction, U points |
---|
434 | idir3 = (/ jpwe, jpsw, jpnw /) |
---|
435 | llsend2(idir3) = llsend2(idir3) .OR. lsend_bdyext(jbdy,2,idir3,ir) ! nei might search point towards its ea bdy |
---|
436 | llrecv2( : ) = llrecv2( : ) .OR. lrecv_bdyint(jbdy,2, : ,ir) ! possibly every direction, U points |
---|
437 | idir3 = (/ jpea, jpse, jpne /) |
---|
438 | llrecv2(idir3) = llrecv2(idir3) .OR. lrecv_bdyext(jbdy,2,idir3,ir) ! might search point towards east bdy |
---|
439 | END IF |
---|
440 | END DO |
---|
441 | IF( ANY(llsend2) .OR. ANY(llrecv2) ) THEN ! if need to send/recv in at least one direction |
---|
442 | CALL lbc_lnk( 'bdyice', u_ice, 'U', -1.0_wp, kfillmode=jpfillnothing ,lsend=llsend2, lrecv=llrecv2 ) |
---|
443 | END IF |
---|
444 | CASE ( 'V' ) |
---|
445 | IF( nn_hls > 1 .AND. ir == 1 ) CYCLE ! at least 2 halos will be corrected -> no need to correct rim 1 before rim 0 |
---|
446 | IF( nn_hls == 1 ) THEN ; llsend3(:) = .false. ; llrecv3(:) = .false. ; END IF |
---|
447 | DO jbdy = 1, nb_bdy |
---|
448 | IF( cn_ice(jbdy) == 'frs' .AND. nn_ice_dta(jbdy) /= 0 ) THEN |
---|
449 | llsend3( : ) = llsend3( : ) .OR. lsend_bdyint(jbdy,3, : ,ir) ! possibly every direction, V points |
---|
450 | idir3 = (/ jpso, jpsw, jpse /) |
---|
451 | llsend3(idir3) = llsend3(idir3) .OR. lsend_bdyext(jbdy,3,idir3,ir) ! nei might search point towards its no bdy |
---|
452 | llrecv3( : ) = llrecv3( : ) .OR. lrecv_bdyint(jbdy,3, : ,ir) ! possibly every direction, V points |
---|
453 | idir3 = (/ jpno, jpnw, jpne /) |
---|
454 | llrecv3(idir3) = llrecv3(idir3) .OR. lrecv_bdyext(jbdy,3,idir3,ir) ! might search point towards north bdy |
---|
455 | END IF |
---|
456 | END DO |
---|
457 | IF( ANY(llsend3) .OR. ANY(llrecv3) ) THEN ! if need to send/recv in at least one direction |
---|
458 | CALL lbc_lnk( 'bdyice', v_ice, 'V', -1.0_wp, kfillmode=jpfillnothing ,lsend=llsend3, lrecv=llrecv3 ) |
---|
459 | END IF |
---|
460 | END SELECT |
---|
461 | END DO ! ir |
---|
462 | ! |
---|
463 | IF( ln_timing ) CALL timing_stop('bdy_ice_dyn') |
---|
464 | ! |
---|
465 | END SUBROUTINE bdy_ice_dyn |
---|
466 | |
---|
467 | #else |
---|
468 | !!--------------------------------------------------------------------------------- |
---|
469 | !! Default option Empty module |
---|
470 | !!--------------------------------------------------------------------------------- |
---|
471 | CONTAINS |
---|
472 | SUBROUTINE bdy_ice( kt ) ! Empty routine |
---|
473 | IMPLICIT NONE |
---|
474 | INTEGER, INTENT( in ) :: kt |
---|
475 | WRITE(*,*) 'bdy_ice: You should not have seen this print! error?', kt |
---|
476 | END SUBROUTINE bdy_ice |
---|
477 | #endif |
---|
478 | |
---|
479 | !!================================================================================= |
---|
480 | END MODULE bdyice |
---|