New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
bdyice_lim.F90 in trunk/NEMOGCM/NEMO/OPA_SRC/BDY – NEMO

source: trunk/NEMOGCM/NEMO/OPA_SRC/BDY/bdyice_lim.F90 @ 4767

Last change on this file since 4767 was 4767, checked in by jamesharle, 10 years ago

Addition of import module dom_ice_2 to allow compilation of BDY and LIM2 (see ticket #1383)

File size: 20.7 KB
Line 
1MODULE 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 par_ice
29   USE ice             ! LIM_3 ice variables
30   USE dom_ice         ! sea-ice domain
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   REAL(wp) ::   epsi20 = 1.e-20_wp  ! module constants
48   REAL(wp) ::   epsi10 = 1.e-10_wp  ! min area allowed by ice model
49   !!----------------------------------------------------------------------
50   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
51   !! $Id: bdyice.F90 2715 2011-03-30 15:58:35Z rblod $
52   !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
53   !!----------------------------------------------------------------------
54CONTAINS
55
56   SUBROUTINE bdy_ice_lim( kt )
57      !!----------------------------------------------------------------------
58      !!                  ***  SUBROUTINE bdy_ice_lim  ***
59      !!
60      !! ** Purpose : - Apply open boundary conditions for ice (LIM2 and LIM3)
61      !!
62      !!----------------------------------------------------------------------
63      INTEGER, INTENT( in ) :: kt     ! Main time step counter
64      !!
65      INTEGER               :: ib_bdy ! Loop index
66      DO ib_bdy=1, nb_bdy
67
68         SELECT CASE( cn_ice_lim(ib_bdy) )
69         CASE('none')
70            CYCLE
71         CASE('frs')
72            CALL bdy_ice_frs( idx_bdy(ib_bdy), dta_bdy(ib_bdy), kt, ib_bdy )
73         CASE DEFAULT
74            CALL ctl_stop( 'bdy_ice_lim : unrecognised option for open boundaries for ice fields' )
75         END SELECT
76      ENDDO
77
78   END SUBROUTINE bdy_ice_lim
79
80   SUBROUTINE bdy_ice_frs( idx, dta, kt, ib_bdy )
81      !!------------------------------------------------------------------------------
82      !!                 ***  SUBROUTINE bdy_ice_frs  ***
83      !!                   
84      !! ** Purpose : Apply the Flow Relaxation Scheme for sea-ice fields in the case
85      !!              of unstructured open boundaries.
86      !!
87      !! Reference : Engedahl H., 1995: Use of the flow relaxation scheme in a three-
88      !!             dimensional baroclinic ocean model with realistic topography. Tellus, 365-382.
89      !!------------------------------------------------------------------------------
90      TYPE(OBC_INDEX), INTENT(in) ::   idx  ! OBC indices
91      TYPE(OBC_DATA),  INTENT(in) ::   dta  ! OBC external data
92      INTEGER,         INTENT(in) ::   kt   ! main time-step counter
93      INTEGER,         INTENT(in) ::   ib_bdy  ! BDY set index      !!
94
95      INTEGER  ::   jpbound            ! 0 = incoming ice
96                                       ! 1 = outgoing ice
97      INTEGER  ::   jb, jk, jgrd, jl   ! dummy loop indices
98      INTEGER  ::   ji, jj, ii, ij     ! local scalar
99      REAL(wp) ::   zwgt, zwgt1        ! local scalar
100      REAL(wp) ::   zinda, ztmelts, zdh
101
102      !!------------------------------------------------------------------------------
103      !
104      IF( nn_timing == 1 ) CALL timing_start('bdy_ice_frs')
105      !
106      jgrd = 1      ! Everything is at T-points here
107      !
108#if defined key_lim2
109      DO jb = 1, idx%nblen(jgrd)
110         ji    = idx%nbi(jb,jgrd)
111         jj    = idx%nbj(jb,jgrd)
112         zwgt  = idx%nbw(jb,jgrd)
113         zwgt1 = 1.e0 - idx%nbw(jb,jgrd)
114         frld (ji,jj) = ( frld (ji,jj) * zwgt1 + dta%frld (jb) * zwgt ) * tmask(ji,jj,1)     ! Leads fraction
115         hicif(ji,jj) = ( hicif(ji,jj) * zwgt1 + dta%hicif(jb) * zwgt ) * tmask(ji,jj,1)     ! Ice depth
116         hsnif(ji,jj) = ( hsnif(ji,jj) * zwgt1 + dta%hsnif(jb) * zwgt ) * tmask(ji,jj,1)     ! Snow depth
117
118!         zinda = 1.0 - MAX( 0.0_wp , SIGN ( 1.0_wp , - frld(ji,jj) ) ) ! 0 if no ice
119!         !------------------------------
120!         ! Sea ice surface temperature
121!         !------------------------------
122!         sist(ji,jj)   = zinda * 270.0 + ( 1.0 - zinda ) * tfu(ji,jj)
123!         !-----------------------------------------------
124!         ! Ice/snow temperatures and energy stored in brines
125!         !-----------------------------------------------
126!         !!! TO BE CONTIUNED (as LIM3 below) !!!
127!            zindhe = MAX( 0.e0, SIGN( 1.e0, fcor(1,jj) ) )              ! = 0 for SH, =1 for NH
128!
129!               ! Recover in situ values.
130!               zindb         = MAX( rzero, SIGN( rone, zs0a(ji,jj) - epsi06 ) )
131!               zacrith       = 1.0 - ( zindhe * acrit(1) + ( 1.0 - zindhe ) * acrit(2) )
132!               zs0a (ji,jj)  = zindb * MIN( zs0a(ji,jj), zacrith )
133!               hsnif(ji,jj)  = zindb * ( zs0sn(ji,jj) /MAX( zs0a(ji,jj), epsi16 ) )
134!               hicif(ji,jj)  = zindb * ( zs0ice(ji,jj)/MAX( zs0a(ji,jj), epsi16 ) )
135!               zindsn        = MAX( rzero, SIGN( rone, hsnif(ji,jj) - epsi06 ) )
136!               zindic        = MAX( rzero, SIGN( rone, hicif(ji,jj) - epsi03 ) )
137!               zindb         = MAX( zindsn, zindic )
138!               zs0a (ji,jj)  = zindb * zs0a(ji,jj)
139!               frld (ji,jj)  = 1.0 - zs0a(ji,jj)
140!               hsnif(ji,jj)  = zindsn * hsnif(ji,jj)
141!               hicif(ji,jj)  = zindic * hicif(ji,jj)
142!               zusvosn       = 1.0/MAX( hsnif(ji,jj) * zs0a(ji,jj), epsi16 )
143!               zusvoic       = 1.0/MAX( hicif(ji,jj) * zs0a(ji,jj), epsi16 )
144!               zignm         = MAX( rzero,  SIGN( rone, hsndif - hsnif(ji,jj) ) )
145!               zrtt          = 173.15 * rone
146!               ztsn          =          zignm   * tbif(ji,jj,1)  &
147!                              + ( 1.0 - zignm ) * MIN( MAX( zrtt, rt0_snow * zusvosn * zs0c0(ji,jj)) , tfu(ji,jj) )
148!               ztic1          = MIN( MAX( zrtt, rt0_ice * zusvoic * zs0c1(ji,jj) ) , tfu(ji,jj) )
149!               ztic2          = MIN( MAX( zrtt, rt0_ice * zusvoic * zs0c2(ji,jj) ) , tfu(ji,jj) )
150!
151!               tbif(ji,jj,1) = zindsn * ztsn  + ( 1.0 - zindsn ) * tfu(ji,jj)               
152!               tbif(ji,jj,2) = zindic * ztic1 + ( 1.0 - zindic ) * tfu(ji,jj)
153!               tbif(ji,jj,3) = zindic * ztic2 + ( 1.0 - zindic ) * tfu(ji,jj)
154!               qstoif(ji,jj) = zindb  * xlic * zs0st(ji,jj) /  MAX( zs0a(ji,jj), epsi16 )
155
156      END DO
157
158      CALL lbc_bdy_lnk( frld,  'T', 1., ib_bdy )                                         ! lateral boundary conditions
159      CALL lbc_bdy_lnk( hicif, 'T', 1., ib_bdy )
160      CALL lbc_bdy_lnk( hsnif, 'T', 1., ib_bdy )
161
162      vt_i(:,:) = hicif(:,:) * frld(:,:)
163      vt_s(:,:) = hsnif(:,:) * frld(:,:)
164      !
165#elif defined key_lim3
166
167      DO jl = 1, jpl
168         DO jb = 1, idx%nblen(jgrd)
169            ji    = idx%nbi(jb,jgrd)
170            jj    = idx%nbj(jb,jgrd)
171            zwgt  = idx%nbw(jb,jgrd)
172            zwgt1 = 1.e0 - idx%nbw(jb,jgrd)
173            a_i (ji,jj,jl) = ( a_i (ji,jj,jl) * zwgt1 + dta%a_i (jb,jl) * zwgt ) * tmask(ji,jj,1)  ! Leads fraction
174            ht_i(ji,jj,jl) = ( ht_i(ji,jj,jl) * zwgt1 + dta%ht_i(jb,jl) * zwgt ) * tmask(ji,jj,1)  ! Ice depth
175            ht_s(ji,jj,jl) = ( ht_s(ji,jj,jl) * zwgt1 + dta%ht_s(jb,jl) * zwgt ) * tmask(ji,jj,1)  ! Snow depth
176
177            ! -----------------
178            ! Pathological case
179            ! -----------------
180            ! In case a) snow load would be in excess or b) ice is coming into a warmer environment that would lead to
181            ! very large transformation from snow to ice (see limthd_dh.F90)
182
183            ! Then, a) transfer the snow excess into the ice (different from limthd_dh)
184            zdh = MAX( 0._wp, ( rhosn * ht_s(ji,jj,jl) + ( rhoic - rau0 ) * ht_i(ji,jj,jl) ) * r1_rau0 )
185            ! Or, b) transfer all the snow into ice (if incoming ice is likely to melt as it comes into a warmer environment)
186            !zdh = MAX( 0._wp, ht_s(ji,jj,jl) * rhosn / rhoic )
187
188            ! recompute ht_i, ht_s
189            ht_i(ji,jj,jl) = MIN( hi_max(jl), ht_i(ji,jj,jl) + zdh )
190            ht_s(ji,jj,jl) = MAX( 0._wp, ht_s(ji,jj,jl) - zdh * rhoic / rhosn ) 
191
192         ENDDO
193         CALL lbc_bdy_lnk(  a_i(:,:,jl), 'T', 1., ib_bdy )
194         CALL lbc_bdy_lnk( ht_i(:,:,jl), 'T', 1., ib_bdy )
195         CALL lbc_bdy_lnk( ht_s(:,:,jl), 'T', 1., ib_bdy )
196      ENDDO
197      ! retrieve at_i
198      at_i(:,:) = 0._wp
199      DO jl = 1, jpl
200         at_i(:,:) = a_i(:,:,jl) + at_i(:,:)
201      END DO
202
203      DO jl = 1, jpl
204         DO jb = 1, idx%nblen(jgrd)
205            ji    = idx%nbi(jb,jgrd)
206            jj    = idx%nbj(jb,jgrd)
207
208            ! condition on ice thickness depends on the ice velocity
209            ! if velocity is outward (strictly), then ice thickness, volume... must be equal to adjacent values
210            jpbound = 0; ii = ji; ij = jj;
211
212            IF ( u_ice(ji+1,jj  ) < 0. .AND. umask(ji-1,jj  ,1) == 0. ) jpbound = 1; ii = ji+1; ij = jj
213            IF ( u_ice(ji-1,jj  ) > 0. .AND. umask(ji+1,jj  ,1) == 0. ) jpbound = 1; ii = ji-1; ij = jj
214            IF ( v_ice(ji  ,jj+1) < 0. .AND. vmask(ji  ,jj-1,1) == 0. ) jpbound = 1; ii = ji  ; ij = jj+1
215            IF ( v_ice(ji  ,jj-1) > 0. .AND. vmask(ji  ,jj+1,1) == 0. ) jpbound = 1; ii = ji  ; ij = jj-1
216
217            zinda = 1.0 - MAX( 0.0_wp , SIGN ( 1.0_wp , - at_i(ii,ij) + 0.01 ) ) ! 0 if no ice
218
219            ! concentration and thickness
220            a_i (ji,jj,jl) = a_i (ii,ij,jl) * zinda
221            ht_i(ji,jj,jl) = ht_i(ii,ij,jl) * zinda
222            ht_s(ji,jj,jl) = ht_s(ii,ij,jl) * zinda
223
224            ! Ice and snow volumes
225            v_i(ji,jj,jl) = ht_i(ji,jj,jl) * a_i(ji,jj,jl)
226            v_s(ji,jj,jl) = ht_s(ji,jj,jl) * a_i(ji,jj,jl)
227
228            SELECT CASE( jpbound )
229
230            CASE( 0 ) ! velocity is inward
231
232               ! Ice salinity, age, temperature
233               sm_i(ji,jj,jl)   = zinda * rn_ice_sal(ib_bdy)  + ( 1.0 - zinda ) * s_i_min
234               o_i(ji,jj,jl)    = zinda * rn_ice_age(ib_bdy)  + ( 1.0 - zinda )
235               t_su(ji,jj,jl)   = zinda * rn_ice_tem(ib_bdy)  + ( 1.0 - zinda ) * rn_ice_tem(ib_bdy)
236               DO jk = 1, nlay_s
237                  t_s(ji,jj,jk,jl) = zinda * rn_ice_tem(ib_bdy) + ( 1.0 - zinda ) * rtt
238               END DO
239               DO jk = 1, nlay_i
240                  t_i(ji,jj,jk,jl) = zinda * rn_ice_tem(ib_bdy) + ( 1.0 - zinda ) * rtt 
241                  s_i(ji,jj,jk,jl) = zinda * rn_ice_sal(ib_bdy) + ( 1.0 - zinda ) * s_i_min
242               END DO
243               
244            CASE( 1 ) ! velocity is outward
245 
246               ! Ice salinity, age, temperature
247               sm_i(ji,jj,jl)   = zinda * sm_i(ii,ij,jl)  + ( 1.0 - zinda ) * s_i_min
248               o_i(ji,jj,jl)    = zinda * o_i(ii,ij,jl)   + ( 1.0 - zinda )
249               t_su(ji,jj,jl)   = zinda * t_su(ii,ij,jl)  + ( 1.0 - zinda ) * rtt
250               DO jk = 1, nlay_s
251                  t_s(ji,jj,jk,jl) = zinda * t_s(ii,ij,jk,jl) + ( 1.0 - zinda ) * rtt
252               END DO
253               DO jk = 1, nlay_i
254                  t_i(ji,jj,jk,jl) = zinda * t_i(ii,ij,jk,jl) + ( 1.0 - zinda ) * rtt
255                  s_i(ji,jj,jk,jl) = zinda * s_i(ii,ij,jk,jl) + ( 1.0 - zinda ) * s_i_min
256               END DO
257
258            END SELECT
259
260            ! if salinity is constant, then overwrite rn_ice_sal
261            IF( num_sal == 1 ) THEN
262               sm_i(ji,jj,jl)   = bulk_sal
263               s_i (ji,jj,:,jl) = bulk_sal
264            ENDIF
265
266            ! contents
267            smv_i(ji,jj,jl)  = MIN( sm_i(ji,jj,jl) , sss_m(ji,jj) ) * v_i(ji,jj,jl)
268            oa_i(ji,jj,jl)   = o_i(ji,jj,jl) * a_i(ji,jj,jl)
269            DO jk = 1, nlay_s
270               ! Snow energy of melting
271               e_s(ji,jj,jk,jl) = zinda * rhosn * ( cpic * ( rtt - t_s(ji,jj,jk,jl) ) + lfus )
272               ! Change dimensions
273               e_s(ji,jj,jk,jl) = e_s(ji,jj,jk,jl) / unit_fac
274               ! Multiply by volume, so that heat content in 10^9 Joules
275               e_s(ji,jj,jk,jl) = e_s(ji,jj,jk,jl) * area(ji,jj) * v_s(ji,jj,jl) / nlay_s
276            END DO
277            DO jk = 1, nlay_i
278               ztmelts          = - tmut * s_i(ji,jj,jk,jl) + rtt !Melting temperature in K                 
279               ! heat content per unit volume
280               e_i(ji,jj,jk,jl) = zinda * rhoic * &
281                  (   cpic    * ( ztmelts - t_i(ji,jj,jk,jl) ) &
282                  +   lfus    * ( 1.0 - (ztmelts-rtt) / MIN((t_i(ji,jj,jk,jl)-rtt),-epsi20) ) &
283                  - rcp      * ( ztmelts - rtt ) )
284               ! Correct dimensions to avoid big values
285               e_i(ji,jj,jk,jl) = e_i(ji,jj,jk,jl) / unit_fac 
286               ! Mutliply by ice volume, and divide by number of layers to get heat content in 10^9 J
287               e_i(ji,jj,jk,jl) = e_i(ji,jj,jk,jl) * area(ji,jj) * a_i(ji,jj,jl) * ht_i(ji,jj,jl) / nlay_i
288            END DO
289
290
291         END DO !jb
292 
293         CALL lbc_bdy_lnk(  a_i(:,:,jl), 'T', 1., ib_bdy )                                         ! lateral boundary conditions
294         CALL lbc_bdy_lnk( ht_i(:,:,jl), 'T', 1., ib_bdy )
295         CALL lbc_bdy_lnk( ht_s(:,:,jl), 'T', 1., ib_bdy )
296         CALL lbc_bdy_lnk(  v_i(:,:,jl), 'T', 1., ib_bdy )
297         CALL lbc_bdy_lnk(  v_s(:,:,jl), 'T', 1., ib_bdy )
298 
299         CALL lbc_bdy_lnk( smv_i(:,:,jl), 'T', 1., ib_bdy )
300         CALL lbc_bdy_lnk(  sm_i(:,:,jl), 'T', 1., ib_bdy )
301         CALL lbc_bdy_lnk(  oa_i(:,:,jl), 'T', 1., ib_bdy )
302         CALL lbc_bdy_lnk(   o_i(:,:,jl), 'T', 1., ib_bdy )
303         CALL lbc_bdy_lnk(  t_su(:,:,jl), 'T', 1., ib_bdy )
304         DO jk = 1, nlay_s
305            CALL lbc_bdy_lnk(t_s(:,:,jk,jl), 'T', 1., ib_bdy )
306            CALL lbc_bdy_lnk(e_s(:,:,jk,jl), 'T', 1., ib_bdy )
307         END DO
308         DO jk = 1, nlay_i
309            CALL lbc_bdy_lnk(t_i(:,:,jk,jl), 'T', 1., ib_bdy )
310            CALL lbc_bdy_lnk(e_i(:,:,jk,jl), 'T', 1., ib_bdy )
311         END DO
312
313      END DO !jl
314   
315#endif
316      !     
317      IF( nn_timing == 1 ) CALL timing_stop('bdy_ice_frs')
318      !
319   END SUBROUTINE bdy_ice_frs
320
321
322   SUBROUTINE bdy_ice_lim_dyn( cd_type )
323      !!------------------------------------------------------------------------------
324      !!                 ***  SUBROUTINE bdy_ice_lim_dyn  ***
325      !!                   
326      !! ** Purpose : Apply dynamics boundary conditions for sea-ice in the cas of unstructured open boundaries.
327      !!              u_ice and v_ice are equal to the value of the adjacent grid point if this latter is not ice free
328      !!              if adjacent grid point is ice free, then u_ice and v_ice are equal to ocean velocities
329      !!
330      !! 2013-06 : C. Rousset
331      !!------------------------------------------------------------------------------
332      !!
333      CHARACTER(len=1), INTENT(in)  ::   cd_type   ! nature of velocity grid-points
334      INTEGER  ::   jb, jgrd   ! dummy loop indices
335      INTEGER  ::   ji, jj             ! local scalar
336      INTEGER  ::   ib_bdy ! Loop index
337      REAL(wp) ::   zmsk1, zmsk2, zflag, zinda 
338     !!------------------------------------------------------------------------------
339      !
340      IF( nn_timing == 1 ) CALL timing_start('bdy_ice_lim_dyn')
341      !
342      DO ib_bdy=1, nb_bdy
343         !
344         SELECT CASE( cn_ice_lim(ib_bdy) )
345
346         CASE('none')
347
348            CYCLE
349           
350         CASE('frs')
351           
352
353            SELECT CASE ( cd_type )
354
355            CASE ( 'U' )
356               
357               jgrd = 2      ! u velocity
358               DO jb = 1, idx_bdy(ib_bdy)%nblen(jgrd)
359                  ji    = idx_bdy(ib_bdy)%nbi(jb,jgrd)
360                  jj    = idx_bdy(ib_bdy)%nbj(jb,jgrd)
361                  zflag = idx_bdy(ib_bdy)%flagu(jb,jgrd)
362                 
363                  IF ( ABS( zflag ) == 1. ) THEN  ! eastern and western boundaries
364                     ! one of the two zmsk is always 0 (because of zflag)
365                     zmsk1 = 1._wp - MAX( 0.0_wp, SIGN ( 1.0_wp , - vt_i(ji+1,jj) ) ) ! 0 if no ice
366                     zmsk2 = 1._wp - MAX( 0.0_wp, SIGN ( 1.0_wp , - vt_i(ji-1,jj) ) ) ! 0 if no ice
367                     
368                     ! u_ice = u_ice of the adjacent grid point except if this grid point is ice-free (then u_ice = u_oce)
369                     u_ice (ji,jj) = u_ice(ji+1,jj) * 0.5 * ABS( zflag + 1._wp ) * zmsk1 + &
370                        &            u_ice(ji-1,jj) * 0.5 * ABS( zflag - 1._wp ) * zmsk2 + &
371                        &            u_oce(ji  ,jj) * ( 1._wp - MIN( 1._wp, zmsk1 + zmsk2 ) )
372                  ELSE                             ! everywhere else
373                     !u_ice(ji,jj) = u_oce(ji,jj)
374                     u_ice(ji,jj) = 0._wp
375                  ENDIF
376                  ! mask ice velocities
377                  zinda = 1.0 - MAX( 0.0_wp , SIGN ( 1.0_wp , - at_i(ji,jj) + 0.01 ) ) ! 0 if no ice
378                  u_ice(ji,jj) = zinda * u_ice(ji,jj)
379                 
380               ENDDO
381
382               CALL lbc_bdy_lnk( u_ice(:,:), 'U', -1., ib_bdy )
383               
384            CASE ( 'V' )
385               
386               jgrd = 3      ! v velocity
387               DO jb = 1, idx_bdy(ib_bdy)%nblen(jgrd)
388                  ji    = idx_bdy(ib_bdy)%nbi(jb,jgrd)
389                  jj    = idx_bdy(ib_bdy)%nbj(jb,jgrd)
390                  zflag = idx_bdy(ib_bdy)%flagv(jb,jgrd)
391                 
392                  IF ( ABS( zflag ) == 1. ) THEN  ! northern and southern boundaries
393                     ! one of the two zmsk is always 0 (because of zflag)
394                     zmsk1 = 1._wp - MAX( 0.0_wp, SIGN ( 1.0_wp , - vt_i(ji,jj+1) ) ) ! 0 if no ice
395                     zmsk2 = 1._wp - MAX( 0.0_wp, SIGN ( 1.0_wp , - vt_i(ji,jj-1) ) ) ! 0 if no ice
396                     
397                     ! u_ice = u_ice of the adjacent grid point except if this grid point is ice-free (then u_ice = u_oce)
398                     v_ice (ji,jj) = v_ice(ji,jj+1) * 0.5 * ABS( zflag + 1._wp ) * zmsk1 + &
399                        &            v_ice(ji,jj-1) * 0.5 * ABS( zflag - 1._wp ) * zmsk2 + &
400                        &            v_oce(ji,jj  ) * ( 1._wp - MIN( 1._wp, zmsk1 + zmsk2 ) )
401                  ELSE                             ! everywhere else
402                     !v_ice(ji,jj) = v_oce(ji,jj)
403                     v_ice(ji,jj) = 0._wp
404                  ENDIF
405                  ! mask ice velocities
406                  zinda = 1.0 - MAX( 0.0_wp , SIGN ( 1.0_wp , - at_i(ji,jj) + 0.01 ) ) ! 0 if no ice
407                  v_ice(ji,jj) = zinda * v_ice(ji,jj)
408                 
409               ENDDO
410               
411               CALL lbc_bdy_lnk( v_ice(:,:), 'V', -1., ib_bdy )
412               
413            END SELECT
414           
415         CASE DEFAULT
416            CALL ctl_stop( 'bdy_ice_lim_dyn : unrecognised option for open boundaries for ice fields' )
417         END SELECT
418         
419      ENDDO
420
421      IF( nn_timing == 1 ) CALL timing_stop('bdy_ice_lim_dyn')
422     
423    END SUBROUTINE bdy_ice_lim_dyn
424
425#else
426   !!---------------------------------------------------------------------------------
427   !!   Default option                                                    Empty module
428   !!---------------------------------------------------------------------------------
429CONTAINS
430   SUBROUTINE bdy_ice_lim( kt )      ! Empty routine
431      WRITE(*,*) 'bdy_ice_lim: You should not have seen this print! error?', kt
432   END SUBROUTINE bdy_ice_lim
433#endif
434
435   !!=================================================================================
436END MODULE bdyice_lim
Note: See TracBrowser for help on using the repository browser.