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 @ 5985

Last change on this file since 5985 was 5985, checked in by timgraham, 8 years ago

Reinstate keywords before upgrading to head of trunk

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