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 branches/UKMO/dev_r5518_AMM15_package/NEMOGCM/NEMO/OPA_SRC/BDY – NEMO

source: branches/UKMO/dev_r5518_AMM15_package/NEMOGCM/NEMO/OPA_SRC/BDY/bdyice_lim.F90 @ 10248

Last change on this file since 10248 was 10248, checked in by kingr, 6 years ago

Merged 2015/nemo_v3_6_STABLE@6232

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