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/2016/dev_v3_6_STABLE_r6506_AGRIF_LIM3/NEMOGCM/NEMO/OPA_SRC/BDY – NEMO

source: branches/2016/dev_v3_6_STABLE_r6506_AGRIF_LIM3/NEMOGCM/NEMO/OPA_SRC/BDY/bdyice_lim.F90 @ 6746

Last change on this file since 6746 was 6746, checked in by clem, 8 years ago

landfast ice parameterization + update from trunk + removing useless dom_ice.F90 and limmsh.F90 and limwri_dimg.h90

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