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 branches/UKMO/dev_r8126_LIM3_couple/NEMOGCM/NEMO/OPA_SRC/BDY – NEMO

source: branches/UKMO/dev_r8126_LIM3_couple/NEMOGCM/NEMO/OPA_SRC/BDY/bdyice.F90 @ 8879

Last change on this file since 8879 was 8879, checked in by frrh, 6 years ago

Merge in http://fcm3/projects/NEMO.xm/log/branches/UKMO/dev_r8183_ICEMODEL_svn_removed
revisions 8738:8847 inclusive.

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