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 @ 4689

Last change on this file since 4689 was 4689, checked in by clem, 10 years ago

new version of LIM3 with perfect conservation of heat, see ticket #1352

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