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 licence (./LICENSE) |
---|
45 | !!---------------------------------------------------------------------- |
---|
46 | CONTAINS |
---|
47 | |
---|
48 | SUBROUTINE bdy_ice( kt ) |
---|
49 | !!---------------------------------------------------------------------- |
---|
50 | !! *** SUBROUTINE bdy_ice *** |
---|
51 | !! |
---|
52 | !! ** Purpose : - Apply open boundary conditions for ice (SI3) |
---|
53 | !! |
---|
54 | !!---------------------------------------------------------------------- |
---|
55 | INTEGER, INTENT(in) :: kt ! Main time step counter |
---|
56 | ! |
---|
57 | INTEGER :: ib_bdy ! Loop index |
---|
58 | !!---------------------------------------------------------------------- |
---|
59 | ! |
---|
60 | IF( ln_timing ) CALL timing_start('bdy_ice') |
---|
61 | ! |
---|
62 | CALL ice_var_glo2eqv |
---|
63 | ! |
---|
64 | DO ib_bdy = 1, nb_bdy |
---|
65 | ! |
---|
66 | SELECT CASE( cn_ice(ib_bdy) ) |
---|
67 | CASE('none') ; CYCLE |
---|
68 | CASE('frs' ) ; CALL bdy_ice_frs( idx_bdy(ib_bdy), dta_bdy(ib_bdy), kt, ib_bdy ) |
---|
69 | CASE DEFAULT |
---|
70 | CALL ctl_stop( 'bdy_ice : unrecognised option for open boundaries for ice fields' ) |
---|
71 | END SELECT |
---|
72 | ! |
---|
73 | END DO |
---|
74 | ! |
---|
75 | CALL ice_cor( kt , 0 ) ! -- In case categories are out of bounds, do a remapping |
---|
76 | ! ! i.e. inputs have not the same ice thickness distribution (set by rn_himean) |
---|
77 | ! ! than the regional simulation |
---|
78 | CALL ice_var_agg(1) |
---|
79 | ! |
---|
80 | IF( ln_icectl ) CALL ice_prt( kt, iiceprt, jiceprt, 1, ' - ice thermo bdy - ' ) |
---|
81 | IF( ln_timing ) CALL timing_stop('bdy_ice') |
---|
82 | ! |
---|
83 | END SUBROUTINE bdy_ice |
---|
84 | |
---|
85 | |
---|
86 | SUBROUTINE bdy_ice_frs( idx, dta, kt, ib_bdy ) |
---|
87 | !!------------------------------------------------------------------------------ |
---|
88 | !! *** SUBROUTINE bdy_ice_frs *** |
---|
89 | !! |
---|
90 | !! ** Purpose : Apply the Flow Relaxation Scheme for sea-ice fields in the case |
---|
91 | !! of unstructured open boundaries. |
---|
92 | !! |
---|
93 | !! Reference : Engedahl H., 1995: Use of the flow relaxation scheme in a three- |
---|
94 | !! dimensional baroclinic ocean model with realistic topography. Tellus, 365-382. |
---|
95 | !!------------------------------------------------------------------------------ |
---|
96 | TYPE(OBC_INDEX), INTENT(in) :: idx ! OBC indices |
---|
97 | TYPE(OBC_DATA), INTENT(in) :: dta ! OBC external data |
---|
98 | INTEGER, INTENT(in) :: kt ! main time-step counter |
---|
99 | INTEGER, INTENT(in) :: ib_bdy ! BDY set index |
---|
100 | ! |
---|
101 | INTEGER :: jpbound ! 0 = incoming ice |
---|
102 | ! ! 1 = outgoing ice |
---|
103 | INTEGER :: jb, jk, jgrd, jl ! dummy loop indices |
---|
104 | INTEGER :: ji, jj, ii, ij ! local scalar |
---|
105 | REAL(wp) :: zwgt, zwgt1 ! local scalar |
---|
106 | REAL(wp) :: ztmelts, zdh |
---|
107 | !!------------------------------------------------------------------------------ |
---|
108 | ! |
---|
109 | jgrd = 1 ! Everything is at T-points here |
---|
110 | ! |
---|
111 | DO jl = 1, jpl |
---|
112 | DO jb = 1, idx%nblenrim(jgrd) |
---|
113 | ji = idx%nbi(jb,jgrd) |
---|
114 | jj = idx%nbj(jb,jgrd) |
---|
115 | zwgt = idx%nbw(jb,jgrd) |
---|
116 | zwgt1 = 1.e0 - idx%nbw(jb,jgrd) |
---|
117 | a_i(ji,jj,jl) = ( a_i(ji,jj,jl) * zwgt1 + dta%a_i(jb,jl) * zwgt ) * tmask(ji,jj,1) ! Leads fraction |
---|
118 | h_i(ji,jj,jl) = ( h_i(ji,jj,jl) * zwgt1 + dta%h_i(jb,jl) * zwgt ) * tmask(ji,jj,1) ! Ice depth |
---|
119 | h_s(ji,jj,jl) = ( h_s(ji,jj,jl) * zwgt1 + dta%h_s(jb,jl) * zwgt ) * tmask(ji,jj,1) ! Snow depth |
---|
120 | |
---|
121 | ! ----------------- |
---|
122 | ! Pathological case |
---|
123 | ! ----------------- |
---|
124 | ! In case a) snow load would be in excess or b) ice is coming into a warmer environment that would lead to |
---|
125 | ! very large transformation from snow to ice (see icethd_dh.F90) |
---|
126 | |
---|
127 | ! Then, a) transfer the snow excess into the ice (different from icethd_dh) |
---|
128 | zdh = MAX( 0._wp, ( rhosn * h_s(ji,jj,jl) + ( rhoic - rau0 ) * h_i(ji,jj,jl) ) * r1_rau0 ) |
---|
129 | ! Or, b) transfer all the snow into ice (if incoming ice is likely to melt as it comes into a warmer environment) |
---|
130 | !zdh = MAX( 0._wp, h_s(ji,jj,jl) * rhosn / rhoic ) |
---|
131 | |
---|
132 | ! recompute h_i, h_s |
---|
133 | h_i(ji,jj,jl) = MIN( hi_max(jl), h_i(ji,jj,jl) + zdh ) |
---|
134 | h_s(ji,jj,jl) = MAX( 0._wp, h_s(ji,jj,jl) - zdh * rhoic / rhosn ) |
---|
135 | |
---|
136 | ENDDO |
---|
137 | CALL lbc_bdy_lnk( a_i(:,:,jl), 'T', 1., ib_bdy ) |
---|
138 | CALL lbc_bdy_lnk( h_i(:,:,jl), 'T', 1., ib_bdy ) |
---|
139 | CALL lbc_bdy_lnk( h_s(:,:,jl), 'T', 1., ib_bdy ) |
---|
140 | ENDDO |
---|
141 | |
---|
142 | DO jl = 1, jpl |
---|
143 | DO jb = 1, idx%nblenrim(jgrd) |
---|
144 | ji = idx%nbi(jb,jgrd) |
---|
145 | jj = idx%nbj(jb,jgrd) |
---|
146 | |
---|
147 | ! condition on ice thickness depends on the ice velocity |
---|
148 | ! if velocity is outward (strictly), then ice thickness, volume... must be equal to adjacent values |
---|
149 | jpbound = 0 ; ii = ji ; ij = jj |
---|
150 | ! |
---|
151 | IF( u_ice(ji+1,jj ) < 0. .AND. umask(ji-1,jj ,1) == 0. ) jpbound = 1; ii = ji+1; ij = jj |
---|
152 | IF( u_ice(ji-1,jj ) > 0. .AND. umask(ji+1,jj ,1) == 0. ) jpbound = 1; ii = ji-1; ij = jj |
---|
153 | IF( v_ice(ji ,jj+1) < 0. .AND. vmask(ji ,jj-1,1) == 0. ) jpbound = 1; ii = ji ; ij = jj+1 |
---|
154 | IF( v_ice(ji ,jj-1) > 0. .AND. vmask(ji ,jj+1,1) == 0. ) jpbound = 1; ii = ji ; ij = jj-1 |
---|
155 | ! |
---|
156 | IF( nn_ice_dta(ib_bdy) == 0 ) jpbound = 0; ii = ji; ij = jj ! case ice boundaries = initial conditions |
---|
157 | ! ! do not make state variables dependent on velocity |
---|
158 | ! |
---|
159 | IF( a_i(ii,ij,jl) > 0._wp ) THEN ! there is ice at the boundary |
---|
160 | ! |
---|
161 | a_i(ji,jj,jl) = a_i(ii,ij,jl) ! concentration |
---|
162 | h_i(ji,jj,jl) = h_i(ii,ij,jl) ! thickness ice |
---|
163 | h_s(ji,jj,jl) = h_s(ii,ij,jl) ! thickness snw |
---|
164 | ! |
---|
165 | SELECT CASE( jpbound ) |
---|
166 | ! |
---|
167 | CASE( 0 ) ! velocity is inward |
---|
168 | ! |
---|
169 | oa_i(ji,jj, jl) = rn_ice_age(ib_bdy) * a_i(ji,jj,jl) ! age |
---|
170 | a_ip(ji,jj, jl) = 0._wp ! pond concentration |
---|
171 | v_ip(ji,jj, jl) = 0._wp ! pond volume |
---|
172 | t_su(ji,jj, jl) = rn_ice_tem(ib_bdy) ! temperature surface |
---|
173 | t_s (ji,jj,:,jl) = rn_ice_tem(ib_bdy) ! temperature snw |
---|
174 | t_i (ji,jj,:,jl) = rn_ice_tem(ib_bdy) ! temperature ice |
---|
175 | s_i (ji,jj, jl) = rn_ice_sal(ib_bdy) ! salinity |
---|
176 | sz_i(ji,jj,:,jl) = rn_ice_sal(ib_bdy) ! salinity profile |
---|
177 | ! |
---|
178 | CASE( 1 ) ! velocity is outward |
---|
179 | ! |
---|
180 | oa_i(ji,jj, jl) = oa_i(ii,ij, jl) ! age |
---|
181 | a_ip(ji,jj, jl) = a_ip(ii,ij, jl) ! pond concentration |
---|
182 | v_ip(ji,jj, jl) = v_ip(ii,ij, jl) ! pond volume |
---|
183 | t_su(ji,jj, jl) = t_su(ii,ij, jl) ! temperature surface |
---|
184 | t_s (ji,jj,:,jl) = t_s (ii,ij,:,jl) ! temperature snw |
---|
185 | t_i (ji,jj,:,jl) = t_i (ii,ij,:,jl) ! temperature ice |
---|
186 | s_i (ji,jj, jl) = s_i (ii,ij, jl) ! salinity |
---|
187 | sz_i(ji,jj,:,jl) = sz_i(ii,ij,:,jl) ! salinity profile |
---|
188 | ! |
---|
189 | END SELECT |
---|
190 | ! |
---|
191 | IF( nn_icesal == 1 ) THEN ! if constant salinity |
---|
192 | s_i (ji,jj ,jl) = rn_icesal |
---|
193 | sz_i(ji,jj,:,jl) = rn_icesal |
---|
194 | ENDIF |
---|
195 | ! |
---|
196 | ! global fields |
---|
197 | v_i (ji,jj,jl) = h_i(ji,jj,jl) * a_i(ji,jj,jl) ! volume ice |
---|
198 | v_s (ji,jj,jl) = h_s(ji,jj,jl) * a_i(ji,jj,jl) ! volume snw |
---|
199 | sv_i(ji,jj,jl) = MIN( s_i(ji,jj,jl) , sss_m(ji,jj) ) * v_i(ji,jj,jl) ! salt content |
---|
200 | DO jk = 1, nlay_s |
---|
201 | e_s(ji,jj,jk,jl) = rhosn * ( cpic * ( rt0 - t_s(ji,jj,jk,jl) ) + lfus ) ! enthalpy in J/m3 |
---|
202 | e_s(ji,jj,jk,jl) = e_s(ji,jj,jk,jl) * v_s(ji,jj,jl) * r1_nlay_s ! enthalpy in J/m2 |
---|
203 | END DO |
---|
204 | DO jk = 1, nlay_i |
---|
205 | ztmelts = - tmut * sz_i(ji,jj,jk,jl) ! Melting temperature in C |
---|
206 | 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 |
---|
207 | ! |
---|
208 | e_i(ji,jj,jk,jl) = rhoic * ( cpic * ( ztmelts - ( t_i(ji,jj,jk,jl) - rt0 ) ) & ! enthalpy in J/m3 |
---|
209 | & + lfus * ( 1._wp - ztmelts / ( t_i(ji,jj,jk,jl) - rt0 ) ) & |
---|
210 | & - rcp * ztmelts ) |
---|
211 | e_i(ji,jj,jk,jl) = e_i(ji,jj,jk,jl) * v_i(ji,jj,jl) * r1_nlay_i ! enthalpy in J/m2 |
---|
212 | END DO |
---|
213 | ! |
---|
214 | ELSE ! no ice at the boundary |
---|
215 | ! |
---|
216 | a_i (ji,jj, jl) = 0._wp |
---|
217 | h_i (ji,jj, jl) = 0._wp |
---|
218 | h_s (ji,jj, jl) = 0._wp |
---|
219 | oa_i(ji,jj, jl) = 0._wp |
---|
220 | a_ip(ji,jj, jl) = 0._wp |
---|
221 | v_ip(ji,jj, jl) = 0._wp |
---|
222 | t_su(ji,jj, jl) = rt0 |
---|
223 | t_s (ji,jj,:,jl) = rt0 |
---|
224 | t_i (ji,jj,:,jl) = rt0 |
---|
225 | |
---|
226 | IF( nn_icesal == 1 ) THEN ! if constant salinity |
---|
227 | s_i (ji,jj ,jl) = rn_icesal |
---|
228 | sz_i(ji,jj,:,jl) = rn_icesal |
---|
229 | ELSE ! if variable salinity |
---|
230 | s_i (ji,jj,jl) = rn_simin |
---|
231 | sz_i(ji,jj,:,jl) = rn_simin |
---|
232 | ENDIF |
---|
233 | ! |
---|
234 | ! global fields |
---|
235 | v_i (ji,jj, jl) = 0._wp |
---|
236 | v_s (ji,jj, jl) = 0._wp |
---|
237 | sv_i(ji,jj, jl) = 0._wp |
---|
238 | e_s (ji,jj,:,jl) = 0._wp |
---|
239 | e_i (ji,jj,:,jl) = 0._wp |
---|
240 | |
---|
241 | ENDIF |
---|
242 | |
---|
243 | END DO |
---|
244 | ! |
---|
245 | CALL lbc_bdy_lnk( a_i (:,:,jl), 'T', 1., ib_bdy ) |
---|
246 | CALL lbc_bdy_lnk( h_i (:,:,jl), 'T', 1., ib_bdy ) |
---|
247 | CALL lbc_bdy_lnk( h_s (:,:,jl), 'T', 1., ib_bdy ) |
---|
248 | CALL lbc_bdy_lnk( oa_i(:,:,jl), 'T', 1., ib_bdy ) |
---|
249 | CALL lbc_bdy_lnk( a_ip(:,:,jl), 'T', 1., ib_bdy ) |
---|
250 | CALL lbc_bdy_lnk( v_ip(:,:,jl), 'T', 1., ib_bdy ) |
---|
251 | CALL lbc_bdy_lnk( s_i (:,:,jl), 'T', 1., ib_bdy ) |
---|
252 | CALL lbc_bdy_lnk( t_su(:,:,jl), 'T', 1., ib_bdy ) |
---|
253 | CALL lbc_bdy_lnk( v_i (:,:,jl), 'T', 1., ib_bdy ) |
---|
254 | CALL lbc_bdy_lnk( v_s (:,:,jl), 'T', 1., ib_bdy ) |
---|
255 | CALL lbc_bdy_lnk( sv_i(:,:,jl), 'T', 1., ib_bdy ) |
---|
256 | DO jk = 1, nlay_s |
---|
257 | CALL lbc_bdy_lnk(t_s(:,:,jk,jl), 'T', 1., ib_bdy ) |
---|
258 | CALL lbc_bdy_lnk(e_s(:,:,jk,jl), 'T', 1., ib_bdy ) |
---|
259 | END DO |
---|
260 | DO jk = 1, nlay_i |
---|
261 | CALL lbc_bdy_lnk(t_i(:,:,jk,jl), 'T', 1., ib_bdy ) |
---|
262 | CALL lbc_bdy_lnk(e_i(:,:,jk,jl), 'T', 1., ib_bdy ) |
---|
263 | END DO |
---|
264 | ! |
---|
265 | END DO ! jl |
---|
266 | ! |
---|
267 | END SUBROUTINE bdy_ice_frs |
---|
268 | |
---|
269 | |
---|
270 | SUBROUTINE bdy_ice_dyn( cd_type ) |
---|
271 | !!------------------------------------------------------------------------------ |
---|
272 | !! *** SUBROUTINE bdy_ice_dyn *** |
---|
273 | !! |
---|
274 | !! ** Purpose : Apply dynamics boundary conditions for sea-ice in the cas of unstructured open boundaries. |
---|
275 | !! u_ice and v_ice are equal to the value of the adjacent grid point if this latter is not ice free |
---|
276 | !! if adjacent grid point is ice free, then u_ice and v_ice are equal to ocean velocities |
---|
277 | !! |
---|
278 | !! 2013-06 : C. Rousset |
---|
279 | !!------------------------------------------------------------------------------ |
---|
280 | CHARACTER(len=1), INTENT(in) :: cd_type ! nature of velocity grid-points |
---|
281 | ! |
---|
282 | INTEGER :: jb, jgrd ! dummy loop indices |
---|
283 | INTEGER :: ji, jj ! local scalar |
---|
284 | INTEGER :: ib_bdy ! Loop index |
---|
285 | REAL(wp) :: zmsk1, zmsk2, zflag |
---|
286 | !!------------------------------------------------------------------------------ |
---|
287 | ! |
---|
288 | DO ib_bdy=1, nb_bdy |
---|
289 | ! |
---|
290 | SELECT CASE( cn_ice(ib_bdy) ) |
---|
291 | ! |
---|
292 | CASE('none') |
---|
293 | CYCLE |
---|
294 | ! |
---|
295 | CASE('frs') |
---|
296 | ! |
---|
297 | IF( nn_ice_dta(ib_bdy) == 0 ) CYCLE ! case ice boundaries = initial conditions |
---|
298 | ! ! do not change ice velocity (it is only computed by rheology) |
---|
299 | SELECT CASE ( cd_type ) |
---|
300 | ! |
---|
301 | CASE ( 'U' ) |
---|
302 | jgrd = 2 ! u velocity |
---|
303 | DO jb = 1, idx_bdy(ib_bdy)%nblenrim(jgrd) |
---|
304 | ji = idx_bdy(ib_bdy)%nbi(jb,jgrd) |
---|
305 | jj = idx_bdy(ib_bdy)%nbj(jb,jgrd) |
---|
306 | zflag = idx_bdy(ib_bdy)%flagu(jb,jgrd) |
---|
307 | ! |
---|
308 | IF ( ABS( zflag ) == 1. ) THEN ! eastern and western boundaries |
---|
309 | ! one of the two zmsk is always 0 (because of zflag) |
---|
310 | zmsk1 = 1._wp - MAX( 0.0_wp, SIGN ( 1.0_wp , - vt_i(ji+1,jj) ) ) ! 0 if no ice |
---|
311 | zmsk2 = 1._wp - MAX( 0.0_wp, SIGN ( 1.0_wp , - vt_i(ji-1,jj) ) ) ! 0 if no ice |
---|
312 | ! |
---|
313 | ! u_ice = u_ice of the adjacent grid point except if this grid point is ice-free (then do not change u_ice) |
---|
314 | u_ice (ji,jj) = u_ice(ji+1,jj) * 0.5_wp * ABS( zflag + 1._wp ) * zmsk1 + & |
---|
315 | & u_ice(ji-1,jj) * 0.5_wp * ABS( zflag - 1._wp ) * zmsk2 + & |
---|
316 | & u_ice(ji ,jj) * ( 1._wp - MIN( 1._wp, zmsk1 + zmsk2 ) ) |
---|
317 | ELSE ! everywhere else |
---|
318 | u_ice(ji,jj) = 0._wp |
---|
319 | ENDIF |
---|
320 | ! |
---|
321 | END DO |
---|
322 | CALL lbc_bdy_lnk( u_ice(:,:), 'U', -1., ib_bdy ) |
---|
323 | ! |
---|
324 | CASE ( 'V' ) |
---|
325 | jgrd = 3 ! v velocity |
---|
326 | DO jb = 1, idx_bdy(ib_bdy)%nblenrim(jgrd) |
---|
327 | ji = idx_bdy(ib_bdy)%nbi(jb,jgrd) |
---|
328 | jj = idx_bdy(ib_bdy)%nbj(jb,jgrd) |
---|
329 | zflag = idx_bdy(ib_bdy)%flagv(jb,jgrd) |
---|
330 | ! |
---|
331 | IF ( ABS( zflag ) == 1. ) THEN ! northern and southern boundaries |
---|
332 | ! one of the two zmsk is always 0 (because of zflag) |
---|
333 | zmsk1 = 1._wp - MAX( 0.0_wp, SIGN ( 1.0_wp , - vt_i(ji,jj+1) ) ) ! 0 if no ice |
---|
334 | zmsk2 = 1._wp - MAX( 0.0_wp, SIGN ( 1.0_wp , - vt_i(ji,jj-1) ) ) ! 0 if no ice |
---|
335 | ! |
---|
336 | ! v_ice = v_ice of the adjacent grid point except if this grid point is ice-free (then do not change v_ice) |
---|
337 | v_ice (ji,jj) = v_ice(ji,jj+1) * 0.5_wp * ABS( zflag + 1._wp ) * zmsk1 + & |
---|
338 | & v_ice(ji,jj-1) * 0.5_wp * ABS( zflag - 1._wp ) * zmsk2 + & |
---|
339 | & v_ice(ji,jj ) * ( 1._wp - MIN( 1._wp, zmsk1 + zmsk2 ) ) |
---|
340 | ELSE ! everywhere else |
---|
341 | v_ice(ji,jj) = 0._wp |
---|
342 | ENDIF |
---|
343 | ! |
---|
344 | END DO |
---|
345 | CALL lbc_bdy_lnk( v_ice(:,:), 'V', -1., ib_bdy ) |
---|
346 | ! |
---|
347 | END SELECT |
---|
348 | ! |
---|
349 | CASE DEFAULT |
---|
350 | CALL ctl_stop( 'bdy_ice_dyn : unrecognised option for open boundaries for ice fields' ) |
---|
351 | END SELECT |
---|
352 | ! |
---|
353 | END DO |
---|
354 | ! |
---|
355 | END SUBROUTINE bdy_ice_dyn |
---|
356 | |
---|
357 | #else |
---|
358 | !!--------------------------------------------------------------------------------- |
---|
359 | !! Default option Empty module |
---|
360 | !!--------------------------------------------------------------------------------- |
---|
361 | CONTAINS |
---|
362 | SUBROUTINE bdy_ice( kt ) ! Empty routine |
---|
363 | WRITE(*,*) 'bdy_ice: You should not have seen this print! error?', kt |
---|
364 | END SUBROUTINE bdy_ice |
---|
365 | #endif |
---|
366 | |
---|
367 | !!================================================================================= |
---|
368 | END MODULE bdyice |
---|