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/2020/dev_r13296_HPC-07_mocavero_mpi3/src/OCE/BDY – NEMO

source: NEMO/branches/2020/dev_r13296_HPC-07_mocavero_mpi3/src/OCE/BDY/bdyice.F90 @ 13630

Last change on this file since 13630 was 13630, checked in by mocavero, 4 years ago

Add neighborhood collectives calls in the NEMO src - ticket #2496

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