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

source: branches/2016/dev_r6522_SIMPLIF_3/NEMOGCM/NEMO/OPA_SRC/BDY/bdyice_lim.F90 @ 6862

Last change on this file since 6862 was 6862, checked in by lovato, 8 years ago

#1729 - trunk: removed key_bdy from the code and set usage of ln_bdy. Tested with SETTE.

  • Property svn:keywords set to Id
File size: 18.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_lim2 || defined key_lim3
11   !!----------------------------------------------------------------------
12   !!   'key_lim2'                                                 LIM-2 sea ice model
13   !!   'key_lim3'                                                 LIM-3 sea ice model
14   !!----------------------------------------------------------------------
15   !!   bdy_ice_lim        : Application of open boundaries to ice
16   !!   bdy_ice_frs        : Application of Flow Relaxation Scheme
17   !!----------------------------------------------------------------------
18   USE timing          ! Timing
19   USE phycst          ! physical constant
20   USE eosbn2          ! equation of state
21   USE oce             ! ocean dynamics and tracers variables
22#if defined key_lim2
23   USE par_ice_2
24   USE ice_2           ! LIM_2 ice variables
25   USE dom_ice_2       ! sea-ice domain
26#elif defined key_lim3
27   USE ice             ! LIM_3 ice variables
28   USE dom_ice         ! sea-ice domain
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      !
62      INTEGER ::   ib_bdy   ! Loop index
63      !!----------------------------------------------------------------------
64      !
65#if defined key_lim3
66      CALL lim_var_glo2eqv
67#endif
68      !
69      DO ib_bdy=1, nb_bdy
70         !
71         SELECT CASE( cn_ice_lim(ib_bdy) )
72         CASE('none')
73            CYCLE
74         CASE('frs')
75            CALL bdy_ice_frs( idx_bdy(ib_bdy), dta_bdy(ib_bdy), kt, ib_bdy )
76         CASE DEFAULT
77            CALL ctl_stop( 'bdy_ice_lim : unrecognised option for open boundaries for ice fields' )
78         END SELECT
79         !
80      END DO
81      !
82#if defined key_lim3
83      CALL lim_var_zapsmall
84      CALL lim_var_agg(1)
85#endif
86      !
87   END SUBROUTINE bdy_ice_lim
88
89
90   SUBROUTINE bdy_ice_frs( idx, dta, kt, ib_bdy )
91      !!------------------------------------------------------------------------------
92      !!                 ***  SUBROUTINE bdy_ice_frs  ***
93      !!                   
94      !! ** Purpose : Apply the Flow Relaxation Scheme for sea-ice fields in the case
95      !!              of unstructured open boundaries.
96      !!
97      !! Reference : Engedahl H., 1995: Use of the flow relaxation scheme in a three-
98      !!             dimensional baroclinic ocean model with realistic topography. Tellus, 365-382.
99      !!------------------------------------------------------------------------------
100      TYPE(OBC_INDEX), INTENT(in) ::   idx     ! OBC indices
101      TYPE(OBC_DATA),  INTENT(in) ::   dta     ! OBC external data
102      INTEGER,         INTENT(in) ::   kt      ! main time-step counter
103      INTEGER,         INTENT(in) ::   ib_bdy  ! BDY set index
104      !
105      INTEGER  ::   jpbound            ! 0 = incoming ice
106      !                                ! 1 = outgoing ice
107      INTEGER  ::   jb, jk, jgrd, jl   ! dummy loop indices
108      INTEGER  ::   ji, jj, ii, ij     ! local scalar
109      REAL(wp) ::   zwgt, zwgt1        ! local scalar
110      REAL(wp) ::   ztmelts, zdh
111#if  defined key_lim2 && ! defined key_lim2_vp && defined key_agrif
112     USE ice_2, vt_s => hsnm
113     USE ice_2, vt_i => hicm
114#endif
115      !!------------------------------------------------------------------------------
116      !
117      IF( nn_timing == 1 )   CALL timing_start('bdy_ice_frs')
118      !
119      jgrd = 1      ! Everything is at T-points here
120      !
121#if defined key_lim2
122      DO jb = 1, idx%nblen(jgrd)
123         ji    = idx%nbi(jb,jgrd)
124         jj    = idx%nbj(jb,jgrd)
125         zwgt  = idx%nbw(jb,jgrd)
126         zwgt1 = 1.e0 - idx%nbw(jb,jgrd)
127         frld (ji,jj) = ( frld (ji,jj) * zwgt1 + dta%frld (jb) * zwgt ) * tmask(ji,jj,1)     ! Leads fraction
128         hicif(ji,jj) = ( hicif(ji,jj) * zwgt1 + dta%hicif(jb) * zwgt ) * tmask(ji,jj,1)     ! Ice depth
129         hsnif(ji,jj) = ( hsnif(ji,jj) * zwgt1 + dta%hsnif(jb) * zwgt ) * tmask(ji,jj,1)     ! Snow depth
130      END DO
131
132      CALL lbc_bdy_lnk( frld,  'T', 1., ib_bdy )                                         ! lateral boundary conditions
133      CALL lbc_bdy_lnk( hicif, 'T', 1., ib_bdy )
134      CALL lbc_bdy_lnk( hsnif, 'T', 1., ib_bdy )
135
136      vt_i(:,:) = hicif(:,:) * frld(:,:)
137      vt_s(:,:) = hsnif(:,:) * frld(:,:)
138      !
139#elif defined key_lim3
140
141      DO jl = 1, jpl
142         DO jb = 1, idx%nblen(jgrd)
143            ji    = idx%nbi(jb,jgrd)
144            jj    = idx%nbj(jb,jgrd)
145            zwgt  = idx%nbw(jb,jgrd)
146            zwgt1 = 1.e0 - idx%nbw(jb,jgrd)
147            a_i (ji,jj,jl) = ( a_i (ji,jj,jl) * zwgt1 + dta%a_i (jb,jl) * zwgt ) * tmask(ji,jj,1)  ! Leads fraction
148            ht_i(ji,jj,jl) = ( ht_i(ji,jj,jl) * zwgt1 + dta%ht_i(jb,jl) * zwgt ) * tmask(ji,jj,1)  ! Ice depth
149            ht_s(ji,jj,jl) = ( ht_s(ji,jj,jl) * zwgt1 + dta%ht_s(jb,jl) * zwgt ) * tmask(ji,jj,1)  ! Snow depth
150
151            ! -----------------
152            ! Pathological case
153            ! -----------------
154            ! In case a) snow load would be in excess or b) ice is coming into a warmer environment that would lead to
155            ! very large transformation from snow to ice (see limthd_dh.F90)
156
157            ! Then, a) transfer the snow excess into the ice (different from limthd_dh)
158            zdh = MAX( 0._wp, ( rhosn * ht_s(ji,jj,jl) + ( rhoic - rau0 ) * ht_i(ji,jj,jl) ) * r1_rau0 )
159            ! Or, b) transfer all the snow into ice (if incoming ice is likely to melt as it comes into a warmer environment)
160            !zdh = MAX( 0._wp, ht_s(ji,jj,jl) * rhosn / rhoic )
161
162            ! recompute ht_i, ht_s
163            ht_i(ji,jj,jl) = MIN( hi_max(jl), ht_i(ji,jj,jl) + zdh )
164            ht_s(ji,jj,jl) = MAX( 0._wp, ht_s(ji,jj,jl) - zdh * rhoic / rhosn ) 
165
166         ENDDO
167         CALL lbc_bdy_lnk(  a_i(:,:,jl), 'T', 1., ib_bdy )
168         CALL lbc_bdy_lnk( ht_i(:,:,jl), 'T', 1., ib_bdy )
169         CALL lbc_bdy_lnk( ht_s(:,:,jl), 'T', 1., ib_bdy )
170      ENDDO
171      ! retrieve at_i
172      at_i(:,:) = 0._wp
173      DO jl = 1, jpl
174         at_i(:,:) = a_i(:,:,jl) + at_i(:,:)
175      END DO
176
177      DO jl = 1, jpl
178         DO jb = 1, idx%nblen(jgrd)
179            ji    = idx%nbi(jb,jgrd)
180            jj    = idx%nbj(jb,jgrd)
181
182            ! condition on ice thickness depends on the ice velocity
183            ! if velocity is outward (strictly), then ice thickness, volume... must be equal to adjacent values
184            jpbound = 0   ;   ii = ji   ;   ij = jj
185            !
186            IF( u_ice(ji+1,jj  ) < 0. .AND. umask(ji-1,jj  ,1) == 0. ) jpbound = 1; ii = ji+1; ij = jj
187            IF( u_ice(ji-1,jj  ) > 0. .AND. umask(ji+1,jj  ,1) == 0. ) jpbound = 1; ii = ji-1; ij = jj
188            IF( v_ice(ji  ,jj+1) < 0. .AND. vmask(ji  ,jj-1,1) == 0. ) jpbound = 1; ii = ji  ; ij = jj+1
189            IF( v_ice(ji  ,jj-1) > 0. .AND. vmask(ji  ,jj+1,1) == 0. ) jpbound = 1; ii = ji  ; ij = jj-1
190            !
191            IF( nn_ice_lim_dta(ib_bdy) == 0 ) jpbound = 0; ii = ji; ij = jj   ! case ice boundaries = initial conditions
192            !                                                                 !      do not make state variables dependent on velocity
193            !
194            rswitch = MAX( 0.0_wp , SIGN ( 1.0_wp , at_i(ii,ij) - 0.01 ) ) ! 0 if no ice
195            !
196            ! concentration and thickness
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
200            !
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            !
205            SELECT CASE( jpbound )
206            !
207            CASE( 0 )   ! velocity is inward
208               !
209               ! Ice salinity, age, temperature
210               sm_i(ji,jj,jl)   = rswitch * rn_ice_sal(ib_bdy)  + ( 1.0 - rswitch ) * rn_simin
211               oa_i(ji,jj,jl)   = rswitch * rn_ice_age(ib_bdy) * a_i(ji,jj,jl)
212               t_su(ji,jj,jl)   = rswitch * rn_ice_tem(ib_bdy)  + ( 1.0 - rswitch ) * rn_ice_tem(ib_bdy)
213               DO jk = 1, nlay_s
214                  t_s(ji,jj,jk,jl) = rswitch * rn_ice_tem(ib_bdy) + ( 1.0 - rswitch ) * rt0
215               END DO
216               DO jk = 1, nlay_i
217                  t_i(ji,jj,jk,jl) = rswitch * rn_ice_tem(ib_bdy) + ( 1.0 - rswitch ) * rt0 
218                  s_i(ji,jj,jk,jl) = rswitch * rn_ice_sal(ib_bdy) + ( 1.0 - rswitch ) * rn_simin
219               END DO
220               !
221            CASE( 1 )   ! velocity is outward
222               !
223               ! Ice salinity, age, temperature
224               sm_i(ji,jj,jl)   = rswitch * sm_i(ii,ij,jl)  + ( 1.0 - rswitch ) * rn_simin
225               oa_i(ji,jj,jl)   = rswitch * oa_i(ii,ij,jl)
226               t_su(ji,jj,jl)   = rswitch * t_su(ii,ij,jl)  + ( 1.0 - rswitch ) * rt0
227               DO jk = 1, nlay_s
228                  t_s(ji,jj,jk,jl) = rswitch * t_s(ii,ij,jk,jl) + ( 1.0 - rswitch ) * rt0
229               END DO
230               DO jk = 1, nlay_i
231                  t_i(ji,jj,jk,jl) = rswitch * t_i(ii,ij,jk,jl) + ( 1.0 - rswitch ) * rt0
232                  s_i(ji,jj,jk,jl) = rswitch * s_i(ii,ij,jk,jl) + ( 1.0 - rswitch ) * rn_simin
233               END DO
234               !
235            END SELECT
236            !
237            IF( nn_icesal == 1 ) THEN     ! constant salinity : overwrite rn_ice_sal
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      CHARACTER(len=1), INTENT(in)  ::   cd_type   ! nature of velocity grid-points
302      !
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            CYCLE
317            !
318         CASE('frs')
319            !
320            IF( nn_ice_lim_dta(ib_bdy) == 0 ) CYCLE            ! case ice boundaries = initial conditions
321            !                                                  !      do not change ice velocity (it is only computed by rheology)
322            SELECT CASE ( cd_type )
323            !     
324            CASE ( 'U' ) 
325               jgrd = 2      ! u velocity
326               DO jb = 1, idx_bdy(ib_bdy)%nblen(jgrd)
327                  ji    = idx_bdy(ib_bdy)%nbi(jb,jgrd)
328                  jj    = idx_bdy(ib_bdy)%nbj(jb,jgrd)
329                  zflag = idx_bdy(ib_bdy)%flagu(jb,jgrd)
330                  !
331                  IF ( ABS( zflag ) == 1. ) THEN  ! eastern and western boundaries
332                     ! one of the two zmsk is always 0 (because of zflag)
333                     zmsk1 = 1._wp - MAX( 0.0_wp, SIGN ( 1.0_wp , - vt_i(ji+1,jj) ) ) ! 0 if no ice
334                     zmsk2 = 1._wp - MAX( 0.0_wp, SIGN ( 1.0_wp , - vt_i(ji-1,jj) ) ) ! 0 if no ice
335                     
336                     ! u_ice = u_ice of the adjacent grid point except if this grid point is ice-free (then u_ice = u_oce)
337                     u_ice (ji,jj) = u_ice(ji+1,jj) * 0.5_wp * ABS( zflag + 1._wp ) * zmsk1 + &
338                        &            u_ice(ji-1,jj) * 0.5_wp * ABS( zflag - 1._wp ) * zmsk2 + &
339                        &            u_oce(ji  ,jj) * ( 1._wp - MIN( 1._wp, zmsk1 + zmsk2 ) )
340                  ELSE                             ! everywhere else
341                     !u_ice(ji,jj) = u_oce(ji,jj)
342                     u_ice(ji,jj) = 0._wp
343                  ENDIF
344                  ! mask ice velocities
345                  rswitch = MAX( 0.0_wp , SIGN ( 1.0_wp , at_i(ji,jj) - 0.01_wp ) ) ! 0 if no ice
346                  u_ice(ji,jj) = rswitch * u_ice(ji,jj)
347                  !
348               END DO
349               CALL lbc_bdy_lnk( u_ice(:,:), 'U', -1., ib_bdy )
350               !
351            CASE ( 'V' )
352               jgrd = 3      ! v velocity
353               DO jb = 1, idx_bdy(ib_bdy)%nblen(jgrd)
354                  ji    = idx_bdy(ib_bdy)%nbi(jb,jgrd)
355                  jj    = idx_bdy(ib_bdy)%nbj(jb,jgrd)
356                  zflag = idx_bdy(ib_bdy)%flagv(jb,jgrd)
357                  !
358                  IF ( ABS( zflag ) == 1. ) THEN  ! northern and southern boundaries
359                     ! one of the two zmsk is always 0 (because of zflag)
360                     zmsk1 = 1._wp - MAX( 0.0_wp, SIGN ( 1.0_wp , - vt_i(ji,jj+1) ) ) ! 0 if no ice
361                     zmsk2 = 1._wp - MAX( 0.0_wp, SIGN ( 1.0_wp , - vt_i(ji,jj-1) ) ) ! 0 if no ice
362                     
363                     ! u_ice = u_ice of the adjacent grid point except if this grid point is ice-free (then u_ice = u_oce)
364                     v_ice (ji,jj) = v_ice(ji,jj+1) * 0.5_wp * ABS( zflag + 1._wp ) * zmsk1 + &
365                        &            v_ice(ji,jj-1) * 0.5_wp * ABS( zflag - 1._wp ) * zmsk2 + &
366                        &            v_oce(ji,jj  ) * ( 1._wp - MIN( 1._wp, zmsk1 + zmsk2 ) )
367                  ELSE                             ! everywhere else
368                     !v_ice(ji,jj) = v_oce(ji,jj)
369                     v_ice(ji,jj) = 0._wp
370                  ENDIF
371                  ! mask ice velocities
372                  rswitch = MAX( 0.0_wp , SIGN ( 1.0_wp , at_i(ji,jj) - 0.01 ) ) ! 0 if no ice
373                  v_ice(ji,jj) = rswitch * v_ice(ji,jj)
374                  !
375               END DO
376               CALL lbc_bdy_lnk( v_ice(:,:), 'V', -1., ib_bdy )
377               !
378            END SELECT
379            !
380         CASE DEFAULT
381            CALL ctl_stop( 'bdy_ice_lim_dyn : unrecognised option for open boundaries for ice fields' )
382         END SELECT
383         !
384      END DO
385      !
386      IF( nn_timing == 1 ) CALL timing_stop('bdy_ice_lim_dyn')
387      !
388    END SUBROUTINE bdy_ice_lim_dyn
389
390#else
391   !!---------------------------------------------------------------------------------
392   !!   Default option                                                    Empty module
393   !!---------------------------------------------------------------------------------
394CONTAINS
395   SUBROUTINE bdy_ice_lim( kt )      ! Empty routine
396      WRITE(*,*) 'bdy_ice_lim: You should not have seen this print! error?', kt
397   END SUBROUTINE bdy_ice_lim
398#endif
399
400   !!=================================================================================
401END MODULE bdyice_lim
Note: See TracBrowser for help on using the repository browser.