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 NEMO/branches/UKMO/NEMO_4.0.1_fix_cpl/src/OCE/BDY – NEMO

source: NEMO/branches/UKMO/NEMO_4.0.1_fix_cpl/src/OCE/BDY/bdyice.F90 @ 12402

Last change on this file since 12402 was 12402, checked in by dancopsey, 4 years ago

Add melt pond lid code.

File size: 25.6 KB
RevLine 
[8586]1MODULE bdyice
2   !!======================================================================
3   !!                       ***  MODULE  bdyice  ***
[9656]4   !! Unstructured Open Boundary Cond. :  Open boundary conditions for sea-ice (SI3)
[8586]5   !!======================================================================
6   !!  History :  3.3  !  2010-09 (D. Storkey)  Original code
[9656]7   !!             3.4  !  2012-01 (C. Rousset)  add new sea ice model
8   !!             4.0  !  2018    (C. Rousset)  SI3 compatibility
[8586]9   !!----------------------------------------------------------------------
[9570]10#if defined key_si3
[8586]11   !!----------------------------------------------------------------------
[9656]12   !!   'key_si3'                                          SI3 sea ice model
[8586]13   !!----------------------------------------------------------------------
14   !!   bdy_ice        : Application of open boundaries to ice
15   !!   bdy_ice_frs    : Application of Flow Relaxation Scheme
16   !!----------------------------------------------------------------------
17   USE oce             ! ocean dynamics and tracers variables
[8637]18   USE ice             ! sea-ice: variables
19   USE icevar          ! sea-ice: operations
[9880]20   USE icecor          ! sea-ice: corrections
[8637]21   USE icectl          ! sea-ice: control prints
[8586]22   USE phycst          ! physical constant
23   USE eosbn2          ! equation of state
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   !
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   USE timing          ! Timing
34
35   IMPLICIT NONE
36   PRIVATE
37
38   PUBLIC   bdy_ice     ! routine called in sbcmod
[9656]39   PUBLIC   bdy_ice_dyn ! routine called in icedyn_rhg_evp
[8586]40
41   !!----------------------------------------------------------------------
[9598]42   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
[10069]43   !! $Id$
[10068]44   !! Software governed by the CeCILL license (see ./LICENSE)
[8586]45   !!----------------------------------------------------------------------
46CONTAINS
47
48   SUBROUTINE bdy_ice( kt )
49      !!----------------------------------------------------------------------
50      !!                  ***  SUBROUTINE bdy_ice  ***
51      !!
[9890]52      !! ** Purpose : Apply open boundary conditions for sea ice
[8586]53      !!
54      !!----------------------------------------------------------------------
55      INTEGER, INTENT(in) ::   kt   ! Main time step counter
56      !
[11536]57      INTEGER ::   jbdy, ir                             ! BDY set index, rim index
58      INTEGER ::   ibeg, iend                           ! length of rim to be treated (rim 0 or rim 1)
59      LOGICAL ::   llrim0                               ! indicate if rim 0 is treated
60      LOGICAL, DIMENSION(4)  :: llsend1, llrecv1        ! indicate how communications are to be carried out
[8586]61      !!----------------------------------------------------------------------
[11041]62      ! controls
63      IF( ln_timing    )   CALL timing_start('bdy_ice_thd')                                                            ! timing
64      IF( ln_icediachk )   CALL ice_cons_hsm(0,'bdy_ice_thd', rdiag_v, rdiag_s, rdiag_t, rdiag_fv, rdiag_fs, rdiag_ft) ! conservation
[11536]65      IF( ln_icediachk )   CALL ice_cons2D  (0,'bdy_ice_thd',  diag_v,  diag_s,  diag_t,  diag_fv,  diag_fs,  diag_ft) ! conservation
[8586]66      !
67      CALL ice_var_glo2eqv
68      !
[11536]69      llsend1(:) = .false.   ;   llrecv1(:) = .false.
70      DO ir = 1, 0, -1   ! treat rim 1 before rim 0
71         IF( ir == 0 ) THEN   ;   llrim0 = .TRUE.
72         ELSE                 ;   llrim0 = .FALSE.
73         END IF
74         DO jbdy = 1, nb_bdy
75            !
76            SELECT CASE( cn_ice(jbdy) )
77            CASE('none')   ;   CYCLE
78            CASE('frs' )   ;   CALL bdy_ice_frs( idx_bdy(jbdy), dta_bdy(jbdy), kt, jbdy, llrim0 )
79            CASE DEFAULT
80               CALL ctl_stop( 'bdy_ice : unrecognised option for open boundaries for ice fields' )
81            END SELECT
82            !
83         END DO
[8586]84         !
[11536]85         ! Update bdy points       
86         IF( nn_hls > 1 .AND. ir == 1 ) CYCLE   ! at least 2 halos will be corrected -> no need to correct rim 1 before rim 0
87         IF( nn_hls == 1 ) THEN   ;   llsend1(:) = .false.   ;   llrecv1(:) = .false.   ;   END IF
88         DO jbdy = 1, nb_bdy
89            IF( cn_ice(jbdy) == 'frs' ) THEN
90               llsend1(:) = llsend1(:) .OR. lsend_bdyint(jbdy,1,:,ir)   ! possibly every direction, T points
91               llrecv1(:) = llrecv1(:) .OR. lrecv_bdyint(jbdy,1,:,ir)   ! possibly every direction, T points
92            END IF
93         END DO   ! jbdy
94         IF( ANY(llsend1) .OR. ANY(llrecv1) ) THEN   ! if need to send/recv in at least one direction
95            ! exchange 3d arrays
96            CALL lbc_lnk_multi( 'bdyice', a_i , 'T', 1., h_i , 'T', 1., h_s , 'T', 1., oa_i, 'T', 1. &
[12402]97                 &                      , a_ip, 'T', 1., v_ip, 'T', 1., lh_ip,'T', 1., s_i , 'T', 1., t_su, 'T', 1. &
[11536]98                 &                      , v_i , 'T', 1., v_s , 'T', 1., sv_i, 'T', 1.                &
99                 &                      , kfillmode=jpfillnothing ,lsend=llsend1, lrecv=llrecv1      )
100            ! exchange 4d arrays :   third dimension = 1   and then   third dimension = jpk
101            CALL lbc_lnk_multi( 'bdyice', t_s , 'T', 1., e_s , 'T', 1., kfillmode=jpfillnothing ,lsend=llsend1, lrecv=llrecv1 )
102            CALL lbc_lnk_multi( 'bdyice', t_i , 'T', 1., e_i , 'T', 1., kfillmode=jpfillnothing ,lsend=llsend1, lrecv=llrecv1 )
103         END IF
104      END DO   ! ir
[8586]105      !
[9880]106      CALL ice_cor( kt , 0 )      ! -- In case categories are out of bounds, do a remapping
[9885]107      !                           !    i.e. inputs have not the same ice thickness distribution (set by rn_himean)
108      !                           !         than the regional simulation
[9124]109      CALL ice_var_agg(1)
110      !
[11041]111      ! controls
[11536]112      IF( ln_icectl    )   CALL ice_prt     ( kt, iiceprt, jiceprt, 1, ' - ice thermo bdy - ' )                        ! prints
[11041]113      IF( ln_icediachk )   CALL ice_cons_hsm(1,'bdy_ice_thd', rdiag_v, rdiag_s, rdiag_t, rdiag_fv, rdiag_fs, rdiag_ft) ! conservation
[11536]114      IF( ln_icediachk )   CALL ice_cons2D  (1,'bdy_ice_thd',  diag_v,  diag_s,  diag_t,  diag_fv,  diag_fs,  diag_ft) ! conservation
[11041]115      IF( ln_timing    )   CALL timing_stop ('bdy_ice_thd')                                                            ! timing
[8586]116      !
117   END SUBROUTINE bdy_ice
118
119
[11536]120   SUBROUTINE bdy_ice_frs( idx, dta, kt, jbdy, llrim0 )
[8586]121      !!------------------------------------------------------------------------------
122      !!                 ***  SUBROUTINE bdy_ice_frs  ***
123      !!                   
[9890]124      !! ** Purpose : Apply the Flow Relaxation Scheme for sea-ice fields
[8586]125      !!
126      !! Reference : Engedahl H., 1995: Use of the flow relaxation scheme in a three-
127      !!             dimensional baroclinic ocean model with realistic topography. Tellus, 365-382.
128      !!------------------------------------------------------------------------------
[11536]129      TYPE(OBC_INDEX), INTENT(in) ::   idx      ! OBC indices
130      TYPE(OBC_DATA),  INTENT(in) ::   dta      ! OBC external data
131      INTEGER,         INTENT(in) ::   kt       ! main time-step counter
132      INTEGER,         INTENT(in) ::   jbdy     ! BDY set index
133      LOGICAL,         INTENT(in) ::   llrim0   ! indicate if rim 0 is treated
[8586]134      !
135      INTEGER  ::   jpbound            ! 0 = incoming ice
136      !                                ! 1 = outgoing ice
[11536]137      INTEGER  ::   ibeg, iend         ! length of rim to be treated (rim 0 or rim 1)
[9890]138      INTEGER  ::   i_bdy, jgrd        ! dummy loop indices
139      INTEGER  ::   ji, jj, jk, jl, ib, jb
[8586]140      REAL(wp) ::   zwgt, zwgt1        ! local scalar
141      REAL(wp) ::   ztmelts, zdh
[11536]142      REAL(wp), POINTER  :: flagu, flagv              ! short cuts
[8586]143      !!------------------------------------------------------------------------------
144      !
145      jgrd = 1      ! Everything is at T-points here
[11536]146      IF( llrim0 ) THEN   ;   ibeg = 1                       ;   iend = idx%nblenrim0(jgrd)
147      ELSE                ;   ibeg = idx%nblenrim0(jgrd)+1   ;   iend = idx%nblenrim(jgrd)
148      END IF
[8586]149      !
150      DO jl = 1, jpl
[11536]151         DO i_bdy = ibeg, iend
[9890]152            ji    = idx%nbi(i_bdy,jgrd)
153            jj    = idx%nbj(i_bdy,jgrd)
154            zwgt  = idx%nbw(i_bdy,jgrd)
155            zwgt1 = 1.e0 - idx%nbw(i_bdy,jgrd)
[11536]156            a_i (ji,jj,  jl) = ( a_i (ji,jj,  jl) * zwgt1 + dta%a_i(i_bdy,jl) * zwgt ) * tmask(ji,jj,1)  ! Ice  concentration
157            h_i (ji,jj,  jl) = ( h_i (ji,jj,  jl) * zwgt1 + dta%h_i(i_bdy,jl) * zwgt ) * tmask(ji,jj,1)  ! Ice  depth
158            h_s (ji,jj,  jl) = ( h_s (ji,jj,  jl) * zwgt1 + dta%h_s(i_bdy,jl) * zwgt ) * tmask(ji,jj,1)  ! Snow depth
159            t_i (ji,jj,:,jl) = ( t_i (ji,jj,:,jl) * zwgt1 + dta%t_i(i_bdy,jl) * zwgt ) * tmask(ji,jj,1)  ! Ice  temperature
160            t_s (ji,jj,:,jl) = ( t_s (ji,jj,:,jl) * zwgt1 + dta%t_s(i_bdy,jl) * zwgt ) * tmask(ji,jj,1)  ! Snow temperature
161            t_su(ji,jj,  jl) = ( t_su(ji,jj,  jl) * zwgt1 + dta%tsu(i_bdy,jl) * zwgt ) * tmask(ji,jj,1)  ! Surf temperature
162            s_i (ji,jj,  jl) = ( s_i (ji,jj,  jl) * zwgt1 + dta%s_i(i_bdy,jl) * zwgt ) * tmask(ji,jj,1)  ! Ice  salinity
163            a_ip(ji,jj,  jl) = ( a_ip(ji,jj,  jl) * zwgt1 + dta%aip(i_bdy,jl) * zwgt ) * tmask(ji,jj,1)  ! Ice  pond concentration
164            h_ip(ji,jj,  jl) = ( h_ip(ji,jj,  jl) * zwgt1 + dta%hip(i_bdy,jl) * zwgt ) * tmask(ji,jj,1)  ! Ice  pond depth
[12402]165            lh_ip(ji,jj, jl) = ( lh_ip(ji,jj, jl) * zwgt1 + dta%hip(i_bdy,jl) * zwgt ) * tmask(ji,jj,1)  ! Ice  pond lid thickness
[11536]166            !
167            sz_i(ji,jj,:,jl) = s_i(ji,jj,jl)
168            !
169            ! make sure ponds = 0 if no ponds scheme
170            IF( .NOT.ln_pnd ) THEN
171               a_ip(ji,jj,jl) = 0._wp
172               h_ip(ji,jj,jl) = 0._wp
[12402]173               lh_ip(ji,jj,jl) = 0._wp
[11536]174            ENDIF
175            !
[8586]176            ! -----------------
177            ! Pathological case
178            ! -----------------
179            ! In case a) snow load would be in excess or b) ice is coming into a warmer environment that would lead to
180            ! very large transformation from snow to ice (see icethd_dh.F90)
181
182            ! Then, a) transfer the snow excess into the ice (different from icethd_dh)
[9935]183            zdh = MAX( 0._wp, ( rhos * h_s(ji,jj,jl) + ( rhoi - rau0 ) * h_i(ji,jj,jl) ) * r1_rau0 )
[8586]184            ! Or, b) transfer all the snow into ice (if incoming ice is likely to melt as it comes into a warmer environment)
[9935]185            !zdh = MAX( 0._wp, h_s(ji,jj,jl) * rhos / rhoi )
[8586]186
187            ! recompute h_i, h_s
188            h_i(ji,jj,jl) = MIN( hi_max(jl), h_i(ji,jj,jl) + zdh )
[9935]189            h_s(ji,jj,jl) = MAX( 0._wp, h_s(ji,jj,jl) - zdh * rhoi / rhos ) 
[11536]190            !
[8586]191         ENDDO
192      ENDDO
193
194      DO jl = 1, jpl
[11536]195         DO i_bdy = ibeg, iend
[9890]196            ji = idx%nbi(i_bdy,jgrd)
197            jj = idx%nbj(i_bdy,jgrd)
[11536]198            flagu => idx%flagu(i_bdy,jgrd)
199            flagv => idx%flagv(i_bdy,jgrd)
[8586]200            ! condition on ice thickness depends on the ice velocity
201            ! if velocity is outward (strictly), then ice thickness, volume... must be equal to adjacent values
[9890]202            jpbound = 0   ;   ib = ji   ;   jb = jj
[8586]203            !
[11536]204            IF( flagu ==  1. )   THEN
205               IF( ji+1 > jpi )   CYCLE
206               IF( u_ice(ji  ,jj  ) < 0. )   jpbound = 1 ; ib = ji+1
207            END IF
208            IF( flagu == -1. )   THEN
209               IF( ji-1 < 1   )   CYCLE
210               IF( u_ice(ji-1,jj  ) < 0. )   jpbound = 1 ; ib = ji-1
211            END IF
212            IF( flagv ==  1. )   THEN
213               IF( jj+1 > jpj )   CYCLE
214               IF( v_ice(ji  ,jj  ) < 0. )   jpbound = 1 ; jb = jj+1
215            END IF
216            IF( flagv == -1. )   THEN
217               IF( jj-1 < 1   )   CYCLE
218               IF( v_ice(ji  ,jj-1) < 0. )   jpbound = 1 ; jb = jj-1
219            END IF
[8586]220            !
[9890]221            IF( nn_ice_dta(jbdy) == 0 )   jpbound = 0 ; ib = ji ; jb = jj   ! case ice boundaries = initial conditions
222            !                                                               !      do not make state variables dependent on velocity
[8586]223            !
[9890]224            IF( a_i(ib,jb,jl) > 0._wp ) THEN   ! there is ice at the boundary
[8586]225               !
[11536]226               a_i (ji,jj,  jl) = a_i (ib,jb,  jl)
227               h_i (ji,jj,  jl) = h_i (ib,jb,  jl)
228               h_s (ji,jj,  jl) = h_s (ib,jb,  jl)
229               t_i (ji,jj,:,jl) = t_i (ib,jb,:,jl)
230               t_s (ji,jj,:,jl) = t_s (ib,jb,:,jl)
231               t_su(ji,jj,  jl) = t_su(ib,jb,  jl)
232               s_i (ji,jj,  jl) = s_i (ib,jb,  jl)
233               a_ip(ji,jj,  jl) = a_ip(ib,jb,  jl)
234               h_ip(ji,jj,  jl) = h_ip(ib,jb,  jl)
[12402]235               lh_ip(ji,jj, jl) = lh_ip(ib,jb, jl)
[8586]236               !
[11536]237               sz_i(ji,jj,:,jl) = sz_i(ib,jb,:,jl)
[9885]238               !
[11536]239               ! ice age
240               IF    ( jpbound == 0 ) THEN  ! velocity is inward
241                  oa_i(ji,jj,jl) = rice_age(jbdy) * a_i(ji,jj,jl)
242               ELSEIF( jpbound == 1 ) THEN  ! velocity is outward
243                  oa_i(ji,jj,jl) = oa_i(ib,jb,jl)
244               ENDIF
245               !
[9888]246               IF( nn_icesal == 1 ) THEN     ! if constant salinity
247                  s_i (ji,jj  ,jl) = rn_icesal
248                  sz_i(ji,jj,:,jl) = rn_icesal
249               ENDIF
[8586]250               !
[9888]251               ! global fields
252               v_i (ji,jj,jl) = h_i(ji,jj,jl) * a_i(ji,jj,jl)                       ! volume ice
253               v_s (ji,jj,jl) = h_s(ji,jj,jl) * a_i(ji,jj,jl)                       ! volume snw
254               sv_i(ji,jj,jl) = MIN( s_i(ji,jj,jl) , sss_m(ji,jj) ) * v_i(ji,jj,jl) ! salt content
[8586]255               DO jk = 1, nlay_s
[9935]256                  e_s(ji,jj,jk,jl) = rhos * ( rcpi * ( rt0 - t_s(ji,jj,jk,jl) ) + rLfus )   ! enthalpy in J/m3
[9888]257                  e_s(ji,jj,jk,jl) = e_s(ji,jj,jk,jl) * v_s(ji,jj,jl) * r1_nlay_s           ! enthalpy in J/m2
258               END DO               
[8586]259               DO jk = 1, nlay_i
[9935]260                  ztmelts          = - rTmlt  * sz_i(ji,jj,jk,jl)             ! Melting temperature in C
[9888]261                  t_i(ji,jj,jk,jl) = MIN( t_i(ji,jj,jk,jl), ztmelts + rt0 )   ! Force t_i to be lower than melting point => likely conservation issue
262                  !
[9935]263                  e_i(ji,jj,jk,jl) = rhoi * ( rcpi  * ( ztmelts - ( t_i(ji,jj,jk,jl) - rt0 ) )           &   ! enthalpy in J/m3
264                     &                      + rLfus * ( 1._wp - ztmelts / ( t_i(ji,jj,jk,jl) - rt0 ) )   &
265                     &                      - rcp   *   ztmelts )                 
[9888]266                  e_i(ji,jj,jk,jl) = e_i(ji,jj,jk,jl) * v_i(ji,jj,jl) * r1_nlay_i                            ! enthalpy in J/m2
[8586]267               END DO
268               !
[11536]269               ! melt ponds
270               IF( a_i(ji,jj,jl) > epsi10 ) THEN
271                  a_ip_frac(ji,jj,jl) = a_ip(ji,jj,jl) / a_i (ji,jj,jl)
272               ELSE
273                  a_ip_frac(ji,jj,jl) = 0._wp
274               ENDIF
275               v_ip(ji,jj,jl) = h_ip(ji,jj,jl) * a_ip(ji,jj,jl)
276               !
[9888]277            ELSE   ! no ice at the boundary
[9885]278               !
[9888]279               a_i (ji,jj,  jl) = 0._wp
280               h_i (ji,jj,  jl) = 0._wp
281               h_s (ji,jj,  jl) = 0._wp
282               oa_i(ji,jj,  jl) = 0._wp
283               a_ip(ji,jj,  jl) = 0._wp
284               v_ip(ji,jj,  jl) = 0._wp
[12402]285               lh_ip(ji,jj, jl) = 0._wp
[9888]286               t_su(ji,jj,  jl) = rt0
287               t_s (ji,jj,:,jl) = rt0
288               t_i (ji,jj,:,jl) = rt0 
[11536]289
290               a_ip_frac(ji,jj,jl) = 0._wp
291               h_ip     (ji,jj,jl) = 0._wp
292               a_ip     (ji,jj,jl) = 0._wp
293               v_ip     (ji,jj,jl) = 0._wp
[9888]294               
295               IF( nn_icesal == 1 ) THEN     ! if constant salinity
296                  s_i (ji,jj  ,jl) = rn_icesal
297                  sz_i(ji,jj,:,jl) = rn_icesal
298               ELSE                          ! if variable salinity
299                  s_i (ji,jj,jl)   = rn_simin
300                  sz_i(ji,jj,:,jl) = rn_simin
301               ENDIF
302               !
303               ! global fields
304               v_i (ji,jj,  jl) = 0._wp
305               v_s (ji,jj,  jl) = 0._wp
306               sv_i(ji,jj,  jl) = 0._wp
307               e_s (ji,jj,:,jl) = 0._wp
308               e_i (ji,jj,:,jl) = 0._wp
309
[8586]310            ENDIF
[9888]311                       
[8586]312         END DO
313         !
[9888]314      END DO ! jl
[8586]315      !     
316   END SUBROUTINE bdy_ice_frs
317
318
319   SUBROUTINE bdy_ice_dyn( cd_type )
320      !!------------------------------------------------------------------------------
321      !!                 ***  SUBROUTINE bdy_ice_dyn  ***
322      !!                   
[9890]323      !! ** Purpose : Apply dynamics boundary conditions for sea-ice.
[8586]324      !!
[9890]325      !! ** Method :  if this adjacent grid point is not ice free, then u_ice and v_ice take its value
326      !!              if                          is     ice free, then u_ice and v_ice are unchanged by BDY
327      !!                                                           they keep values calculated in rheology
328      !!
[8586]329      !!------------------------------------------------------------------------------
330      CHARACTER(len=1), INTENT(in)  ::   cd_type   ! nature of velocity grid-points
331      !
[11536]332      INTEGER  ::   i_bdy, jgrd       ! dummy loop indices
333      INTEGER  ::   ji, jj            ! local scalar
334      INTEGER  ::   jbdy, ir     ! BDY set index, rim index
335      INTEGER  ::   ibeg, iend   ! length of rim to be treated (rim 0 or rim 1)
[8586]336      REAL(wp) ::   zmsk1, zmsk2, zflag
[11536]337      LOGICAL, DIMENSION(4) :: llsend2, llrecv2, llsend3, llrecv3  ! indicate how communications are to be carried out
[8586]338      !!------------------------------------------------------------------------------
[9905]339      IF( ln_timing )   CALL timing_start('bdy_ice_dyn')
[8586]340      !
[11536]341      llsend2(:) = .false.   ;   llrecv2(:) = .false.
342      llsend3(:) = .false.   ;   llrecv3(:) = .false.
343      DO ir = 1, 0, -1
344         DO jbdy = 1, nb_bdy
[8586]345            !
[11536]346            SELECT CASE( cn_ice(jbdy) )
[8586]347               !
[11536]348            CASE('none')
349               CYCLE
350               !
351            CASE('frs')
352               !
353               IF( nn_ice_dta(jbdy) == 0 ) CYCLE            ! case ice boundaries = initial conditions
354               !                                            !      do not change ice velocity (it is only computed by rheology)
355               SELECT CASE ( cd_type )
356                  !     
357               CASE ( 'U' ) 
358                  jgrd = 2      ! u velocity
359                  IF( ir == 0 ) THEN   ;   ibeg = 1                                 ;   iend = idx_bdy(jbdy)%nblenrim0(jgrd)
360                  ELSE                 ;   ibeg = idx_bdy(jbdy)%nblenrim0(jgrd)+1   ;   iend = idx_bdy(jbdy)%nblenrim(jgrd)
361                  END IF
362                  DO i_bdy = ibeg, iend
363                     ji    = idx_bdy(jbdy)%nbi(i_bdy,jgrd)
364                     jj    = idx_bdy(jbdy)%nbj(i_bdy,jgrd)
365                     zflag = idx_bdy(jbdy)%flagu(i_bdy,jgrd)
366                     !     i-1  i   i    |  !        i  i i+1 |  !          i  i i+1 |
367                     !      >  ice  >    |  !        o  > ice |  !          o  >  o  |     
368                     ! => set at u_ice(i-1) !  => set to O       !  => unchanged
369                     IF( zflag == -1. .AND. ji > 1 .AND. ji < jpi )   THEN 
370                        IF    ( vt_i(ji  ,jj) > 0. )   THEN   ;   u_ice(ji,jj) = u_ice(ji-1,jj) 
371                        ELSEIF( vt_i(ji+1,jj) > 0. )   THEN   ;   u_ice(ji,jj) = 0._wp
372                        END IF
373                     END IF
374                     ! |    i  i+1 i+1        !  |  i   i i+1        !  | i  i i+1
375                     ! |    >  ice  >         !  | ice  >  o         !  | o  >  o   
376                     ! => set at u_ice(i+1)   !     => set to O      !     =>  unchanged
377                     IF( zflag ==  1. .AND. ji+1 < jpi+1 )   THEN
378                        IF    ( vt_i(ji+1,jj) > 0. )   THEN   ;   u_ice(ji,jj) = u_ice(ji+1,jj)
379                        ELSEIF( vt_i(ji  ,jj) > 0. )   THEN   ;   u_ice(ji,jj) = 0._wp
380                        END IF
381                     END IF
382                     !
383                     IF( zflag ==  0. )   u_ice(ji,jj) = 0._wp   ! u_ice = 0  if north/south bdy 
384                     !
385                  END DO
[8586]386                  !
[11536]387               CASE ( 'V' )
388                  jgrd = 3      ! v velocity
389                  IF( ir == 0 ) THEN   ;   ibeg = 1                                 ;   iend = idx_bdy(jbdy)%nblenrim0(jgrd)
390                  ELSE                 ;   ibeg = idx_bdy(jbdy)%nblenrim0(jgrd)+1   ;   iend = idx_bdy(jbdy)%nblenrim(jgrd)
391                  END IF
392                  DO i_bdy = ibeg, iend
393                     ji    = idx_bdy(jbdy)%nbi(i_bdy,jgrd)
394                     jj    = idx_bdy(jbdy)%nbj(i_bdy,jgrd)
395                     zflag = idx_bdy(jbdy)%flagv(i_bdy,jgrd)
396                     !                         !      ice   (jj+1)       !       o    (jj+1)
397                     !       ^    (jj  )       !       ^    (jj  )       !       ^    (jj  )       
398                     !      ice   (jj  )       !       o    (jj  )       !       o    (jj  )       
399                     !       ^    (jj-1)       !                         !
400                     ! => set to u_ice(jj-1)   !  =>   set to 0          !   => unchanged       
401                     IF( zflag == -1. .AND. jj > 1 .AND. jj < jpj )   THEN                 
402                        IF    ( vt_i(ji,jj  ) > 0. )   THEN   ;   v_ice(ji,jj) = v_ice(ji,jj-1)
403                        ELSEIF( vt_i(ji,jj+1) > 0. )   THEN   ;   v_ice(ji,jj) = 0._wp
404                        END IF
405                     END IF
406                     !       ^    (jj+1)       !                         !             
407                     !      ice   (jj+1)       !       o    (jj+1)       !       o    (jj+1)       
408                     !       ^    (jj  )       !       ^    (jj  )       !       ^    (jj  )
409                     !   ________________      !  ____ice___(jj  )_      !  _____o____(jj  )
410                     ! => set to u_ice(jj+1)   !    => set to 0          !    => unchanged 
411                     IF( zflag ==  1. .AND. jj < jpj )   THEN             
412                        IF    ( vt_i(ji,jj+1) > 0. )   THEN   ;   v_ice(ji,jj) = v_ice(ji,jj+1)
413                        ELSEIF( vt_i(ji,jj  ) > 0. )   THEN   ;   v_ice(ji,jj) = 0._wp
414                        END IF
415                     END IF
416                     !
417                     IF( zflag ==  0. )   v_ice(ji,jj) = 0._wp   ! v_ice = 0  if west/east bdy 
418                     !
419                  END DO
[8586]420                  !
[11536]421               END SELECT
[8586]422               !
[11536]423            CASE DEFAULT
424               CALL ctl_stop( 'bdy_ice_dyn : unrecognised option for open boundaries for ice fields' )
[8586]425            END SELECT
426            !
[11536]427         END DO    ! jbdy
428         !
429         SELECT CASE ( cd_type )       
430         CASE ( 'U' ) 
431         IF( nn_hls > 1 .AND. ir == 1 ) CYCLE   ! at least 2 halos will be corrected -> no need to correct rim 1 before rim 0
432         IF( nn_hls == 1 ) THEN   ;   llsend2(:) = .false.   ;   llrecv2(:) = .false.   ;   END IF
433            DO jbdy = 1, nb_bdy
434               IF( cn_ice(jbdy) == 'frs' .AND. nn_ice_dta(jbdy) /= 0 ) THEN
435                  llsend2(:) = llsend2(:) .OR. lsend_bdyint(jbdy,2,:,ir)   ! possibly every direction, U points
436                  llsend2(1) = llsend2(1) .OR. lsend_bdyext(jbdy,2,1,ir)   ! neighbour might search point towards its west bdy
437                  llrecv2(:) = llrecv2(:) .OR. lrecv_bdyint(jbdy,2,:,ir)   ! possibly every direction, U points
438                  llrecv2(2) = llrecv2(2) .OR. lrecv_bdyext(jbdy,2,2,ir)   ! might search point towards east bdy
439               END IF
440            END DO
441            IF( ANY(llsend2) .OR. ANY(llrecv2) ) THEN   ! if need to send/recv in at least one direction
442               CALL lbc_lnk( 'bdyice', u_ice, 'U', -1., kfillmode=jpfillnothing ,lsend=llsend2, lrecv=llrecv2 )
443            END IF
444         CASE ( 'V' )
445         IF( nn_hls > 1 .AND. ir == 1 ) CYCLE   ! at least 2 halos will be corrected -> no need to correct rim 1 before rim 0
446         IF( nn_hls == 1 ) THEN   ;   llsend3(:) = .false.   ;   llrecv3(:) = .false.   ;   END IF
447            DO jbdy = 1, nb_bdy
448               IF( cn_ice(jbdy) == 'frs' .AND. nn_ice_dta(jbdy) /= 0 ) THEN
449                  llsend3(:) = llsend3(:) .OR. lsend_bdyint(jbdy,3,:,ir)   ! possibly every direction, V points
450                  llsend3(3) = llsend3(3) .OR. lsend_bdyext(jbdy,3,3,ir)   ! neighbour might search point towards its south bdy
451                  llrecv3(:) = llrecv3(:) .OR. lrecv_bdyint(jbdy,3,:,ir)   ! possibly every direction, V points
452                  llrecv3(4) = llrecv3(4) .OR. lrecv_bdyext(jbdy,3,4,ir)   ! might search point towards north bdy
453               END IF
454            END DO
455            IF( ANY(llsend3) .OR. ANY(llrecv3) ) THEN   ! if need to send/recv in at least one direction
456               CALL lbc_lnk( 'bdyice', v_ice, 'V', -1., kfillmode=jpfillnothing ,lsend=llsend3, lrecv=llrecv3 )
457            END IF
[8586]458         END SELECT
[11536]459      END DO   ! ir
[8586]460      !
[9905]461      IF( ln_timing )   CALL timing_stop('bdy_ice_dyn')
462      !
[8586]463    END SUBROUTINE bdy_ice_dyn
464
465#else
466   !!---------------------------------------------------------------------------------
467   !!   Default option                                                    Empty module
468   !!---------------------------------------------------------------------------------
469CONTAINS
470   SUBROUTINE bdy_ice( kt )      ! Empty routine
[9927]471      IMPLICIT NONE
472      INTEGER, INTENT( in ) :: kt
[8586]473      WRITE(*,*) 'bdy_ice: You should not have seen this print! error?', kt
474   END SUBROUTINE bdy_ice
475#endif
476
477   !!=================================================================================
478END MODULE bdyice
Note: See TracBrowser for help on using the repository browser.