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

source: trunk/NEMOGCM/NEMO/OPA_SRC/BDY/bdyice_lim.F90 @ 5143

Last change on this file since 5143 was 5143, checked in by clem, 9 years ago

LIM3: BDY conditions equal to initial state is now working

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#endif
31   USE par_oce         ! ocean parameters
32   USE dom_oce         ! ocean space and time domain variables
33   USE sbc_oce         ! Surface boundary condition: ocean fields
34   USE bdy_oce         ! ocean open boundary conditions
35   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
36   USE in_out_manager  ! write to numout file
37   USE lib_mpp         ! distributed memory computing
38   USE lib_fortran     ! to use key_nosignedzero
39
40   IMPLICIT NONE
41   PRIVATE
42
43   PUBLIC   bdy_ice_lim     ! routine called in sbcmod
44   PUBLIC   bdy_ice_lim_dyn ! routine called in limrhg
45
46   !!----------------------------------------------------------------------
47   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
48   !! $Id: bdyice.F90 2715 2011-03-30 15:58:35Z rblod $
49   !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
50   !!----------------------------------------------------------------------
51CONTAINS
52
53   SUBROUTINE bdy_ice_lim( kt )
54      !!----------------------------------------------------------------------
55      !!                  ***  SUBROUTINE bdy_ice_lim  ***
56      !!
57      !! ** Purpose : - Apply open boundary conditions for ice (LIM2 and LIM3)
58      !!
59      !!----------------------------------------------------------------------
60      INTEGER, INTENT( in ) :: kt     ! Main time step counter
61      INTEGER               :: ib_bdy ! Loop index
62
63      DO ib_bdy=1, nb_bdy
64
65         SELECT CASE( cn_ice_lim(ib_bdy) )
66         CASE('none')
67            CYCLE
68         CASE('frs')
69            CALL bdy_ice_frs( idx_bdy(ib_bdy), dta_bdy(ib_bdy), kt, ib_bdy )
70         CASE DEFAULT
71            CALL ctl_stop( 'bdy_ice_lim : unrecognised option for open boundaries for ice fields' )
72         END SELECT
73
74      END DO
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            IF( nn_ice_lim_dta(ib_bdy) == 0 ) jpbound = 0; ii = ji; ij = jj   ! case ice boundaries = initial conditions
177                                                                              !      do not make state variables dependent on velocity
178               
179
180            rswitch = 1.0 - MAX( 0.0_wp , SIGN ( 1.0_wp , - at_i(ii,ij) + 0.01 ) ) ! 0 if no ice
181
182            ! concentration and thickness
183            a_i (ji,jj,jl) = a_i (ii,ij,jl) * rswitch
184            ht_i(ji,jj,jl) = ht_i(ii,ij,jl) * rswitch
185            ht_s(ji,jj,jl) = ht_s(ii,ij,jl) * rswitch
186
187            ! Ice and snow volumes
188            v_i(ji,jj,jl) = ht_i(ji,jj,jl) * a_i(ji,jj,jl)
189            v_s(ji,jj,jl) = ht_s(ji,jj,jl) * a_i(ji,jj,jl)
190
191            SELECT CASE( jpbound )
192
193            CASE( 0 ) ! velocity is inward
194
195               ! Ice salinity, age, temperature
196               sm_i(ji,jj,jl)   = rswitch * rn_ice_sal(ib_bdy)  + ( 1.0 - rswitch ) * rn_simin
197               o_i(ji,jj,jl)    = rswitch * rn_ice_age(ib_bdy)  + ( 1.0 - rswitch )
198               t_su(ji,jj,jl)   = rswitch * rn_ice_tem(ib_bdy)  + ( 1.0 - rswitch ) * rn_ice_tem(ib_bdy)
199               DO jk = 1, nlay_s
200                  t_s(ji,jj,jk,jl) = rswitch * rn_ice_tem(ib_bdy) + ( 1.0 - rswitch ) * rt0
201               END DO
202               DO jk = 1, nlay_i
203                  t_i(ji,jj,jk,jl) = rswitch * rn_ice_tem(ib_bdy) + ( 1.0 - rswitch ) * rt0 
204                  s_i(ji,jj,jk,jl) = rswitch * rn_ice_sal(ib_bdy) + ( 1.0 - rswitch ) * rn_simin
205               END DO
206               
207            CASE( 1 ) ! velocity is outward
208 
209               ! Ice salinity, age, temperature
210               sm_i(ji,jj,jl)   = rswitch * sm_i(ii,ij,jl)  + ( 1.0 - rswitch ) * rn_simin
211               o_i(ji,jj,jl)    = rswitch * o_i(ii,ij,jl)   + ( 1.0 - rswitch )
212               t_su(ji,jj,jl)   = rswitch * t_su(ii,ij,jl)  + ( 1.0 - rswitch ) * rt0
213               DO jk = 1, nlay_s
214                  t_s(ji,jj,jk,jl) = rswitch * t_s(ii,ij,jk,jl) + ( 1.0 - rswitch ) * rt0
215               END DO
216               DO jk = 1, nlay_i
217                  t_i(ji,jj,jk,jl) = rswitch * t_i(ii,ij,jk,jl) + ( 1.0 - rswitch ) * rt0
218                  s_i(ji,jj,jk,jl) = rswitch * s_i(ii,ij,jk,jl) + ( 1.0 - rswitch ) * rn_simin
219               END DO
220
221            END SELECT
222
223            ! if salinity is constant, then overwrite rn_ice_sal
224            IF( nn_icesal == 1 ) THEN
225               sm_i(ji,jj,jl)   = rn_icesal
226               s_i (ji,jj,:,jl) = rn_icesal
227            ENDIF
228
229            ! contents
230            smv_i(ji,jj,jl)  = MIN( sm_i(ji,jj,jl) , sss_m(ji,jj) ) * v_i(ji,jj,jl)
231            oa_i(ji,jj,jl)   = o_i(ji,jj,jl) * a_i(ji,jj,jl)
232            DO jk = 1, nlay_s
233               ! Snow energy of melting
234               e_s(ji,jj,jk,jl) = rswitch * rhosn * ( cpic * ( rt0 - t_s(ji,jj,jk,jl) ) + lfus )
235               ! Multiply by volume, so that heat content in J/m2
236               e_s(ji,jj,jk,jl) = e_s(ji,jj,jk,jl) * v_s(ji,jj,jl) * r1_nlay_s
237            END DO
238            DO jk = 1, nlay_i
239               ztmelts          = - tmut * s_i(ji,jj,jk,jl) + rt0 !Melting temperature in K                 
240               ! heat content per unit volume
241               e_i(ji,jj,jk,jl) = rswitch * rhoic * &
242                  (   cpic    * ( ztmelts - t_i(ji,jj,jk,jl) ) &
243                  +   lfus    * ( 1.0 - (ztmelts-rt0) / MIN((t_i(ji,jj,jk,jl)-rt0),-epsi20) ) &
244                  - rcp      * ( ztmelts - rt0 ) )
245               ! Mutliply by ice volume, and divide by number of layers to get heat content in J/m2
246               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
247            END DO
248
249         END DO
250 
251         CALL lbc_bdy_lnk(  a_i(:,:,jl), 'T', 1., ib_bdy )
252         CALL lbc_bdy_lnk( ht_i(:,:,jl), 'T', 1., ib_bdy )
253         CALL lbc_bdy_lnk( ht_s(:,:,jl), 'T', 1., ib_bdy )
254         CALL lbc_bdy_lnk(  v_i(:,:,jl), 'T', 1., ib_bdy )
255         CALL lbc_bdy_lnk(  v_s(:,:,jl), 'T', 1., ib_bdy )
256 
257         CALL lbc_bdy_lnk( smv_i(:,:,jl), 'T', 1., ib_bdy )
258         CALL lbc_bdy_lnk(  sm_i(:,:,jl), 'T', 1., ib_bdy )
259         CALL lbc_bdy_lnk(  oa_i(:,:,jl), 'T', 1., ib_bdy )
260         CALL lbc_bdy_lnk(   o_i(:,:,jl), 'T', 1., ib_bdy )
261         CALL lbc_bdy_lnk(  t_su(:,:,jl), 'T', 1., ib_bdy )
262         DO jk = 1, nlay_s
263            CALL lbc_bdy_lnk(t_s(:,:,jk,jl), 'T', 1., ib_bdy )
264            CALL lbc_bdy_lnk(e_s(:,:,jk,jl), 'T', 1., ib_bdy )
265         END DO
266         DO jk = 1, nlay_i
267            CALL lbc_bdy_lnk(t_i(:,:,jk,jl), 'T', 1., ib_bdy )
268            CALL lbc_bdy_lnk(e_i(:,:,jk,jl), 'T', 1., ib_bdy )
269         END DO
270
271      END DO !jl
272   
273#endif
274      !     
275      IF( nn_timing == 1 ) CALL timing_stop('bdy_ice_frs')
276      !
277   END SUBROUTINE bdy_ice_frs
278
279
280   SUBROUTINE bdy_ice_lim_dyn( cd_type )
281      !!------------------------------------------------------------------------------
282      !!                 ***  SUBROUTINE bdy_ice_lim_dyn  ***
283      !!                   
284      !! ** Purpose : Apply dynamics boundary conditions for sea-ice in the cas of unstructured open boundaries.
285      !!              u_ice and v_ice are equal to the value of the adjacent grid point if this latter is not ice free
286      !!              if adjacent grid point is ice free, then u_ice and v_ice are equal to ocean velocities
287      !!
288      !! 2013-06 : C. Rousset
289      !!------------------------------------------------------------------------------
290      !!
291      CHARACTER(len=1), INTENT(in)  ::   cd_type   ! nature of velocity grid-points
292      INTEGER  ::   jb, jgrd           ! dummy loop indices
293      INTEGER  ::   ji, jj             ! local scalar
294      INTEGER  ::   ib_bdy             ! Loop index
295      REAL(wp) ::   zmsk1, zmsk2, zflag
296     !!------------------------------------------------------------------------------
297      !
298      IF( nn_timing == 1 ) CALL timing_start('bdy_ice_lim_dyn')
299      !
300      DO ib_bdy=1, nb_bdy
301         !
302         SELECT CASE( cn_ice_lim(ib_bdy) )
303
304         CASE('none')
305
306            CYCLE
307           
308         CASE('frs')
309           
310            IF( nn_ice_lim_dta(ib_bdy) == 0 ) CYCLE            ! case ice boundaries = initial conditions
311                                                               !      do not change ice velocity (it is only computed by rheology)
312 
313            SELECT CASE ( cd_type )
314               
315            CASE ( 'U' )
316               
317               jgrd = 2      ! u velocity
318               DO jb = 1, idx_bdy(ib_bdy)%nblen(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)%flagu(jb,jgrd)
322                 
323                  IF ( ABS( zflag ) == 1. ) THEN  ! eastern and western 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+1,jj) ) ) ! 0 if no ice
326                     zmsk2 = 1._wp - MAX( 0.0_wp, SIGN ( 1.0_wp , - vt_i(ji-1,jj) ) ) ! 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                     u_ice (ji,jj) = u_ice(ji+1,jj) * 0.5_wp * ABS( zflag + 1._wp ) * zmsk1 + &
330                        &            u_ice(ji-1,jj) * 0.5_wp * ABS( zflag - 1._wp ) * zmsk2 + &
331                        &            u_oce(ji  ,jj) * ( 1._wp - MIN( 1._wp, zmsk1 + zmsk2 ) )
332                  ELSE                             ! everywhere else
333                     !u_ice(ji,jj) = u_oce(ji,jj)
334                     u_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_wp ) ) ! 0 if no ice
338                  u_ice(ji,jj) = rswitch * u_ice(ji,jj)
339                 
340               ENDDO
341               
342               CALL lbc_bdy_lnk( u_ice(:,:), 'U', -1., ib_bdy )
343               
344            CASE ( 'V' )
345               
346               jgrd = 3      ! v velocity
347               DO jb = 1, idx_bdy(ib_bdy)%nblen(jgrd)
348                  ji    = idx_bdy(ib_bdy)%nbi(jb,jgrd)
349                  jj    = idx_bdy(ib_bdy)%nbj(jb,jgrd)
350                  zflag = idx_bdy(ib_bdy)%flagv(jb,jgrd)
351                 
352                  IF ( ABS( zflag ) == 1. ) THEN  ! northern and southern boundaries
353                     ! one of the two zmsk is always 0 (because of zflag)
354                     zmsk1 = 1._wp - MAX( 0.0_wp, SIGN ( 1.0_wp , - vt_i(ji,jj+1) ) ) ! 0 if no ice
355                     zmsk2 = 1._wp - MAX( 0.0_wp, SIGN ( 1.0_wp , - vt_i(ji,jj-1) ) ) ! 0 if no ice
356                     
357                     ! u_ice = u_ice of the adjacent grid point except if this grid point is ice-free (then u_ice = u_oce)
358                     v_ice (ji,jj) = v_ice(ji,jj+1) * 0.5_wp * ABS( zflag + 1._wp ) * zmsk1 + &
359                        &            v_ice(ji,jj-1) * 0.5_wp * ABS( zflag - 1._wp ) * zmsk2 + &
360                        &            v_oce(ji,jj  ) * ( 1._wp - MIN( 1._wp, zmsk1 + zmsk2 ) )
361                  ELSE                             ! everywhere else
362                     !v_ice(ji,jj) = v_oce(ji,jj)
363                     v_ice(ji,jj) = 0._wp
364                  ENDIF
365                  ! mask ice velocities
366                  rswitch = MAX( 0.0_wp , SIGN ( 1.0_wp , at_i(ji,jj) - 0.01 ) ) ! 0 if no ice
367                  v_ice(ji,jj) = rswitch * v_ice(ji,jj)
368                 
369               ENDDO
370               
371               CALL lbc_bdy_lnk( v_ice(:,:), 'V', -1., ib_bdy )
372                 
373            END SELECT
374           
375         CASE DEFAULT
376            CALL ctl_stop( 'bdy_ice_lim_dyn : unrecognised option for open boundaries for ice fields' )
377         END SELECT
378         
379      ENDDO
380
381      IF( nn_timing == 1 ) CALL timing_stop('bdy_ice_lim_dyn')
382     
383    END SUBROUTINE bdy_ice_lim_dyn
384
385#else
386   !!---------------------------------------------------------------------------------
387   !!   Default option                                                    Empty module
388   !!---------------------------------------------------------------------------------
389CONTAINS
390   SUBROUTINE bdy_ice_lim( kt )      ! Empty routine
391      WRITE(*,*) 'bdy_ice_lim: You should not have seen this print! error?', kt
392   END SUBROUTINE bdy_ice_lim
393#endif
394
395   !!=================================================================================
396END MODULE bdyice_lim
Note: See TracBrowser for help on using the repository browser.