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

Last change on this file since 11510 was 11510, checked in by clem, 5 years ago

add ponds for bdy

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