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/NERC/dev_r5840_BDY_MSK/NEMOGCM/NEMO/OPA_SRC/BDY – NEMO

source: branches/NERC/dev_r5840_BDY_MSK/NEMOGCM/NEMO/OPA_SRC/BDY/bdyini.F90 @ 6948

Last change on this file since 6948 was 6948, checked in by jamesharle, 8 years ago

updated masking code in bdyini.F90 and dommsk.F90 to properly account for bdy_msk

  • Property svn:keywords set to Id
File size: 82.1 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  ::   iwe, ies, iso, ino, 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
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     ENDIF
385
386      ! -------------------------------------------------
387      ! Initialise indices arrays for open boundaries
388      ! -------------------------------------------------
389
390      ! Work out global dimensions of boundary data
391      ! ---------------------------------------------
392      REWIND( numnam_cfg )     
393
394      !!----------------------------------------------------------------------
395
396             
397               
398      nblendta(:,:) = 0
399      nbdysege = 0
400      nbdysegw = 0
401      nbdysegn = 0
402      nbdysegs = 0
403      icount   = 0 ! count user defined segments
404      ! Dimensions below are used to allocate arrays to read external data
405      jpbdtas = 1 ! Maximum size of boundary data (structured case)
406      jpbdtau = 1 ! Maximum size of boundary data (unstructured case)
407
408      DO ib_bdy = 1, nb_bdy
409
410         IF( .NOT. ln_coords_file(ib_bdy) ) THEN ! Work out size of global arrays from namelist parameters
411 
412            icount = icount + 1
413            ! No REWIND here because may need to read more than one nambdy_index namelist.
414            ! Read only namelist_cfg to avoid unseccessfull overwrite
415!!          REWIND( numnam_ref )              ! Namelist nambdy_index in reference namelist : Open boundaries indexes
416!!          READ  ( numnam_ref, namrun, IOSTAT = ios, ERR = 903)
417!!903       IF( ios /= 0 ) CALL ctl_nam ( ios , 'nambdy_index in reference namelist', lwp )
418
419!!          REWIND( numnam_cfg )              ! Namelist nambdy_index in configuration namelist : Open boundaries indexes
420            READ  ( numnam_cfg, nambdy_index, IOSTAT = ios, ERR = 904 )
421904         IF( ios /= 0 ) CALL ctl_nam ( ios , 'nambdy_index in configuration namelist', lwp )
422            IF(lwm) WRITE ( numond, nambdy_index )
423
424            SELECT CASE ( TRIM(ctypebdy) )
425              CASE( 'N' )
426                 IF( nbdyind == -1 ) THEN  ! Automatic boundary definition: if nbdysegX = -1
427                    nbdyind  = jpjglo - 2  ! set boundary to whole side of model domain.
428                    nbdybeg  = 2
429                    nbdyend  = jpiglo - 1
430                 ENDIF
431                 nbdysegn = nbdysegn + 1
432                 npckgn(nbdysegn) = ib_bdy ! Save bdy package number
433                 jpjnob(nbdysegn) = nbdyind
434                 jpindt(nbdysegn) = nbdybeg
435                 jpinft(nbdysegn) = nbdyend
436                 !
437              CASE( 'S' )
438                 IF( nbdyind == -1 ) THEN  ! Automatic boundary definition: if nbdysegX = -1
439                    nbdyind  = 2           ! set boundary to whole side of model domain.
440                    nbdybeg  = 2
441                    nbdyend  = jpiglo - 1
442                 ENDIF
443                 nbdysegs = nbdysegs + 1
444                 npckgs(nbdysegs) = ib_bdy ! Save bdy package number
445                 jpjsob(nbdysegs) = nbdyind
446                 jpisdt(nbdysegs) = nbdybeg
447                 jpisft(nbdysegs) = nbdyend
448                 !
449              CASE( 'E' )
450                 IF( nbdyind == -1 ) THEN  ! Automatic boundary definition: if nbdysegX = -1
451                    nbdyind  = jpiglo - 2  ! set boundary to whole side of model domain.
452                    nbdybeg  = 2
453                    nbdyend  = jpjglo - 1
454                 ENDIF
455                 nbdysege = nbdysege + 1 
456                 npckge(nbdysege) = ib_bdy ! Save bdy package number
457                 jpieob(nbdysege) = nbdyind
458                 jpjedt(nbdysege) = nbdybeg
459                 jpjeft(nbdysege) = nbdyend
460                 !
461              CASE( 'W' )
462                 IF( nbdyind == -1 ) THEN  ! Automatic boundary definition: if nbdysegX = -1
463                    nbdyind  = 2           ! set boundary to whole side of model domain.
464                    nbdybeg  = 2
465                    nbdyend  = jpjglo - 1
466                 ENDIF
467                 nbdysegw = nbdysegw + 1
468                 npckgw(nbdysegw) = ib_bdy ! Save bdy package number
469                 jpiwob(nbdysegw) = nbdyind
470                 jpjwdt(nbdysegw) = nbdybeg
471                 jpjwft(nbdysegw) = nbdyend
472                 !
473              CASE DEFAULT   ;   CALL ctl_stop( 'ctypebdy must be N, S, E or W' )
474            END SELECT
475
476            ! For simplicity we assume that in case of straight bdy, arrays have the same length
477            ! (even if it is true that last tangential velocity points
478            ! are useless). This simplifies a little bit boundary data format (and agrees with format
479            ! used so far in obc package)
480
481            nblendta(1:jpbgrd,ib_bdy) =  (nbdyend - nbdybeg + 1) * nn_rimwidth(ib_bdy)
482            jpbdtas = MAX(jpbdtas, (nbdyend - nbdybeg + 1))
483            IF (lwp.and.(nn_rimwidth(ib_bdy)>nrimmax)) &
484            & CALL ctl_stop( 'rimwidth must be lower than nrimmax' )
485
486         ELSE            ! Read size of arrays in boundary coordinates file.
487            CALL iom_open( cn_coords_file(ib_bdy), inum )
488            DO igrd = 1, jpbgrd
489               id_dummy = iom_varid( inum, 'nbi'//cgrid(igrd), kdimsz=kdimsz ) 
490               !clem nblendta(igrd,ib_bdy) = kdimsz(1)
491               !clem jpbdtau = MAX(jpbdtau, kdimsz(1))
492               nblendta(igrd,ib_bdy) = MAXVAL(kdimsz)
493               jpbdtau = MAX(jpbdtau, MAXVAL(kdimsz))
494            ENDDO
495            CALL iom_close( inum )
496
497         ENDIF
498
499      ENDDO ! ib_bdy
500
501      IF (nb_bdy>0) THEN
502         jpbdta = MAXVAL(nblendta(1:jpbgrd,1:nb_bdy))
503
504         ! Allocate arrays
505         !---------------
506         ALLOCATE( nbidta(jpbdta, jpbgrd, nb_bdy), nbjdta(jpbdta, jpbgrd, nb_bdy),    &
507            &      nbrdta(jpbdta, jpbgrd, nb_bdy) )
508
509         ALLOCATE( dta_global(jpbdtau, 1, jpk) )
510         IF ( icount>0 ) ALLOCATE( dta_global2(jpbdtas, nrimmax, jpk) )
511         !
512      ENDIF
513
514      ! Now look for crossings in user (namelist) defined open boundary segments:
515      !--------------------------------------------------------------------------
516      IF ( icount>0 ) CALL bdy_ctl_seg
517
518      ! Calculate global boundary index arrays or read in from file
519      !------------------------------------------------------------               
520      ! 1. Read global index arrays from boundary coordinates file.
521      DO ib_bdy = 1, nb_bdy
522
523         IF( ln_coords_file(ib_bdy) ) THEN
524
525            CALL iom_open( cn_coords_file(ib_bdy), inum )
526            DO igrd = 1, jpbgrd
527               CALL iom_get( inum, jpdom_unknown, 'nbi'//cgrid(igrd), dta_global(1:nblendta(igrd,ib_bdy),:,1) )
528               DO ii = 1,nblendta(igrd,ib_bdy)
529                  nbidta(ii,igrd,ib_bdy) = INT( dta_global(ii,1,1) )
530               END DO
531               CALL iom_get( inum, jpdom_unknown, 'nbj'//cgrid(igrd), dta_global(1:nblendta(igrd,ib_bdy),:,1) )
532               DO ii = 1,nblendta(igrd,ib_bdy)
533                  nbjdta(ii,igrd,ib_bdy) = INT( dta_global(ii,1,1) )
534               END DO
535               CALL iom_get( inum, jpdom_unknown, 'nbr'//cgrid(igrd), dta_global(1:nblendta(igrd,ib_bdy),:,1) )
536               DO ii = 1,nblendta(igrd,ib_bdy)
537                  nbrdta(ii,igrd,ib_bdy) = INT( dta_global(ii,1,1) )
538               END DO
539
540               ibr_max = MAXVAL( nbrdta(:,igrd,ib_bdy) )
541               IF(lwp) WRITE(numout,*)
542               IF(lwp) WRITE(numout,*) ' Maximum rimwidth in file is ', ibr_max
543               IF(lwp) WRITE(numout,*) ' nn_rimwidth from namelist is ', nn_rimwidth(ib_bdy)
544               IF (ibr_max < nn_rimwidth(ib_bdy))   &
545                     CALL ctl_stop( 'nn_rimwidth is larger than maximum rimwidth in file',cn_coords_file(ib_bdy) )
546            END DO
547            CALL iom_close( inum )
548
549         ENDIF
550
551      ENDDO     
552   
553      ! 2. Now fill indices corresponding to straight open boundary arrays:
554      ! East
555      !-----
556      DO iseg = 1, nbdysege
557         ib_bdy = npckge(iseg)
558         !
559         ! ------------ T points -------------
560         igrd=1
561         icount=0
562         DO ir = 1, nn_rimwidth(ib_bdy)
563            DO ij = jpjedt(iseg), jpjeft(iseg)
564               icount = icount + 1
565               nbidta(icount, igrd, ib_bdy) = jpieob(iseg) + 2 - ir
566               nbjdta(icount, igrd, ib_bdy) = ij
567               nbrdta(icount, igrd, ib_bdy) = ir
568            ENDDO
569         ENDDO
570         !
571         ! ------------ U points -------------
572         igrd=2
573         icount=0
574         DO ir = 1, nn_rimwidth(ib_bdy)
575            DO ij = jpjedt(iseg), jpjeft(iseg)
576               icount = icount + 1
577               nbidta(icount, igrd, ib_bdy) = jpieob(iseg) + 1 - ir
578               nbjdta(icount, igrd, ib_bdy) = ij
579               nbrdta(icount, igrd, ib_bdy) = ir
580            ENDDO
581         ENDDO
582         !
583         ! ------------ V points -------------
584         igrd=3
585         icount=0
586         DO ir = 1, nn_rimwidth(ib_bdy)
587!            DO ij = jpjedt(iseg), jpjeft(iseg) - 1
588            DO ij = jpjedt(iseg), jpjeft(iseg)
589               icount = icount + 1
590               nbidta(icount, igrd, ib_bdy) = jpieob(iseg) + 2 - ir
591               nbjdta(icount, igrd, ib_bdy) = ij
592               nbrdta(icount, igrd, ib_bdy) = ir
593            ENDDO
594            nbidta(icount, igrd, ib_bdy) = -ib_bdy ! Discount this point
595            nbjdta(icount, igrd, ib_bdy) = -ib_bdy ! Discount this point
596         ENDDO
597      ENDDO
598      !
599      ! West
600      !-----
601      DO iseg = 1, nbdysegw
602         ib_bdy = npckgw(iseg)
603         !
604         ! ------------ T points -------------
605         igrd=1
606         icount=0
607         DO ir = 1, nn_rimwidth(ib_bdy)
608            DO ij = jpjwdt(iseg), jpjwft(iseg)
609               icount = icount + 1
610               nbidta(icount, igrd, ib_bdy) = jpiwob(iseg) + ir - 1
611               nbjdta(icount, igrd, ib_bdy) = ij
612               nbrdta(icount, igrd, ib_bdy) = ir
613            ENDDO
614         ENDDO
615         !
616         ! ------------ U points -------------
617         igrd=2
618         icount=0
619         DO ir = 1, nn_rimwidth(ib_bdy)
620            DO ij = jpjwdt(iseg), jpjwft(iseg)
621               icount = icount + 1
622               nbidta(icount, igrd, ib_bdy) = jpiwob(iseg) + ir - 1
623               nbjdta(icount, igrd, ib_bdy) = ij
624               nbrdta(icount, igrd, ib_bdy) = ir
625            ENDDO
626         ENDDO
627         !
628         ! ------------ V points -------------
629         igrd=3
630         icount=0
631         DO ir = 1, nn_rimwidth(ib_bdy)
632!            DO ij = jpjwdt(iseg), jpjwft(iseg) - 1
633            DO ij = jpjwdt(iseg), jpjwft(iseg)
634               icount = icount + 1
635               nbidta(icount, igrd, ib_bdy) = jpiwob(iseg) + ir - 1
636               nbjdta(icount, igrd, ib_bdy) = ij
637               nbrdta(icount, igrd, ib_bdy) = ir
638            ENDDO
639            nbidta(icount, igrd, ib_bdy) = -ib_bdy ! Discount this point
640            nbjdta(icount, igrd, ib_bdy) = -ib_bdy ! Discount this point
641         ENDDO
642      ENDDO
643      !
644      ! North
645      !-----
646      DO iseg = 1, nbdysegn
647         ib_bdy = npckgn(iseg)
648         !
649         ! ------------ T points -------------
650         igrd=1
651         icount=0
652         DO ir = 1, nn_rimwidth(ib_bdy)
653            DO ii = jpindt(iseg), jpinft(iseg)
654               icount = icount + 1
655               nbidta(icount, igrd, ib_bdy) = ii
656               nbjdta(icount, igrd, ib_bdy) = jpjnob(iseg) + 2 - ir 
657               nbrdta(icount, igrd, ib_bdy) = ir
658            ENDDO
659         ENDDO
660         !
661         ! ------------ U points -------------
662         igrd=2
663         icount=0
664         DO ir = 1, nn_rimwidth(ib_bdy)
665!            DO ii = jpindt(iseg), jpinft(iseg) - 1
666            DO ii = jpindt(iseg), jpinft(iseg)
667               icount = icount + 1
668               nbidta(icount, igrd, ib_bdy) = ii
669               nbjdta(icount, igrd, ib_bdy) = jpjnob(iseg) + 2 - ir
670               nbrdta(icount, igrd, ib_bdy) = ir
671            ENDDO
672            nbidta(icount, igrd, ib_bdy) = -ib_bdy ! Discount this point
673            nbjdta(icount, igrd, ib_bdy) = -ib_bdy ! Discount this point
674         ENDDO
675         !
676         ! ------------ V points -------------
677         igrd=3
678         icount=0
679         DO ir = 1, nn_rimwidth(ib_bdy)
680            DO ii = jpindt(iseg), jpinft(iseg)
681               icount = icount + 1
682               nbidta(icount, igrd, ib_bdy) = ii
683               nbjdta(icount, igrd, ib_bdy) = jpjnob(iseg) + 1 - ir
684               nbrdta(icount, igrd, ib_bdy) = ir
685            ENDDO
686         ENDDO
687      ENDDO
688      !
689      ! South
690      !-----
691      DO iseg = 1, nbdysegs
692         ib_bdy = npckgs(iseg)
693         !
694         ! ------------ T points -------------
695         igrd=1
696         icount=0
697         DO ir = 1, nn_rimwidth(ib_bdy)
698            DO ii = jpisdt(iseg), jpisft(iseg)
699               icount = icount + 1
700               nbidta(icount, igrd, ib_bdy) = ii
701               nbjdta(icount, igrd, ib_bdy) = jpjsob(iseg) + ir - 1
702               nbrdta(icount, igrd, ib_bdy) = ir
703            ENDDO
704         ENDDO
705         !
706         ! ------------ U points -------------
707         igrd=2
708         icount=0
709         DO ir = 1, nn_rimwidth(ib_bdy)
710!            DO ii = jpisdt(iseg), jpisft(iseg) - 1
711            DO ii = jpisdt(iseg), jpisft(iseg)
712               icount = icount + 1
713               nbidta(icount, igrd, ib_bdy) = ii
714               nbjdta(icount, igrd, ib_bdy) = jpjsob(iseg) + ir - 1
715               nbrdta(icount, igrd, ib_bdy) = ir
716            ENDDO
717            nbidta(icount, igrd, ib_bdy) = -ib_bdy ! Discount this point
718            nbjdta(icount, igrd, ib_bdy) = -ib_bdy ! Discount this point
719         ENDDO
720         !
721         ! ------------ V points -------------
722         igrd=3
723         icount=0
724         DO ir = 1, nn_rimwidth(ib_bdy)
725            DO ii = jpisdt(iseg), jpisft(iseg)
726               icount = icount + 1
727               nbidta(icount, igrd, ib_bdy) = ii
728               nbjdta(icount, igrd, ib_bdy) = jpjsob(iseg) + ir - 1
729               nbrdta(icount, igrd, ib_bdy) = ir
730            ENDDO
731         ENDDO
732      ENDDO
733
734      !  Deal with duplicated points
735      !-----------------------------
736      ! We assign negative indices to duplicated points (to remove them from bdy points to be updated)
737      ! if their distance to the bdy is greater than the other
738      ! If their distance are the same, just keep only one to avoid updating a point twice
739      DO igrd = 1, jpbgrd
740         DO ib_bdy1 = 1, nb_bdy
741            DO ib_bdy2 = 1, nb_bdy
742               IF (ib_bdy1/=ib_bdy2) THEN
743                  DO ib1 = 1, nblendta(igrd,ib_bdy1)
744                     DO ib2 = 1, nblendta(igrd,ib_bdy2)
745                        IF ((nbidta(ib1, igrd, ib_bdy1)==nbidta(ib2, igrd, ib_bdy2)).AND. &
746                        &   (nbjdta(ib1, igrd, ib_bdy1)==nbjdta(ib2, igrd, ib_bdy2))) THEN
747!                           IF ((lwp).AND.(igrd==1)) WRITE(numout,*) ' found coincident point ji, jj:', &
748!                                                       &              nbidta(ib1, igrd, ib_bdy1),      &
749!                                                       &              nbjdta(ib2, igrd, ib_bdy2)
750                           ! keep only points with the lowest distance to boundary:
751                           IF (nbrdta(ib1, igrd, ib_bdy1)<nbrdta(ib2, igrd, ib_bdy2)) THEN
752                             nbidta(ib2, igrd, ib_bdy2) =-ib_bdy2
753                             nbjdta(ib2, igrd, ib_bdy2) =-ib_bdy2
754                           ELSEIF (nbrdta(ib1, igrd, ib_bdy1)>nbrdta(ib2, igrd, ib_bdy2)) THEN
755                             nbidta(ib1, igrd, ib_bdy1) =-ib_bdy1
756                             nbjdta(ib1, igrd, ib_bdy1) =-ib_bdy1
757                           ! Arbitrary choice if distances are the same:
758                           ELSE
759                             nbidta(ib1, igrd, ib_bdy1) =-ib_bdy1
760                             nbjdta(ib1, igrd, ib_bdy1) =-ib_bdy1
761                           ENDIF
762                        END IF
763                     END DO
764                  END DO
765               ENDIF
766            END DO
767         END DO
768      END DO
769
770      ! Work out dimensions of boundary data on each processor
771      ! ------------------------------------------------------
772
773      ! Rather assume that boundary data indices are given on global domain
774      ! TO BE DISCUSSED ?
775!      iw = mig(1) + 1            ! if monotasking and no zoom, iw=2
776!      ie = mig(1) + nlci-1 - 1   ! if monotasking and no zoom, ie=jpim1
777!      is = mjg(1) + 1            ! if monotasking and no zoom, is=2
778!      in = mjg(1) + nlcj-1 - 1   ! if monotasking and no zoom, in=jpjm1     
779      iwe = mig(1) - jpizoom + 2         ! if monotasking and no zoom, iw=2
780      ies = mig(1) + nlci - jpizoom - 1  ! if monotasking and no zoom, ie=jpim1
781      iso = mjg(1) - jpjzoom + 2         ! if monotasking and no zoom, is=2
782      ino = mjg(1) + nlcj - jpjzoom - 1  ! if monotasking and no zoom, in=jpjm1
783
784      ALLOCATE( nbondi_bdy(nb_bdy))
785      ALLOCATE( nbondj_bdy(nb_bdy))
786      nbondi_bdy(:)=2
787      nbondj_bdy(:)=2
788      ALLOCATE( nbondi_bdy_b(nb_bdy))
789      ALLOCATE( nbondj_bdy_b(nb_bdy))
790      nbondi_bdy_b(:)=2
791      nbondj_bdy_b(:)=2
792
793      ! Work out dimensions of boundary data on each neighbour process
794      IF(nbondi .eq. 0) THEN
795         iw_b(1) = jpizoom + nimppt(nowe+1)
796         ie_b(1) = jpizoom + nimppt(nowe+1)+nlcit(nowe+1)-3
797         is_b(1) = jpjzoom + njmppt(nowe+1)
798         in_b(1) = jpjzoom + njmppt(nowe+1)+nlcjt(nowe+1)-3
799
800         iw_b(2) = jpizoom + nimppt(noea+1)
801         ie_b(2) = jpizoom + nimppt(noea+1)+nlcit(noea+1)-3
802         is_b(2) = jpjzoom + njmppt(noea+1)
803         in_b(2) = jpjzoom + njmppt(noea+1)+nlcjt(noea+1)-3
804      ELSEIF(nbondi .eq. 1) THEN
805         iw_b(1) = jpizoom + nimppt(nowe+1)
806         ie_b(1) = jpizoom + nimppt(nowe+1)+nlcit(nowe+1)-3
807         is_b(1) = jpjzoom + njmppt(nowe+1)
808         in_b(1) = jpjzoom + njmppt(nowe+1)+nlcjt(nowe+1)-3
809      ELSEIF(nbondi .eq. -1) THEN
810         iw_b(2) = jpizoom + nimppt(noea+1)
811         ie_b(2) = jpizoom + nimppt(noea+1)+nlcit(noea+1)-3
812         is_b(2) = jpjzoom + njmppt(noea+1)
813         in_b(2) = jpjzoom + njmppt(noea+1)+nlcjt(noea+1)-3
814      ENDIF
815
816      IF(nbondj .eq. 0) THEN
817         iw_b(3) = jpizoom + nimppt(noso+1)
818         ie_b(3) = jpizoom + nimppt(noso+1)+nlcit(noso+1)-3
819         is_b(3) = jpjzoom + njmppt(noso+1)
820         in_b(3) = jpjzoom + njmppt(noso+1)+nlcjt(noso+1)-3
821
822         iw_b(4) = jpizoom + nimppt(nono+1)
823         ie_b(4) = jpizoom + nimppt(nono+1)+nlcit(nono+1)-3
824         is_b(4) = jpjzoom + njmppt(nono+1)
825         in_b(4) = jpjzoom + njmppt(nono+1)+nlcjt(nono+1)-3
826      ELSEIF(nbondj .eq. 1) THEN
827         iw_b(3) = jpizoom + nimppt(noso+1)
828         ie_b(3) = jpizoom + nimppt(noso+1)+nlcit(noso+1)-3
829         is_b(3) = jpjzoom + njmppt(noso+1)
830         in_b(3) = jpjzoom + njmppt(noso+1)+nlcjt(noso+1)-3
831      ELSEIF(nbondj .eq. -1) THEN
832         iw_b(4) = jpizoom + nimppt(nono+1)
833         ie_b(4) = jpizoom + nimppt(nono+1)+nlcit(nono+1)-3
834         is_b(4) = jpjzoom + njmppt(nono+1)
835         in_b(4) = jpjzoom + njmppt(nono+1)+nlcjt(nono+1)-3
836      ENDIF
837
838      DO ib_bdy = 1, nb_bdy
839         DO igrd = 1, jpbgrd
840            icount  = 0
841            icountr = 0
842            idx_bdy(ib_bdy)%nblen(igrd)    = 0
843            idx_bdy(ib_bdy)%nblenrim(igrd) = 0
844            DO ib = 1, nblendta(igrd,ib_bdy)
845               ! check that data is in correct order in file
846               ibm1 = MAX(1,ib-1)
847               IF(lwp) THEN         ! Since all procs read global data only need to do this check on one proc...
848                  IF( nbrdta(ib,igrd,ib_bdy) < nbrdta(ibm1,igrd,ib_bdy) ) THEN
849                     CALL ctl_stop('bdy_init : ERROR : boundary data in file must be defined ', &
850                          &        ' in order of distance from edge nbr A utility for re-ordering ', &
851                          &        ' boundary coordinates and data files exists in the TOOLS/OBC directory')
852                  ENDIF   
853               ENDIF
854               ! check if point is in local domain
855               IF(  nbidta(ib,igrd,ib_bdy) >= iwe .AND. nbidta(ib,igrd,ib_bdy) <= ies .AND.   &
856                  & nbjdta(ib,igrd,ib_bdy) >= iso .AND. nbjdta(ib,igrd,ib_bdy) <= ino      ) THEN
857                  !
858                  icount = icount  + 1
859                  !
860                  IF( nbrdta(ib,igrd,ib_bdy) == 1 )   icountr = icountr+1
861               ENDIF
862            ENDDO
863            idx_bdy(ib_bdy)%nblenrim(igrd) = icountr !: length of rim boundary data on each proc
864            idx_bdy(ib_bdy)%nblen   (igrd) = icount  !: length of boundary data on each proc       
865         ENDDO  ! igrd
866
867         ! Allocate index arrays for this boundary set
868         !--------------------------------------------
869         ilen1 = MAXVAL(idx_bdy(ib_bdy)%nblen(:))
870         ALLOCATE( idx_bdy(ib_bdy)%nbi(ilen1,jpbgrd) )
871         ALLOCATE( idx_bdy(ib_bdy)%nbj(ilen1,jpbgrd) )
872         ALLOCATE( idx_bdy(ib_bdy)%nbr(ilen1,jpbgrd) )
873         ALLOCATE( idx_bdy(ib_bdy)%nbd(ilen1,jpbgrd) )
874         ALLOCATE( idx_bdy(ib_bdy)%nbdout(ilen1,jpbgrd) )
875         ALLOCATE( idx_bdy(ib_bdy)%nbmap(ilen1,jpbgrd) )
876         ALLOCATE( idx_bdy(ib_bdy)%nbw(ilen1,jpbgrd) )
877         ALLOCATE( idx_bdy(ib_bdy)%flagu(ilen1,jpbgrd) )
878         ALLOCATE( idx_bdy(ib_bdy)%flagv(ilen1,jpbgrd) )
879
880         ! Dispatch mapping indices and discrete distances on each processor
881         ! -----------------------------------------------------------------
882
883         com_east = 0
884         com_west = 0
885         com_south = 0
886         com_north = 0
887
888         com_east_b = 0
889         com_west_b = 0
890         com_south_b = 0
891         com_north_b = 0
892
893         DO igrd = 1, jpbgrd
894            icount  = 0
895            ! Loop on rimwidth to ensure outermost points come first in the local arrays.
896            DO ir=1, nn_rimwidth(ib_bdy)
897               DO ib = 1, nblendta(igrd,ib_bdy)
898                  ! check if point is in local domain and equals ir
899                  IF(  nbidta(ib,igrd,ib_bdy) >= iwe .AND. nbidta(ib,igrd,ib_bdy) <= ies .AND.   &
900                     & nbjdta(ib,igrd,ib_bdy) >= iso .AND. nbjdta(ib,igrd,ib_bdy) <= ino .AND.   &
901                     & nbrdta(ib,igrd,ib_bdy) == ir  ) THEN
902                     !
903                     icount = icount  + 1
904
905                     ! Rather assume that boundary data indices are given on global domain
906                     ! TO BE DISCUSSED ?
907!                     idx_bdy(ib_bdy)%nbi(icount,igrd)   = nbidta(ib,igrd,ib_bdy)- mig(1)+1
908!                     idx_bdy(ib_bdy)%nbj(icount,igrd)   = nbjdta(ib,igrd,ib_bdy)- mjg(1)+1
909                     idx_bdy(ib_bdy)%nbi(icount,igrd)   = nbidta(ib,igrd,ib_bdy)- mig(1)+jpizoom
910                     idx_bdy(ib_bdy)%nbj(icount,igrd)   = nbjdta(ib,igrd,ib_bdy)- mjg(1)+jpjzoom
911                     ! check if point has to be sent
912                     ii = idx_bdy(ib_bdy)%nbi(icount,igrd)
913                     ij = idx_bdy(ib_bdy)%nbj(icount,igrd)
914                     if((com_east .ne. 1) .and. (ii .eq. (nlci-1)) .and. (nbondi .le. 0)) then
915                        com_east = 1
916                     elseif((com_west .ne. 1) .and. (ii .eq. 2) .and. (nbondi .ge. 0) .and. (nbondi .ne. 2)) then
917                        com_west = 1
918                     endif
919                     if((com_south .ne. 1) .and. (ij .eq. 2) .and. (nbondj .ge. 0) .and. (nbondj .ne. 2)) then
920                        com_south = 1
921                     elseif((com_north .ne. 1) .and. (ij .eq. (nlcj-1)) .and. (nbondj .le. 0)) then
922                        com_north = 1
923                     endif
924                     idx_bdy(ib_bdy)%nbr(icount,igrd)   = nbrdta(ib,igrd,ib_bdy)
925                     idx_bdy(ib_bdy)%nbmap(icount,igrd) = ib
926                  ENDIF
927                  ! check if point has to be received from a neighbour
928                  IF(nbondi .eq. 0) THEN
929                     IF( nbidta(ib,igrd,ib_bdy) >= iw_b(1) .AND. nbidta(ib,igrd,ib_bdy) <= ie_b(1) .AND.   &
930                       & nbjdta(ib,igrd,ib_bdy) >= is_b(1) .AND. nbjdta(ib,igrd,ib_bdy) <= in_b(1) .AND.   &
931                       & nbrdta(ib,igrd,ib_bdy) == ir  ) THEN
932                       ii = nbidta(ib,igrd,ib_bdy)- iw_b(1)+2
933                       if((com_west_b .ne. 1) .and. (ii .eq. (nlcit(nowe+1)-1))) then
934                          ij = nbjdta(ib,igrd,ib_bdy) - is_b(1)+2
935                          if((ij .eq. 2) .and. (nbondj .eq. 0 .or. nbondj .eq. 1)) then
936                            com_south = 1
937                          elseif((ij .eq. nlcjt(nowe+1)-1) .and. (nbondj .eq. 0 .or. nbondj .eq. -1)) then
938                            com_north = 1
939                          endif
940                          com_west_b = 1
941                       endif
942                     ENDIF
943                     IF( nbidta(ib,igrd,ib_bdy) >= iw_b(2) .AND. nbidta(ib,igrd,ib_bdy) <= ie_b(2) .AND.   &
944                       & nbjdta(ib,igrd,ib_bdy) >= is_b(2) .AND. nbjdta(ib,igrd,ib_bdy) <= in_b(2) .AND.   &
945                       & nbrdta(ib,igrd,ib_bdy) == ir  ) THEN
946                       ii = nbidta(ib,igrd,ib_bdy)- iw_b(2)+2
947                       if((com_east_b .ne. 1) .and. (ii .eq. 2)) then
948                          ij = nbjdta(ib,igrd,ib_bdy) - is_b(2)+2
949                          if((ij .eq. 2) .and. (nbondj .eq. 0 .or. nbondj .eq. 1)) then
950                            com_south = 1
951                          elseif((ij .eq. nlcjt(noea+1)-1) .and. (nbondj .eq. 0 .or. nbondj .eq. -1)) then
952                            com_north = 1
953                          endif
954                          com_east_b = 1
955                       endif
956                     ENDIF
957                  ELSEIF(nbondi .eq. 1) THEN
958                     IF( nbidta(ib,igrd,ib_bdy) >= iw_b(1) .AND. nbidta(ib,igrd,ib_bdy) <= ie_b(1) .AND.   &
959                       & nbjdta(ib,igrd,ib_bdy) >= is_b(1) .AND. nbjdta(ib,igrd,ib_bdy) <= in_b(1) .AND.   &
960                       & nbrdta(ib,igrd,ib_bdy) == ir  ) THEN
961                       ii = nbidta(ib,igrd,ib_bdy)- iw_b(1)+2
962                       if((com_west_b .ne. 1) .and. (ii .eq. (nlcit(nowe+1)-1))) then
963                          ij = nbjdta(ib,igrd,ib_bdy) - is_b(1)+2
964                          if((ij .eq. 2) .and. (nbondj .eq. 0 .or. nbondj .eq. 1)) then
965                            com_south = 1
966                          elseif((ij .eq. nlcjt(nowe+1)-1) .and. (nbondj .eq. 0 .or. nbondj .eq. -1)) then
967                            com_north = 1
968                          endif
969                          com_west_b = 1
970                       endif
971                     ENDIF
972                  ELSEIF(nbondi .eq. -1) THEN
973                     IF( nbidta(ib,igrd,ib_bdy) >= iw_b(2) .AND. nbidta(ib,igrd,ib_bdy) <= ie_b(2) .AND.   &
974                       & nbjdta(ib,igrd,ib_bdy) >= is_b(2) .AND. nbjdta(ib,igrd,ib_bdy) <= in_b(2) .AND.   &
975                       & nbrdta(ib,igrd,ib_bdy) == ir  ) THEN
976                       ii = nbidta(ib,igrd,ib_bdy)- iw_b(2)+2
977                       if((com_east_b .ne. 1) .and. (ii .eq. 2)) then
978                          ij = nbjdta(ib,igrd,ib_bdy) - is_b(2)+2
979                          if((ij .eq. 2) .and. (nbondj .eq. 0 .or. nbondj .eq. 1)) then
980                            com_south = 1
981                          elseif((ij .eq. nlcjt(noea+1)-1) .and. (nbondj .eq. 0 .or. nbondj .eq. -1)) then
982                            com_north = 1
983                          endif
984                          com_east_b = 1
985                       endif
986                     ENDIF
987                  ENDIF
988                  IF(nbondj .eq. 0) THEN
989                     IF(com_north_b .ne. 1 .AND. (nbidta(ib,igrd,ib_bdy) == iw_b(4)-1  &
990                       & .OR. nbidta(ib,igrd,ib_bdy) == ie_b(4)+1) .AND. &
991                       & nbjdta(ib,igrd,ib_bdy) == is_b(4) .AND. nbrdta(ib,igrd,ib_bdy) == ir) THEN
992                       com_north_b = 1 
993                     ENDIF
994                     IF(com_south_b .ne. 1 .AND. (nbidta(ib,igrd,ib_bdy) == iw_b(3)-1  &
995                       &.OR. nbidta(ib,igrd,ib_bdy) == ie_b(3)+1) .AND. &
996                       & nbjdta(ib,igrd,ib_bdy) == in_b(3) .AND. nbrdta(ib,igrd,ib_bdy) == ir) THEN
997                       com_south_b = 1 
998                     ENDIF
999                     IF( nbidta(ib,igrd,ib_bdy) >= iw_b(3) .AND. nbidta(ib,igrd,ib_bdy) <= ie_b(3) .AND.   &
1000                       & nbjdta(ib,igrd,ib_bdy) >= is_b(3) .AND. nbjdta(ib,igrd,ib_bdy) <= in_b(3) .AND.   &
1001                       & nbrdta(ib,igrd,ib_bdy) == ir  ) THEN
1002                       ij = nbjdta(ib,igrd,ib_bdy)- is_b(3)+2
1003                       if((com_south_b .ne. 1) .and. (ij .eq. (nlcjt(noso+1)-1))) then
1004                          com_south_b = 1
1005                       endif
1006                     ENDIF
1007                     IF( nbidta(ib,igrd,ib_bdy) >= iw_b(4) .AND. nbidta(ib,igrd,ib_bdy) <= ie_b(4) .AND.   &
1008                       & nbjdta(ib,igrd,ib_bdy) >= is_b(4) .AND. nbjdta(ib,igrd,ib_bdy) <= in_b(4) .AND.   &
1009                       & nbrdta(ib,igrd,ib_bdy) == ir  ) THEN
1010                       ij = nbjdta(ib,igrd,ib_bdy)- is_b(4)+2
1011                       if((com_north_b .ne. 1) .and. (ij .eq. 2)) then
1012                          com_north_b = 1
1013                       endif
1014                     ENDIF
1015                  ELSEIF(nbondj .eq. 1) THEN
1016                     IF( com_south_b .ne. 1 .AND. (nbidta(ib,igrd,ib_bdy) == iw_b(3)-1 .OR. &
1017                       & nbidta(ib,igrd,ib_bdy) == ie_b(3)+1) .AND. &
1018                       & nbjdta(ib,igrd,ib_bdy) == in_b(3) .AND. nbrdta(ib,igrd,ib_bdy) == ir) THEN
1019                       com_south_b = 1 
1020                     ENDIF
1021                     IF( nbidta(ib,igrd,ib_bdy) >= iw_b(3) .AND. nbidta(ib,igrd,ib_bdy) <= ie_b(3) .AND.   &
1022                       & nbjdta(ib,igrd,ib_bdy) >= is_b(3) .AND. nbjdta(ib,igrd,ib_bdy) <= in_b(3) .AND.   &
1023                       & nbrdta(ib,igrd,ib_bdy) == ir  ) THEN
1024                       ij = nbjdta(ib,igrd,ib_bdy)- is_b(3)+2
1025                       if((com_south_b .ne. 1) .and. (ij .eq. (nlcjt(noso+1)-1))) then
1026                          com_south_b = 1
1027                       endif
1028                     ENDIF
1029                  ELSEIF(nbondj .eq. -1) THEN
1030                     IF(com_north_b .ne. 1 .AND. (nbidta(ib,igrd,ib_bdy) == iw_b(4)-1  &
1031                       & .OR. nbidta(ib,igrd,ib_bdy) == ie_b(4)+1) .AND. &
1032                       & nbjdta(ib,igrd,ib_bdy) == is_b(4) .AND. nbrdta(ib,igrd,ib_bdy) == ir) THEN
1033                       com_north_b = 1 
1034                     ENDIF
1035                     IF( nbidta(ib,igrd,ib_bdy) >= iw_b(4) .AND. nbidta(ib,igrd,ib_bdy) <= ie_b(4) .AND.   &
1036                       & nbjdta(ib,igrd,ib_bdy) >= is_b(4) .AND. nbjdta(ib,igrd,ib_bdy) <= in_b(4) .AND.   &
1037                       & nbrdta(ib,igrd,ib_bdy) == ir  ) THEN
1038                       ij = nbjdta(ib,igrd,ib_bdy)- is_b(4)+2
1039                       if((com_north_b .ne. 1) .and. (ij .eq. 2)) then
1040                          com_north_b = 1
1041                       endif
1042                     ENDIF
1043                  ENDIF
1044               ENDDO
1045            ENDDO
1046         ENDDO 
1047
1048         ! definition of the i- and j- direction local boundaries arrays
1049         ! used for sending the boudaries
1050         IF((com_east .eq. 1) .and. (com_west .eq. 1)) THEN
1051            nbondi_bdy(ib_bdy) = 0
1052         ELSEIF ((com_east .eq. 1) .and. (com_west .eq. 0)) THEN
1053            nbondi_bdy(ib_bdy) = -1
1054         ELSEIF ((com_east .eq. 0) .and. (com_west .eq. 1)) THEN
1055            nbondi_bdy(ib_bdy) = 1
1056         ENDIF
1057
1058         IF((com_north .eq. 1) .and. (com_south .eq. 1)) THEN
1059            nbondj_bdy(ib_bdy) = 0
1060         ELSEIF ((com_north .eq. 1) .and. (com_south .eq. 0)) THEN
1061            nbondj_bdy(ib_bdy) = -1
1062         ELSEIF ((com_north .eq. 0) .and. (com_south .eq. 1)) THEN
1063            nbondj_bdy(ib_bdy) = 1
1064         ENDIF
1065
1066         ! definition of the i- and j- direction local boundaries arrays
1067         ! used for receiving the boudaries
1068         IF((com_east_b .eq. 1) .and. (com_west_b .eq. 1)) THEN
1069            nbondi_bdy_b(ib_bdy) = 0
1070         ELSEIF ((com_east_b .eq. 1) .and. (com_west_b .eq. 0)) THEN
1071            nbondi_bdy_b(ib_bdy) = -1
1072         ELSEIF ((com_east_b .eq. 0) .and. (com_west_b .eq. 1)) THEN
1073            nbondi_bdy_b(ib_bdy) = 1
1074         ENDIF
1075
1076         IF((com_north_b .eq. 1) .and. (com_south_b .eq. 1)) THEN
1077            nbondj_bdy_b(ib_bdy) = 0
1078         ELSEIF ((com_north_b .eq. 1) .and. (com_south_b .eq. 0)) THEN
1079            nbondj_bdy_b(ib_bdy) = -1
1080         ELSEIF ((com_north_b .eq. 0) .and. (com_south_b .eq. 1)) THEN
1081            nbondj_bdy_b(ib_bdy) = 1
1082         ENDIF
1083
1084         ! Compute rim weights for FRS scheme
1085         ! ----------------------------------
1086         DO igrd = 1, jpbgrd
1087            DO ib = 1, idx_bdy(ib_bdy)%nblen(igrd)
1088               nbr => idx_bdy(ib_bdy)%nbr(ib,igrd)
1089               idx_bdy(ib_bdy)%nbw(ib,igrd) = 1.- TANH( FLOAT( nbr - 1 ) *0.5 )      ! tanh formulation
1090!               idx_bdy(ib_bdy)%nbw(ib,igrd) = (FLOAT(nn_rimwidth(ib_bdy)+1-nbr)/FLOAT(nn_rimwidth(ib_bdy)))**2.  ! quadratic
1091!               idx_bdy(ib_bdy)%nbw(ib,igrd) =  FLOAT(nn_rimwidth(ib_bdy)+1-nbr)/FLOAT(nn_rimwidth(ib_bdy))       ! linear
1092            END DO
1093         END DO 
1094
1095         ! Compute damping coefficients
1096         ! ----------------------------
1097         DO igrd = 1, jpbgrd
1098            DO ib = 1, idx_bdy(ib_bdy)%nblen(igrd)
1099               nbr => idx_bdy(ib_bdy)%nbr(ib,igrd)
1100               idx_bdy(ib_bdy)%nbd(ib,igrd) = 1. / ( rn_time_dmp(ib_bdy) * rday ) & 
1101               & *(FLOAT(nn_rimwidth(ib_bdy)+1-nbr)/FLOAT(nn_rimwidth(ib_bdy)))**2.   ! quadratic
1102               idx_bdy(ib_bdy)%nbdout(ib,igrd) = 1. / ( rn_time_dmp_out(ib_bdy) * rday ) & 
1103               & *(FLOAT(nn_rimwidth(ib_bdy)+1-nbr)/FLOAT(nn_rimwidth(ib_bdy)))**2.   ! quadratic
1104            END DO
1105         END DO
1106
1107      ENDDO
1108
1109      ! ------------------------------------------------------
1110      ! Initialise masks and find normal/tangential directions
1111      ! ------------------------------------------------------
1112
1113      ! Read global 2D mask at T-points: bdytmask
1114      ! -----------------------------------------
1115      ! bdytmask = 1  on the computational domain AND on open boundaries
1116      !          = 0  elsewhere   
1117 
1118      bdytmask(:,:) = ssmask(:,:)
1119
1120      IF( ln_mask_file ) THEN
1121         CALL iom_open( cn_mask_file, inum )
1122         CALL iom_get ( inum, jpdom_data, 'bdy_msk', bdytmask(:,:) )
1123         CALL iom_close( inum )
1124
1125         ! Derive mask on U and V grid from mask on T grid
1126         bdyumask(:,:) = 0.e0
1127         bdyvmask(:,:) = 0.e0
1128         DO ij=1, jpjm1
1129            DO ii=1, jpim1
1130               bdyumask(ii,ij)=bdytmask(ii,ij)*bdytmask(ii+1, ij )
1131               bdyvmask(ii,ij)=bdytmask(ii,ij)*bdytmask(ii  ,ij+1) 
1132            END DO
1133         END DO
1134         CALL lbc_lnk( bdyumask(:,:), 'U', 1. )   ;   CALL lbc_lnk( bdyvmask(:,:), 'V', 1. )      ! Lateral boundary cond.
1135
1136      ENDIF ! ln_mask_file=.TRUE.
1137     
1138      IF( .not. ln_mask_file ) THEN
1139         ! If .not. ln_mask_file then we need to derive mask on U and V grid
1140         ! from mask on T grid here.
1141         bdyumask(:,:) = 0.e0
1142         bdyvmask(:,:) = 0.e0
1143         DO ij=1, jpjm1
1144            DO ii=1, jpim1
1145               bdyumask(ii,ij)=bdytmask(ii,ij)*bdytmask(ii+1, ij )
1146               bdyvmask(ii,ij)=bdytmask(ii,ij)*bdytmask(ii  ,ij+1) 
1147            END DO
1148         END DO
1149         CALL lbc_lnk( bdyumask(:,:), 'U', 1. )   ;   CALL lbc_lnk( bdyvmask(:,:), 'V', 1. )      ! Lateral boundary cond.
1150      ENDIF
1151
1152      ! bdy masks and bmask are now set to zero on boundary points:
1153      igrd = 1       ! In the free surface case, bmask is at T-points
1154      DO ib_bdy = 1, nb_bdy
1155        DO ib = 1, idx_bdy(ib_bdy)%nblenrim(igrd)     
1156          bmask(idx_bdy(ib_bdy)%nbi(ib,igrd), idx_bdy(ib_bdy)%nbj(ib,igrd)) = 0.e0
1157        ENDDO
1158      ENDDO
1159      !
1160      igrd = 1
1161      DO ib_bdy = 1, nb_bdy
1162        DO ib = 1, idx_bdy(ib_bdy)%nblenrim(igrd)     
1163          bdytmask(idx_bdy(ib_bdy)%nbi(ib,igrd), idx_bdy(ib_bdy)%nbj(ib,igrd)) = 0.e0
1164        ENDDO
1165      ENDDO
1166      !
1167      igrd = 2
1168      DO ib_bdy = 1, nb_bdy
1169        DO ib = 1, idx_bdy(ib_bdy)%nblenrim(igrd)
1170          bdyumask(idx_bdy(ib_bdy)%nbi(ib,igrd), idx_bdy(ib_bdy)%nbj(ib,igrd)) = 0.e0
1171        ENDDO
1172      ENDDO
1173      !
1174      igrd = 3
1175      DO ib_bdy = 1, nb_bdy
1176        DO ib = 1, idx_bdy(ib_bdy)%nblenrim(igrd)
1177          bdyvmask(idx_bdy(ib_bdy)%nbi(ib,igrd), idx_bdy(ib_bdy)%nbj(ib,igrd)) = 0.e0
1178        ENDDO
1179      ENDDO
1180
1181      ! For the flagu/flagv calculation below we require a version of fmask without
1182      ! the land boundary condition (shlat) included:
1183      CALL wrk_alloc(jpi,jpj,zfmask) 
1184      DO ij = 2, jpjm1
1185         DO ii = 2, jpim1
1186            zfmask(ii,ij) = tmask(ii,ij  ,1) * tmask(ii+1,ij  ,1)   &
1187           &              * tmask(ii,ij+1,1) * tmask(ii+1,ij+1,1)
1188         END DO     
1189      END DO
1190
1191      ! Lateral boundary conditions
1192      CALL lbc_lnk( zfmask       , 'F', 1. )
1193      CALL lbc_lnk( fmask        , 'F', 1. )   ;   CALL lbc_lnk( bdytmask(:,:), 'T', 1. )
1194      CALL lbc_lnk( bdyumask(:,:), 'U', 1. )   ;   CALL lbc_lnk( bdyvmask(:,:), 'V', 1. )
1195
1196      DO ib_bdy = 1, nb_bdy       ! Indices and directions of rim velocity components
1197
1198         idx_bdy(ib_bdy)%flagu(:,:) = 0.e0
1199         idx_bdy(ib_bdy)%flagv(:,:) = 0.e0
1200         icount = 0 
1201
1202         ! Calculate relationship of U direction to the local orientation of the boundary
1203         ! flagu = -1 : u component is normal to the dynamical boundary and its direction is outward
1204         ! flagu =  0 : u is tangential
1205         ! flagu =  1 : u is normal to the boundary and is direction is inward
1206 
1207         DO igrd = 1,jpbgrd 
1208            SELECT CASE( igrd )
1209               CASE( 1 )
1210                  pmask => umask(:,:,1)
1211                  i_offset = 0
1212               CASE( 2 ) 
1213                  pmask => bdytmask
1214                  i_offset = 1
1215               CASE( 3 ) 
1216                  pmask => zfmask(:,:)
1217                  i_offset = 0
1218            END SELECT
1219            icount = 0
1220            DO ib = 1, idx_bdy(ib_bdy)%nblenrim(igrd) 
1221               nbi => idx_bdy(ib_bdy)%nbi(ib,igrd)
1222               nbj => idx_bdy(ib_bdy)%nbj(ib,igrd)
1223               zefl = pmask(nbi+i_offset-1,nbj)
1224               zwfl = pmask(nbi+i_offset,nbj)
1225               ! This error check only works if you are using the bdyXmask arrays
1226               IF( i_offset == 1 .and. zefl + zwfl == 2 ) THEN
1227                  icount = icount + 1
1228                  IF(lwp) WRITE(numout,*) 'Problem with igrd = ',igrd,' at (global) nbi, nbj : ',mig(nbi),mjg(nbj)
1229               ELSE
1230                  idx_bdy(ib_bdy)%flagu(ib,igrd) = -zefl + zwfl
1231               ENDIF
1232            END DO
1233            IF( icount /= 0 ) THEN
1234               IF(lwp) WRITE(numout,*)
1235               IF(lwp) WRITE(numout,*) ' E R R O R : Some ',cgrid(igrd),' grid points,',   &
1236                  ' are not boundary points (flagu calculation). Check nbi, nbj, indices for boundary set ',ib_bdy
1237               IF(lwp) WRITE(numout,*) ' ========== '
1238               IF(lwp) WRITE(numout,*)
1239               nstop = nstop + 1
1240            ENDIF
1241         END DO
1242
1243         ! Calculate relationship of V direction to the local orientation of the boundary
1244         ! flagv = -1 : v component is normal to the dynamical boundary but its direction is outward
1245         ! flagv =  0 : v is tangential
1246         ! flagv =  1 : v is normal to the boundary and is direction is inward
1247
1248         DO igrd = 1,jpbgrd 
1249            SELECT CASE( igrd )
1250               CASE( 1 )
1251                  pmask => vmask(:,:,1)
1252                  j_offset = 0
1253               CASE( 2 )
1254                  pmask => zfmask(:,:)
1255                  j_offset = 0
1256               CASE( 3 )
1257                  pmask => bdytmask
1258                  j_offset = 1
1259            END SELECT
1260            icount = 0
1261            DO ib = 1, idx_bdy(ib_bdy)%nblenrim(igrd) 
1262               nbi => idx_bdy(ib_bdy)%nbi(ib,igrd)
1263               nbj => idx_bdy(ib_bdy)%nbj(ib,igrd)
1264               znfl = pmask(nbi,nbj+j_offset-1  )
1265               zsfl = pmask(nbi,nbj+j_offset)
1266               ! This error check only works if you are using the bdyXmask arrays
1267               IF( j_offset == 1 .and. znfl + zsfl == 2 ) THEN
1268                  IF(lwp) WRITE(numout,*) 'Problem with igrd = ',igrd,' at (global) nbi, nbj : ',mig(nbi),mjg(nbj)
1269                  icount = icount + 1
1270               ELSE
1271                  idx_bdy(ib_bdy)%flagv(ib,igrd) = -znfl + zsfl
1272               END IF
1273            END DO
1274            IF( icount /= 0 ) THEN
1275               IF(lwp) WRITE(numout,*)
1276               IF(lwp) WRITE(numout,*) ' E R R O R : Some ',cgrid(igrd),' grid points,',   &
1277                  ' are not boundary points (flagv calculation). Check nbi, nbj, indices for boundary set ',ib_bdy
1278               IF(lwp) WRITE(numout,*) ' ========== '
1279               IF(lwp) WRITE(numout,*)
1280               nstop = nstop + 1
1281            ENDIF
1282         END DO
1283
1284      END DO
1285
1286      ! Compute total lateral surface for volume correction:
1287      ! ----------------------------------------------------
1288      ! JC: this must be done at each time step with key_vvl
1289      bdysurftot = 0.e0 
1290      IF( ln_vol ) THEN 
1291         igrd = 2      ! Lateral surface at U-points
1292         DO ib_bdy = 1, nb_bdy
1293            DO ib = 1, idx_bdy(ib_bdy)%nblenrim(igrd)
1294               nbi => idx_bdy(ib_bdy)%nbi(ib,igrd)
1295               nbj => idx_bdy(ib_bdy)%nbj(ib,igrd)
1296               flagu => idx_bdy(ib_bdy)%flagu(ib,igrd)
1297               bdysurftot = bdysurftot + hu     (nbi  , nbj)                           &
1298                  &                    * e2u    (nbi  , nbj) * ABS( flagu ) &
1299                  &                    * tmask_i(nbi  , nbj)                           &
1300                  &                    * tmask_i(nbi+1, nbj)                   
1301            ENDDO
1302         ENDDO
1303
1304         igrd=3 ! Add lateral surface at V-points
1305         DO ib_bdy = 1, nb_bdy
1306            DO ib = 1, idx_bdy(ib_bdy)%nblenrim(igrd)
1307               nbi => idx_bdy(ib_bdy)%nbi(ib,igrd)
1308               nbj => idx_bdy(ib_bdy)%nbj(ib,igrd)
1309               flagv => idx_bdy(ib_bdy)%flagv(ib,igrd)
1310               bdysurftot = bdysurftot + hv     (nbi, nbj  )                           &
1311                  &                    * e1v    (nbi, nbj  ) * ABS( flagv ) &
1312                  &                    * tmask_i(nbi, nbj  )                           &
1313                  &                    * tmask_i(nbi, nbj+1)
1314            ENDDO
1315         ENDDO
1316         !
1317         IF( lk_mpp )   CALL mpp_sum( bdysurftot )      ! sum over the global domain
1318      END IF   
1319      !
1320      ! Tidy up
1321      !--------
1322      IF (nb_bdy>0) THEN
1323         DEALLOCATE(nbidta, nbjdta, nbrdta)
1324      ENDIF
1325
1326      CALL wrk_dealloc(jpi,jpj,zfmask) 
1327
1328      IF( nn_timing == 1 ) CALL timing_stop('bdy_init')
1329
1330   END SUBROUTINE bdy_init
1331
1332   SUBROUTINE bdy_ctl_seg
1333      !!----------------------------------------------------------------------
1334      !!                 ***  ROUTINE bdy_ctl_seg  ***
1335      !!
1336      !! ** Purpose :   Check straight open boundary segments location
1337      !!
1338      !! ** Method  :   - Look for open boundary corners
1339      !!                - Check that segments start or end on land
1340      !!----------------------------------------------------------------------
1341      INTEGER  ::   ib, ib1, ib2, ji ,jj, itest 
1342      INTEGER, DIMENSION(jp_nseg,2) :: icorne, icornw, icornn, icorns 
1343      REAL(wp), DIMENSION(2) ::   ztestmask
1344      !!----------------------------------------------------------------------
1345      !
1346      IF (lwp) WRITE(numout,*) ' '
1347      IF (lwp) WRITE(numout,*) 'bdy_ctl_seg: Check analytical segments'
1348      IF (lwp) WRITE(numout,*) '~~~~~~~~~~~~'
1349      !
1350      IF(lwp) WRITE(numout,*) 'Number of east  segments     : ', nbdysege
1351      IF(lwp) WRITE(numout,*) 'Number of west  segments     : ', nbdysegw
1352      IF(lwp) WRITE(numout,*) 'Number of north segments     : ', nbdysegn
1353      IF(lwp) WRITE(numout,*) 'Number of south segments     : ', nbdysegs
1354      ! 1. Check bounds
1355      !----------------
1356      DO ib = 1, nbdysegn
1357         IF (lwp) WRITE(numout,*) '**check north seg bounds pckg: ', npckgn(ib)
1358         IF ((jpjnob(ib).ge.jpjglo-1).or.& 
1359            &(jpjnob(ib).le.1))        CALL ctl_stop( 'nbdyind out of domain' )
1360         IF (jpindt(ib).ge.jpinft(ib)) CALL ctl_stop( 'Bdy start index is greater than end index' )
1361         IF (jpindt(ib).le.1     )     CALL ctl_stop( 'Start index out of domain' )
1362         IF (jpinft(ib).ge.jpiglo)     CALL ctl_stop( 'End index out of domain' )
1363      END DO
1364      !
1365      DO ib = 1, nbdysegs
1366         IF (lwp) WRITE(numout,*) '**check south seg bounds pckg: ', npckgs(ib)
1367         IF ((jpjsob(ib).ge.jpjglo-1).or.& 
1368            &(jpjsob(ib).le.1))        CALL ctl_stop( 'nbdyind out of domain' )
1369         IF (jpisdt(ib).ge.jpisft(ib)) CALL ctl_stop( 'Bdy start index is greater than end index' )
1370         IF (jpisdt(ib).le.1     )     CALL ctl_stop( 'Start index out of domain' )
1371         IF (jpisft(ib).ge.jpiglo)     CALL ctl_stop( 'End index out of domain' )
1372      END DO
1373      !
1374      DO ib = 1, nbdysege
1375         IF (lwp) WRITE(numout,*) '**check east  seg bounds pckg: ', npckge(ib)
1376         IF ((jpieob(ib).ge.jpiglo-1).or.& 
1377            &(jpieob(ib).le.1))        CALL ctl_stop( 'nbdyind out of domain' )
1378         IF (jpjedt(ib).ge.jpjeft(ib)) CALL ctl_stop( 'Bdy start index is greater than end index' )
1379         IF (jpjedt(ib).le.1     )     CALL ctl_stop( 'Start index out of domain' )
1380         IF (jpjeft(ib).ge.jpjglo)     CALL ctl_stop( 'End index out of domain' )
1381      END DO
1382      !
1383      DO ib = 1, nbdysegw
1384         IF (lwp) WRITE(numout,*) '**check west  seg bounds pckg: ', npckgw(ib)
1385         IF ((jpiwob(ib).ge.jpiglo-1).or.& 
1386            &(jpiwob(ib).le.1))        CALL ctl_stop( 'nbdyind out of domain' )
1387         IF (jpjwdt(ib).ge.jpjwft(ib)) CALL ctl_stop( 'Bdy start index is greater than end index' )
1388         IF (jpjwdt(ib).le.1     )     CALL ctl_stop( 'Start index out of domain' )
1389         IF (jpjwft(ib).ge.jpjglo)     CALL ctl_stop( 'End index out of domain' )
1390      ENDDO
1391      !
1392      !     
1393      ! 2. Look for segment crossings
1394      !------------------------------
1395      IF (lwp) WRITE(numout,*) '**Look for segments corners  :'
1396      !
1397      itest = 0 ! corner number
1398      !
1399      ! flag to detect if start or end of open boundary belongs to a corner
1400      ! if not (=0), it must be on land.
1401      ! if a corner is detected, save bdy package number for further tests
1402      icorne(:,:)=0. ; icornw(:,:)=0. ; icornn(:,:)=0. ; icorns(:,:)=0.
1403      ! South/West crossings
1404      IF ((nbdysegw > 0).AND.(nbdysegs > 0)) THEN
1405         DO ib1 = 1, nbdysegw       
1406            DO ib2 = 1, nbdysegs
1407               IF (( jpisdt(ib2)<=jpiwob(ib1)).AND. &
1408                &  ( jpisft(ib2)>=jpiwob(ib1)).AND. &
1409                &  ( jpjwdt(ib1)<=jpjsob(ib2)).AND. &
1410                &  ( jpjwft(ib1)>=jpjsob(ib2))) THEN
1411                  IF ((jpjwdt(ib1)==jpjsob(ib2)).AND.(jpisdt(ib2)==jpiwob(ib1))) THEN 
1412                     ! We have a possible South-West corner                     
1413!                     WRITE(numout,*) ' Found a South-West corner at (i,j): ', jpisdt(ib2), jpjwdt(ib1)
1414!                     WRITE(numout,*) ' between segments: ', npckgw(ib1), npckgs(ib2)
1415                     icornw(ib1,1) = npckgs(ib2)
1416                     icorns(ib2,1) = npckgw(ib1)
1417                  ELSEIF ((jpisft(ib2)==jpiwob(ib1)).AND.(jpjwft(ib1)==jpjsob(ib2))) THEN
1418                     IF(lwp) WRITE(numout,*)
1419                     IF(lwp) WRITE(numout,*) ' E R R O R : Found an acute open boundary corner at point (i,j)= ', &
1420                     &                                     jpisft(ib2), jpjwft(ib1)
1421                     IF(lwp) WRITE(numout,*) ' ==========  Not allowed yet'
1422                     IF(lwp) WRITE(numout,*) '             Crossing problem with West segment: ',npckgw(ib1), & 
1423                     &                                                    ' and South segment: ',npckgs(ib2)
1424                     IF(lwp) WRITE(numout,*)
1425                     nstop = nstop + 1
1426                  ELSE
1427                     IF(lwp) WRITE(numout,*)
1428                     IF(lwp) WRITE(numout,*) ' E R R O R : Check South and West Open boundary indices'
1429                     IF(lwp) WRITE(numout,*) ' ==========  Crossing problem with West segment: ',npckgw(ib1) , &
1430                     &                                                    ' and South segment: ',npckgs(ib2)
1431                     IF(lwp) WRITE(numout,*)
1432                     nstop = nstop+1
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                     IF(lwp) WRITE(numout,*)
1455                     IF(lwp) WRITE(numout,*) ' E R R O R : Found an acute open boundary corner at point (i,j)= ', &
1456                     &                                     jpisdt(ib2), jpjeft(ib1)
1457                     IF(lwp) WRITE(numout,*) ' ==========  Not allowed yet'
1458                     IF(lwp) WRITE(numout,*) '             Crossing problem with East segment: ',npckge(ib1), &
1459                     &                                                    ' and South segment: ',npckgs(ib2)
1460                     IF(lwp) WRITE(numout,*)
1461                     nstop = nstop + 1
1462                  ELSE
1463                     IF(lwp) WRITE(numout,*)
1464                     IF(lwp) WRITE(numout,*) ' E R R O R : Check South and East Open boundary indices'
1465                     IF(lwp) WRITE(numout,*) ' ==========  Crossing problem with East segment: ',npckge(ib1), &
1466                     &                                                    ' and South segment: ',npckgs(ib2)
1467                     IF(lwp) WRITE(numout,*)
1468                     nstop = nstop + 1
1469                  END IF
1470               END IF
1471            END DO
1472         END DO
1473      END IF
1474      !
1475      ! North/West crossings
1476      IF ((nbdysegn > 0).AND.(nbdysegw > 0)) THEN
1477         DO ib1 = 1, nbdysegw       
1478            DO ib2 = 1, nbdysegn
1479               IF (( jpindt(ib2)<=jpiwob(ib1)  ).AND. &
1480                &  ( jpinft(ib2)>=jpiwob(ib1)  ).AND. &
1481                &  ( jpjwdt(ib1)<=jpjnob(ib2)+1).AND. &
1482                &  ( jpjwft(ib1)>=jpjnob(ib2)+1)) THEN
1483                  IF ((jpjwft(ib1)==jpjnob(ib2)+1).AND.(jpindt(ib2)==jpiwob(ib1))) THEN
1484                     ! We have a possible North-West corner
1485!                     WRITE(numout,*) ' Found a North-West corner at (i,j): ', jpindt(ib2), jpjwft(ib1)
1486!                     WRITE(numout,*) ' between segments: ', npckgw(ib1), npckgn(ib2)
1487                     icornw(ib1,2) = npckgn(ib2)
1488                     icornn(ib2,1) = npckgw(ib1)
1489                  ELSEIF ((jpjwdt(ib1)==jpjnob(ib2)+1).AND.(jpinft(ib2)==jpiwob(ib1))) THEN
1490                     IF(lwp) WRITE(numout,*)
1491                     IF(lwp) WRITE(numout,*) ' E R R O R : Found an acute open boundary corner at point (i,j)= ', &
1492                     &                                     jpinft(ib2), jpjwdt(ib1)
1493                     IF(lwp) WRITE(numout,*) ' ==========  Not allowed yet'
1494                     IF(lwp) WRITE(numout,*) '             Crossing problem with West segment: ',npckgw(ib1), &
1495                     &                                                    ' and North segment: ',npckgn(ib2)
1496                     IF(lwp) WRITE(numout,*)
1497                     nstop = nstop + 1
1498                  ELSE
1499                     IF(lwp) WRITE(numout,*)
1500                     IF(lwp) WRITE(numout,*) ' E R R O R : Check North and West Open boundary indices'
1501                     IF(lwp) WRITE(numout,*) ' ==========  Crossing problem with West segment: ',npckgw(ib1), &
1502                     &                                                    ' and North segment: ',npckgn(ib2)
1503                     IF(lwp) WRITE(numout,*)
1504                     nstop = nstop + 1
1505                  END IF
1506               END IF
1507            END DO
1508         END DO
1509      END IF
1510      !
1511      ! North/East crossings
1512      IF ((nbdysegn > 0).AND.(nbdysege > 0)) THEN
1513         DO ib1 = 1, nbdysege       
1514            DO ib2 = 1, nbdysegn
1515               IF (( jpindt(ib2)<=jpieob(ib1)+1).AND. &
1516                &  ( jpinft(ib2)>=jpieob(ib1)+1).AND. &
1517                &  ( jpjedt(ib1)<=jpjnob(ib2)+1).AND. &
1518                &  ( jpjeft(ib1)>=jpjnob(ib2)+1)) THEN
1519                  IF ((jpjeft(ib1)==jpjnob(ib2)+1).AND.(jpinft(ib2)==jpieob(ib1)+1)) THEN
1520                     ! We have a possible North-East corner
1521!                     WRITE(numout,*) ' Found a North-East corner at (i,j): ', jpinft(ib2), jpjeft(ib1)
1522!                     WRITE(numout,*) ' between segments: ', npckge(ib1), npckgn(ib2)
1523                     icorne(ib1,2) = npckgn(ib2)
1524                     icornn(ib2,2) = npckge(ib1)
1525                  ELSEIF ((jpjedt(ib1)==jpjnob(ib2)+1).AND.(jpindt(ib2)==jpieob(ib1)+1)) THEN
1526                     IF(lwp) WRITE(numout,*)
1527                     IF(lwp) WRITE(numout,*) ' E R R O R : Found an acute open boundary corner at point (i,j)= ', &
1528                     &                                     jpindt(ib2), jpjedt(ib1)
1529                     IF(lwp) WRITE(numout,*) ' ==========  Not allowed yet'
1530                     IF(lwp) WRITE(numout,*) '             Crossing problem with East segment: ',npckge(ib1), &
1531                     &                                                    ' and North segment: ',npckgn(ib2)
1532                     IF(lwp) WRITE(numout,*)
1533                     nstop = nstop + 1
1534                  ELSE
1535                     IF(lwp) WRITE(numout,*)
1536                     IF(lwp) WRITE(numout,*) ' E R R O R : Check North and East Open boundary indices'
1537                     IF(lwp) WRITE(numout,*) ' ==========  Crossing problem with East segment: ',npckge(ib1), &
1538                     &                                                    ' and North segment: ',npckgn(ib2)
1539                     IF(lwp) WRITE(numout,*)
1540                     nstop = nstop + 1
1541                  END IF
1542               END IF
1543            END DO
1544         END DO
1545      END IF
1546      !
1547      ! 3. Check if segment extremities are on land
1548      !--------------------------------------------
1549      !
1550      ! West segments
1551      DO ib = 1, nbdysegw
1552         ! get mask at boundary extremities:
1553         ztestmask(1:2)=0.
1554         DO ji = 1, jpi
1555            DO jj = 1, jpj             
1556              IF (((ji + nimpp - 1) == jpiwob(ib)).AND. & 
1557               &  ((jj + njmpp - 1) == jpjwdt(ib))) ztestmask(1)=tmask(ji,jj,1)
1558              IF (((ji + nimpp - 1) == jpiwob(ib)).AND. & 
1559               &  ((jj + njmpp - 1) == jpjwft(ib))) ztestmask(2)=tmask(ji,jj,1) 
1560            END DO
1561         END DO
1562         IF( lk_mpp )   CALL mpp_sum( ztestmask, 2 )   ! sum over the global domain
1563
1564         IF (ztestmask(1)==1) THEN
1565            IF (icornw(ib,1)==0) THEN
1566               IF(lwp) WRITE(numout,*)
1567               IF(lwp) WRITE(numout,*) ' E R R O R : Open boundary segment ', npckgw(ib)
1568               IF(lwp) WRITE(numout,*) ' ==========  does not start on land or on a corner'                                                 
1569               IF(lwp) WRITE(numout,*)
1570               nstop = nstop + 1
1571            ELSE
1572               ! This is a corner
1573               IF(lwp) WRITE(numout,*) 'Found a South-West corner at (i,j): ', jpiwob(ib), jpjwdt(ib)
1574               CALL bdy_ctl_corn(npckgw(ib), icornw(ib,1))
1575               itest=itest+1
1576            ENDIF
1577         ENDIF
1578         IF (ztestmask(2)==1) THEN
1579            IF (icornw(ib,2)==0) THEN
1580               IF(lwp) WRITE(numout,*)
1581               IF(lwp) WRITE(numout,*) ' E R R O R : Open boundary segment ', npckgw(ib)
1582               IF(lwp) WRITE(numout,*) ' ==========  does not end on land or on a corner'                                                 
1583               IF(lwp) WRITE(numout,*)
1584               nstop = nstop + 1
1585            ELSE
1586               ! This is a corner
1587               IF(lwp) WRITE(numout,*) 'Found a North-West corner at (i,j): ', jpiwob(ib), jpjwft(ib)
1588               CALL bdy_ctl_corn(npckgw(ib), icornw(ib,2))
1589               itest=itest+1
1590            ENDIF
1591         ENDIF
1592      END DO
1593      !
1594      ! East segments
1595      DO ib = 1, nbdysege
1596         ! get mask at boundary extremities:
1597         ztestmask(1:2)=0.
1598         DO ji = 1, jpi
1599            DO jj = 1, jpj             
1600              IF (((ji + nimpp - 1) == jpieob(ib)+1).AND. & 
1601               &  ((jj + njmpp - 1) == jpjedt(ib))) ztestmask(1)=tmask(ji,jj,1)
1602              IF (((ji + nimpp - 1) == jpieob(ib)+1).AND. & 
1603               &  ((jj + njmpp - 1) == jpjeft(ib))) ztestmask(2)=tmask(ji,jj,1) 
1604            END DO
1605         END DO
1606         IF( lk_mpp )   CALL mpp_sum( ztestmask, 2 )   ! sum over the global domain
1607
1608         IF (ztestmask(1)==1) THEN
1609            IF (icorne(ib,1)==0) THEN
1610               IF(lwp) WRITE(numout,*)
1611               IF(lwp) WRITE(numout,*) ' E R R O R : Open boundary segment ', npckge(ib)
1612               IF(lwp) WRITE(numout,*) ' ==========  does not start on land or on a corner'                                                 
1613               IF(lwp) WRITE(numout,*)
1614               nstop = nstop + 1 
1615            ELSE
1616               ! This is a corner
1617               IF(lwp) WRITE(numout,*) 'Found a South-East corner at (i,j): ', jpieob(ib)+1, jpjedt(ib)
1618               CALL bdy_ctl_corn(npckge(ib), icorne(ib,1))
1619               itest=itest+1
1620            ENDIF
1621         ENDIF
1622         IF (ztestmask(2)==1) THEN
1623            IF (icorne(ib,2)==0) THEN
1624               IF(lwp) WRITE(numout,*)
1625               IF(lwp) WRITE(numout,*) ' E R R O R : Open boundary segment ', npckge(ib)
1626               IF(lwp) WRITE(numout,*) ' ==========  does not end on land or on a corner'                                                 
1627               IF(lwp) WRITE(numout,*)
1628               nstop = nstop + 1
1629            ELSE
1630               ! This is a corner
1631               IF(lwp) WRITE(numout,*) 'Found a North-East corner at (i,j): ', jpieob(ib)+1, jpjeft(ib)
1632               CALL bdy_ctl_corn(npckge(ib), icorne(ib,2))
1633               itest=itest+1
1634            ENDIF
1635         ENDIF
1636      END DO
1637      !
1638      ! South segments
1639      DO ib = 1, nbdysegs
1640         ! get mask at boundary extremities:
1641         ztestmask(1:2)=0.
1642         DO ji = 1, jpi
1643            DO jj = 1, jpj             
1644              IF (((jj + njmpp - 1) == jpjsob(ib)).AND. & 
1645               &  ((ji + nimpp - 1) == jpisdt(ib))) ztestmask(1)=tmask(ji,jj,1)
1646              IF (((jj + njmpp - 1) == jpjsob(ib)).AND. & 
1647               &  ((ji + nimpp - 1) == jpisft(ib))) ztestmask(2)=tmask(ji,jj,1) 
1648            END DO
1649         END DO
1650         IF( lk_mpp )   CALL mpp_sum( ztestmask, 2 )   ! sum over the global domain
1651
1652         IF ((ztestmask(1)==1).AND.(icorns(ib,1)==0)) THEN
1653            IF(lwp) WRITE(numout,*)
1654            IF(lwp) WRITE(numout,*) ' E R R O R : Open boundary segment ', npckgs(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         ENDIF
1659         IF ((ztestmask(2)==1).AND.(icorns(ib,2)==0)) THEN
1660            IF(lwp) WRITE(numout,*)
1661            IF(lwp) WRITE(numout,*) ' E R R O R : Open boundary segment ', npckgs(ib)
1662            IF(lwp) WRITE(numout,*) ' ==========  does not end on land or on a corner'                                                 
1663            IF(lwp) WRITE(numout,*)
1664            nstop = nstop + 1
1665         ENDIF
1666      END DO
1667      !
1668      ! North segments
1669      DO ib = 1, nbdysegn
1670         ! get mask at boundary extremities:
1671         ztestmask(1:2)=0.
1672         DO ji = 1, jpi
1673            DO jj = 1, jpj             
1674              IF (((jj + njmpp - 1) == jpjnob(ib)+1).AND. & 
1675               &  ((ji + nimpp - 1) == jpindt(ib))) ztestmask(1)=tmask(ji,jj,1)
1676              IF (((jj + njmpp - 1) == jpjnob(ib)+1).AND. & 
1677               &  ((ji + nimpp - 1) == jpinft(ib))) ztestmask(2)=tmask(ji,jj,1) 
1678            END DO
1679         END DO
1680         IF( lk_mpp )   CALL mpp_sum( ztestmask, 2 )   ! sum over the global domain
1681
1682         IF ((ztestmask(1)==1).AND.(icornn(ib,1)==0)) THEN
1683            IF(lwp) WRITE(numout,*)
1684            IF(lwp) WRITE(numout,*) ' E R R O R : Open boundary segment ', npckgn(ib)
1685            IF(lwp) WRITE(numout,*) ' ==========  does not start on land'                                                 
1686            IF(lwp) WRITE(numout,*)
1687            nstop = nstop + 1
1688         ENDIF
1689         IF ((ztestmask(2)==1).AND.(icornn(ib,2)==0)) THEN
1690            IF(lwp) WRITE(numout,*)
1691            IF(lwp) WRITE(numout,*) ' E R R O R : Open boundary segment ', npckgn(ib)
1692            IF(lwp) WRITE(numout,*) ' ==========  does not end on land'                                                 
1693            IF(lwp) WRITE(numout,*)
1694            nstop = nstop + 1
1695         ENDIF
1696      END DO
1697      !
1698      IF ((itest==0).AND.(lwp)) WRITE(numout,*) 'NO open boundary corner found'
1699      !
1700      ! Other tests TBD:
1701      ! segments completly on land
1702      ! optimized open boundary array length according to landmask
1703      ! Nudging layers that overlap with interior domain
1704      !
1705   END SUBROUTINE bdy_ctl_seg
1706
1707   SUBROUTINE bdy_ctl_corn( ib1, ib2 )
1708      !!----------------------------------------------------------------------
1709      !!                 ***  ROUTINE bdy_ctl_corn  ***
1710      !!
1711      !! ** Purpose :   Check numerical schemes consistency between
1712      !!                segments having a common corner
1713      !!
1714      !! ** Method  :   
1715      !!----------------------------------------------------------------------
1716      INTEGER, INTENT(in)  ::   ib1, ib2
1717      INTEGER :: itest
1718      !!----------------------------------------------------------------------
1719      itest = 0
1720
1721      IF (cn_dyn2d(ib1)/=cn_dyn2d(ib2)) itest = itest + 1
1722      IF (cn_dyn3d(ib1)/=cn_dyn3d(ib2)) itest = itest + 1
1723      IF (cn_tra(ib1)/=cn_tra(ib2)) itest = itest + 1
1724      !
1725      IF (nn_dyn2d_dta(ib1)/=nn_dyn2d_dta(ib2)) itest = itest + 1
1726      IF (nn_dyn3d_dta(ib1)/=nn_dyn3d_dta(ib2)) itest = itest + 1
1727      IF (nn_tra_dta(ib1)/=nn_tra_dta(ib2)) itest = itest + 1
1728      !
1729      IF (nn_rimwidth(ib1)/=nn_rimwidth(ib2)) itest = itest + 1   
1730      !
1731      IF ( itest>0 ) THEN
1732         IF(lwp) WRITE(numout,*) ' E R R O R : Segments ', ib1, 'and ', ib2
1733         IF(lwp) WRITE(numout,*) ' ==========  have different open bdy schemes'                                                 
1734         IF(lwp) WRITE(numout,*)
1735         nstop = nstop + 1
1736      ENDIF
1737      !
1738   END SUBROUTINE bdy_ctl_corn
1739
1740#else
1741   !!---------------------------------------------------------------------------------
1742   !!   Dummy module                                   NO open boundaries
1743   !!---------------------------------------------------------------------------------
1744CONTAINS
1745   SUBROUTINE bdy_init      ! Dummy routine
1746   END SUBROUTINE bdy_init
1747#endif
1748
1749   !!=================================================================================
1750END MODULE bdyini
Note: See TracBrowser for help on using the repository browser.