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.F90 in NEMO/branches/2019/dev_r10984_HPC-13_IRRMANN_BDY_optimization/src/OCE/BDY – NEMO

source: NEMO/branches/2019/dev_r10984_HPC-13_IRRMANN_BDY_optimization/src/OCE/BDY/bdyice.F90 @ 11049

Last change on this file since 11049 was 11049, checked in by girrmann, 5 years ago

dev_r10984_HPC-13 : CYCLE instruction is not systematic anymore, computation is done on the halo whenever possible and overwritten by lbc_bdy instruction, see #2285

  • Property svn:keywords set to Id
File size: 20.0 KB
Line 
1MODULE 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   !!----------------------------------------------------------------------
46CONTAINS
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   ! BDY set index
58      !!----------------------------------------------------------------------
59      ! controls
60      IF( ln_timing    )   CALL timing_start('bdy_ice_thd')                                                            ! timing
61      IF( ln_icediachk )   CALL ice_cons_hsm(0,'bdy_ice_thd', rdiag_v, rdiag_s, rdiag_t, rdiag_fv, rdiag_fs, rdiag_ft) ! conservation
62      !
63      CALL ice_var_glo2eqv
64      !
65      DO jbdy = 1, nb_bdy
66         !
67         SELECT CASE( cn_ice(jbdy) )
68         CASE('none')   ;   CYCLE
69         CASE('frs' )   ;   CALL bdy_ice_frs( idx_bdy(jbdy), dta_bdy(jbdy), kt, jbdy )
70         CASE DEFAULT
71            CALL ctl_stop( 'bdy_ice : unrecognised option for open boundaries for ice fields' )
72         END SELECT
73         !
74      END DO
75      !
76      CALL ice_cor( kt , 0 )      ! -- In case categories are out of bounds, do a remapping
77      !                           !    i.e. inputs have not the same ice thickness distribution (set by rn_himean)
78      !                           !         than the regional simulation
79      CALL ice_var_agg(1)
80      !
81      ! controls
82      IF( ln_icediachk )   CALL ice_cons_hsm(1,'bdy_ice_thd', rdiag_v, rdiag_s, rdiag_t, rdiag_fv, rdiag_fs, rdiag_ft) ! conservation
83      IF( ln_icectl    )   CALL ice_prt     ( kt, iiceprt, jiceprt, 1, ' - ice thermo bdy - ' )                        ! prints
84      IF( ln_timing    )   CALL timing_stop ('bdy_ice_thd')                                                            ! timing
85      !
86   END SUBROUTINE bdy_ice
87
88
89   SUBROUTINE bdy_ice_frs( idx, dta, kt, jbdy )
90      !!------------------------------------------------------------------------------
91      !!                 ***  SUBROUTINE bdy_ice_frs  ***
92      !!                   
93      !! ** Purpose : Apply the Flow Relaxation Scheme for sea-ice fields
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
101      INTEGER,         INTENT(in) ::   jbdy    ! BDY set index
102      !
103      INTEGER  ::   jpbound            ! 0 = incoming ice
104      !                                ! 1 = outgoing ice
105      INTEGER  ::   i_bdy, jgrd        ! dummy loop indices
106      INTEGER  ::   ji, jj, jk, jl, ib, jb
107      REAL(wp) ::   zwgt, zwgt1        ! local scalar
108      REAL(wp) ::   ztmelts, zdh
109      REAL(wp), POINTER  :: flagu, flagv              ! short cuts
110      !!------------------------------------------------------------------------------
111      !
112      jgrd = 1      ! Everything is at T-points here
113      !
114      DO jl = 1, jpl
115         DO i_bdy = 1, idx%nblenrim(jgrd)
116            ji    = idx%nbi(i_bdy,jgrd)
117            jj    = idx%nbj(i_bdy,jgrd)
118            zwgt  = idx%nbw(i_bdy,jgrd)
119            zwgt1 = 1.e0 - idx%nbw(i_bdy,jgrd)
120            a_i(ji,jj,jl) = ( a_i(ji,jj,jl) * zwgt1 + dta%a_i(i_bdy,jl) * zwgt ) * tmask(ji,jj,1)  ! Leads fraction
121            h_i(ji,jj,jl) = ( h_i(ji,jj,jl) * zwgt1 + dta%h_i(i_bdy,jl) * zwgt ) * tmask(ji,jj,1)  ! Ice depth
122            h_s(ji,jj,jl) = ( h_s(ji,jj,jl) * zwgt1 + dta%h_s(i_bdy,jl) * zwgt ) * tmask(ji,jj,1)  ! Snow depth
123
124            ! -----------------
125            ! Pathological case
126            ! -----------------
127            ! In case a) snow load would be in excess or b) ice is coming into a warmer environment that would lead to
128            ! very large transformation from snow to ice (see icethd_dh.F90)
129
130            ! Then, a) transfer the snow excess into the ice (different from icethd_dh)
131            zdh = MAX( 0._wp, ( rhos * h_s(ji,jj,jl) + ( rhoi - rau0 ) * h_i(ji,jj,jl) ) * r1_rau0 )
132            ! Or, b) transfer all the snow into ice (if incoming ice is likely to melt as it comes into a warmer environment)
133            !zdh = MAX( 0._wp, h_s(ji,jj,jl) * rhos / rhoi )
134
135            ! recompute h_i, h_s
136            h_i(ji,jj,jl) = MIN( hi_max(jl), h_i(ji,jj,jl) + zdh )
137            h_s(ji,jj,jl) = MAX( 0._wp, h_s(ji,jj,jl) - zdh * rhoi / rhos ) 
138
139         ENDDO
140      ENDDO
141      CALL lbc_bdy_lnk( 'bdyice', a_i(:,:,:), 'T', 1., jbdy )
142      CALL lbc_bdy_lnk( 'bdyice', h_i(:,:,:), 'T', 1., jbdy )
143      CALL lbc_bdy_lnk( 'bdyice', h_s(:,:,:), 'T', 1., jbdy )
144
145      DO jl = 1, jpl
146         DO i_bdy = 1, idx%nblenrim(jgrd)
147            ji = idx%nbi(i_bdy,jgrd)
148            jj = idx%nbj(i_bdy,jgrd)
149            flagu => idx%flagu(i_bdy,jgrd)
150            flagv => idx%flagv(i_bdy,jgrd)
151            ! condition on ice thickness depends on the ice velocity
152            ! if velocity is outward (strictly), then ice thickness, volume... must be equal to adjacent values
153            jpbound = 0   ;   ib = ji   ;   jb = jj
154            !
155            IF( flagu ==  1. )   THEN
156               IF( ji+1 > jpi  )   CYCLE
157               IF( u_ice(ji  ,jj  ) < 0. )   jpbound = 1 ; ib = ji+1
158            END IF
159            IF( flagu == -1. )   THEN
160               IF( ji-1 < 1    )   CYCLE
161               IF( u_ice(ji-1,jj  ) < 0. )   jpbound = 1 ; ib = ji-1
162            END IF
163            IF( flagv ==  1. )   THEN
164               IF( ji+1 > jpj )   CYCLE
165               IF( v_ice(ji  ,jj  ) < 0. )   jpbound = 1 ; jb = jj+1
166            END IF
167            IF( flagv == -1. )   THEN
168               IF( jj-1 < 1   )   CYCLE
169               IF( v_ice(ji  ,jj-1) < 0. )   jpbound = 1 ; jb = jj-1
170            END IF
171            !
172            IF( nn_ice_dta(jbdy) == 0 )   jpbound = 0 ; ib = ji ; jb = jj   ! case ice boundaries = initial conditions
173            !                                                               !      do not make state variables dependent on velocity
174            !
175            IF( a_i(ib,jb,jl) > 0._wp ) THEN   ! there is ice at the boundary
176               !
177               a_i(ji,jj,jl) = a_i(ib,jb,jl) ! concentration
178               h_i(ji,jj,jl) = h_i(ib,jb,jl) ! thickness ice
179               h_s(ji,jj,jl) = h_s(ib,jb,jl) ! thickness snw
180               !
181               SELECT CASE( jpbound )
182                  !
183               CASE( 0 )   ! velocity is inward
184                  !
185                  oa_i(ji,jj,  jl) = rn_ice_age(jbdy) * a_i(ji,jj,jl) ! age
186                  a_ip(ji,jj,  jl) = 0._wp                            ! pond concentration
187                  v_ip(ji,jj,  jl) = 0._wp                            ! pond volume
188                  t_su(ji,jj,  jl) = rn_ice_tem(jbdy)                 ! temperature surface
189                  t_s (ji,jj,:,jl) = rn_ice_tem(jbdy)                 ! temperature snw
190                  t_i (ji,jj,:,jl) = rn_ice_tem(jbdy)                 ! temperature ice
191                  s_i (ji,jj,  jl) = rn_ice_sal(jbdy)                 ! salinity
192                  sz_i(ji,jj,:,jl) = rn_ice_sal(jbdy)                 ! salinity profile
193                  !
194               CASE( 1 )   ! velocity is outward
195                  !
196                  oa_i(ji,jj,  jl) = oa_i(ib,jb,  jl) ! age
197                  a_ip(ji,jj,  jl) = a_ip(ib,jb,  jl) ! pond concentration
198                  v_ip(ji,jj,  jl) = v_ip(ib,jb,  jl) ! pond volume
199                  t_su(ji,jj,  jl) = t_su(ib,jb,  jl) ! temperature surface
200                  t_s (ji,jj,:,jl) = t_s (ib,jb,:,jl) ! temperature snw
201                  t_i (ji,jj,:,jl) = t_i (ib,jb,:,jl) ! temperature ice
202                  s_i (ji,jj,  jl) = s_i (ib,jb,  jl) ! salinity
203                  sz_i(ji,jj,:,jl) = sz_i(ib,jb,:,jl) ! salinity profile
204                  !
205               END SELECT
206               !
207               IF( nn_icesal == 1 ) THEN     ! if constant salinity
208                  s_i (ji,jj  ,jl) = rn_icesal
209                  sz_i(ji,jj,:,jl) = rn_icesal
210               ENDIF
211               !
212               ! global fields
213               v_i (ji,jj,jl) = h_i(ji,jj,jl) * a_i(ji,jj,jl)                       ! volume ice
214               v_s (ji,jj,jl) = h_s(ji,jj,jl) * a_i(ji,jj,jl)                       ! volume snw
215               sv_i(ji,jj,jl) = MIN( s_i(ji,jj,jl) , sss_m(ji,jj) ) * v_i(ji,jj,jl) ! salt content
216               DO jk = 1, nlay_s
217                  e_s(ji,jj,jk,jl) = rhos * ( rcpi * ( rt0 - t_s(ji,jj,jk,jl) ) + rLfus )   ! enthalpy in J/m3
218                  e_s(ji,jj,jk,jl) = e_s(ji,jj,jk,jl) * v_s(ji,jj,jl) * r1_nlay_s           ! enthalpy in J/m2
219               END DO               
220               DO jk = 1, nlay_i
221                  ztmelts          = - rTmlt  * sz_i(ji,jj,jk,jl)             ! Melting temperature in C
222                  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
223                  !
224                  e_i(ji,jj,jk,jl) = rhoi * ( rcpi  * ( ztmelts - ( t_i(ji,jj,jk,jl) - rt0 ) )           &   ! enthalpy in J/m3
225                     &                      + rLfus * ( 1._wp - ztmelts / ( t_i(ji,jj,jk,jl) - rt0 ) )   &
226                     &                      - rcp   *   ztmelts )                 
227                  e_i(ji,jj,jk,jl) = e_i(ji,jj,jk,jl) * v_i(ji,jj,jl) * r1_nlay_i                            ! enthalpy in J/m2
228               END DO
229               !
230            ELSE   ! no ice at the boundary
231               !
232               a_i (ji,jj,  jl) = 0._wp
233               h_i (ji,jj,  jl) = 0._wp
234               h_s (ji,jj,  jl) = 0._wp
235               oa_i(ji,jj,  jl) = 0._wp
236               a_ip(ji,jj,  jl) = 0._wp
237               v_ip(ji,jj,  jl) = 0._wp
238               t_su(ji,jj,  jl) = rt0
239               t_s (ji,jj,:,jl) = rt0
240               t_i (ji,jj,:,jl) = rt0 
241               
242               IF( nn_icesal == 1 ) THEN     ! if constant salinity
243                  s_i (ji,jj  ,jl) = rn_icesal
244                  sz_i(ji,jj,:,jl) = rn_icesal
245               ELSE                          ! if variable salinity
246                  s_i (ji,jj,jl)   = rn_simin
247                  sz_i(ji,jj,:,jl) = rn_simin
248               ENDIF
249               !
250               ! global fields
251               v_i (ji,jj,  jl) = 0._wp
252               v_s (ji,jj,  jl) = 0._wp
253               sv_i(ji,jj,  jl) = 0._wp
254               e_s (ji,jj,:,jl) = 0._wp
255               e_i (ji,jj,:,jl) = 0._wp
256
257            ENDIF
258                       
259         END DO
260         !
261      END DO ! jl
262
263      CALL lbc_bdy_lnk( 'bdyice', a_i (:,:,:)  , 'T', 1., jbdy )
264      CALL lbc_bdy_lnk( 'bdyice', h_i (:,:,:)  , 'T', 1., jbdy )
265      CALL lbc_bdy_lnk( 'bdyice', h_s (:,:,:)  , 'T', 1., jbdy )
266      CALL lbc_bdy_lnk( 'bdyice', oa_i(:,:,:)  , 'T', 1., jbdy )
267      CALL lbc_bdy_lnk( 'bdyice', a_ip(:,:,:)  , 'T', 1., jbdy )
268      CALL lbc_bdy_lnk( 'bdyice', v_ip(:,:,:)  , 'T', 1., jbdy )
269      CALL lbc_bdy_lnk( 'bdyice', s_i (:,:,:)  , 'T', 1., jbdy )
270      CALL lbc_bdy_lnk( 'bdyice', t_su(:,:,:)  , 'T', 1., jbdy )
271      CALL lbc_bdy_lnk( 'bdyice', v_i (:,:,:)  , 'T', 1., jbdy )
272      CALL lbc_bdy_lnk( 'bdyice', v_s (:,:,:)  , 'T', 1., jbdy )
273      CALL lbc_bdy_lnk( 'bdyice', sv_i(:,:,:)  , 'T', 1., jbdy )
274      CALL lbc_bdy_lnk( 'bdyice', t_s (:,:,:,:), 'T', 1., jbdy )
275      CALL lbc_bdy_lnk( 'bdyice', e_s (:,:,:,:), 'T', 1., jbdy )
276      CALL lbc_bdy_lnk( 'bdyice', t_i (:,:,:,:), 'T', 1., jbdy )
277      CALL lbc_bdy_lnk( 'bdyice', e_i (:,:,:,:), 'T', 1., jbdy )
278      !     
279   END SUBROUTINE bdy_ice_frs
280
281
282   SUBROUTINE bdy_ice_dyn( cd_type )
283      !!------------------------------------------------------------------------------
284      !!                 ***  SUBROUTINE bdy_ice_dyn  ***
285      !!                   
286      !! ** Purpose : Apply dynamics boundary conditions for sea-ice.
287      !!
288      !! ** Method :  if this adjacent grid point is not ice free, then u_ice and v_ice take its value
289      !!              if                          is     ice free, then u_ice and v_ice are unchanged by BDY
290      !!                                                           they keep values calculated in rheology
291      !!
292      !!------------------------------------------------------------------------------
293      CHARACTER(len=1), INTENT(in)  ::   cd_type   ! nature of velocity grid-points
294      !
295      INTEGER  ::   i_bdy, jgrd      ! dummy loop indices
296      INTEGER  ::   ji, jj           ! local scalar
297      INTEGER  ::   jbdy             ! BDY set index
298      REAL(wp) ::   zmsk1, zmsk2, zflag
299      !!------------------------------------------------------------------------------
300      IF( ln_timing )   CALL timing_start('bdy_ice_dyn')
301      !
302      DO jbdy=1, nb_bdy
303         !
304         SELECT CASE( cn_ice(jbdy) )
305         !
306         CASE('none')
307            CYCLE
308            !
309         CASE('frs')
310            !
311            IF( nn_ice_dta(jbdy) == 0 ) CYCLE            ! case ice boundaries = initial conditions
312            !                                            !      do not change ice velocity (it is only computed by rheology)
313            SELECT CASE ( cd_type )
314            !     
315            CASE ( 'U' ) 
316               jgrd = 2      ! u velocity
317               DO i_bdy = 1, idx_bdy(jbdy)%nblenrim(jgrd)
318                  ji    = idx_bdy(jbdy)%nbi(i_bdy,jgrd)
319                  jj    = idx_bdy(jbdy)%nbj(i_bdy,jgrd)
320                  zflag = idx_bdy(jbdy)%flagu(i_bdy,jgrd)
321                  !     i-1  i   i    |  !        i  i i+1 |  !          i  i i+1 |
322                  !      >  ice  >    |  !        o  > ice |  !          o  >  o  |     
323                  ! => set at u_ice(i-1) !  => set to O       !  => unchanged
324                  IF( zflag == -1. .AND. ji > 1 .AND. ji < jpi )   THEN 
325                     IF    ( vt_i(ji  ,jj) > 0. )   THEN   ;   u_ice(ji,jj) = u_ice(ji-1,jj) 
326                     ELSEIF( vt_i(ji+1,jj) > 0. )   THEN   ;   u_ice(ji,jj) = 0._wp
327                     END IF
328                  END IF
329                  ! |    i  i+1 i+1        !  |  i   i i+1        !  | i  i i+1
330                  ! |    >  ice  >         !  | ice  >  o         !  | o  >  o   
331                  ! => set at u_ice(i+1)   !     => set to O      !     =>  unchanged
332                  IF( zflag ==  1. .AND. ji+1 < jpi+1 )   THEN
333                     IF    ( vt_i(ji+1,jj) > 0. )   THEN   ;   u_ice(ji,jj) = u_ice(ji+1,jj)
334                     ELSEIF( vt_i(ji  ,jj) > 0. )   THEN   ;   u_ice(ji,jj) = 0._wp
335                     END IF
336                  END IF
337                  !
338                  IF( zflag ==  0. )   u_ice(ji,jj) = 0._wp   ! u_ice = 0  if north/south bdy 
339                  !
340               END DO
341               CALL lbc_bdy_lnk( 'bdyice', u_ice(:,:), 'U', -1., jbdy )
342               !
343            CASE ( 'V' )
344               jgrd = 3      ! v velocity
345               DO i_bdy = 1, idx_bdy(jbdy)%nblenrim(jgrd)
346                  ji    = idx_bdy(jbdy)%nbi(i_bdy,jgrd)
347                  jj    = idx_bdy(jbdy)%nbj(i_bdy,jgrd)
348                  zflag = idx_bdy(jbdy)%flagv(i_bdy,jgrd)
349                  !    ¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨¨     !  ¨¨¨¨ïce¨¨¨(jj+1)¨¨     ! ¨¨¨¨¨¨ö¨¨¨¨(jj+1)       
350                  !       ^    (jj  )       !       ^    (jj  )       !       ^    (jj  )       
351                  !      ice   (jj  )       !       o    (jj  )       !       o    (jj  )       
352                  !       ^    (jj-1)       !                         !
353                  ! => set to u_ice(jj-1)   !  =>   set to 0          !   => unchanged       
354                  IF( zflag == -1. .AND. jj > 1 .AND. jj < jpj )   THEN                 
355                    IF    ( vt_i(ji,jj  ) > 0. )   THEN   ;   v_ice(ji,jj) = v_ice(ji,jj-1)
356                    ELSEIF( vt_i(ji,jj+1) > 0. )   THEN   ;   v_ice(ji,jj) = 0._wp
357                    END IF
358                  END IF 
359                  !       ^    (jj+1)       !                         !             
360                  !      ice   (jj+1)       !       o    (jj+1)       !       o    (jj+1)       
361                  !       ^    (jj  )       !       ^    (jj  )       !       ^    (jj  )
362                  !   ________________      !  ____ice___(jj  )_      !  _____o____(jj  )
363                  ! => set to u_ice(jj+1)   !    => set to 0          !    => unchanged 
364                  IF( zflag ==  1. .AND. jj < jpj )   THEN             
365                     IF    ( vt_i(ji,jj+1) > 0. )   THEN   ;   v_ice(ji,jj) = v_ice(ji,jj+1)
366                     ELSEIF( vt_i(ji,jj  ) > 0. )   THEN   ;   v_ice(ji,jj) = 0._wp
367                     END IF
368                  END IF                                         
369                  !
370                  IF( zflag ==  0. )   v_ice(ji,jj) = 0._wp   ! v_ice = 0  if west/east bdy 
371                  !
372               END DO
373               CALL lbc_bdy_lnk( 'bdyice', v_ice(:,:), 'V', -1., jbdy )
374               !
375            END SELECT
376            !
377         CASE DEFAULT
378            CALL ctl_stop( 'bdy_ice_dyn : unrecognised option for open boundaries for ice fields' )
379         END SELECT
380         !
381      END DO
382      !
383      IF( ln_timing )   CALL timing_stop('bdy_ice_dyn')
384      !
385    END SUBROUTINE bdy_ice_dyn
386
387#else
388   !!---------------------------------------------------------------------------------
389   !!   Default option                                                    Empty module
390   !!---------------------------------------------------------------------------------
391CONTAINS
392   SUBROUTINE bdy_ice( kt )      ! Empty routine
393      IMPLICIT NONE
394      INTEGER, INTENT( in ) :: kt
395      WRITE(*,*) 'bdy_ice: You should not have seen this print! error?', kt
396   END SUBROUTINE bdy_ice
397#endif
398
399   !!=================================================================================
400END MODULE bdyice
Note: See TracBrowser for help on using the repository browser.