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

source: branches/2014/dev_r4650_UKMO7_STARTHOUR/NEMOGCM/NEMO/OPA_SRC/BDY/bdyice_lim.F90 @ 5075

Last change on this file since 5075 was 5075, checked in by timgraham, 9 years ago

Upgraded branch to current head of trunk (r5072) so it can be used with the trunk

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