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

source: branches/UKMO/dev_r8183_ICEMODEL_svn_removed/NEMOGCM/NEMO/OPA_SRC/BDY/bdyice.F90 @ 8738

Last change on this file since 8738 was 8738, checked in by dancopsey, 6 years ago

Merged in main ICEMODEL branch (branches/2017/dev_r8183_ICEMODEL) up to revision 8588

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