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

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

Fix compile errors.

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