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

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

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

dev_r10984_HPC-13 : Step 1, boundary is now detected all over the local domain, this does not change the result. Improve bdy treatment for bdy_rnf in bdytra.F90, this changes the result when keyword runoff is specified in namelist

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