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.
bdyini.F90 in branches/2014/dev_r4621_NOC4_BDY_VERT_INTERP/NEMOGCM/NEMO/OPA_SRC/BDY – NEMO

source: branches/2014/dev_r4621_NOC4_BDY_VERT_INTERP/NEMOGCM/NEMO/OPA_SRC/BDY/bdyini.F90 @ 5620

Last change on this file since 5620 was 5620, checked in by jamesharle, 9 years ago

Merge with r5619 of trunk, update to unstructured BDY interpolation in
fldread.F90. Structured BDY interpolation incomplete.

  • Property svn:keywords set to Id
File size: 83.8 KB
Line 
1MODULE bdyini
2   !!======================================================================
3   !!                       ***  MODULE  bdyini  ***
4   !! Unstructured open boundaries : initialisation
5   !!======================================================================
6   !! History :  1.0  !  2005-01  (J. Chanut, A. Sellar)  Original code
7   !!             -   !  2007-01  (D. Storkey) Update to use IOM module
8   !!             -   !  2007-01  (D. Storkey) Tidal forcing
9   !!            3.0  !  2008-04  (NEMO team)  add in the reference version
10   !!            3.3  !  2010-09  (E.O'Dea) updates for Shelf configurations
11   !!            3.3  !  2010-09  (D.Storkey) add ice boundary conditions
12   !!            3.4  !  2011     (D. Storkey) rewrite in preparation for OBC-BDY merge
13   !!            3.4  !  2012     (J. Chanut) straight open boundary case update
14   !!            3.5  !  2012     (S. Mocavero, I. Epicoco) Updates for the
15   !!                             optimization of BDY communications
16   !!----------------------------------------------------------------------
17#if defined key_bdy
18   !!----------------------------------------------------------------------
19   !!   'key_bdy'                     Unstructured Open Boundary Conditions
20   !!----------------------------------------------------------------------
21   !!   bdy_init       : Initialization of unstructured open boundaries
22   !!----------------------------------------------------------------------
23   USE wrk_nemo        ! Memory Allocation
24   USE timing          ! Timing
25   USE oce             ! ocean dynamics and tracers variables
26   USE dom_oce         ! ocean space and time domain
27   USE bdy_oce         ! unstructured open boundary conditions
28   USE in_out_manager  ! I/O units
29   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
30   USE lib_mpp         ! for mpp_sum 
31   USE iom             ! I/O
32   USE sbctide, ONLY: lk_tide ! Tidal forcing or not
33   USE phycst, ONLY: rday
34
35   IMPLICIT NONE
36   PRIVATE
37
38   PUBLIC   bdy_init   ! routine called in nemo_init
39
40   INTEGER, PARAMETER          :: jp_nseg = 100
41   INTEGER, PARAMETER          :: nrimmax = 20 ! maximum rimwidth in structured
42                                               ! open boundary data files
43   ! Straight open boundary segment parameters:
44   INTEGER  :: nbdysege, nbdysegw, nbdysegn, nbdysegs 
45   INTEGER, DIMENSION(jp_nseg) :: jpieob, jpjedt, jpjeft, npckge
46   INTEGER, DIMENSION(jp_nseg) :: jpiwob, jpjwdt, jpjwft, npckgw
47   INTEGER, DIMENSION(jp_nseg) :: jpjnob, jpindt, jpinft, npckgn
48   INTEGER, DIMENSION(jp_nseg) :: jpjsob, jpisdt, jpisft, npckgs
49   !!----------------------------------------------------------------------
50   !! NEMO/OPA 4.0 , NEMO Consortium (2011)
51   !! $Id$
52   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
53   !!----------------------------------------------------------------------
54CONTAINS
55   
56   SUBROUTINE bdy_init
57      !!----------------------------------------------------------------------
58      !!                 ***  ROUTINE bdy_init  ***
59      !!         
60      !! ** Purpose :   Initialization of the dynamics and tracer fields with
61      !!              unstructured open boundaries.
62      !!
63      !! ** Method  :   Read initialization arrays (mask, indices) to identify
64      !!              an unstructured open boundary
65      !!
66      !! ** Input   :  bdy_init.nc, input file for unstructured open boundaries
67      !!----------------------------------------------------------------------     
68      ! namelist variables
69      !-------------------
70      CHARACTER(LEN=80),DIMENSION(jpbgrd)  ::   clfile
71      CHARACTER(LEN=1)   ::   ctypebdy
72      INTEGER :: nbdyind, nbdybeg, nbdyend
73
74      ! local variables
75      !-------------------
76      INTEGER  ::   ib_bdy, ii, ij, ik, igrd, ib, ir, iseg ! dummy loop indices
77      INTEGER  ::   icount, icountr, ibr_max, ilen1, ibm1  ! local integers
78      INTEGER  ::   iw, ie, is, in, inum, id_dummy         !   -       -
79      INTEGER  ::   igrd_start, igrd_end, jpbdta           !   -       -
80      INTEGER  ::   jpbdtau, jpbdtas                       !   -       -
81      INTEGER  ::   ib_bdy1, ib_bdy2, ib1, ib2             !   -       -
82      INTEGER  ::   i_offset, j_offset                     !   -       -
83      INTEGER, POINTER  ::  nbi, nbj, nbr                  ! short cuts
84      REAL(wp), POINTER  ::  flagu, flagv                  !    -   -
85      REAL(wp), POINTER, DIMENSION(:,:)       ::   pmask    ! pointer to 2D mask fields
86      REAL(wp) ::   zefl, zwfl, znfl, zsfl                 ! local scalars
87      INTEGER, DIMENSION (2)                  ::   kdimsz
88      INTEGER, DIMENSION(jpbgrd,jp_bdy)       ::   nblendta         ! Length of index arrays
89      INTEGER, ALLOCATABLE, DIMENSION(:,:,:)  ::   nbidta, nbjdta   ! Index arrays: i and j indices of bdy dta
90      INTEGER, ALLOCATABLE, DIMENSION(:,:,:)  ::   nbrdta           ! Discrete distance from rim points
91      CHARACTER(LEN=1),DIMENSION(jpbgrd)      ::   cgrid
92      INTEGER :: com_east, com_west, com_south, com_north          ! Flags for boundaries sending
93      INTEGER :: com_east_b, com_west_b, com_south_b, com_north_b  ! Flags for boundaries receiving
94      INTEGER :: iw_b(4), ie_b(4), is_b(4), in_b(4)                ! Arrays for neighbours coordinates
95      REAL(wp), POINTER, DIMENSION(:,:)       ::   zfmask  ! temporary fmask array excluding coastal boundary condition (shlat)
96
97      !!
98      NAMELIST/nambdy/ nb_bdy, ln_coords_file, cn_coords_file,                 &
99         &             ln_mask_file, cn_mask_file, cn_dyn2d, nn_dyn2d_dta,     &
100         &             cn_dyn3d, nn_dyn3d_dta, cn_tra, nn_tra_dta,             & 
101         &             ln_tra_dmp, ln_dyn3d_dmp, rn_time_dmp, rn_time_dmp_out, &
102         &             cn_ice_lim, nn_ice_lim_dta,                           &
103         &             rn_ice_tem, rn_ice_sal, rn_ice_age,                 &
104         &             ln_vol, nn_volctl, nn_rimwidth, nb_jpk_bdy
105      !!
106      NAMELIST/nambdy_index/ ctypebdy, nbdyind, nbdybeg, nbdyend
107      INTEGER  ::   ios                 ! Local integer output status for namelist read
108      !!----------------------------------------------------------------------
109
110      IF( nn_timing == 1 ) CALL timing_start('bdy_init')
111
112      IF(lwp) WRITE(numout,*)
113      IF(lwp) WRITE(numout,*) 'bdy_init : initialization of open boundaries'
114      IF(lwp) WRITE(numout,*) '~~~~~~~~'
115      !
116
117      IF( jperio /= 0 )   CALL ctl_stop( 'Cyclic or symmetric,',   &
118         &                               ' and general open boundary condition are not compatible' )
119
120      cgrid= (/'t','u','v'/)
121     
122      ! ------------------------
123      ! Read namelist parameters
124      ! ------------------------
125
126      REWIND( numnam_ref )              ! Namelist nambdy in reference namelist :Unstructured open boundaries 
127      READ  ( numnam_ref, nambdy, IOSTAT = ios, ERR = 901)
128901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'nambdy in reference namelist', lwp )
129
130      REWIND( numnam_cfg )              ! Namelist nambdy in configuration namelist :Unstructured open boundaries
131      READ  ( numnam_cfg, nambdy, IOSTAT = ios, ERR = 902 )
132902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'nambdy in configuration namelist', lwp )
133      IF(lwm) WRITE ( numond, nambdy )
134
135      ! -----------------------------------------
136      ! Check and write out namelist parameters
137      ! -----------------------------------------
138      !                                   ! control prints
139      IF(lwp) WRITE(numout,*) '         nambdy'
140
141      IF( nb_bdy .eq. 0 ) THEN
142        IF(lwp) WRITE(numout,*) 'nb_bdy = 0, NO OPEN BOUNDARIES APPLIED.'
143      ELSE
144        IF(lwp) WRITE(numout,*) 'Number of open boundary sets : ',nb_bdy
145      ENDIF
146
147      DO ib_bdy = 1,nb_bdy
148        IF(lwp) WRITE(numout,*) ' ' 
149        IF(lwp) WRITE(numout,*) '------ Open boundary data set ',ib_bdy,'------' 
150
151        IF( ln_coords_file(ib_bdy) ) THEN
152           IF(lwp) WRITE(numout,*) 'Boundary definition read from file '//TRIM(cn_coords_file(ib_bdy))
153        ELSE
154           IF(lwp) WRITE(numout,*) 'Boundary defined in namelist.'
155        ENDIF
156        IF(lwp) WRITE(numout,*)
157       
158        IF(lwp) WRITE(numout,*) 'Boundary conditions for barotropic solution:  '
159        SELECT CASE( cn_dyn2d(ib_bdy) )                 
160          CASE('none')         
161             IF(lwp) WRITE(numout,*) '      no open boundary condition'       
162             dta_bdy(ib_bdy)%ll_ssh = .false.
163             dta_bdy(ib_bdy)%ll_u2d = .false.
164             dta_bdy(ib_bdy)%ll_v2d = .false.
165          CASE('frs')         
166             IF(lwp) WRITE(numout,*) '      Flow Relaxation Scheme'
167             dta_bdy(ib_bdy)%ll_ssh = .false.
168             dta_bdy(ib_bdy)%ll_u2d = .true.
169             dta_bdy(ib_bdy)%ll_v2d = .true.
170          CASE('flather')     
171             IF(lwp) WRITE(numout,*) '      Flather radiation condition'
172             dta_bdy(ib_bdy)%ll_ssh = .true.
173             dta_bdy(ib_bdy)%ll_u2d = .true.
174             dta_bdy(ib_bdy)%ll_v2d = .true.
175          CASE('orlanski')     
176             IF(lwp) WRITE(numout,*) '      Orlanski (fully oblique) radiation condition with adaptive nudging'
177             dta_bdy(ib_bdy)%ll_ssh = .false.
178             dta_bdy(ib_bdy)%ll_u2d = .true.
179             dta_bdy(ib_bdy)%ll_v2d = .true.
180          CASE('orlanski_npo') 
181             IF(lwp) WRITE(numout,*) '      Orlanski (NPO) radiation condition with adaptive nudging'
182             dta_bdy(ib_bdy)%ll_ssh = .false.
183             dta_bdy(ib_bdy)%ll_u2d = .true.
184             dta_bdy(ib_bdy)%ll_v2d = .true.
185          CASE DEFAULT   ;   CALL ctl_stop( 'unrecognised value for cn_dyn2d' )
186        END SELECT
187        IF( cn_dyn2d(ib_bdy) /= 'none' ) THEN
188           SELECT CASE( nn_dyn2d_dta(ib_bdy) )                   !
189              CASE( 0 )      ;   IF(lwp) WRITE(numout,*) '      initial state used for bdy data'       
190              CASE( 1 )      ;   IF(lwp) WRITE(numout,*) '      boundary data taken from file'
191              CASE( 2 )      ;   IF(lwp) WRITE(numout,*) '      tidal harmonic forcing taken from file'
192              CASE( 3 )      ;   IF(lwp) WRITE(numout,*) '      boundary data AND tidal harmonic forcing taken from files'
193              CASE DEFAULT   ;   CALL ctl_stop( 'nn_dyn2d_dta must be between 0 and 3' )
194           END SELECT
195           IF (( nn_dyn2d_dta(ib_bdy) .ge. 2 ).AND.(.NOT.lk_tide)) THEN
196             CALL ctl_stop( 'You must activate key_tide to add tidal forcing at open boundaries' )
197           ENDIF
198        ENDIF
199        IF(lwp) WRITE(numout,*)
200
201        IF(lwp) WRITE(numout,*) 'Boundary conditions for baroclinic velocities:  '
202        SELECT CASE( cn_dyn3d(ib_bdy) )                 
203          CASE('none')
204             IF(lwp) WRITE(numout,*) '      no open boundary condition'       
205             dta_bdy(ib_bdy)%ll_u3d = .false.
206             dta_bdy(ib_bdy)%ll_v3d = .false.
207          CASE('frs')       
208             IF(lwp) WRITE(numout,*) '      Flow Relaxation Scheme'
209             dta_bdy(ib_bdy)%ll_u3d = .true.
210             dta_bdy(ib_bdy)%ll_v3d = .true.
211          CASE('specified')
212             IF(lwp) WRITE(numout,*) '      Specified value'
213             dta_bdy(ib_bdy)%ll_u3d = .true.
214             dta_bdy(ib_bdy)%ll_v3d = .true.
215          CASE('zero')
216             IF(lwp) WRITE(numout,*) '      Zero baroclinic velocities (runoff case)'
217             dta_bdy(ib_bdy)%ll_u3d = .false.
218             dta_bdy(ib_bdy)%ll_v3d = .false.
219          CASE('orlanski')
220             IF(lwp) WRITE(numout,*) '      Orlanski (fully oblique) radiation condition with adaptive nudging'
221             dta_bdy(ib_bdy)%ll_u3d = .true.
222             dta_bdy(ib_bdy)%ll_v3d = .true.
223          CASE('orlanski_npo')
224             IF(lwp) WRITE(numout,*) '      Orlanski (NPO) radiation condition with adaptive nudging'
225             dta_bdy(ib_bdy)%ll_u3d = .true.
226             dta_bdy(ib_bdy)%ll_v3d = .true.
227          CASE DEFAULT   ;   CALL ctl_stop( 'unrecognised value for cn_dyn3d' )
228        END SELECT
229        IF( cn_dyn3d(ib_bdy) /= 'none' ) THEN
230           SELECT CASE( nn_dyn3d_dta(ib_bdy) )                   !
231              CASE( 0 )      ;   IF(lwp) WRITE(numout,*) '      initial state used for bdy data'       
232              CASE( 1 )      ;   IF(lwp) WRITE(numout,*) '      boundary data taken from file'
233              CASE DEFAULT   ;   CALL ctl_stop( 'nn_dyn3d_dta must be 0 or 1' )
234           END SELECT
235        ENDIF
236
237        IF ( ln_dyn3d_dmp(ib_bdy) ) THEN
238           IF ( cn_dyn3d(ib_bdy) == 'none' ) THEN
239              IF(lwp) WRITE(numout,*) 'No open boundary condition for baroclinic velocities: ln_dyn3d_dmp is set to .false.'
240              ln_dyn3d_dmp(ib_bdy)=.false.
241           ELSEIF ( cn_dyn3d(ib_bdy) == 'frs' ) THEN
242              CALL ctl_stop( 'Use FRS OR relaxation' )
243           ELSE
244              IF(lwp) WRITE(numout,*) '      + baroclinic velocities relaxation zone'
245              IF(lwp) WRITE(numout,*) '      Damping time scale: ',rn_time_dmp(ib_bdy),' days'
246              IF((lwp).AND.rn_time_dmp(ib_bdy)<0) CALL ctl_stop( 'Time scale must be positive' )
247              dta_bdy(ib_bdy)%ll_u3d = .true.
248              dta_bdy(ib_bdy)%ll_v3d = .true.
249           ENDIF
250        ELSE
251           IF(lwp) WRITE(numout,*) '      NO relaxation on baroclinic velocities'
252        ENDIF
253        IF(lwp) WRITE(numout,*)
254
255        IF(lwp) WRITE(numout,*) 'Boundary conditions for temperature and salinity:  '
256        SELECT CASE( cn_tra(ib_bdy) )                 
257          CASE('none')
258             IF(lwp) WRITE(numout,*) '      no open boundary condition'       
259             dta_bdy(ib_bdy)%ll_tem = .false.
260             dta_bdy(ib_bdy)%ll_sal = .false.
261          CASE('frs')
262             IF(lwp) WRITE(numout,*) '      Flow Relaxation Scheme'
263             dta_bdy(ib_bdy)%ll_tem = .true.
264             dta_bdy(ib_bdy)%ll_sal = .true.
265          CASE('specified')
266             IF(lwp) WRITE(numout,*) '      Specified value'
267             dta_bdy(ib_bdy)%ll_tem = .true.
268             dta_bdy(ib_bdy)%ll_sal = .true.
269          CASE('neumann')
270             IF(lwp) WRITE(numout,*) '      Neumann conditions'
271             dta_bdy(ib_bdy)%ll_tem = .false.
272             dta_bdy(ib_bdy)%ll_sal = .false.
273          CASE('runoff')
274             IF(lwp) WRITE(numout,*) '      Runoff conditions : Neumann for T and specified to 0.1 for salinity'
275             dta_bdy(ib_bdy)%ll_tem = .false.
276             dta_bdy(ib_bdy)%ll_sal = .false.
277          CASE('orlanski')
278             IF(lwp) WRITE(numout,*) '      Orlanski (fully oblique) radiation condition with adaptive nudging'
279             dta_bdy(ib_bdy)%ll_tem = .true.
280             dta_bdy(ib_bdy)%ll_sal = .true.
281          CASE('orlanski_npo')
282             IF(lwp) WRITE(numout,*) '      Orlanski (NPO) radiation condition with adaptive nudging'
283             dta_bdy(ib_bdy)%ll_tem = .true.
284             dta_bdy(ib_bdy)%ll_sal = .true.
285          CASE DEFAULT   ;   CALL ctl_stop( 'unrecognised value for cn_tra' )
286        END SELECT
287        IF( cn_tra(ib_bdy) /= 'none' ) THEN
288           SELECT CASE( nn_tra_dta(ib_bdy) )                   !
289              CASE( 0 )      ;   IF(lwp) WRITE(numout,*) '      initial state used for bdy data'       
290              CASE( 1 )      ;   IF(lwp) WRITE(numout,*) '      boundary data taken from file'
291              CASE DEFAULT   ;   CALL ctl_stop( 'nn_tra_dta must be 0 or 1' )
292           END SELECT
293        ENDIF
294
295        IF ( ln_tra_dmp(ib_bdy) ) THEN
296           IF ( cn_tra(ib_bdy) == 'none' ) THEN
297              IF(lwp) WRITE(numout,*) 'No open boundary condition for tracers: ln_tra_dmp is set to .false.'
298              ln_tra_dmp(ib_bdy)=.false.
299           ELSEIF ( cn_tra(ib_bdy) == 'frs' ) THEN
300              CALL ctl_stop( 'Use FRS OR relaxation' )
301           ELSE
302              IF(lwp) WRITE(numout,*) '      + T/S relaxation zone'
303              IF(lwp) WRITE(numout,*) '      Damping time scale: ',rn_time_dmp(ib_bdy),' days'
304              IF(lwp) WRITE(numout,*) '      Outflow damping time scale: ',rn_time_dmp_out(ib_bdy),' days'
305              IF((lwp).AND.rn_time_dmp(ib_bdy)<0) CALL ctl_stop( 'Time scale must be positive' )
306              dta_bdy(ib_bdy)%ll_tem = .true.
307              dta_bdy(ib_bdy)%ll_sal = .true.
308           ENDIF
309        ELSE
310           IF(lwp) WRITE(numout,*) '      NO T/S relaxation'
311        ENDIF
312        IF(lwp) WRITE(numout,*)
313
314#if defined key_lim2
315        IF(lwp) WRITE(numout,*) 'Boundary conditions for sea ice:  '
316        SELECT CASE( cn_ice_lim(ib_bdy) )                 
317          CASE('none')
318             IF(lwp) WRITE(numout,*) '      no open boundary condition'       
319             dta_bdy(ib_bdy)%ll_frld  = .false.
320             dta_bdy(ib_bdy)%ll_hicif = .false.
321             dta_bdy(ib_bdy)%ll_hsnif = .false.
322          CASE('frs')
323             IF(lwp) WRITE(numout,*) '      Flow Relaxation Scheme'
324             dta_bdy(ib_bdy)%ll_frld  = .true.
325             dta_bdy(ib_bdy)%ll_hicif = .true.
326             dta_bdy(ib_bdy)%ll_hsnif = .true.
327          CASE DEFAULT   ;   CALL ctl_stop( 'unrecognised value for cn_ice_lim' )
328        END SELECT
329        IF( cn_ice_lim(ib_bdy) /= 'none' ) THEN
330           SELECT CASE( nn_ice_lim_dta(ib_bdy) )                   !
331              CASE( 0 )      ;   IF(lwp) WRITE(numout,*) '      initial state used for bdy data'       
332              CASE( 1 )      ;   IF(lwp) WRITE(numout,*) '      boundary data taken from file'
333              CASE DEFAULT   ;   CALL ctl_stop( 'nn_ice_lim_dta must be 0 or 1' )
334           END SELECT
335        ENDIF
336        IF(lwp) WRITE(numout,*)
337#elif defined key_lim3
338        IF(lwp) WRITE(numout,*) 'Boundary conditions for sea ice:  '
339        SELECT CASE( cn_ice_lim(ib_bdy) )                 
340          CASE('none')
341             IF(lwp) WRITE(numout,*) '      no open boundary condition'       
342             dta_bdy(ib_bdy)%ll_a_i  = .false.
343             dta_bdy(ib_bdy)%ll_ht_i = .false.
344             dta_bdy(ib_bdy)%ll_ht_s = .false.
345          CASE('frs')
346             IF(lwp) WRITE(numout,*) '      Flow Relaxation Scheme'
347             dta_bdy(ib_bdy)%ll_a_i  = .true.
348             dta_bdy(ib_bdy)%ll_ht_i = .true.
349             dta_bdy(ib_bdy)%ll_ht_s = .true.
350          CASE DEFAULT   ;   CALL ctl_stop( 'unrecognised value for cn_ice_lim' )
351        END SELECT
352        IF( cn_ice_lim(ib_bdy) /= 'none' ) THEN
353           SELECT CASE( nn_ice_lim_dta(ib_bdy) )                   !
354              CASE( 0 )      ;   IF(lwp) WRITE(numout,*) '      initial state used for bdy data'       
355              CASE( 1 )      ;   IF(lwp) WRITE(numout,*) '      boundary data taken from file'
356              CASE DEFAULT   ;   CALL ctl_stop( 'nn_ice_lim_dta must be 0 or 1' )
357           END SELECT
358        ENDIF
359        IF(lwp) WRITE(numout,*)
360        IF(lwp) WRITE(numout,*) '      tem of bdy sea-ice = ', rn_ice_tem(ib_bdy)         
361        IF(lwp) WRITE(numout,*) '      sal of bdy sea-ice = ', rn_ice_sal(ib_bdy)         
362        IF(lwp) WRITE(numout,*) '      age of bdy sea-ice = ', rn_ice_age(ib_bdy)         
363#endif
364
365        IF(lwp) WRITE(numout,*) '      Width of relaxation zone = ', nn_rimwidth(ib_bdy)
366        IF(lwp) WRITE(numout,*)
367
368      ENDDO
369
370     IF (nb_bdy .gt. 0) THEN
371        IF( ln_vol ) THEN                     ! check volume conservation (nn_volctl value)
372          IF(lwp) WRITE(numout,*) 'Volume correction applied at open boundaries'
373          IF(lwp) WRITE(numout,*)
374          SELECT CASE ( nn_volctl )
375            CASE( 1 )      ;   IF(lwp) WRITE(numout,*) '      The total volume will be constant'
376            CASE( 0 )      ;   IF(lwp) WRITE(numout,*) '      The total volume will vary according to the surface E-P flux'
377            CASE DEFAULT   ;   CALL ctl_stop( 'nn_volctl must be 0 or 1' )
378          END SELECT
379          IF(lwp) WRITE(numout,*)
380        ELSE
381          IF(lwp) WRITE(numout,*) 'No volume correction applied at open boundaries'
382          IF(lwp) WRITE(numout,*)
383        ENDIF
384        IF( nb_jpk_bdy > 0 ) THEN
385           IF(lwp) WRITE(numout,*) '*** open boundary will be interpolate in the vertical onto the native grid ***'
386        ELSE
387           IF(lwp) WRITE(numout,*) '*** open boundary will be read straight onto the native grid without vertical interpolation ***'
388        ENDIF
389     ENDIF
390
391      ! -------------------------------------------------
392      ! Initialise indices arrays for open boundaries
393      ! -------------------------------------------------
394
395      ! Work out global dimensions of boundary data
396      ! ---------------------------------------------
397      REWIND( numnam_cfg )     
398
399      !!----------------------------------------------------------------------
400
401             
402               
403      nblendta(:,:) = 0
404      nbdysege = 0
405      nbdysegw = 0
406      nbdysegn = 0
407      nbdysegs = 0
408      icount   = 0 ! count user defined segments
409      ! Dimensions below are used to allocate arrays to read external data
410      jpbdtas = 1 ! Maximum size of boundary data (structured case)
411      jpbdtau = 1 ! Maximum size of boundary data (unstructured case)
412
413      DO ib_bdy = 1, nb_bdy
414
415         IF( .NOT. ln_coords_file(ib_bdy) ) THEN ! Work out size of global arrays from namelist parameters
416 
417            icount = icount + 1
418            ! No REWIND here because may need to read more than one nambdy_index namelist.
419            ! Read only namelist_cfg to avoid unseccessfull overwrite
420!!          REWIND( numnam_ref )              ! Namelist nambdy_index in reference namelist : Open boundaries indexes
421!!          READ  ( numnam_ref, namrun, IOSTAT = ios, ERR = 903)
422!!903       IF( ios /= 0 ) CALL ctl_nam ( ios , 'nambdy_index in reference namelist', lwp )
423
424!!          REWIND( numnam_cfg )              ! Namelist nambdy_index in configuration namelist : Open boundaries indexes
425            READ  ( numnam_cfg, nambdy_index, IOSTAT = ios, ERR = 904 )
426904         IF( ios /= 0 ) CALL ctl_nam ( ios , 'nambdy_index in configuration namelist', lwp )
427            IF(lwm) WRITE ( numond, nambdy_index )
428
429            SELECT CASE ( TRIM(ctypebdy) )
430              CASE( 'N' )
431                 IF( nbdyind == -1 ) THEN  ! Automatic boundary definition: if nbdysegX = -1
432                    nbdyind  = jpjglo - 2  ! set boundary to whole side of model domain.
433                    nbdybeg  = 2
434                    nbdyend  = jpiglo - 1
435                 ENDIF
436                 nbdysegn = nbdysegn + 1
437                 npckgn(nbdysegn) = ib_bdy ! Save bdy package number
438                 jpjnob(nbdysegn) = nbdyind
439                 jpindt(nbdysegn) = nbdybeg
440                 jpinft(nbdysegn) = nbdyend
441                 !
442              CASE( 'S' )
443                 IF( nbdyind == -1 ) THEN  ! Automatic boundary definition: if nbdysegX = -1
444                    nbdyind  = 2           ! set boundary to whole side of model domain.
445                    nbdybeg  = 2
446                    nbdyend  = jpiglo - 1
447                 ENDIF
448                 nbdysegs = nbdysegs + 1
449                 npckgs(nbdysegs) = ib_bdy ! Save bdy package number
450                 jpjsob(nbdysegs) = nbdyind
451                 jpisdt(nbdysegs) = nbdybeg
452                 jpisft(nbdysegs) = nbdyend
453                 !
454              CASE( 'E' )
455                 IF( nbdyind == -1 ) THEN  ! Automatic boundary definition: if nbdysegX = -1
456                    nbdyind  = jpiglo - 2  ! set boundary to whole side of model domain.
457                    nbdybeg  = 2
458                    nbdyend  = jpjglo - 1
459                 ENDIF
460                 nbdysege = nbdysege + 1 
461                 npckge(nbdysege) = ib_bdy ! Save bdy package number
462                 jpieob(nbdysege) = nbdyind
463                 jpjedt(nbdysege) = nbdybeg
464                 jpjeft(nbdysege) = nbdyend
465                 !
466              CASE( 'W' )
467                 IF( nbdyind == -1 ) THEN  ! Automatic boundary definition: if nbdysegX = -1
468                    nbdyind  = 2           ! set boundary to whole side of model domain.
469                    nbdybeg  = 2
470                    nbdyend  = jpjglo - 1
471                 ENDIF
472                 nbdysegw = nbdysegw + 1
473                 npckgw(nbdysegw) = ib_bdy ! Save bdy package number
474                 jpiwob(nbdysegw) = nbdyind
475                 jpjwdt(nbdysegw) = nbdybeg
476                 jpjwft(nbdysegw) = nbdyend
477                 !
478              CASE DEFAULT   ;   CALL ctl_stop( 'ctypebdy must be N, S, E or W' )
479            END SELECT
480
481            ! For simplicity we assume that in case of straight bdy, arrays have the same length
482            ! (even if it is true that last tangential velocity points
483            ! are useless). This simplifies a little bit boundary data format (and agrees with format
484            ! used so far in obc package)
485
486            nblendta(1:jpbgrd,ib_bdy) =  (nbdyend - nbdybeg + 1) * nn_rimwidth(ib_bdy)
487            jpbdtas = MAX(jpbdtas, (nbdyend - nbdybeg + 1))
488            IF (lwp.and.(nn_rimwidth(ib_bdy)>nrimmax)) &
489            & CALL ctl_stop( 'rimwidth must be lower than nrimmax' )
490
491         ELSE            ! Read size of arrays in boundary coordinates file.
492            CALL iom_open( cn_coords_file(ib_bdy), inum )
493            DO igrd = 1, jpbgrd
494               id_dummy = iom_varid( inum, 'nbi'//cgrid(igrd), kdimsz=kdimsz ) 
495               !clem nblendta(igrd,ib_bdy) = kdimsz(1)
496               !clem jpbdtau = MAX(jpbdtau, kdimsz(1))
497               nblendta(igrd,ib_bdy) = MAXVAL(kdimsz)
498               jpbdtau = MAX(jpbdtau, MAXVAL(kdimsz))
499            ENDDO
500            CALL iom_close( inum )
501
502         ENDIF
503
504      ENDDO ! ib_bdy
505
506      IF (nb_bdy>0) THEN
507         jpbdta = MAXVAL(nblendta(1:jpbgrd,1:nb_bdy))
508
509         ! Allocate arrays
510         !---------------
511         ALLOCATE( nbidta(jpbdta, jpbgrd, nb_bdy), nbjdta(jpbdta, jpbgrd, nb_bdy),    &
512            &      nbrdta(jpbdta, jpbgrd, nb_bdy) )
513
514         IF( nb_jpk_bdy>0 ) THEN
515            ALLOCATE( dta_global(jpbdtau, 1, nb_jpk_bdy) )
516            ALLOCATE( dta_global_z(jpbdtau, 1, nb_jpk_bdy) )
517         ELSE
518            ALLOCATE( dta_global(jpbdtau, 1, jpk) )
519            ALLOCATE( dta_global_z(jpbdtau, 1, jpk) )
520         ENDIF
521
522         IF ( icount>0 ) THEN
523            IF( nb_jpk_bdy>0 ) THEN
524               ALLOCATE( dta_global2(jpbdtas, nrimmax, nb_jpk_bdy) )
525               ALLOCATE( dta_global2_z(jpbdtas, nrimmax, nb_jpk_bdy) )
526            ELSE
527               ALLOCATE( dta_global2(jpbdtas, nrimmax, jpk) )
528               ALLOCATE( dta_global2_z(jpbdtas, nrimmax, jpk) )
529            ENDIF
530         ENDIF
531         !
532      ENDIF
533
534      ! Now look for crossings in user (namelist) defined open boundary segments:
535      !--------------------------------------------------------------------------
536      IF ( icount>0 ) CALL bdy_ctl_seg
537
538      ! Calculate global boundary index arrays or read in from file
539      !------------------------------------------------------------               
540      ! 1. Read global index arrays from boundary coordinates file.
541      DO ib_bdy = 1, nb_bdy
542
543         IF( ln_coords_file(ib_bdy) ) THEN
544
545            CALL iom_open( cn_coords_file(ib_bdy), inum )
546            DO igrd = 1, jpbgrd
547               CALL iom_get( inum, jpdom_unknown, 'nbi'//cgrid(igrd), dta_global(1:nblendta(igrd,ib_bdy),:,1) )
548               DO ii = 1,nblendta(igrd,ib_bdy)
549                  nbidta(ii,igrd,ib_bdy) = INT( dta_global(ii,1,1) )
550               END DO
551               CALL iom_get( inum, jpdom_unknown, 'nbj'//cgrid(igrd), dta_global(1:nblendta(igrd,ib_bdy),:,1) )
552               DO ii = 1,nblendta(igrd,ib_bdy)
553                  nbjdta(ii,igrd,ib_bdy) = INT( dta_global(ii,1,1) )
554               END DO
555               CALL iom_get( inum, jpdom_unknown, 'nbr'//cgrid(igrd), dta_global(1:nblendta(igrd,ib_bdy),:,1) )
556               DO ii = 1,nblendta(igrd,ib_bdy)
557                  nbrdta(ii,igrd,ib_bdy) = INT( dta_global(ii,1,1) )
558               END DO
559
560               ibr_max = MAXVAL( nbrdta(:,igrd,ib_bdy) )
561               IF(lwp) WRITE(numout,*)
562               IF(lwp) WRITE(numout,*) ' Maximum rimwidth in file is ', ibr_max
563               IF(lwp) WRITE(numout,*) ' nn_rimwidth from namelist is ', nn_rimwidth(ib_bdy)
564               IF (ibr_max < nn_rimwidth(ib_bdy))   &
565                     CALL ctl_stop( 'nn_rimwidth is larger than maximum rimwidth in file',cn_coords_file(ib_bdy) )
566            END DO
567            CALL iom_close( inum )
568
569         ENDIF
570
571      ENDDO     
572   
573      ! 2. Now fill indices corresponding to straight open boundary arrays:
574      ! East
575      !-----
576      DO iseg = 1, nbdysege
577         ib_bdy = npckge(iseg)
578         !
579         ! ------------ T points -------------
580         igrd=1
581         icount=0
582         DO ir = 1, nn_rimwidth(ib_bdy)
583            DO ij = jpjedt(iseg), jpjeft(iseg)
584               icount = icount + 1
585               nbidta(icount, igrd, ib_bdy) = jpieob(iseg) + 2 - ir
586               nbjdta(icount, igrd, ib_bdy) = ij
587               nbrdta(icount, igrd, ib_bdy) = ir
588            ENDDO
589         ENDDO
590         !
591         ! ------------ U points -------------
592         igrd=2
593         icount=0
594         DO ir = 1, nn_rimwidth(ib_bdy)
595            DO ij = jpjedt(iseg), jpjeft(iseg)
596               icount = icount + 1
597               nbidta(icount, igrd, ib_bdy) = jpieob(iseg) + 1 - ir
598               nbjdta(icount, igrd, ib_bdy) = ij
599               nbrdta(icount, igrd, ib_bdy) = ir
600            ENDDO
601         ENDDO
602         !
603         ! ------------ V points -------------
604         igrd=3
605         icount=0
606         DO ir = 1, nn_rimwidth(ib_bdy)
607!            DO ij = jpjedt(iseg), jpjeft(iseg) - 1
608            DO ij = jpjedt(iseg), jpjeft(iseg)
609               icount = icount + 1
610               nbidta(icount, igrd, ib_bdy) = jpieob(iseg) + 2 - ir
611               nbjdta(icount, igrd, ib_bdy) = ij
612               nbrdta(icount, igrd, ib_bdy) = ir
613            ENDDO
614            nbidta(icount, igrd, ib_bdy) = -ib_bdy ! Discount this point
615            nbjdta(icount, igrd, ib_bdy) = -ib_bdy ! Discount this point
616         ENDDO
617      ENDDO
618      !
619      ! West
620      !-----
621      DO iseg = 1, nbdysegw
622         ib_bdy = npckgw(iseg)
623         !
624         ! ------------ T points -------------
625         igrd=1
626         icount=0
627         DO ir = 1, nn_rimwidth(ib_bdy)
628            DO ij = jpjwdt(iseg), jpjwft(iseg)
629               icount = icount + 1
630               nbidta(icount, igrd, ib_bdy) = jpiwob(iseg) + ir - 1
631               nbjdta(icount, igrd, ib_bdy) = ij
632               nbrdta(icount, igrd, ib_bdy) = ir
633            ENDDO
634         ENDDO
635         !
636         ! ------------ U points -------------
637         igrd=2
638         icount=0
639         DO ir = 1, nn_rimwidth(ib_bdy)
640            DO ij = jpjwdt(iseg), jpjwft(iseg)
641               icount = icount + 1
642               nbidta(icount, igrd, ib_bdy) = jpiwob(iseg) + ir - 1
643               nbjdta(icount, igrd, ib_bdy) = ij
644               nbrdta(icount, igrd, ib_bdy) = ir
645            ENDDO
646         ENDDO
647         !
648         ! ------------ V points -------------
649         igrd=3
650         icount=0
651         DO ir = 1, nn_rimwidth(ib_bdy)
652!            DO ij = jpjwdt(iseg), jpjwft(iseg) - 1
653            DO ij = jpjwdt(iseg), jpjwft(iseg)
654               icount = icount + 1
655               nbidta(icount, igrd, ib_bdy) = jpiwob(iseg) + ir - 1
656               nbjdta(icount, igrd, ib_bdy) = ij
657               nbrdta(icount, igrd, ib_bdy) = ir
658            ENDDO
659            nbidta(icount, igrd, ib_bdy) = -ib_bdy ! Discount this point
660            nbjdta(icount, igrd, ib_bdy) = -ib_bdy ! Discount this point
661         ENDDO
662      ENDDO
663      !
664      ! North
665      !-----
666      DO iseg = 1, nbdysegn
667         ib_bdy = npckgn(iseg)
668         !
669         ! ------------ T points -------------
670         igrd=1
671         icount=0
672         DO ir = 1, nn_rimwidth(ib_bdy)
673            DO ii = jpindt(iseg), jpinft(iseg)
674               icount = icount + 1
675               nbidta(icount, igrd, ib_bdy) = ii
676               nbjdta(icount, igrd, ib_bdy) = jpjnob(iseg) + 2 - ir 
677               nbrdta(icount, igrd, ib_bdy) = ir
678            ENDDO
679         ENDDO
680         !
681         ! ------------ U points -------------
682         igrd=2
683         icount=0
684         DO ir = 1, nn_rimwidth(ib_bdy)
685!            DO ii = jpindt(iseg), jpinft(iseg) - 1
686            DO ii = jpindt(iseg), jpinft(iseg)
687               icount = icount + 1
688               nbidta(icount, igrd, ib_bdy) = ii
689               nbjdta(icount, igrd, ib_bdy) = jpjnob(iseg) + 2 - ir
690               nbrdta(icount, igrd, ib_bdy) = ir
691            ENDDO
692            nbidta(icount, igrd, ib_bdy) = -ib_bdy ! Discount this point
693            nbjdta(icount, igrd, ib_bdy) = -ib_bdy ! Discount this point
694         ENDDO
695         !
696         ! ------------ V points -------------
697         igrd=3
698         icount=0
699         DO ir = 1, nn_rimwidth(ib_bdy)
700            DO ii = jpindt(iseg), jpinft(iseg)
701               icount = icount + 1
702               nbidta(icount, igrd, ib_bdy) = ii
703               nbjdta(icount, igrd, ib_bdy) = jpjnob(iseg) + 1 - ir
704               nbrdta(icount, igrd, ib_bdy) = ir
705            ENDDO
706         ENDDO
707      ENDDO
708      !
709      ! South
710      !-----
711      DO iseg = 1, nbdysegs
712         ib_bdy = npckgs(iseg)
713         !
714         ! ------------ T points -------------
715         igrd=1
716         icount=0
717         DO ir = 1, nn_rimwidth(ib_bdy)
718            DO ii = jpisdt(iseg), jpisft(iseg)
719               icount = icount + 1
720               nbidta(icount, igrd, ib_bdy) = ii
721               nbjdta(icount, igrd, ib_bdy) = jpjsob(iseg) + ir - 1
722               nbrdta(icount, igrd, ib_bdy) = ir
723            ENDDO
724         ENDDO
725         !
726         ! ------------ U points -------------
727         igrd=2
728         icount=0
729         DO ir = 1, nn_rimwidth(ib_bdy)
730!            DO ii = jpisdt(iseg), jpisft(iseg) - 1
731            DO ii = jpisdt(iseg), jpisft(iseg)
732               icount = icount + 1
733               nbidta(icount, igrd, ib_bdy) = ii
734               nbjdta(icount, igrd, ib_bdy) = jpjsob(iseg) + ir - 1
735               nbrdta(icount, igrd, ib_bdy) = ir
736            ENDDO
737            nbidta(icount, igrd, ib_bdy) = -ib_bdy ! Discount this point
738            nbjdta(icount, igrd, ib_bdy) = -ib_bdy ! Discount this point
739         ENDDO
740         !
741         ! ------------ V points -------------
742         igrd=3
743         icount=0
744         DO ir = 1, nn_rimwidth(ib_bdy)
745            DO ii = jpisdt(iseg), jpisft(iseg)
746               icount = icount + 1
747               nbidta(icount, igrd, ib_bdy) = ii
748               nbjdta(icount, igrd, ib_bdy) = jpjsob(iseg) + ir - 1
749               nbrdta(icount, igrd, ib_bdy) = ir
750            ENDDO
751         ENDDO
752      ENDDO
753
754      !  Deal with duplicated points
755      !-----------------------------
756      ! We assign negative indices to duplicated points (to remove them from bdy points to be updated)
757      ! if their distance to the bdy is greater than the other
758      ! If their distance are the same, just keep only one to avoid updating a point twice
759      DO igrd = 1, jpbgrd
760         DO ib_bdy1 = 1, nb_bdy
761            DO ib_bdy2 = 1, nb_bdy
762               IF (ib_bdy1/=ib_bdy2) THEN
763                  DO ib1 = 1, nblendta(igrd,ib_bdy1)
764                     DO ib2 = 1, nblendta(igrd,ib_bdy2)
765                        IF ((nbidta(ib1, igrd, ib_bdy1)==nbidta(ib2, igrd, ib_bdy2)).AND. &
766                        &   (nbjdta(ib1, igrd, ib_bdy1)==nbjdta(ib2, igrd, ib_bdy2))) THEN
767!                           IF ((lwp).AND.(igrd==1)) WRITE(numout,*) ' found coincident point ji, jj:', &
768!                                                       &              nbidta(ib1, igrd, ib_bdy1),      &
769!                                                       &              nbjdta(ib2, igrd, ib_bdy2)
770                           ! keep only points with the lowest distance to boundary:
771                           IF (nbrdta(ib1, igrd, ib_bdy1)<nbrdta(ib2, igrd, ib_bdy2)) THEN
772                             nbidta(ib2, igrd, ib_bdy2) =-ib_bdy2
773                             nbjdta(ib2, igrd, ib_bdy2) =-ib_bdy2
774                           ELSEIF (nbrdta(ib1, igrd, ib_bdy1)>nbrdta(ib2, igrd, ib_bdy2)) THEN
775                             nbidta(ib1, igrd, ib_bdy1) =-ib_bdy1
776                             nbjdta(ib1, igrd, ib_bdy1) =-ib_bdy1
777                           ! Arbitrary choice if distances are the same:
778                           ELSE
779                             nbidta(ib1, igrd, ib_bdy1) =-ib_bdy1
780                             nbjdta(ib1, igrd, ib_bdy1) =-ib_bdy1
781                           ENDIF
782                        END IF
783                     END DO
784                  END DO
785               ENDIF
786            END DO
787         END DO
788      END DO
789
790      ! Work out dimensions of boundary data on each processor
791      ! ------------------------------------------------------
792
793      ! Rather assume that boundary data indices are given on global domain
794      ! TO BE DISCUSSED ?
795!      iw = mig(1) + 1            ! if monotasking and no zoom, iw=2
796!      ie = mig(1) + nlci-1 - 1   ! if monotasking and no zoom, ie=jpim1
797!      is = mjg(1) + 1            ! if monotasking and no zoom, is=2
798!      in = mjg(1) + nlcj-1 - 1   ! if monotasking and no zoom, in=jpjm1     
799      iw = mig(1) - jpizoom + 2         ! if monotasking and no zoom, iw=2
800      ie = mig(1) + nlci - jpizoom - 1  ! if monotasking and no zoom, ie=jpim1
801      is = mjg(1) - jpjzoom + 2         ! if monotasking and no zoom, is=2
802      in = mjg(1) + nlcj - jpjzoom - 1  ! if monotasking and no zoom, in=jpjm1
803
804      ALLOCATE( nbondi_bdy(nb_bdy))
805      ALLOCATE( nbondj_bdy(nb_bdy))
806      nbondi_bdy(:)=2
807      nbondj_bdy(:)=2
808      ALLOCATE( nbondi_bdy_b(nb_bdy))
809      ALLOCATE( nbondj_bdy_b(nb_bdy))
810      nbondi_bdy_b(:)=2
811      nbondj_bdy_b(:)=2
812
813      ! Work out dimensions of boundary data on each neighbour process
814      IF(nbondi .eq. 0) THEN
815         iw_b(1) = jpizoom + nimppt(nowe+1)
816         ie_b(1) = jpizoom + nimppt(nowe+1)+nlcit(nowe+1)-3
817         is_b(1) = jpjzoom + njmppt(nowe+1)
818         in_b(1) = jpjzoom + njmppt(nowe+1)+nlcjt(nowe+1)-3
819
820         iw_b(2) = jpizoom + nimppt(noea+1)
821         ie_b(2) = jpizoom + nimppt(noea+1)+nlcit(noea+1)-3
822         is_b(2) = jpjzoom + njmppt(noea+1)
823         in_b(2) = jpjzoom + njmppt(noea+1)+nlcjt(noea+1)-3
824      ELSEIF(nbondi .eq. 1) THEN
825         iw_b(1) = jpizoom + nimppt(nowe+1)
826         ie_b(1) = jpizoom + nimppt(nowe+1)+nlcit(nowe+1)-3
827         is_b(1) = jpjzoom + njmppt(nowe+1)
828         in_b(1) = jpjzoom + njmppt(nowe+1)+nlcjt(nowe+1)-3
829      ELSEIF(nbondi .eq. -1) THEN
830         iw_b(2) = jpizoom + nimppt(noea+1)
831         ie_b(2) = jpizoom + nimppt(noea+1)+nlcit(noea+1)-3
832         is_b(2) = jpjzoom + njmppt(noea+1)
833         in_b(2) = jpjzoom + njmppt(noea+1)+nlcjt(noea+1)-3
834      ENDIF
835
836      IF(nbondj .eq. 0) THEN
837         iw_b(3) = jpizoom + nimppt(noso+1)
838         ie_b(3) = jpizoom + nimppt(noso+1)+nlcit(noso+1)-3
839         is_b(3) = jpjzoom + njmppt(noso+1)
840         in_b(3) = jpjzoom + njmppt(noso+1)+nlcjt(noso+1)-3
841
842         iw_b(4) = jpizoom + nimppt(nono+1)
843         ie_b(4) = jpizoom + nimppt(nono+1)+nlcit(nono+1)-3
844         is_b(4) = jpjzoom + njmppt(nono+1)
845         in_b(4) = jpjzoom + njmppt(nono+1)+nlcjt(nono+1)-3
846      ELSEIF(nbondj .eq. 1) THEN
847         iw_b(3) = jpizoom + nimppt(noso+1)
848         ie_b(3) = jpizoom + nimppt(noso+1)+nlcit(noso+1)-3
849         is_b(3) = jpjzoom + njmppt(noso+1)
850         in_b(3) = jpjzoom + njmppt(noso+1)+nlcjt(noso+1)-3
851      ELSEIF(nbondj .eq. -1) THEN
852         iw_b(4) = jpizoom + nimppt(nono+1)
853         ie_b(4) = jpizoom + nimppt(nono+1)+nlcit(nono+1)-3
854         is_b(4) = jpjzoom + njmppt(nono+1)
855         in_b(4) = jpjzoom + njmppt(nono+1)+nlcjt(nono+1)-3
856      ENDIF
857
858      DO ib_bdy = 1, nb_bdy
859         DO igrd = 1, jpbgrd
860            icount  = 0
861            icountr = 0
862            idx_bdy(ib_bdy)%nblen(igrd)    = 0
863            idx_bdy(ib_bdy)%nblenrim(igrd) = 0
864            DO ib = 1, nblendta(igrd,ib_bdy)
865               ! check that data is in correct order in file
866               ibm1 = MAX(1,ib-1)
867               IF(lwp) THEN         ! Since all procs read global data only need to do this check on one proc...
868                  IF( nbrdta(ib,igrd,ib_bdy) < nbrdta(ibm1,igrd,ib_bdy) ) THEN
869                     CALL ctl_stop('bdy_init : ERROR : boundary data in file must be defined ', &
870                          &        ' in order of distance from edge nbr A utility for re-ordering ', &
871                          &        ' boundary coordinates and data files exists in the TOOLS/OBC directory')
872                  ENDIF   
873               ENDIF
874               ! check if point is in local domain
875               IF(  nbidta(ib,igrd,ib_bdy) >= iw .AND. nbidta(ib,igrd,ib_bdy) <= ie .AND.   &
876                  & nbjdta(ib,igrd,ib_bdy) >= is .AND. nbjdta(ib,igrd,ib_bdy) <= in       ) THEN
877                  !
878                  icount = icount  + 1
879                  !
880                  IF( nbrdta(ib,igrd,ib_bdy) == 1 )   icountr = icountr+1
881               ENDIF
882            ENDDO
883            idx_bdy(ib_bdy)%nblenrim(igrd) = icountr !: length of rim boundary data on each proc
884            idx_bdy(ib_bdy)%nblen   (igrd) = icount  !: length of boundary data on each proc       
885         ENDDO  ! igrd
886
887         ! Allocate index arrays for this boundary set
888         !--------------------------------------------
889         ilen1 = MAXVAL(idx_bdy(ib_bdy)%nblen(:))
890         ALLOCATE( idx_bdy(ib_bdy)%nbi(ilen1,jpbgrd) )
891         ALLOCATE( idx_bdy(ib_bdy)%nbj(ilen1,jpbgrd) )
892         ALLOCATE( idx_bdy(ib_bdy)%nbr(ilen1,jpbgrd) )
893         ALLOCATE( idx_bdy(ib_bdy)%nbd(ilen1,jpbgrd) )
894         ALLOCATE( idx_bdy(ib_bdy)%nbdout(ilen1,jpbgrd) )
895         ALLOCATE( idx_bdy(ib_bdy)%nbmap(ilen1,jpbgrd) )
896         ALLOCATE( idx_bdy(ib_bdy)%nbw(ilen1,jpbgrd) )
897         ALLOCATE( idx_bdy(ib_bdy)%flagu(ilen1,jpbgrd) )
898         ALLOCATE( idx_bdy(ib_bdy)%flagv(ilen1,jpbgrd) )
899
900         ! Dispatch mapping indices and discrete distances on each processor
901         ! -----------------------------------------------------------------
902
903         com_east = 0
904         com_west = 0
905         com_south = 0
906         com_north = 0
907
908         com_east_b = 0
909         com_west_b = 0
910         com_south_b = 0
911         com_north_b = 0
912         DO igrd = 1, jpbgrd
913            icount  = 0
914            ! Loop on rimwidth to ensure outermost points come first in the local arrays.
915            DO ir=1, nn_rimwidth(ib_bdy)
916               DO ib = 1, nblendta(igrd,ib_bdy)
917                  ! check if point is in local domain and equals ir
918                  IF(  nbidta(ib,igrd,ib_bdy) >= iw .AND. nbidta(ib,igrd,ib_bdy) <= ie .AND.   &
919                     & nbjdta(ib,igrd,ib_bdy) >= is .AND. nbjdta(ib,igrd,ib_bdy) <= in .AND.   &
920                     & nbrdta(ib,igrd,ib_bdy) == ir  ) THEN
921                     !
922                     icount = icount  + 1
923
924                     ! Rather assume that boundary data indices are given on global domain
925                     ! TO BE DISCUSSED ?
926!                     idx_bdy(ib_bdy)%nbi(icount,igrd)   = nbidta(ib,igrd,ib_bdy)- mig(1)+1
927!                     idx_bdy(ib_bdy)%nbj(icount,igrd)   = nbjdta(ib,igrd,ib_bdy)- mjg(1)+1
928                     idx_bdy(ib_bdy)%nbi(icount,igrd)   = nbidta(ib,igrd,ib_bdy)- mig(1)+jpizoom
929                     idx_bdy(ib_bdy)%nbj(icount,igrd)   = nbjdta(ib,igrd,ib_bdy)- mjg(1)+jpjzoom
930                     ! check if point has to be sent
931                     ii = idx_bdy(ib_bdy)%nbi(icount,igrd)
932                     ij = idx_bdy(ib_bdy)%nbj(icount,igrd)
933                     if((com_east .ne. 1) .and. (ii .eq. (nlci-1)) .and. (nbondi .le. 0)) then
934                        com_east = 1
935                     elseif((com_west .ne. 1) .and. (ii .eq. 2) .and. (nbondi .ge. 0) .and. (nbondi .ne. 2)) then
936                        com_west = 1
937                     endif
938                     if((com_south .ne. 1) .and. (ij .eq. 2) .and. (nbondj .ge. 0) .and. (nbondj .ne. 2)) then
939                        com_south = 1
940                     elseif((com_north .ne. 1) .and. (ij .eq. (nlcj-1)) .and. (nbondj .le. 0)) then
941                        com_north = 1
942                     endif
943                     idx_bdy(ib_bdy)%nbr(icount,igrd)   = nbrdta(ib,igrd,ib_bdy)
944                     idx_bdy(ib_bdy)%nbmap(icount,igrd) = ib
945                  ENDIF
946                  ! check if point has to be received from a neighbour
947                  IF(nbondi .eq. 0) THEN
948                     IF( nbidta(ib,igrd,ib_bdy) >= iw_b(1) .AND. nbidta(ib,igrd,ib_bdy) <= ie_b(1) .AND.   &
949                       & nbjdta(ib,igrd,ib_bdy) >= is_b(1) .AND. nbjdta(ib,igrd,ib_bdy) <= in_b(1) .AND.   &
950                       & nbrdta(ib,igrd,ib_bdy) == ir  ) THEN
951                       ii = nbidta(ib,igrd,ib_bdy)- iw_b(1)+2
952                       if((com_west_b .ne. 1) .and. (ii .eq. (nlcit(nowe+1)-1))) then
953                          ij = nbjdta(ib,igrd,ib_bdy) - is_b(1)+2
954                          if((ij .eq. 2) .and. (nbondj .eq. 0 .or. nbondj .eq. 1)) then
955                            com_south = 1
956                          elseif((ij .eq. nlcjt(nowe+1)-1) .and. (nbondj .eq. 0 .or. nbondj .eq. -1)) then
957                            com_north = 1
958                          endif
959                          com_west_b = 1
960                       endif
961                     ENDIF
962                     IF( nbidta(ib,igrd,ib_bdy) >= iw_b(2) .AND. nbidta(ib,igrd,ib_bdy) <= ie_b(2) .AND.   &
963                       & nbjdta(ib,igrd,ib_bdy) >= is_b(2) .AND. nbjdta(ib,igrd,ib_bdy) <= in_b(2) .AND.   &
964                       & nbrdta(ib,igrd,ib_bdy) == ir  ) THEN
965                       ii = nbidta(ib,igrd,ib_bdy)- iw_b(2)+2
966                       if((com_east_b .ne. 1) .and. (ii .eq. 2)) then
967                          ij = nbjdta(ib,igrd,ib_bdy) - is_b(2)+2
968                          if((ij .eq. 2) .and. (nbondj .eq. 0 .or. nbondj .eq. 1)) then
969                            com_south = 1
970                          elseif((ij .eq. nlcjt(noea+1)-1) .and. (nbondj .eq. 0 .or. nbondj .eq. -1)) then
971                            com_north = 1
972                          endif
973                          com_east_b = 1
974                       endif
975                     ENDIF
976                  ELSEIF(nbondi .eq. 1) THEN
977                     IF( nbidta(ib,igrd,ib_bdy) >= iw_b(1) .AND. nbidta(ib,igrd,ib_bdy) <= ie_b(1) .AND.   &
978                       & nbjdta(ib,igrd,ib_bdy) >= is_b(1) .AND. nbjdta(ib,igrd,ib_bdy) <= in_b(1) .AND.   &
979                       & nbrdta(ib,igrd,ib_bdy) == ir  ) THEN
980                       ii = nbidta(ib,igrd,ib_bdy)- iw_b(1)+2
981                       if((com_west_b .ne. 1) .and. (ii .eq. (nlcit(nowe+1)-1))) then
982                          ij = nbjdta(ib,igrd,ib_bdy) - is_b(1)+2
983                          if((ij .eq. 2) .and. (nbondj .eq. 0 .or. nbondj .eq. 1)) then
984                            com_south = 1
985                          elseif((ij .eq. nlcjt(nowe+1)-1) .and. (nbondj .eq. 0 .or. nbondj .eq. -1)) then
986                            com_north = 1
987                          endif
988                          com_west_b = 1
989                       endif
990                     ENDIF
991                  ELSEIF(nbondi .eq. -1) THEN
992                     IF( nbidta(ib,igrd,ib_bdy) >= iw_b(2) .AND. nbidta(ib,igrd,ib_bdy) <= ie_b(2) .AND.   &
993                       & nbjdta(ib,igrd,ib_bdy) >= is_b(2) .AND. nbjdta(ib,igrd,ib_bdy) <= in_b(2) .AND.   &
994                       & nbrdta(ib,igrd,ib_bdy) == ir  ) THEN
995                       ii = nbidta(ib,igrd,ib_bdy)- iw_b(2)+2
996                       if((com_east_b .ne. 1) .and. (ii .eq. 2)) then
997                          ij = nbjdta(ib,igrd,ib_bdy) - is_b(2)+2
998                          if((ij .eq. 2) .and. (nbondj .eq. 0 .or. nbondj .eq. 1)) then
999                            com_south = 1
1000                          elseif((ij .eq. nlcjt(noea+1)-1) .and. (nbondj .eq. 0 .or. nbondj .eq. -1)) then
1001                            com_north = 1
1002                          endif
1003                          com_east_b = 1
1004                       endif
1005                     ENDIF
1006                  ENDIF
1007                  IF(nbondj .eq. 0) THEN
1008                     IF(com_north_b .ne. 1 .AND. (nbidta(ib,igrd,ib_bdy) == iw_b(4)-1  &
1009                       & .OR. nbidta(ib,igrd,ib_bdy) == ie_b(4)+1) .AND. &
1010                       & nbjdta(ib,igrd,ib_bdy) == is_b(4) .AND. nbrdta(ib,igrd,ib_bdy) == ir) THEN
1011                       com_north_b = 1 
1012                     ENDIF
1013                     IF(com_south_b .ne. 1 .AND. (nbidta(ib,igrd,ib_bdy) == iw_b(3)-1  &
1014                       &.OR. nbidta(ib,igrd,ib_bdy) == ie_b(3)+1) .AND. &
1015                       & nbjdta(ib,igrd,ib_bdy) == in_b(3) .AND. nbrdta(ib,igrd,ib_bdy) == ir) THEN
1016                       com_south_b = 1 
1017                     ENDIF
1018                     IF( nbidta(ib,igrd,ib_bdy) >= iw_b(3) .AND. nbidta(ib,igrd,ib_bdy) <= ie_b(3) .AND.   &
1019                       & nbjdta(ib,igrd,ib_bdy) >= is_b(3) .AND. nbjdta(ib,igrd,ib_bdy) <= in_b(3) .AND.   &
1020                       & nbrdta(ib,igrd,ib_bdy) == ir  ) THEN
1021                       ij = nbjdta(ib,igrd,ib_bdy)- is_b(3)+2
1022                       if((com_south_b .ne. 1) .and. (ij .eq. (nlcjt(noso+1)-1))) then
1023                          com_south_b = 1
1024                       endif
1025                     ENDIF
1026                     IF( nbidta(ib,igrd,ib_bdy) >= iw_b(4) .AND. nbidta(ib,igrd,ib_bdy) <= ie_b(4) .AND.   &
1027                       & nbjdta(ib,igrd,ib_bdy) >= is_b(4) .AND. nbjdta(ib,igrd,ib_bdy) <= in_b(4) .AND.   &
1028                       & nbrdta(ib,igrd,ib_bdy) == ir  ) THEN
1029                       ij = nbjdta(ib,igrd,ib_bdy)- is_b(4)+2
1030                       if((com_north_b .ne. 1) .and. (ij .eq. 2)) then
1031                          com_north_b = 1
1032                       endif
1033                     ENDIF
1034                  ELSEIF(nbondj .eq. 1) THEN
1035                     IF( com_south_b .ne. 1 .AND. (nbidta(ib,igrd,ib_bdy) == iw_b(3)-1 .OR. &
1036                       & nbidta(ib,igrd,ib_bdy) == ie_b(3)+1) .AND. &
1037                       & nbjdta(ib,igrd,ib_bdy) == in_b(3) .AND. nbrdta(ib,igrd,ib_bdy) == ir) THEN
1038                       com_south_b = 1 
1039                     ENDIF
1040                     IF( nbidta(ib,igrd,ib_bdy) >= iw_b(3) .AND. nbidta(ib,igrd,ib_bdy) <= ie_b(3) .AND.   &
1041                       & nbjdta(ib,igrd,ib_bdy) >= is_b(3) .AND. nbjdta(ib,igrd,ib_bdy) <= in_b(3) .AND.   &
1042                       & nbrdta(ib,igrd,ib_bdy) == ir  ) THEN
1043                       ij = nbjdta(ib,igrd,ib_bdy)- is_b(3)+2
1044                       if((com_south_b .ne. 1) .and. (ij .eq. (nlcjt(noso+1)-1))) then
1045                          com_south_b = 1
1046                       endif
1047                     ENDIF
1048                  ELSEIF(nbondj .eq. -1) THEN
1049                     IF(com_north_b .ne. 1 .AND. (nbidta(ib,igrd,ib_bdy) == iw_b(4)-1  &
1050                       & .OR. nbidta(ib,igrd,ib_bdy) == ie_b(4)+1) .AND. &
1051                       & nbjdta(ib,igrd,ib_bdy) == is_b(4) .AND. nbrdta(ib,igrd,ib_bdy) == ir) THEN
1052                       com_north_b = 1 
1053                     ENDIF
1054                     IF( nbidta(ib,igrd,ib_bdy) >= iw_b(4) .AND. nbidta(ib,igrd,ib_bdy) <= ie_b(4) .AND.   &
1055                       & nbjdta(ib,igrd,ib_bdy) >= is_b(4) .AND. nbjdta(ib,igrd,ib_bdy) <= in_b(4) .AND.   &
1056                       & nbrdta(ib,igrd,ib_bdy) == ir  ) THEN
1057                       ij = nbjdta(ib,igrd,ib_bdy)- is_b(4)+2
1058                       if((com_north_b .ne. 1) .and. (ij .eq. 2)) then
1059                          com_north_b = 1
1060                       endif
1061                     ENDIF
1062                  ENDIF
1063               ENDDO
1064            ENDDO
1065         ENDDO 
1066
1067         ! definition of the i- and j- direction local boundaries arrays
1068         ! used for sending the boudaries
1069         IF((com_east .eq. 1) .and. (com_west .eq. 1)) THEN
1070            nbondi_bdy(ib_bdy) = 0
1071         ELSEIF ((com_east .eq. 1) .and. (com_west .eq. 0)) THEN
1072            nbondi_bdy(ib_bdy) = -1
1073         ELSEIF ((com_east .eq. 0) .and. (com_west .eq. 1)) THEN
1074            nbondi_bdy(ib_bdy) = 1
1075         ENDIF
1076
1077         IF((com_north .eq. 1) .and. (com_south .eq. 1)) THEN
1078            nbondj_bdy(ib_bdy) = 0
1079         ELSEIF ((com_north .eq. 1) .and. (com_south .eq. 0)) THEN
1080            nbondj_bdy(ib_bdy) = -1
1081         ELSEIF ((com_north .eq. 0) .and. (com_south .eq. 1)) THEN
1082            nbondj_bdy(ib_bdy) = 1
1083         ENDIF
1084
1085         ! definition of the i- and j- direction local boundaries arrays
1086         ! used for receiving the boudaries
1087         IF((com_east_b .eq. 1) .and. (com_west_b .eq. 1)) THEN
1088            nbondi_bdy_b(ib_bdy) = 0
1089         ELSEIF ((com_east_b .eq. 1) .and. (com_west_b .eq. 0)) THEN
1090            nbondi_bdy_b(ib_bdy) = -1
1091         ELSEIF ((com_east_b .eq. 0) .and. (com_west_b .eq. 1)) THEN
1092            nbondi_bdy_b(ib_bdy) = 1
1093         ENDIF
1094
1095         IF((com_north_b .eq. 1) .and. (com_south_b .eq. 1)) THEN
1096            nbondj_bdy_b(ib_bdy) = 0
1097         ELSEIF ((com_north_b .eq. 1) .and. (com_south_b .eq. 0)) THEN
1098            nbondj_bdy_b(ib_bdy) = -1
1099         ELSEIF ((com_north_b .eq. 0) .and. (com_south_b .eq. 1)) THEN
1100            nbondj_bdy_b(ib_bdy) = 1
1101         ENDIF
1102
1103         ! Compute rim weights for FRS scheme
1104         ! ----------------------------------
1105         DO igrd = 1, jpbgrd
1106            DO ib = 1, idx_bdy(ib_bdy)%nblen(igrd)
1107               nbr => idx_bdy(ib_bdy)%nbr(ib,igrd)
1108               idx_bdy(ib_bdy)%nbw(ib,igrd) = 1.- TANH( FLOAT( nbr - 1 ) *0.5 )      ! tanh formulation
1109!               idx_bdy(ib_bdy)%nbw(ib,igrd) = (FLOAT(nn_rimwidth(ib_bdy)+1-nbr)/FLOAT(nn_rimwidth(ib_bdy)))**2.  ! quadratic
1110!               idx_bdy(ib_bdy)%nbw(ib,igrd) =  FLOAT(nn_rimwidth(ib_bdy)+1-nbr)/FLOAT(nn_rimwidth(ib_bdy))       ! linear
1111            END DO
1112         END DO 
1113
1114         ! Compute damping coefficients
1115         ! ----------------------------
1116         DO igrd = 1, jpbgrd
1117            DO ib = 1, idx_bdy(ib_bdy)%nblen(igrd)
1118               nbr => idx_bdy(ib_bdy)%nbr(ib,igrd)
1119               idx_bdy(ib_bdy)%nbd(ib,igrd) = 1. / ( rn_time_dmp(ib_bdy) * rday ) & 
1120               & *(FLOAT(nn_rimwidth(ib_bdy)+1-nbr)/FLOAT(nn_rimwidth(ib_bdy)))**2.   ! quadratic
1121               idx_bdy(ib_bdy)%nbdout(ib,igrd) = 1. / ( rn_time_dmp_out(ib_bdy) * rday ) & 
1122               & *(FLOAT(nn_rimwidth(ib_bdy)+1-nbr)/FLOAT(nn_rimwidth(ib_bdy)))**2.   ! quadratic
1123            END DO
1124         END DO
1125
1126      ENDDO
1127
1128      ! ------------------------------------------------------
1129      ! Initialise masks and find normal/tangential directions
1130      ! ------------------------------------------------------
1131
1132      ! Read global 2D mask at T-points: bdytmask
1133      ! -----------------------------------------
1134      ! bdytmask = 1  on the computational domain AND on open boundaries
1135      !          = 0  elsewhere   
1136 
1137      IF( ln_mask_file ) THEN
1138         CALL iom_open( cn_mask_file, inum )
1139         CALL iom_get ( inum, jpdom_data, 'bdy_msk', bdytmask(:,:) )
1140         CALL iom_close( inum )
1141
1142         ! Derive mask on U and V grid from mask on T grid
1143         bdyumask(:,:) = 0.e0
1144         bdyvmask(:,:) = 0.e0
1145         DO ij=1, jpjm1
1146            DO ii=1, jpim1
1147               bdyumask(ii,ij)=bdytmask(ii,ij)*bdytmask(ii+1, ij )
1148               bdyvmask(ii,ij)=bdytmask(ii,ij)*bdytmask(ii  ,ij+1) 
1149            END DO
1150         END DO
1151         CALL lbc_lnk( bdyumask(:,:), 'U', 1. )   ;   CALL lbc_lnk( bdyvmask(:,:), 'V', 1. )      ! Lateral boundary cond.
1152
1153
1154         ! Mask corrections
1155         ! ----------------
1156         DO ik = 1, jpkm1
1157            DO ij = 1, jpj
1158               DO ii = 1, jpi
1159                  tmask(ii,ij,ik) = tmask(ii,ij,ik) * bdytmask(ii,ij)
1160                  umask(ii,ij,ik) = umask(ii,ij,ik) * bdyumask(ii,ij)
1161                  vmask(ii,ij,ik) = vmask(ii,ij,ik) * bdyvmask(ii,ij)
1162                  bmask(ii,ij)    = bmask(ii,ij)    * bdytmask(ii,ij)
1163               END DO     
1164            END DO
1165         END DO
1166
1167         DO ik = 1, jpkm1
1168            DO ij = 2, jpjm1
1169               DO ii = 2, jpim1
1170                  fmask(ii,ij,ik) = fmask(ii,ij,ik) * bdytmask(ii,ij  ) * bdytmask(ii+1,ij  )   &
1171                     &                              * bdytmask(ii,ij+1) * bdytmask(ii+1,ij+1)
1172               END DO     
1173            END DO
1174         END DO
1175
1176         tmask_i (:,:) = ssmask(:,:) * tmask_i(:,:)
1177
1178      ENDIF ! ln_mask_file=.TRUE.
1179     
1180      bdytmask(:,:) = ssmask(:,:)
1181      IF( .not. ln_mask_file ) THEN
1182         ! If .not. ln_mask_file then we need to derive mask on U and V grid
1183         ! from mask on T grid here.
1184         bdyumask(:,:) = 0.e0
1185         bdyvmask(:,:) = 0.e0
1186         DO ij=1, jpjm1
1187            DO ii=1, jpim1
1188               bdyumask(ii,ij)=bdytmask(ii,ij)*bdytmask(ii+1, ij )
1189               bdyvmask(ii,ij)=bdytmask(ii,ij)*bdytmask(ii  ,ij+1) 
1190            END DO
1191         END DO
1192         CALL lbc_lnk( bdyumask(:,:), 'U', 1. )   ;   CALL lbc_lnk( bdyvmask(:,:), 'V', 1. )      ! Lateral boundary cond.
1193      ENDIF
1194
1195      ! bdy masks and bmask are now set to zero on boundary points:
1196      igrd = 1       ! In the free surface case, bmask is at T-points
1197      DO ib_bdy = 1, nb_bdy
1198        DO ib = 1, idx_bdy(ib_bdy)%nblenrim(igrd)     
1199          bmask(idx_bdy(ib_bdy)%nbi(ib,igrd), idx_bdy(ib_bdy)%nbj(ib,igrd)) = 0.e0
1200        ENDDO
1201      ENDDO
1202      !
1203      igrd = 1
1204      DO ib_bdy = 1, nb_bdy
1205        DO ib = 1, idx_bdy(ib_bdy)%nblenrim(igrd)     
1206          bdytmask(idx_bdy(ib_bdy)%nbi(ib,igrd), idx_bdy(ib_bdy)%nbj(ib,igrd)) = 0.e0
1207        ENDDO
1208      ENDDO
1209      !
1210      igrd = 2
1211      DO ib_bdy = 1, nb_bdy
1212        DO ib = 1, idx_bdy(ib_bdy)%nblenrim(igrd)
1213          bdyumask(idx_bdy(ib_bdy)%nbi(ib,igrd), idx_bdy(ib_bdy)%nbj(ib,igrd)) = 0.e0
1214        ENDDO
1215      ENDDO
1216      !
1217      igrd = 3
1218      DO ib_bdy = 1, nb_bdy
1219        DO ib = 1, idx_bdy(ib_bdy)%nblenrim(igrd)
1220          bdyvmask(idx_bdy(ib_bdy)%nbi(ib,igrd), idx_bdy(ib_bdy)%nbj(ib,igrd)) = 0.e0
1221        ENDDO
1222      ENDDO
1223
1224      ! For the flagu/flagv calculation below we require a version of fmask without
1225      ! the land boundary condition (shlat) included:
1226      CALL wrk_alloc(jpi,jpj,zfmask) 
1227      DO ij = 2, jpjm1
1228         DO ii = 2, jpim1
1229            zfmask(ii,ij) = tmask(ii,ij  ,1) * tmask(ii+1,ij  ,1)   &
1230           &              * tmask(ii,ij+1,1) * tmask(ii+1,ij+1,1)
1231         END DO     
1232      END DO
1233
1234      ! Lateral boundary conditions
1235      CALL lbc_lnk( zfmask       , 'F', 1. )
1236      CALL lbc_lnk( fmask        , 'F', 1. )   ;   CALL lbc_lnk( bdytmask(:,:), 'T', 1. )
1237      CALL lbc_lnk( bdyumask(:,:), 'U', 1. )   ;   CALL lbc_lnk( bdyvmask(:,:), 'V', 1. )
1238
1239      DO ib_bdy = 1, nb_bdy       ! Indices and directions of rim velocity components
1240
1241         idx_bdy(ib_bdy)%flagu(:,:) = 0.e0
1242         idx_bdy(ib_bdy)%flagv(:,:) = 0.e0
1243         icount = 0 
1244
1245         ! Calculate relationship of U direction to the local orientation of the boundary
1246         ! flagu = -1 : u component is normal to the dynamical boundary and its direction is outward
1247         ! flagu =  0 : u is tangential
1248         ! flagu =  1 : u is normal to the boundary and is direction is inward
1249 
1250         DO igrd = 1,jpbgrd 
1251            SELECT CASE( igrd )
1252               CASE( 1 )
1253                  pmask => umask(:,:,1)
1254                  i_offset = 0
1255               CASE( 2 ) 
1256                  pmask => bdytmask
1257                  i_offset = 1
1258               CASE( 3 ) 
1259                  pmask => zfmask(:,:)
1260                  i_offset = 0
1261            END SELECT
1262            icount = 0
1263            DO ib = 1, idx_bdy(ib_bdy)%nblenrim(igrd) 
1264               nbi => idx_bdy(ib_bdy)%nbi(ib,igrd)
1265               nbj => idx_bdy(ib_bdy)%nbj(ib,igrd)
1266               zefl = pmask(nbi+i_offset-1,nbj)
1267               zwfl = pmask(nbi+i_offset,nbj)
1268               ! This error check only works if you are using the bdyXmask arrays
1269               IF( i_offset == 1 .and. zefl + zwfl == 2 ) THEN
1270                  icount = icount + 1
1271                  IF(lwp) WRITE(numout,*) 'Problem with igrd = ',igrd,' at (global) nbi, nbj : ',mig(nbi),mjg(nbj)
1272               ELSE
1273                  idx_bdy(ib_bdy)%flagu(ib,igrd) = -zefl + zwfl
1274               ENDIF
1275            END DO
1276            IF( icount /= 0 ) THEN
1277               IF(lwp) WRITE(numout,*)
1278               IF(lwp) WRITE(numout,*) ' E R R O R : Some ',cgrid(igrd),' grid points,',   &
1279                  ' are not boundary points (flagu calculation). Check nbi, nbj, indices for boundary set ',ib_bdy
1280               IF(lwp) WRITE(numout,*) ' ========== '
1281               IF(lwp) WRITE(numout,*)
1282               nstop = nstop + 1
1283            ENDIF
1284         END DO
1285
1286         ! Calculate relationship of V direction to the local orientation of the boundary
1287         ! flagv = -1 : v component is normal to the dynamical boundary but its direction is outward
1288         ! flagv =  0 : v is tangential
1289         ! flagv =  1 : v is normal to the boundary and is direction is inward
1290
1291         DO igrd = 1,jpbgrd 
1292            SELECT CASE( igrd )
1293               CASE( 1 )
1294                  pmask => vmask(:,:,1)
1295                  j_offset = 0
1296               CASE( 2 )
1297                  pmask => zfmask(:,:)
1298                  j_offset = 0
1299               CASE( 3 )
1300                  pmask => bdytmask
1301                  j_offset = 1
1302            END SELECT
1303            icount = 0
1304            DO ib = 1, idx_bdy(ib_bdy)%nblenrim(igrd) 
1305               nbi => idx_bdy(ib_bdy)%nbi(ib,igrd)
1306               nbj => idx_bdy(ib_bdy)%nbj(ib,igrd)
1307               znfl = pmask(nbi,nbj+j_offset-1  )
1308               zsfl = pmask(nbi,nbj+j_offset)
1309               ! This error check only works if you are using the bdyXmask arrays
1310               IF( j_offset == 1 .and. znfl + zsfl == 2 ) THEN
1311                  IF(lwp) WRITE(numout,*) 'Problem with igrd = ',igrd,' at (global) nbi, nbj : ',mig(nbi),mjg(nbj)
1312                  icount = icount + 1
1313               ELSE
1314                  idx_bdy(ib_bdy)%flagv(ib,igrd) = -znfl + zsfl
1315               END IF
1316            END DO
1317            IF( icount /= 0 ) THEN
1318               IF(lwp) WRITE(numout,*)
1319               IF(lwp) WRITE(numout,*) ' E R R O R : Some ',cgrid(igrd),' grid points,',   &
1320                  ' are not boundary points (flagv calculation). Check nbi, nbj, indices for boundary set ',ib_bdy
1321               IF(lwp) WRITE(numout,*) ' ========== '
1322               IF(lwp) WRITE(numout,*)
1323               nstop = nstop + 1
1324            ENDIF
1325         END DO
1326
1327      END DO
1328
1329      ! Compute total lateral surface for volume correction:
1330      ! ----------------------------------------------------
1331      ! JC: this must be done at each time step with key_vvl
1332      bdysurftot = 0.e0 
1333      IF( ln_vol ) THEN 
1334         igrd = 2      ! Lateral surface at U-points
1335         DO ib_bdy = 1, nb_bdy
1336            DO ib = 1, idx_bdy(ib_bdy)%nblenrim(igrd)
1337               nbi => idx_bdy(ib_bdy)%nbi(ib,igrd)
1338               nbj => idx_bdy(ib_bdy)%nbj(ib,igrd)
1339               flagu => idx_bdy(ib_bdy)%flagu(ib,igrd)
1340               bdysurftot = bdysurftot + hu     (nbi  , nbj)                           &
1341                  &                    * e2u    (nbi  , nbj) * ABS( flagu ) &
1342                  &                    * tmask_i(nbi  , nbj)                           &
1343                  &                    * tmask_i(nbi+1, nbj)                   
1344            ENDDO
1345         ENDDO
1346
1347         igrd=3 ! Add lateral surface at V-points
1348         DO ib_bdy = 1, nb_bdy
1349            DO ib = 1, idx_bdy(ib_bdy)%nblenrim(igrd)
1350               nbi => idx_bdy(ib_bdy)%nbi(ib,igrd)
1351               nbj => idx_bdy(ib_bdy)%nbj(ib,igrd)
1352               flagv => idx_bdy(ib_bdy)%flagv(ib,igrd)
1353               bdysurftot = bdysurftot + hv     (nbi, nbj  )                           &
1354                  &                    * e1v    (nbi, nbj  ) * ABS( flagv ) &
1355                  &                    * tmask_i(nbi, nbj  )                           &
1356                  &                    * tmask_i(nbi, nbj+1)
1357            ENDDO
1358         ENDDO
1359         !
1360         IF( lk_mpp )   CALL mpp_sum( bdysurftot )      ! sum over the global domain
1361      END IF   
1362      !
1363      ! Tidy up
1364      !--------
1365      IF (nb_bdy>0) THEN
1366         DEALLOCATE(nbidta, nbjdta, nbrdta)
1367      ENDIF
1368
1369      CALL wrk_dealloc(jpi,jpj,zfmask) 
1370
1371      IF( nn_timing == 1 ) CALL timing_stop('bdy_init')
1372
1373   END SUBROUTINE bdy_init
1374
1375   SUBROUTINE bdy_ctl_seg
1376      !!----------------------------------------------------------------------
1377      !!                 ***  ROUTINE bdy_ctl_seg  ***
1378      !!
1379      !! ** Purpose :   Check straight open boundary segments location
1380      !!
1381      !! ** Method  :   - Look for open boundary corners
1382      !!                - Check that segments start or end on land
1383      !!----------------------------------------------------------------------
1384      INTEGER  ::   ib, ib1, ib2, ji ,jj, itest 
1385      INTEGER, DIMENSION(jp_nseg,2) :: icorne, icornw, icornn, icorns 
1386      REAL(wp), DIMENSION(2) ::   ztestmask
1387      !!----------------------------------------------------------------------
1388      !
1389      IF (lwp) WRITE(numout,*) ' '
1390      IF (lwp) WRITE(numout,*) 'bdy_ctl_seg: Check analytical segments'
1391      IF (lwp) WRITE(numout,*) '~~~~~~~~~~~~'
1392      !
1393      IF(lwp) WRITE(numout,*) 'Number of east  segments     : ', nbdysege
1394      IF(lwp) WRITE(numout,*) 'Number of west  segments     : ', nbdysegw
1395      IF(lwp) WRITE(numout,*) 'Number of north segments     : ', nbdysegn
1396      IF(lwp) WRITE(numout,*) 'Number of south segments     : ', nbdysegs
1397      ! 1. Check bounds
1398      !----------------
1399      DO ib = 1, nbdysegn
1400         IF (lwp) WRITE(numout,*) '**check north seg bounds pckg: ', npckgn(ib)
1401         IF ((jpjnob(ib).ge.jpjglo-1).or.& 
1402            &(jpjnob(ib).le.1))        CALL ctl_stop( 'nbdyind out of domain' )
1403         IF (jpindt(ib).ge.jpinft(ib)) CALL ctl_stop( 'Bdy start index is greater than end index' )
1404         IF (jpindt(ib).le.1     )     CALL ctl_stop( 'Start index out of domain' )
1405         IF (jpinft(ib).ge.jpiglo)     CALL ctl_stop( 'End index out of domain' )
1406      END DO
1407      !
1408      DO ib = 1, nbdysegs
1409         IF (lwp) WRITE(numout,*) '**check south seg bounds pckg: ', npckgs(ib)
1410         IF ((jpjsob(ib).ge.jpjglo-1).or.& 
1411            &(jpjsob(ib).le.1))        CALL ctl_stop( 'nbdyind out of domain' )
1412         IF (jpisdt(ib).ge.jpisft(ib)) CALL ctl_stop( 'Bdy start index is greater than end index' )
1413         IF (jpisdt(ib).le.1     )     CALL ctl_stop( 'Start index out of domain' )
1414         IF (jpisft(ib).ge.jpiglo)     CALL ctl_stop( 'End index out of domain' )
1415      END DO
1416      !
1417      DO ib = 1, nbdysege
1418         IF (lwp) WRITE(numout,*) '**check east  seg bounds pckg: ', npckge(ib)
1419         IF ((jpieob(ib).ge.jpiglo-1).or.& 
1420            &(jpieob(ib).le.1))        CALL ctl_stop( 'nbdyind out of domain' )
1421         IF (jpjedt(ib).ge.jpjeft(ib)) CALL ctl_stop( 'Bdy start index is greater than end index' )
1422         IF (jpjedt(ib).le.1     )     CALL ctl_stop( 'Start index out of domain' )
1423         IF (jpjeft(ib).ge.jpjglo)     CALL ctl_stop( 'End index out of domain' )
1424      END DO
1425      !
1426      DO ib = 1, nbdysegw
1427         IF (lwp) WRITE(numout,*) '**check west  seg bounds pckg: ', npckgw(ib)
1428         IF ((jpiwob(ib).ge.jpiglo-1).or.& 
1429            &(jpiwob(ib).le.1))        CALL ctl_stop( 'nbdyind out of domain' )
1430         IF (jpjwdt(ib).ge.jpjwft(ib)) CALL ctl_stop( 'Bdy start index is greater than end index' )
1431         IF (jpjwdt(ib).le.1     )     CALL ctl_stop( 'Start index out of domain' )
1432         IF (jpjwft(ib).ge.jpjglo)     CALL ctl_stop( 'End index out of domain' )
1433      ENDDO
1434      !
1435      !     
1436      ! 2. Look for segment crossings
1437      !------------------------------
1438      IF (lwp) WRITE(numout,*) '**Look for segments corners  :'
1439      !
1440      itest = 0 ! corner number
1441      !
1442      ! flag to detect if start or end of open boundary belongs to a corner
1443      ! if not (=0), it must be on land.
1444      ! if a corner is detected, save bdy package number for further tests
1445      icorne(:,:)=0. ; icornw(:,:)=0. ; icornn(:,:)=0. ; icorns(:,:)=0.
1446      ! South/West crossings
1447      IF ((nbdysegw > 0).AND.(nbdysegs > 0)) THEN
1448         DO ib1 = 1, nbdysegw       
1449            DO ib2 = 1, nbdysegs
1450               IF (( jpisdt(ib2)<=jpiwob(ib1)).AND. &
1451                &  ( jpisft(ib2)>=jpiwob(ib1)).AND. &
1452                &  ( jpjwdt(ib1)<=jpjsob(ib2)).AND. &
1453                &  ( jpjwft(ib1)>=jpjsob(ib2))) THEN
1454                  IF ((jpjwdt(ib1)==jpjsob(ib2)).AND.(jpisdt(ib2)==jpiwob(ib1))) THEN 
1455                     ! We have a possible South-West corner                     
1456!                     WRITE(numout,*) ' Found a South-West corner at (i,j): ', jpisdt(ib2), jpjwdt(ib1)
1457!                     WRITE(numout,*) ' between segments: ', npckgw(ib1), npckgs(ib2)
1458                     icornw(ib1,1) = npckgs(ib2)
1459                     icorns(ib2,1) = npckgw(ib1)
1460                  ELSEIF ((jpisft(ib2)==jpiwob(ib1)).AND.(jpjwft(ib1)==jpjsob(ib2))) THEN
1461                     IF(lwp) WRITE(numout,*)
1462                     IF(lwp) WRITE(numout,*) ' E R R O R : Found an acute open boundary corner at point (i,j)= ', &
1463                     &                                     jpisft(ib2), jpjwft(ib1)
1464                     IF(lwp) WRITE(numout,*) ' ==========  Not allowed yet'
1465                     IF(lwp) WRITE(numout,*) '             Crossing problem with West segment: ',npckgw(ib1), & 
1466                     &                                                    ' and South segment: ',npckgs(ib2)
1467                     IF(lwp) WRITE(numout,*)
1468                     nstop = nstop + 1
1469                  ELSE
1470                     IF(lwp) WRITE(numout,*)
1471                     IF(lwp) WRITE(numout,*) ' E R R O R : Check South and West Open boundary indices'
1472                     IF(lwp) WRITE(numout,*) ' ==========  Crossing problem with West segment: ',npckgw(ib1) , &
1473                     &                                                    ' and South segment: ',npckgs(ib2)
1474                     IF(lwp) WRITE(numout,*)
1475                     nstop = nstop+1
1476                  END IF
1477               END IF
1478            END DO
1479         END DO
1480      END IF
1481      !
1482      ! South/East crossings
1483      IF ((nbdysege > 0).AND.(nbdysegs > 0)) THEN
1484         DO ib1 = 1, nbdysege
1485            DO ib2 = 1, nbdysegs
1486               IF (( jpisdt(ib2)<=jpieob(ib1)+1).AND. &
1487                &  ( jpisft(ib2)>=jpieob(ib1)+1).AND. &
1488                &  ( jpjedt(ib1)<=jpjsob(ib2)  ).AND. &
1489                &  ( jpjeft(ib1)>=jpjsob(ib2)  )) THEN
1490                  IF ((jpjedt(ib1)==jpjsob(ib2)).AND.(jpisft(ib2)==jpieob(ib1)+1)) THEN
1491                     ! We have a possible South-East corner
1492!                     WRITE(numout,*) ' Found a South-East corner at (i,j): ', jpisft(ib2), jpjedt(ib1)
1493!                     WRITE(numout,*) ' between segments: ', npckge(ib1), npckgs(ib2)
1494                     icorne(ib1,1) = npckgs(ib2)
1495                     icorns(ib2,2) = npckge(ib1)
1496                  ELSEIF ((jpjeft(ib1)==jpjsob(ib2)).AND.(jpisdt(ib2)==jpieob(ib1)+1)) THEN
1497                     IF(lwp) WRITE(numout,*)
1498                     IF(lwp) WRITE(numout,*) ' E R R O R : Found an acute open boundary corner at point (i,j)= ', &
1499                     &                                     jpisdt(ib2), jpjeft(ib1)
1500                     IF(lwp) WRITE(numout,*) ' ==========  Not allowed yet'
1501                     IF(lwp) WRITE(numout,*) '             Crossing problem with East segment: ',npckge(ib1), &
1502                     &                                                    ' and South segment: ',npckgs(ib2)
1503                     IF(lwp) WRITE(numout,*)
1504                     nstop = nstop + 1
1505                  ELSE
1506                     IF(lwp) WRITE(numout,*)
1507                     IF(lwp) WRITE(numout,*) ' E R R O R : Check South and East Open boundary indices'
1508                     IF(lwp) WRITE(numout,*) ' ==========  Crossing problem with East segment: ',npckge(ib1), &
1509                     &                                                    ' and South segment: ',npckgs(ib2)
1510                     IF(lwp) WRITE(numout,*)
1511                     nstop = nstop + 1
1512                  END IF
1513               END IF
1514            END DO
1515         END DO
1516      END IF
1517      !
1518      ! North/West crossings
1519      IF ((nbdysegn > 0).AND.(nbdysegw > 0)) THEN
1520         DO ib1 = 1, nbdysegw       
1521            DO ib2 = 1, nbdysegn
1522               IF (( jpindt(ib2)<=jpiwob(ib1)  ).AND. &
1523                &  ( jpinft(ib2)>=jpiwob(ib1)  ).AND. &
1524                &  ( jpjwdt(ib1)<=jpjnob(ib2)+1).AND. &
1525                &  ( jpjwft(ib1)>=jpjnob(ib2)+1)) THEN
1526                  IF ((jpjwft(ib1)==jpjnob(ib2)+1).AND.(jpindt(ib2)==jpiwob(ib1))) THEN
1527                     ! We have a possible North-West corner
1528!                     WRITE(numout,*) ' Found a North-West corner at (i,j): ', jpindt(ib2), jpjwft(ib1)
1529!                     WRITE(numout,*) ' between segments: ', npckgw(ib1), npckgn(ib2)
1530                     icornw(ib1,2) = npckgn(ib2)
1531                     icornn(ib2,1) = npckgw(ib1)
1532                  ELSEIF ((jpjwdt(ib1)==jpjnob(ib2)+1).AND.(jpinft(ib2)==jpiwob(ib1))) THEN
1533                     IF(lwp) WRITE(numout,*)
1534                     IF(lwp) WRITE(numout,*) ' E R R O R : Found an acute open boundary corner at point (i,j)= ', &
1535                     &                                     jpinft(ib2), jpjwdt(ib1)
1536                     IF(lwp) WRITE(numout,*) ' ==========  Not allowed yet'
1537                     IF(lwp) WRITE(numout,*) '             Crossing problem with West segment: ',npckgw(ib1), &
1538                     &                                                    ' and North segment: ',npckgn(ib2)
1539                     IF(lwp) WRITE(numout,*)
1540                     nstop = nstop + 1
1541                  ELSE
1542                     IF(lwp) WRITE(numout,*)
1543                     IF(lwp) WRITE(numout,*) ' E R R O R : Check North and West Open boundary indices'
1544                     IF(lwp) WRITE(numout,*) ' ==========  Crossing problem with West segment: ',npckgw(ib1), &
1545                     &                                                    ' and North segment: ',npckgn(ib2)
1546                     IF(lwp) WRITE(numout,*)
1547                     nstop = nstop + 1
1548                  END IF
1549               END IF
1550            END DO
1551         END DO
1552      END IF
1553      !
1554      ! North/East crossings
1555      IF ((nbdysegn > 0).AND.(nbdysege > 0)) THEN
1556         DO ib1 = 1, nbdysege       
1557            DO ib2 = 1, nbdysegn
1558               IF (( jpindt(ib2)<=jpieob(ib1)+1).AND. &
1559                &  ( jpinft(ib2)>=jpieob(ib1)+1).AND. &
1560                &  ( jpjedt(ib1)<=jpjnob(ib2)+1).AND. &
1561                &  ( jpjeft(ib1)>=jpjnob(ib2)+1)) THEN
1562                  IF ((jpjeft(ib1)==jpjnob(ib2)+1).AND.(jpinft(ib2)==jpieob(ib1)+1)) THEN
1563                     ! We have a possible North-East corner
1564!                     WRITE(numout,*) ' Found a North-East corner at (i,j): ', jpinft(ib2), jpjeft(ib1)
1565!                     WRITE(numout,*) ' between segments: ', npckge(ib1), npckgn(ib2)
1566                     icorne(ib1,2) = npckgn(ib2)
1567                     icornn(ib2,2) = npckge(ib1)
1568                  ELSEIF ((jpjedt(ib1)==jpjnob(ib2)+1).AND.(jpindt(ib2)==jpieob(ib1)+1)) THEN
1569                     IF(lwp) WRITE(numout,*)
1570                     IF(lwp) WRITE(numout,*) ' E R R O R : Found an acute open boundary corner at point (i,j)= ', &
1571                     &                                     jpindt(ib2), jpjedt(ib1)
1572                     IF(lwp) WRITE(numout,*) ' ==========  Not allowed yet'
1573                     IF(lwp) WRITE(numout,*) '             Crossing problem with East segment: ',npckge(ib1), &
1574                     &                                                    ' and North segment: ',npckgn(ib2)
1575                     IF(lwp) WRITE(numout,*)
1576                     nstop = nstop + 1
1577                  ELSE
1578                     IF(lwp) WRITE(numout,*)
1579                     IF(lwp) WRITE(numout,*) ' E R R O R : Check North and East Open boundary indices'
1580                     IF(lwp) WRITE(numout,*) ' ==========  Crossing problem with East segment: ',npckge(ib1), &
1581                     &                                                    ' and North segment: ',npckgn(ib2)
1582                     IF(lwp) WRITE(numout,*)
1583                     nstop = nstop + 1
1584                  END IF
1585               END IF
1586            END DO
1587         END DO
1588      END IF
1589      !
1590      ! 3. Check if segment extremities are on land
1591      !--------------------------------------------
1592      !
1593      ! West segments
1594      DO ib = 1, nbdysegw
1595         ! get mask at boundary extremities:
1596         ztestmask(1:2)=0.
1597         DO ji = 1, jpi
1598            DO jj = 1, jpj             
1599              IF (((ji + nimpp - 1) == jpiwob(ib)).AND. & 
1600               &  ((jj + njmpp - 1) == jpjwdt(ib))) ztestmask(1)=tmask(ji,jj,1)
1601              IF (((ji + nimpp - 1) == jpiwob(ib)).AND. & 
1602               &  ((jj + njmpp - 1) == jpjwft(ib))) ztestmask(2)=tmask(ji,jj,1) 
1603            END DO
1604         END DO
1605         IF( lk_mpp )   CALL mpp_sum( ztestmask, 2 )   ! sum over the global domain
1606
1607         IF (ztestmask(1)==1) THEN
1608            IF (icornw(ib,1)==0) THEN
1609               IF(lwp) WRITE(numout,*)
1610               IF(lwp) WRITE(numout,*) ' E R R O R : Open boundary segment ', npckgw(ib)
1611               IF(lwp) WRITE(numout,*) ' ==========  does not start on land or on a corner'                                                 
1612               IF(lwp) WRITE(numout,*)
1613               nstop = nstop + 1
1614            ELSE
1615               ! This is a corner
1616               WRITE(numout,*) 'Found a South-West corner at (i,j): ', jpiwob(ib), jpjwdt(ib)
1617               CALL bdy_ctl_corn(npckgw(ib), icornw(ib,1))
1618               itest=itest+1
1619            ENDIF
1620         ENDIF
1621         IF (ztestmask(2)==1) THEN
1622            IF (icornw(ib,2)==0) THEN
1623               IF(lwp) WRITE(numout,*)
1624               IF(lwp) WRITE(numout,*) ' E R R O R : Open boundary segment ', npckgw(ib)
1625               IF(lwp) WRITE(numout,*) ' ==========  does not end on land or on a corner'                                                 
1626               IF(lwp) WRITE(numout,*)
1627               nstop = nstop + 1
1628            ELSE
1629               ! This is a corner
1630               WRITE(numout,*) 'Found a North-West corner at (i,j): ', jpiwob(ib), jpjwft(ib)
1631               CALL bdy_ctl_corn(npckgw(ib), icornw(ib,2))
1632               itest=itest+1
1633            ENDIF
1634         ENDIF
1635      END DO
1636      !
1637      ! East segments
1638      DO ib = 1, nbdysege
1639         ! get mask at boundary extremities:
1640         ztestmask(1:2)=0.
1641         DO ji = 1, jpi
1642            DO jj = 1, jpj             
1643              IF (((ji + nimpp - 1) == jpieob(ib)+1).AND. & 
1644               &  ((jj + njmpp - 1) == jpjedt(ib))) ztestmask(1)=tmask(ji,jj,1)
1645              IF (((ji + nimpp - 1) == jpieob(ib)+1).AND. & 
1646               &  ((jj + njmpp - 1) == jpjeft(ib))) ztestmask(2)=tmask(ji,jj,1) 
1647            END DO
1648         END DO
1649         IF( lk_mpp )   CALL mpp_sum( ztestmask, 2 )   ! sum over the global domain
1650
1651         IF (ztestmask(1)==1) THEN
1652            IF (icorne(ib,1)==0) THEN
1653               IF(lwp) WRITE(numout,*)
1654               IF(lwp) WRITE(numout,*) ' E R R O R : Open boundary segment ', npckge(ib)
1655               IF(lwp) WRITE(numout,*) ' ==========  does not start on land or on a corner'                                                 
1656               IF(lwp) WRITE(numout,*)
1657               nstop = nstop + 1 
1658            ELSE
1659               ! This is a corner
1660               WRITE(numout,*) 'Found a South-East corner at (i,j): ', jpieob(ib)+1, jpjedt(ib)
1661               CALL bdy_ctl_corn(npckge(ib), icorne(ib,1))
1662               itest=itest+1
1663            ENDIF
1664         ENDIF
1665         IF (ztestmask(2)==1) THEN
1666            IF (icorne(ib,2)==0) THEN
1667               IF(lwp) WRITE(numout,*)
1668               IF(lwp) WRITE(numout,*) ' E R R O R : Open boundary segment ', npckge(ib)
1669               IF(lwp) WRITE(numout,*) ' ==========  does not end on land or on a corner'                                                 
1670               IF(lwp) WRITE(numout,*)
1671               nstop = nstop + 1
1672            ELSE
1673               ! This is a corner
1674               WRITE(numout,*) 'Found a North-East corner at (i,j): ', jpieob(ib)+1, jpjeft(ib)
1675               CALL bdy_ctl_corn(npckge(ib), icorne(ib,2))
1676               itest=itest+1
1677            ENDIF
1678         ENDIF
1679      END DO
1680      !
1681      ! South segments
1682      DO ib = 1, nbdysegs
1683         ! get mask at boundary extremities:
1684         ztestmask(1:2)=0.
1685         DO ji = 1, jpi
1686            DO jj = 1, jpj             
1687              IF (((jj + njmpp - 1) == jpjsob(ib)).AND. & 
1688               &  ((ji + nimpp - 1) == jpisdt(ib))) ztestmask(1)=tmask(ji,jj,1)
1689              IF (((jj + njmpp - 1) == jpjsob(ib)).AND. & 
1690               &  ((ji + nimpp - 1) == jpisft(ib))) ztestmask(2)=tmask(ji,jj,1) 
1691            END DO
1692         END DO
1693         IF( lk_mpp )   CALL mpp_sum( ztestmask, 2 )   ! sum over the global domain
1694
1695         IF ((ztestmask(1)==1).AND.(icorns(ib,1)==0)) THEN
1696            IF(lwp) WRITE(numout,*)
1697            IF(lwp) WRITE(numout,*) ' E R R O R : Open boundary segment ', npckgs(ib)
1698            IF(lwp) WRITE(numout,*) ' ==========  does not start on land or on a corner'                                                 
1699            IF(lwp) WRITE(numout,*)
1700            nstop = nstop + 1
1701         ENDIF
1702         IF ((ztestmask(2)==1).AND.(icorns(ib,2)==0)) THEN
1703            IF(lwp) WRITE(numout,*)
1704            IF(lwp) WRITE(numout,*) ' E R R O R : Open boundary segment ', npckgs(ib)
1705            IF(lwp) WRITE(numout,*) ' ==========  does not end on land or on a corner'                                                 
1706            IF(lwp) WRITE(numout,*)
1707            nstop = nstop + 1
1708         ENDIF
1709      END DO
1710      !
1711      ! North segments
1712      DO ib = 1, nbdysegn
1713         ! get mask at boundary extremities:
1714         ztestmask(1:2)=0.
1715         DO ji = 1, jpi
1716            DO jj = 1, jpj             
1717              IF (((jj + njmpp - 1) == jpjnob(ib)+1).AND. & 
1718               &  ((ji + nimpp - 1) == jpindt(ib))) ztestmask(1)=tmask(ji,jj,1)
1719              IF (((jj + njmpp - 1) == jpjnob(ib)+1).AND. & 
1720               &  ((ji + nimpp - 1) == jpinft(ib))) ztestmask(2)=tmask(ji,jj,1) 
1721            END DO
1722         END DO
1723         IF( lk_mpp )   CALL mpp_sum( ztestmask, 2 )   ! sum over the global domain
1724
1725         IF ((ztestmask(1)==1).AND.(icornn(ib,1)==0)) THEN
1726            IF(lwp) WRITE(numout,*)
1727            IF(lwp) WRITE(numout,*) ' E R R O R : Open boundary segment ', npckgn(ib)
1728            IF(lwp) WRITE(numout,*) ' ==========  does not start on land'                                                 
1729            IF(lwp) WRITE(numout,*)
1730            nstop = nstop + 1
1731         ENDIF
1732         IF ((ztestmask(2)==1).AND.(icornn(ib,2)==0)) THEN
1733            IF(lwp) WRITE(numout,*)
1734            IF(lwp) WRITE(numout,*) ' E R R O R : Open boundary segment ', npckgn(ib)
1735            IF(lwp) WRITE(numout,*) ' ==========  does not end on land'                                                 
1736            IF(lwp) WRITE(numout,*)
1737            nstop = nstop + 1
1738         ENDIF
1739      END DO
1740      !
1741      IF ((itest==0).AND.(lwp)) WRITE(numout,*) 'NO open boundary corner found'
1742      !
1743      ! Other tests TBD:
1744      ! segments completly on land
1745      ! optimized open boundary array length according to landmask
1746      ! Nudging layers that overlap with interior domain
1747      !
1748   END SUBROUTINE bdy_ctl_seg
1749
1750   SUBROUTINE bdy_ctl_corn( ib1, ib2 )
1751      !!----------------------------------------------------------------------
1752      !!                 ***  ROUTINE bdy_ctl_corn  ***
1753      !!
1754      !! ** Purpose :   Check numerical schemes consistency between
1755      !!                segments having a common corner
1756      !!
1757      !! ** Method  :   
1758      !!----------------------------------------------------------------------
1759      INTEGER, INTENT(in)  ::   ib1, ib2
1760      INTEGER :: itest
1761      !!----------------------------------------------------------------------
1762      itest = 0
1763
1764      IF (cn_dyn2d(ib1)/=cn_dyn2d(ib2)) itest = itest + 1
1765      IF (cn_dyn3d(ib1)/=cn_dyn3d(ib2)) itest = itest + 1
1766      IF (cn_tra(ib1)/=cn_tra(ib2)) itest = itest + 1
1767      !
1768      IF (nn_dyn2d_dta(ib1)/=nn_dyn2d_dta(ib2)) itest = itest + 1
1769      IF (nn_dyn3d_dta(ib1)/=nn_dyn3d_dta(ib2)) itest = itest + 1
1770      IF (nn_tra_dta(ib1)/=nn_tra_dta(ib2)) itest = itest + 1
1771      !
1772      IF (nn_rimwidth(ib1)/=nn_rimwidth(ib2)) itest = itest + 1   
1773      !
1774      IF ( itest>0 ) THEN
1775         IF(lwp) WRITE(numout,*) ' E R R O R : Segments ', ib1, 'and ', ib2
1776         IF(lwp) WRITE(numout,*) ' ==========  have different open bdy schemes'                                                 
1777         IF(lwp) WRITE(numout,*)
1778         nstop = nstop + 1
1779      ENDIF
1780      !
1781   END SUBROUTINE bdy_ctl_corn
1782
1783#else
1784   !!---------------------------------------------------------------------------------
1785   !!   Dummy module                                   NO open boundaries
1786   !!---------------------------------------------------------------------------------
1787CONTAINS
1788   SUBROUTINE bdy_init      ! Dummy routine
1789   END SUBROUTINE bdy_init
1790#endif
1791
1792   !!=================================================================================
1793END MODULE bdyini
Note: See TracBrowser for help on using the repository browser.