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

source: branches/2017/dev_r8183_ICEMODEL/NEMOGCM/NEMO/OPA_SRC/BDY/bdyice_lim.F90 @ 8306

Last change on this file since 8306 was 8306, checked in by clem, 7 years ago

step1: remove LIM2 from the code

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