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/2019/dev_r10984_HPC-13_IRRMANN_BDY_optimization/src/OCE/BDY – NEMO

source: NEMO/branches/2019/dev_r10984_HPC-13_IRRMANN_BDY_optimization/src/OCE/BDY/bdyice.F90 @ 11423

Last change on this file since 11423 was 11210, checked in by girrmann, 5 years ago

dev_r10984_HPC-13 : remove some communications in bdy treatment in case nn_hls > 1, see #2285

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