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/trunk/src/OCE/BDY – NEMO

source: NEMO/trunk/src/OCE/BDY/bdyice.F90 @ 9656

Last change on this file since 9656 was 9656, checked in by clem, 6 years ago

remove the remaining references to LIM

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