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

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

dev_r10984_HPC-13 : new implementation of lbc_bdy_lnk in prevision of step 2, regroup communications, see #2285

  • Property svn:keywords set to Id
File size: 76.0 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                                               ! open boundary data files
40   ! Straight open boundary segment parameters:
41   INTEGER  ::   nbdysege, nbdysegw, nbdysegn, nbdysegs 
42   INTEGER, DIMENSION(jp_nseg) ::   jpieob, jpjedt, jpjeft, npckge   !
43   INTEGER, DIMENSION(jp_nseg) ::   jpiwob, jpjwdt, jpjwft, npckgw   !
44   INTEGER, DIMENSION(jp_nseg) ::   jpjnob, jpindt, jpinft, npckgn   !
45   INTEGER, DIMENSION(jp_nseg) ::   jpjsob, jpisdt, jpisft, npckgs   !
46   !!----------------------------------------------------------------------
47   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
48   !! $Id$
49   !! Software governed by the CeCILL license (see ./LICENSE)
50   !!----------------------------------------------------------------------
51CONTAINS
52
53   SUBROUTINE bdy_init
54      !!----------------------------------------------------------------------
55      !!                 ***  ROUTINE bdy_init  ***
56      !!
57      !! ** Purpose :   Initialization of the dynamics and tracer fields with
58      !!              unstructured open boundaries.
59      !!
60      !! ** Method  :   Read initialization arrays (mask, indices) to identify
61      !!              an unstructured open boundary
62      !!
63      !! ** Input   :  bdy_init.nc, input file for unstructured open boundaries
64      !!----------------------------------------------------------------------
65      NAMELIST/nambdy/ ln_bdy, nb_bdy, ln_coords_file, cn_coords_file,         &
66         &             ln_mask_file, cn_mask_file, cn_dyn2d, nn_dyn2d_dta,     &
67         &             cn_dyn3d, nn_dyn3d_dta, cn_tra, nn_tra_dta,             &
68         &             ln_tra_dmp, ln_dyn3d_dmp, rn_time_dmp, rn_time_dmp_out, &
69         &             cn_ice, nn_ice_dta,                                     &
70         &             rn_ice_tem, rn_ice_sal, rn_ice_age,                     &
71         &             ln_vol, nn_volctl, nn_rimwidth, nb_jpk_bdy
72         !
73      INTEGER  ::   ios                 ! Local integer output status for namelist read
74      !!----------------------------------------------------------------------
75
76      ! ------------------------
77      ! Read namelist parameters
78      ! ------------------------
79      REWIND( numnam_ref )              ! Namelist nambdy in reference namelist :Unstructured open boundaries
80      READ  ( numnam_ref, nambdy, IOSTAT = ios, ERR = 901)
81901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'nambdy in reference namelist', lwp )
82      REWIND( numnam_cfg )              ! Namelist nambdy in configuration namelist :Unstructured open boundaries
83      READ  ( numnam_cfg, nambdy, IOSTAT = ios, ERR = 902 )
84902   IF( ios >  0 )   CALL ctl_nam ( ios , 'nambdy in configuration namelist', lwp )
85      IF(lwm) WRITE ( numond, nambdy )
86
87      IF( .NOT. Agrif_Root() ) ln_bdy = .FALSE.   ! forced for Agrif children
88     
89      ! -----------------------------------------
90      ! unstructured open boundaries use control
91      ! -----------------------------------------
92      IF ( ln_bdy ) THEN
93         IF(lwp) WRITE(numout,*)
94         IF(lwp) WRITE(numout,*) 'bdy_init : initialization of open boundaries'
95         IF(lwp) WRITE(numout,*) '~~~~~~~~'
96         !
97         ! Open boundaries definition (arrays and masks)
98         CALL bdy_segs
99         !
100         ! Open boundaries initialisation of external data arrays
101         CALL bdy_dta_init
102         !
103         ! Open boundaries initialisation of tidal harmonic forcing
104         IF( ln_tide ) CALL bdytide_init
105         !
106      ELSE
107         IF(lwp) WRITE(numout,*)
108         IF(lwp) WRITE(numout,*) 'bdy_init : open boundaries not used (ln_bdy = F)'
109         IF(lwp) WRITE(numout,*) '~~~~~~~~'
110         !
111      ENDIF
112      !
113   END SUBROUTINE bdy_init
114
115
116   SUBROUTINE bdy_segs
117      !!----------------------------------------------------------------------
118      !!                 ***  ROUTINE bdy_init  ***
119      !!         
120      !! ** Purpose :   Definition of unstructured open boundaries.
121      !!
122      !! ** Method  :   Read initialization arrays (mask, indices) to identify
123      !!              an unstructured open boundary
124      !!
125      !! ** Input   :  bdy_init.nc, input file for unstructured open boundaries
126      !!----------------------------------------------------------------------     
127      INTEGER  ::   ib_bdy, ii, ij, ik, igrd, ib, ir, iseg ! dummy loop indices
128      INTEGER  ::   icount, icountr, ibr_max, ilen1, ibm1  ! local integers
129      INTEGER  ::   iwe, ies, iso, ino, inum, id_dummy     !   -       -
130      INTEGER  ::   igrd_start, igrd_end, jpbdta           !   -       -
131      INTEGER  ::   jpbdtau, jpbdtas                       !   -       -
132      INTEGER  ::   ib_bdy1, ib_bdy2, ib1, ib2             !   -       -
133      INTEGER  ::   i_offset, j_offset, inbdy, itreat      !   -       -
134      INTEGER , POINTER  ::  nbi, nbj, nbr                 ! short cuts
135      REAL(wp), POINTER, DIMENSION(:,:)       ::   pmask    ! pointer to 2D mask fields
136      REAL(wp) ::   zefl, zwfl, znfl, zsfl                 ! local scalars
137      INTEGER, DIMENSION (2)                  ::   kdimsz
138      INTEGER, DIMENSION(jpbgrd,jp_bdy)       ::   nblendta         ! Length of index arrays
139      INTEGER, ALLOCATABLE, DIMENSION(:,:,:)  ::   nbidta, nbjdta   ! Index arrays: i and j indices of bdy dta
140      INTEGER, ALLOCATABLE, DIMENSION(:,:,:)  ::   nbrdta           ! Discrete distance from rim points
141      CHARACTER(LEN=1),DIMENSION(jpbgrd)      ::   cgrid
142      INTEGER :: jpk_max
143      LOGICAL :: llsend_ea, llsend_we, llsend_so, llsend_no  ! Flags for boundaries sending
144      LOGICAL :: llrecv_ea, llrecv_we, llrecv_so, llrecv_no  ! Flags for boundaries receiving
145      REAL(wp), TARGET, DIMENSION(jpi,jpj) ::   zfmask  ! temporary fmask array excluding coastal boundary condition (shlat)
146      REAL(wp)        , DIMENSION(jpi,jpj) ::   ztmp
147      REAL(wp), POINTER                    :: flagu,     flagv     ! short cuts
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      ! Find lenght of boundaries and rim on local mpi domain
794      !------------------------------------------------------
795      !
796      iwe = mig(1)
797      ies = mig(nlci)
798      iso = mjg(1) 
799      ino = mjg(nlcj) 
800      !
801      DO ib_bdy = 1, nb_bdy
802         DO igrd = 1, jpbgrd
803            icount  = 0   ! initialization of local bdy length
804            icountr = 0   ! initialization of local rim bdy length
805            idx_bdy(ib_bdy)%nblen(igrd)    = 0
806            idx_bdy(ib_bdy)%nblenrim(igrd) = 0
807            DO ib = 1, nblendta(igrd,ib_bdy)
808               ! check that data is in correct order in file
809               ibm1 = MAX(1,ib-1)
810               IF(lwp) THEN         ! Since all procs read global data only need to do this check on one proc...
811                  IF( nbrdta(ib,igrd,ib_bdy) < nbrdta(ibm1,igrd,ib_bdy) ) THEN
812                     CALL ctl_stop('bdy_segs : ERROR : boundary data in file must be defined ', &
813                          &        ' in order of distance from edge nbr A utility for re-ordering ', &
814                          &        ' boundary coordinates and data files exists in the TOOLS/OBC directory')
815                  ENDIF   
816               ENDIF
817               ! check if point is in local domain
818               IF(  nbidta(ib,igrd,ib_bdy) >= iwe .AND. nbidta(ib,igrd,ib_bdy) <= ies .AND.   &
819                  & nbjdta(ib,igrd,ib_bdy) >= iso .AND. nbjdta(ib,igrd,ib_bdy) <= ino      ) THEN
820                  !
821                  icount = icount + 1
822                  IF( nbrdta(ib,igrd,ib_bdy) == 1 )   icountr = icountr + 1
823               ENDIF
824            END DO
825            idx_bdy(ib_bdy)%nblenrim(igrd) = icountr !: length of rim boundary data on each proc
826            idx_bdy(ib_bdy)%nblen   (igrd) = icount  !: length of boundary data on each proc       
827         END DO   ! igrd
828
829         ! Allocate index arrays for this boundary set
830         !--------------------------------------------
831         ilen1 = MAXVAL( idx_bdy(ib_bdy)%nblen(:) )
832         ALLOCATE( idx_bdy(ib_bdy)%nbi   (ilen1,jpbgrd) ,   &
833            &      idx_bdy(ib_bdy)%nbj   (ilen1,jpbgrd) ,   &
834            &      idx_bdy(ib_bdy)%nbr   (ilen1,jpbgrd) ,   &
835            &      idx_bdy(ib_bdy)%nbd   (ilen1,jpbgrd) ,   &
836            &      idx_bdy(ib_bdy)%nbdout(ilen1,jpbgrd) ,   &
837            &      idx_bdy(ib_bdy)%ntreat(ilen1,jpbgrd) ,   &
838            &      idx_bdy(ib_bdy)%nbmap (ilen1,jpbgrd) ,   &
839            &      idx_bdy(ib_bdy)%nbw   (ilen1,jpbgrd) ,   &
840            &      idx_bdy(ib_bdy)%flagu (ilen1,jpbgrd) ,   &
841            &      idx_bdy(ib_bdy)%flagv (ilen1,jpbgrd) )
842
843         ! Dispatch mapping indices and discrete distances on each processor
844         ! -----------------------------------------------------------------
845         DO igrd = 1, jpbgrd
846            icount  = 0
847            ! Outer loop on rimwidth to ensure outermost points come first in the local arrays.
848            DO ir=1, nn_rimwidth(ib_bdy)
849               DO ib = 1, nblendta(igrd,ib_bdy)
850                  ! check if point is in local domain and equals ir
851                  IF(  nbidta(ib,igrd,ib_bdy) >= iwe .AND. nbidta(ib,igrd,ib_bdy) <= ies .AND.   &
852                     & nbjdta(ib,igrd,ib_bdy) >= iso .AND. nbjdta(ib,igrd,ib_bdy) <= ino .AND.   &
853                     & nbrdta(ib,igrd,ib_bdy) == ir  ) THEN
854                     !
855                     icount = icount  + 1
856                     idx_bdy(ib_bdy)%nbi(icount,igrd)   = nbidta(ib,igrd,ib_bdy)- mig(1)+1   ! global to local indexes
857                     idx_bdy(ib_bdy)%nbj(icount,igrd)   = nbjdta(ib,igrd,ib_bdy)- mjg(1)+1   ! global to local indexes
858                     idx_bdy(ib_bdy)%nbr(icount,igrd)   = nbrdta(ib,igrd,ib_bdy)
859                     idx_bdy(ib_bdy)%nbmap(icount,igrd) = ib
860                  ENDIF
861               END DO
862            END DO
863         END DO   ! igrd
864
865      END DO   ! ib_bdy
866
867      ! Initialize array indicating communications in bdy
868      ! -------------------------------------------------
869
870      ! Allocate array indicating if a send    instruction is needed in bdy treatment
871      ALLOCATE( nbondi_bdy(nb_bdy) )
872      ALLOCATE( nbondj_bdy(nb_bdy) )
873      nbondi_bdy(:)=2
874      nbondj_bdy(:)=2
875      ! Allocate array indicating if a receive instruction is needed in bdy treatment
876      ALLOCATE( nbondi_bdy_b(nb_bdy))
877      ALLOCATE( nbondj_bdy_b(nb_bdy))
878      nbondi_bdy_b(:)=2
879      nbondj_bdy_b(:)=2
880
881      DO ib_bdy = 1, nb_bdy
882         ! default : no send
883         llsend_ea = .false.   
884         llsend_we = .false.
885         llsend_so = .false.
886         llsend_no = .false.
887         ! default : no receive
888         llrecv_ea = .false.
889         llrecv_we = .false.
890         llrecv_so = .false.
891         llrecv_no = .false.
892         DO igrd = 1, jpbgrd
893            DO ib = 1, idx_bdy(ib_bdy)%nblenrim(igrd)   ! only the rim triggers communications, see bdy routines
894               ii = idx_bdy(ib_bdy)%nbi(ib,igrd)
895               ij = idx_bdy(ib_bdy)%nbj(ib,igrd)
896               !
897               ! check if point has to be sent     to   a neighbour
898               ! E neighbour and on the inner right side
899               IF( ii == nlci-1 .and. (nbondi == 0 .or. nbondi == -1) )   llsend_ea = .true.
900               ! W neighbour and on the inner left  side
901               IF( ii == 2      .and. (nbondi == 0 .or. nbondi ==  1) )   llsend_we = .true.
902               ! N neighbour and on the inner up   side
903               IF( ij == nlcj-1 .and. (nbondj == 0 .or. nbondj == -1) )   llsend_no = .true.
904               ! S neighbour and on the inner down side
905               IF( ij == 2      .and. (nbondj == 0 .or. nbondj ==  1) )   llsend_so = .true.
906               !
907               ! check if point has to be received from a neighbour
908               ! E neighbour and on the outter right side
909               IF( ii == nlci .and. (nbondi == 0 .or. nbondi == -1) )   llrecv_ea = .true.
910               ! W neighbour and on the outter left  side
911               IF( ii == 1    .and. (nbondi == 0 .or. nbondi ==  1) )   llrecv_we = .true.
912               ! N neighbour and on the outter up   side
913               IF( ij == nlcj .and. (nbondj == 0 .or. nbondj == -1) )   llrecv_no = .true.
914               ! S neighbour and on the outter down side
915               IF( ij == 1    .and. (nbondj == 0 .or. nbondj ==  1) )   llrecv_so = .true.
916               !
917            END DO
918         END DO  ! igrd
919
920         ! definition of the i- and j- direction local boundaries arrays used for sending the boundaries
921         IF(           llsend_ea .and.       llsend_we ) THEN   ;   nbondi_bdy(ib_bdy) =  0
922         ELSEIF(       llsend_ea .and. .not. llsend_we ) THEN   ;   nbondi_bdy(ib_bdy) = -1
923         ELSEIF( .not. llsend_ea .and.       llsend_we ) THEN   ;   nbondi_bdy(ib_bdy) =  1
924         ENDIF
925         IF(           llsend_no .and.       llsend_so ) THEN   ;   nbondj_bdy(ib_bdy) =  0
926         ELSEIF(       llsend_no .and. .not. llsend_so ) THEN   ;   nbondj_bdy(ib_bdy) = -1
927         ELSEIF( .not. llsend_no .and.       llsend_so ) THEN   ;   nbondj_bdy(ib_bdy) =  1
928         ENDIF
929
930         ! definition of the i- and j- direction local boundaries arrays used for receiving the boundaries
931         IF(           llrecv_ea .and.       llrecv_we ) THEN   ;   nbondi_bdy_b(ib_bdy) =  0
932         ELSEIF(       llrecv_ea .and. .not. llrecv_we ) THEN   ;   nbondi_bdy_b(ib_bdy) = -1
933         ELSEIF( .not. llrecv_ea .and.       llrecv_we ) THEN   ;   nbondi_bdy_b(ib_bdy) =  1
934         ENDIF
935         IF(           llrecv_no .and.       llrecv_so ) THEN   ;   nbondj_bdy_b(ib_bdy) =  0
936         ELSEIF(       llrecv_no .and. .not. llrecv_so ) THEN   ;   nbondj_bdy_b(ib_bdy) = -1
937         ELSEIF( .not. llrecv_no .and.       llrecv_so ) THEN   ;   nbondj_bdy_b(ib_bdy) =  1
938         ENDIF
939
940         ! Compute rim weights for FRS scheme
941         ! ----------------------------------
942         DO igrd = 1, jpbgrd
943            DO ib = 1, idx_bdy(ib_bdy)%nblen(igrd)
944               nbr => idx_bdy(ib_bdy)%nbr(ib,igrd)
945               idx_bdy(ib_bdy)%nbw(ib,igrd) = 1.- TANH( REAL( nbr - 1 ) *0.5 )      ! tanh formulation
946!               idx_bdy(ib_bdy)%nbw(ib,igrd) = (REAL(nn_rimwidth(ib_bdy)+1-nbr)/REAL(nn_rimwidth(ib_bdy)))**2.  ! quadratic
947!               idx_bdy(ib_bdy)%nbw(ib,igrd) =  REAL(nn_rimwidth(ib_bdy)+1-nbr)/REAL(nn_rimwidth(ib_bdy))       ! linear
948            END DO
949         END DO 
950
951         ! Compute damping coefficients
952         ! ----------------------------
953         DO igrd = 1, jpbgrd
954            DO ib = 1, idx_bdy(ib_bdy)%nblen(igrd)
955               nbr => idx_bdy(ib_bdy)%nbr(ib,igrd)
956               idx_bdy(ib_bdy)%nbd(ib,igrd) = 1. / ( rn_time_dmp(ib_bdy) * rday ) & 
957               & *(REAL(nn_rimwidth(ib_bdy)+1-nbr)/REAL(nn_rimwidth(ib_bdy)))**2.   ! quadratic
958               idx_bdy(ib_bdy)%nbdout(ib,igrd) = 1. / ( rn_time_dmp_out(ib_bdy) * rday ) & 
959               & *(REAL(nn_rimwidth(ib_bdy)+1-nbr)/REAL(nn_rimwidth(ib_bdy)))**2.   ! quadratic
960            END DO
961         END DO
962
963      END DO ! ib_bdy
964
965      ! ------------------------------------------------------
966      ! Initialise masks and find normal/tangential directions
967      ! ------------------------------------------------------
968
969      ! Read global 2D mask at T-points: bdytmask
970      ! -----------------------------------------
971      ! bdytmask = 1  on the computational domain AND on open boundaries
972      !          = 0  elsewhere   
973 
974      bdytmask(:,:) = ssmask(:,:)
975
976      ! Derive mask on U and V grid from mask on T grid
977
978      bdyumask(:,:) = 0._wp
979      bdyvmask(:,:) = 0._wp
980      DO ij = 1, jpjm1
981         DO ii = 1, jpim1
982            bdyumask(ii,ij) = bdytmask(ii,ij) * bdytmask(ii+1, ij )
983            bdyvmask(ii,ij) = bdytmask(ii,ij) * bdytmask(ii  ,ij+1) 
984         END DO
985      END DO
986      CALL lbc_lnk_multi( 'bdyini', bdyumask, 'U', 1. , bdyvmask, 'V', 1. )   ! Lateral boundary cond.
987
988      ! bdy masks are now set to zero on boundary points:
989      !
990      igrd = 1
991      DO ib_bdy = 1, nb_bdy
992        DO ib = 1, idx_bdy(ib_bdy)%nblenrim(igrd)     
993          bdytmask(idx_bdy(ib_bdy)%nbi(ib,igrd), idx_bdy(ib_bdy)%nbj(ib,igrd)) = 0._wp
994        END DO
995      END DO
996      !
997      igrd = 2
998      DO ib_bdy = 1, nb_bdy
999        DO ib = 1, idx_bdy(ib_bdy)%nblenrim(igrd)
1000          bdyumask(idx_bdy(ib_bdy)%nbi(ib,igrd), idx_bdy(ib_bdy)%nbj(ib,igrd)) = 0._wp
1001        END DO
1002      END DO
1003      !
1004      igrd = 3
1005      DO ib_bdy = 1, nb_bdy
1006        DO ib = 1, idx_bdy(ib_bdy)%nblenrim(igrd)
1007          bdyvmask(idx_bdy(ib_bdy)%nbi(ib,igrd), idx_bdy(ib_bdy)%nbj(ib,igrd)) = 0._wp
1008        END DO
1009      END DO
1010
1011      ! For the flagu/flagv calculation below we require a version of fmask without
1012      ! the land boundary condition (shlat) included:
1013      zfmask(:,:) = 0
1014      DO ij = 2, jpjm1
1015         DO ii = 2, jpim1
1016            zfmask(ii,ij) = tmask(ii,ij  ,1) * tmask(ii+1,ij  ,1)   &
1017           &              * tmask(ii,ij+1,1) * tmask(ii+1,ij+1,1)
1018         END DO     
1019      END DO
1020
1021      ! Lateral boundary conditions
1022      CALL lbc_lnk( 'bdyini', zfmask, 'F', 1. ) 
1023      CALL lbc_lnk_multi( 'bdyini', bdyumask, 'U', 1. , bdyvmask, 'V', 1., bdytmask, 'T', 1. )
1024      DO ib_bdy = 1, nb_bdy       ! Indices and directions of rim velocity components
1025         icount = 0 
1026
1027         ! Calculate relationship of U direction to the local orientation of the boundary
1028         ! flagu = -1 : u component is normal to the dynamical boundary and its direction is outward
1029         ! flagu =  0 : u is tangential
1030         ! flagu =  1 : u is normal to the boundary and is direction is inward
1031 
1032         DO igrd = 1, jpbgrd 
1033            SELECT CASE( igrd )
1034               CASE( 1 )   ;   pmask => umask   (:,:,1)   ;   i_offset = 0
1035               CASE( 2 )   ;   pmask => bdytmask(:,:)     ;   i_offset = 1
1036               CASE( 3 )   ;   pmask => zfmask  (:,:)     ;   i_offset = 0
1037            END SELECT
1038            icount = 0
1039            ztmp(:,:) = 0._wp
1040            DO ib = 1, idx_bdy(ib_bdy)%nblenrim(igrd) 
1041               ii = idx_bdy(ib_bdy)%nbi(ib,igrd)
1042               ij = idx_bdy(ib_bdy)%nbj(ib,igrd)
1043               IF( ii == 1 .OR. ii == jpi .OR. ij == 1 .OR. ij == jpj )  CYCLE
1044               zefl = pmask(ii+i_offset-1,ij)
1045               zwfl = pmask(ii+i_offset  ,ij)
1046               ! This error check only works if you are using the bdyXmask arrays
1047               IF( i_offset == 1 .and. zefl + zwfl == 2. ) THEN
1048                  icount = icount + 1
1049                  IF(lwp) WRITE(numout,*) 'Problem with igrd = ',igrd,' at (global) nbi, nbj : ',mig(ii),mjg(ij)
1050               ELSE
1051                  ztmp(ii,ij) = -zefl + zwfl
1052               ENDIF
1053            END DO
1054            IF( icount /= 0 ) THEN
1055               WRITE(ctmp1,*) ' E R R O R : Some ',cgrid(igrd),' grid points,',   &
1056                  ' are not boundary points (flagu calculation). Check nbi, nbj, indices for boundary set ',ib_bdy
1057               WRITE(ctmp2,*) ' ========== '
1058               CALL ctl_stop( ' ', ctmp1, ctmp2, ' ' )
1059            ENDIF
1060            CALL lbc_lnk( 'bdyini', ztmp, 'T', 1. ) 
1061            DO ib = 1, idx_bdy(ib_bdy)%nblenrim(igrd) 
1062               ii = idx_bdy(ib_bdy)%nbi(ib,igrd)
1063               ij = idx_bdy(ib_bdy)%nbj(ib,igrd)
1064               idx_bdy(ib_bdy)%flagu(ib,igrd) = ztmp(ii,ij)
1065            END DO
1066         END DO
1067
1068         ! Calculate relationship of V direction to the local orientation of the boundary
1069         ! flagv = -1 : v component is normal to the dynamical boundary but its direction is outward
1070         ! flagv =  0 : v is tangential
1071         ! flagv =  1 : v is normal to the boundary and is direction is inward
1072
1073         DO igrd = 1, jpbgrd 
1074            SELECT CASE( igrd )
1075               CASE( 1 )   ;   pmask => vmask (:,:,1)   ;   j_offset = 0
1076               CASE( 2 )   ;   pmask => zfmask(:,:)     ;   j_offset = 0
1077               CASE( 3 )   ;   pmask => bdytmask        ;   j_offset = 1
1078            END SELECT
1079            icount = 0
1080            ztmp(:,:) = 0._wp
1081            DO ib = 1, idx_bdy(ib_bdy)%nblenrim(igrd) 
1082               ii = idx_bdy(ib_bdy)%nbi(ib,igrd)
1083               ij = idx_bdy(ib_bdy)%nbj(ib,igrd)
1084               IF( ii == 1 .OR. ii == jpi .OR. ij == 1 .OR. ij == jpj )  CYCLE
1085               znfl = pmask(ii,ij+j_offset-1)
1086               zsfl = pmask(ii,ij+j_offset  )
1087               ! This error check only works if you are using the bdyXmask arrays
1088               IF( j_offset == 1 .and. znfl + zsfl == 2. ) THEN
1089                  IF(lwp) WRITE(numout,*) 'Problem with igrd = ',igrd,' at (global) nbi, nbj : ',mig(ii),mjg(ij)
1090                  icount = icount + 1
1091               ELSE
1092                  ztmp(ii,ij) = -znfl + zsfl
1093               END IF
1094            END DO
1095            IF( icount /= 0 ) THEN
1096               WRITE(ctmp1,*) ' E R R O R : Some ',cgrid(igrd),' grid points,',   &
1097                  ' are not boundary points (flagv calculation). Check nbi, nbj, indices for boundary set ',ib_bdy
1098               WRITE(ctmp2,*) ' ========== '
1099               CALL ctl_stop( ' ', ctmp1, ctmp2, ' ' )
1100            ENDIF
1101            CALL lbc_lnk( 'bdyini', ztmp, 'T', 1. ) 
1102            DO ib = 1, idx_bdy(ib_bdy)%nblenrim(igrd) 
1103               ii = idx_bdy(ib_bdy)%nbi(ib,igrd)
1104               ij = idx_bdy(ib_bdy)%nbj(ib,igrd)
1105               idx_bdy(ib_bdy)%flagv(ib,igrd) = ztmp(ii,ij)
1106            END DO
1107         END DO
1108         !
1109      END DO
1110     
1111      ! detect corner interior and its orientation index 1 to 4  depending on the orientation
1112      ! detect corner exterior and its orientation index 5 to 8  depending on the orientation
1113      ! detect geometries with 3 neighbours        index 9 to 12 depending on the orientation
1114      ! else                                       index 0
1115      DO ib_bdy = 1, nb_bdy
1116         DO igrd = 1, jpbgrd
1117            SELECT CASE( igrd )
1118               CASE( 1 )   ;   pmask => bdytmask 
1119               CASE( 2 )   ;   pmask => bdyumask 
1120               CASE( 3 )   ;   pmask => bdyvmask 
1121            END SELECT
1122            ztmp(:,:) = 0._wp
1123            DO ib = 1, idx_bdy(ib_bdy)%nblenrim(igrd)
1124               ii = idx_bdy(ib_bdy)%nbi(ib,igrd)
1125               ij = idx_bdy(ib_bdy)%nbj(ib,igrd)
1126               IF( ii == 1 .OR. ii == jpi .OR. ij == 1 .OR. ij == jpj )   CYCLE
1127               llnobdy = pmask(ii  ,ij+1) == 1. 
1128               llsobdy = pmask(ii  ,ij-1) == 1. 
1129               lleabdy = pmask(ii+1,ij  ) == 1. 
1130               llwebdy = pmask(ii-1,ij  ) == 1. 
1131               inbdy  = COUNT( (/ llnobdy, llsobdy, lleabdy, llwebdy /) )
1132               IF( inbdy == 0 )   THEN   ! no neighbours -> interior of a corner
1133                  !               !              !     _____     !     _____     
1134                  !  1 |   o      !  2  o   |    !  3 | x        !  4     x |   
1135                  !    |_x_ _     !    _ _x_|    !    |   o      !      o   |
1136                  IF( pmask(ii+1,ij+1) == 1. )   ztmp(ii,ij) = 1
1137                  IF( pmask(ii-1,ij+1) == 1. )   ztmp(ii,ij) = 2
1138                  IF( pmask(ii+1,ij-1) == 1. )   ztmp(ii,ij) = 3
1139                  IF( pmask(ii-1,ij-1) == 1. )   ztmp(ii,ij) = 4
1140               END IF
1141               IF( inbdy == 1 )   THEN   ! middle of linear bdy
1142                  ztmp(ii,ij) = 0   ! regular treatment with flags
1143               END IF
1144               IF( inbdy == 2 )   THEN   ! exterior of a corner
1145                  !        o      !        o      !    _____|       !       |_____ 
1146                  !  5 ____x o    !  6   o x___   ! 7      x o      !  8   o x     
1147                  !         |     !       |       !        o        !        o
1148                  IF( llnobdy .AND. lleabdy )   ztmp(ii,ij) = 5
1149                  IF( llnobdy .AND. llwebdy )   ztmp(ii,ij) = 6
1150                  IF( llsobdy .AND. lleabdy )   ztmp(ii,ij) = 7
1151                  IF( llsobdy .AND. llwebdy )   ztmp(ii,ij) = 8
1152               END IF
1153               IF( inbdy == 3 )   THEN   ! 3 neighbours __   __
1154                  !    |_  o      !        o  _|  !       |_|     !       o         
1155                  !  9  _| x o    ! 10   o x |_   ! 11   o x o    ! 12  o x o       
1156                  !    |   o      !        o   |  !        o      !    __|¨|__   
1157                  IF( llnobdy .AND. lleabdy .AND. llsobdy )   ztmp(ii,ij) = 9
1158                  IF( llnobdy .AND. llwebdy .AND. llsobdy )   ztmp(ii,ij) = 10
1159                  IF( llwebdy .AND. llsobdy .AND. lleabdy )   ztmp(ii,ij) = 11
1160                  IF( llwebdy .AND. llnobdy .AND. lleabdy )   ztmp(ii,ij) = 12
1161               END IF
1162               IF( inbdy == 4 )   THEN
1163                  WRITE(ctmp1,*) ' E R R O R : Problem with  ',cgrid(igrd) ,' grid points,',   &
1164                       '  some points on boundary set ', ib_bdy, ' have 4 neighbours'
1165                  WRITE(ctmp2,*) ' ========== '
1166                  CALL ctl_stop( ' ', ctmp1, ctmp2, ' ' )
1167               END IF
1168            END DO
1169            CALL lbc_lnk( 'bdyini', ztmp, 'T', 1. )
1170            DO ib = 1, idx_bdy(ib_bdy)%nblenrim(igrd) 
1171               ii = idx_bdy(ib_bdy)%nbi(ib,igrd)
1172               ij = idx_bdy(ib_bdy)%nbj(ib,igrd)
1173               idx_bdy(ib_bdy)%ntreat(ib,igrd) = ztmp(ii,ij)
1174            END DO
1175         END DO
1176      END DO
1177
1178
1179      ALLOCATE( lsend_bdy(nb_bdy,jpbgrd,4), lrecv_bdy(nb_bdy,jpbgrd,4) )
1180      lsend_bdy(:,:,:) = .false.
1181      lrecv_bdy(:,:,:) = .false. 
1182      !
1183      ! Check which boundaries might need communication
1184      DO igrd = 1, jpbgrd
1185         DO ib_bdy = 1, nb_bdy
1186            IF     ( nbondi_bdy  (ib_bdy) ==  0 )   THEN
1187               lsend_bdy(ib_bdy,igrd,1) = .true.
1188               lsend_bdy(ib_bdy,igrd,2) = .true.
1189            ELSE IF( nbondi_bdy  (ib_bdy) ==  1 )   THEN
1190               lsend_bdy(ib_bdy,igrd,1) = .true.
1191            ELSE IF( nbondi_bdy  (ib_bdy) == -1 )   THEN
1192               lsend_bdy(ib_bdy,igrd,2) = .true.
1193            END IF
1194            IF     ( nbondi_bdy_b(ib_bdy) ==  0 )   THEN
1195               lrecv_bdy(ib_bdy,igrd,1) = .true.
1196               lrecv_bdy(ib_bdy,igrd,2) = .true.
1197            ELSE IF( nbondi_bdy_b(ib_bdy) ==  1 )   THEN
1198               lrecv_bdy(ib_bdy,igrd,1) = .true.
1199            ELSE IF( nbondi_bdy_b(ib_bdy) == -1 )   THEN
1200               lrecv_bdy(ib_bdy,igrd,2) = .true.
1201            END IF
1202            IF(      nbondj_bdy  (ib_bdy) ==  0 )   THEN
1203               lsend_bdy(ib_bdy,igrd,3) = .true.
1204               lsend_bdy(ib_bdy,igrd,4) = .true.
1205            ELSE IF( nbondj_bdy  (ib_bdy) ==  1 )   THEN
1206               lsend_bdy(ib_bdy,igrd,3) = .true.
1207            ELSE IF( nbondj_bdy  (ib_bdy) == -1 )   THEN
1208               lsend_bdy(ib_bdy,igrd,4) = .true.
1209            END IF
1210            IF(      nbondj_bdy_b(ib_bdy) ==  0 )   THEN
1211               lrecv_bdy(ib_bdy,igrd,3) = .true.
1212               lrecv_bdy(ib_bdy,igrd,4) = .true.
1213            ELSE IF( nbondj_bdy_b(ib_bdy) ==  1 )   THEN
1214               lrecv_bdy(ib_bdy,igrd,3) = .true.
1215            ELSE IF( nbondj_bdy_b(ib_bdy) == -1 )   THEN
1216               lrecv_bdy(ib_bdy,igrd,4) = .true.
1217            END IF
1218         END DO
1219      END DO
1220      !
1221      ! Tidy up
1222      !--------
1223      IF( nb_bdy>0 )   DEALLOCATE( nbidta, nbjdta, nbrdta )
1224      !
1225   END SUBROUTINE bdy_segs
1226
1227
1228   SUBROUTINE bdy_ctl_seg
1229      !!----------------------------------------------------------------------
1230      !!                 ***  ROUTINE bdy_ctl_seg  ***
1231      !!
1232      !! ** Purpose :   Check straight open boundary segments location
1233      !!
1234      !! ** Method  :   - Look for open boundary corners
1235      !!                - Check that segments start or end on land
1236      !!----------------------------------------------------------------------
1237      INTEGER  ::   ib, ib1, ib2, ji ,jj, itest 
1238      INTEGER, DIMENSION(jp_nseg,2) :: icorne, icornw, icornn, icorns 
1239      REAL(wp), DIMENSION(2) ::   ztestmask
1240      !!----------------------------------------------------------------------
1241      !
1242      IF (lwp) WRITE(numout,*) ' '
1243      IF (lwp) WRITE(numout,*) 'bdy_ctl_seg: Check analytical segments'
1244      IF (lwp) WRITE(numout,*) '~~~~~~~~~~~~'
1245      !
1246      IF(lwp) WRITE(numout,*) 'Number of east  segments     : ', nbdysege
1247      IF(lwp) WRITE(numout,*) 'Number of west  segments     : ', nbdysegw
1248      IF(lwp) WRITE(numout,*) 'Number of north segments     : ', nbdysegn
1249      IF(lwp) WRITE(numout,*) 'Number of south segments     : ', nbdysegs
1250      ! 1. Check bounds
1251      !----------------
1252      DO ib = 1, nbdysegn
1253         IF (lwp) WRITE(numout,*) '**check north seg bounds pckg: ', npckgn(ib)
1254         IF ((jpjnob(ib).ge.jpjglo-1).or.& 
1255            &(jpjnob(ib).le.1))        CALL ctl_stop( 'nbdyind out of domain' )
1256         IF (jpindt(ib).ge.jpinft(ib)) CALL ctl_stop( 'Bdy start index is greater than end index' )
1257         IF (jpindt(ib).le.1     )     CALL ctl_stop( 'Start index out of domain' )
1258         IF (jpinft(ib).ge.jpiglo)     CALL ctl_stop( 'End index out of domain' )
1259      END DO
1260      !
1261      DO ib = 1, nbdysegs
1262         IF (lwp) WRITE(numout,*) '**check south seg bounds pckg: ', npckgs(ib)
1263         IF ((jpjsob(ib).ge.jpjglo-1).or.& 
1264            &(jpjsob(ib).le.1))        CALL ctl_stop( 'nbdyind out of domain' )
1265         IF (jpisdt(ib).ge.jpisft(ib)) CALL ctl_stop( 'Bdy start index is greater than end index' )
1266         IF (jpisdt(ib).le.1     )     CALL ctl_stop( 'Start index out of domain' )
1267         IF (jpisft(ib).ge.jpiglo)     CALL ctl_stop( 'End index out of domain' )
1268      END DO
1269      !
1270      DO ib = 1, nbdysege
1271         IF (lwp) WRITE(numout,*) '**check east  seg bounds pckg: ', npckge(ib)
1272         IF ((jpieob(ib).ge.jpiglo-1).or.& 
1273            &(jpieob(ib).le.1))        CALL ctl_stop( 'nbdyind out of domain' )
1274         IF (jpjedt(ib).ge.jpjeft(ib)) CALL ctl_stop( 'Bdy start index is greater than end index' )
1275         IF (jpjedt(ib).le.1     )     CALL ctl_stop( 'Start index out of domain' )
1276         IF (jpjeft(ib).ge.jpjglo)     CALL ctl_stop( 'End index out of domain' )
1277      END DO
1278      !
1279      DO ib = 1, nbdysegw
1280         IF (lwp) WRITE(numout,*) '**check west  seg bounds pckg: ', npckgw(ib)
1281         IF ((jpiwob(ib).ge.jpiglo-1).or.& 
1282            &(jpiwob(ib).le.1))        CALL ctl_stop( 'nbdyind out of domain' )
1283         IF (jpjwdt(ib).ge.jpjwft(ib)) CALL ctl_stop( 'Bdy start index is greater than end index' )
1284         IF (jpjwdt(ib).le.1     )     CALL ctl_stop( 'Start index out of domain' )
1285         IF (jpjwft(ib).ge.jpjglo)     CALL ctl_stop( 'End index out of domain' )
1286      ENDDO
1287      !
1288      !     
1289      ! 2. Look for segment crossings
1290      !------------------------------
1291      IF (lwp) WRITE(numout,*) '**Look for segments corners  :'
1292      !
1293      itest = 0 ! corner number
1294      !
1295      ! flag to detect if start or end of open boundary belongs to a corner
1296      ! if not (=0), it must be on land.
1297      ! if a corner is detected, save bdy package number for further tests
1298      icorne(:,:)=0. ; icornw(:,:)=0. ; icornn(:,:)=0. ; icorns(:,:)=0.
1299      ! South/West crossings
1300      IF ((nbdysegw > 0).AND.(nbdysegs > 0)) THEN
1301         DO ib1 = 1, nbdysegw       
1302            DO ib2 = 1, nbdysegs
1303               IF (( jpisdt(ib2)<=jpiwob(ib1)).AND. &
1304                &  ( jpisft(ib2)>=jpiwob(ib1)).AND. &
1305                &  ( jpjwdt(ib1)<=jpjsob(ib2)).AND. &
1306                &  ( jpjwft(ib1)>=jpjsob(ib2))) THEN
1307                  IF ((jpjwdt(ib1)==jpjsob(ib2)).AND.(jpisdt(ib2)==jpiwob(ib1))) THEN 
1308                     ! We have a possible South-West corner                     
1309!                     WRITE(numout,*) ' Found a South-West corner at (i,j): ', jpisdt(ib2), jpjwdt(ib1)
1310!                     WRITE(numout,*) ' between segments: ', npckgw(ib1), npckgs(ib2)
1311                     icornw(ib1,1) = npckgs(ib2)
1312                     icorns(ib2,1) = npckgw(ib1)
1313                  ELSEIF ((jpisft(ib2)==jpiwob(ib1)).AND.(jpjwft(ib1)==jpjsob(ib2))) THEN
1314                     WRITE(ctmp1,*) ' E R R O R : Found an acute open boundary corner at point (i,j)= ', &
1315                        &                                     jpisft(ib2), jpjwft(ib1)
1316                     WRITE(ctmp2,*) ' ==========  Not allowed yet'
1317                     WRITE(ctmp3,*) '             Crossing problem with West segment: ',npckgw(ib1), & 
1318                        &                                        ' and South segment: ',npckgs(ib2)
1319                     CALL ctl_stop( ' ', ctmp1, ctmp2, ctmp3, ' ' )
1320                  ELSE
1321                     WRITE(ctmp1,*) ' E R R O R : Check South and West Open boundary indices'
1322                     WRITE(ctmp2,*) ' ==========  Crossing problem with West segment: ',npckgw(ib1) , &
1323                        &                                         ' and South segment: ',npckgs(ib2)
1324                     CALL ctl_stop( ' ', ctmp1, ctmp2, ' ' )
1325                  END IF
1326               END IF
1327            END DO
1328         END DO
1329      END IF
1330      !
1331      ! South/East crossings
1332      IF ((nbdysege > 0).AND.(nbdysegs > 0)) THEN
1333         DO ib1 = 1, nbdysege
1334            DO ib2 = 1, nbdysegs
1335               IF (( jpisdt(ib2)<=jpieob(ib1)+1).AND. &
1336                &  ( jpisft(ib2)>=jpieob(ib1)+1).AND. &
1337                &  ( jpjedt(ib1)<=jpjsob(ib2)  ).AND. &
1338                &  ( jpjeft(ib1)>=jpjsob(ib2)  )) THEN
1339                  IF ((jpjedt(ib1)==jpjsob(ib2)).AND.(jpisft(ib2)==jpieob(ib1)+1)) THEN
1340                     ! We have a possible South-East corner
1341!                     WRITE(numout,*) ' Found a South-East corner at (i,j): ', jpisft(ib2), jpjedt(ib1)
1342!                     WRITE(numout,*) ' between segments: ', npckge(ib1), npckgs(ib2)
1343                     icorne(ib1,1) = npckgs(ib2)
1344                     icorns(ib2,2) = npckge(ib1)
1345                  ELSEIF ((jpjeft(ib1)==jpjsob(ib2)).AND.(jpisdt(ib2)==jpieob(ib1)+1)) THEN
1346                     WRITE(ctmp1,*) ' E R R O R : Found an acute open boundary corner at point (i,j)= ', &
1347                        &                                     jpisdt(ib2), jpjeft(ib1)
1348                     WRITE(ctmp2,*) ' ==========  Not allowed yet'
1349                     WRITE(ctmp3,*) '             Crossing problem with East segment: ',npckge(ib1), &
1350                        &                                        ' and South segment: ',npckgs(ib2)
1351                     CALL ctl_stop( ' ', ctmp1, ctmp2, ctmp3, ' ' )
1352                  ELSE
1353                     WRITE(ctmp1,*) ' E R R O R : Check South and East Open boundary indices'
1354                     WRITE(ctmp2,*) ' ==========  Crossing problem with East segment: ',npckge(ib1), &
1355                     &                                           ' and South segment: ',npckgs(ib2)
1356                     CALL ctl_stop( ' ', ctmp1, ctmp2, ' ' )
1357                  END IF
1358               END IF
1359            END DO
1360         END DO
1361      END IF
1362      !
1363      ! North/West crossings
1364      IF ((nbdysegn > 0).AND.(nbdysegw > 0)) THEN
1365         DO ib1 = 1, nbdysegw       
1366            DO ib2 = 1, nbdysegn
1367               IF (( jpindt(ib2)<=jpiwob(ib1)  ).AND. &
1368                &  ( jpinft(ib2)>=jpiwob(ib1)  ).AND. &
1369                &  ( jpjwdt(ib1)<=jpjnob(ib2)+1).AND. &
1370                &  ( jpjwft(ib1)>=jpjnob(ib2)+1)) THEN
1371                  IF ((jpjwft(ib1)==jpjnob(ib2)+1).AND.(jpindt(ib2)==jpiwob(ib1))) THEN
1372                     ! We have a possible North-West corner
1373!                     WRITE(numout,*) ' Found a North-West corner at (i,j): ', jpindt(ib2), jpjwft(ib1)
1374!                     WRITE(numout,*) ' between segments: ', npckgw(ib1), npckgn(ib2)
1375                     icornw(ib1,2) = npckgn(ib2)
1376                     icornn(ib2,1) = npckgw(ib1)
1377                  ELSEIF ((jpjwdt(ib1)==jpjnob(ib2)+1).AND.(jpinft(ib2)==jpiwob(ib1))) THEN
1378                     WRITE(ctmp1,*) ' E R R O R : Found an acute open boundary corner at point (i,j)= ', &
1379                     &                                     jpinft(ib2), jpjwdt(ib1)
1380                     WRITE(ctmp2,*) ' ==========  Not allowed yet'
1381                     WRITE(ctmp3,*) '             Crossing problem with West segment: ',npckgw(ib1), &
1382                     &                                                    ' and North segment: ',npckgn(ib2)
1383                     CALL ctl_stop( ' ', ctmp1, ctmp2, ctmp3, ' ' )
1384                  ELSE
1385                     WRITE(ctmp1,*) ' E R R O R : Check North and West Open boundary indices'
1386                     WRITE(ctmp2,*) ' ==========  Crossing problem with West segment: ',npckgw(ib1), &
1387                     &                                                    ' and North segment: ',npckgn(ib2)
1388                     CALL ctl_stop( ' ', ctmp1, ctmp2, ' ' )
1389                  END IF
1390               END IF
1391            END DO
1392         END DO
1393      END IF
1394      !
1395      ! North/East crossings
1396      IF ((nbdysegn > 0).AND.(nbdysege > 0)) THEN
1397         DO ib1 = 1, nbdysege       
1398            DO ib2 = 1, nbdysegn
1399               IF (( jpindt(ib2)<=jpieob(ib1)+1).AND. &
1400                &  ( jpinft(ib2)>=jpieob(ib1)+1).AND. &
1401                &  ( jpjedt(ib1)<=jpjnob(ib2)+1).AND. &
1402                &  ( jpjeft(ib1)>=jpjnob(ib2)+1)) THEN
1403                  IF ((jpjeft(ib1)==jpjnob(ib2)+1).AND.(jpinft(ib2)==jpieob(ib1)+1)) THEN
1404                     ! We have a possible North-East corner
1405!                     WRITE(numout,*) ' Found a North-East corner at (i,j): ', jpinft(ib2), jpjeft(ib1)
1406!                     WRITE(numout,*) ' between segments: ', npckge(ib1), npckgn(ib2)
1407                     icorne(ib1,2) = npckgn(ib2)
1408                     icornn(ib2,2) = npckge(ib1)
1409                  ELSEIF ((jpjedt(ib1)==jpjnob(ib2)+1).AND.(jpindt(ib2)==jpieob(ib1)+1)) THEN
1410                     WRITE(ctmp1,*) ' E R R O R : Found an acute open boundary corner at point (i,j)= ', &
1411                     &                                     jpindt(ib2), jpjedt(ib1)
1412                     WRITE(ctmp2,*) ' ==========  Not allowed yet'
1413                     WRITE(ctmp3,*) '             Crossing problem with East segment: ',npckge(ib1), &
1414                     &                                           ' and North segment: ',npckgn(ib2)
1415                     CALL ctl_stop( ' ', ctmp1, ctmp2, ctmp3, ' ' )
1416                  ELSE
1417                     WRITE(ctmp1,*) ' E R R O R : Check North and East Open boundary indices'
1418                     WRITE(ctmp2,*) ' ==========  Crossing problem with East segment: ',npckge(ib1), &
1419                     &                                           ' and North segment: ',npckgn(ib2)
1420                     CALL ctl_stop( ' ', ctmp1, ctmp2, ' ' )
1421                  END IF
1422               END IF
1423            END DO
1424         END DO
1425      END IF
1426      !
1427      ! 3. Check if segment extremities are on land
1428      !--------------------------------------------
1429      !
1430      ! West segments
1431      DO ib = 1, nbdysegw
1432         ! get mask at boundary extremities:
1433         ztestmask(1:2)=0.
1434         DO ji = 1, jpi
1435            DO jj = 1, jpj             
1436              IF (((ji + nimpp - 1) == jpiwob(ib)).AND. & 
1437               &  ((jj + njmpp - 1) == jpjwdt(ib))) ztestmask(1)=tmask(ji,jj,1)
1438              IF (((ji + nimpp - 1) == jpiwob(ib)).AND. & 
1439               &  ((jj + njmpp - 1) == jpjwft(ib))) ztestmask(2)=tmask(ji,jj,1) 
1440            END DO
1441         END DO
1442         CALL mpp_sum( 'bdyini', ztestmask, 2 )   ! sum over the global domain
1443
1444         IF (ztestmask(1)==1) THEN
1445            IF (icornw(ib,1)==0) THEN
1446               WRITE(ctmp1,*) ' E R R O R : Open boundary segment ', npckgw(ib)
1447               WRITE(ctmp2,*) ' ==========  does not start on land or on a corner'                                                 
1448               CALL ctl_stop( ' ', ctmp1, ctmp2, ' ' )
1449            ELSE
1450               ! This is a corner
1451               IF(lwp) WRITE(numout,*) 'Found a South-West corner at (i,j): ', jpiwob(ib), jpjwdt(ib)
1452               CALL bdy_ctl_corn(npckgw(ib), icornw(ib,1))
1453               itest=itest+1
1454            ENDIF
1455         ENDIF
1456         IF (ztestmask(2)==1) THEN
1457            IF (icornw(ib,2)==0) THEN
1458               WRITE(ctmp1,*) ' E R R O R : Open boundary segment ', npckgw(ib)
1459               WRITE(ctmp2,*) ' ==========  does not end on land or on a corner'                                                 
1460               CALL ctl_stop( ' ', ctmp1, ctmp2, ' ' )
1461            ELSE
1462               ! This is a corner
1463               IF(lwp) WRITE(numout,*) 'Found a North-West corner at (i,j): ', jpiwob(ib), jpjwft(ib)
1464               CALL bdy_ctl_corn(npckgw(ib), icornw(ib,2))
1465               itest=itest+1
1466            ENDIF
1467         ENDIF
1468      END DO
1469      !
1470      ! East segments
1471      DO ib = 1, nbdysege
1472         ! get mask at boundary extremities:
1473         ztestmask(1:2)=0.
1474         DO ji = 1, jpi
1475            DO jj = 1, jpj             
1476              IF (((ji + nimpp - 1) == jpieob(ib)+1).AND. & 
1477               &  ((jj + njmpp - 1) == jpjedt(ib))) ztestmask(1)=tmask(ji,jj,1)
1478              IF (((ji + nimpp - 1) == jpieob(ib)+1).AND. & 
1479               &  ((jj + njmpp - 1) == jpjeft(ib))) ztestmask(2)=tmask(ji,jj,1) 
1480            END DO
1481         END DO
1482         CALL mpp_sum( 'bdyini', ztestmask, 2 )   ! sum over the global domain
1483
1484         IF (ztestmask(1)==1) THEN
1485            IF (icorne(ib,1)==0) THEN
1486               WRITE(ctmp1,*) ' E R R O R : Open boundary segment ', npckge(ib)
1487               WRITE(ctmp2,*) ' ==========  does not start on land or on a corner'                                                 
1488               CALL ctl_stop( ' ', ctmp1, ctmp2, ' ' )
1489            ELSE
1490               ! This is a corner
1491               IF(lwp) WRITE(numout,*) 'Found a South-East corner at (i,j): ', jpieob(ib)+1, jpjedt(ib)
1492               CALL bdy_ctl_corn(npckge(ib), icorne(ib,1))
1493               itest=itest+1
1494            ENDIF
1495         ENDIF
1496         IF (ztestmask(2)==1) THEN
1497            IF (icorne(ib,2)==0) THEN
1498               WRITE(ctmp1,*) ' E R R O R : Open boundary segment ', npckge(ib)
1499               WRITE(ctmp2,*) ' ==========  does not end on land or on a corner'                                                 
1500               CALL ctl_stop( ' ', ctmp1, ctmp2, ' ' )
1501            ELSE
1502               ! This is a corner
1503               IF(lwp) WRITE(numout,*) 'Found a North-East corner at (i,j): ', jpieob(ib)+1, jpjeft(ib)
1504               CALL bdy_ctl_corn(npckge(ib), icorne(ib,2))
1505               itest=itest+1
1506            ENDIF
1507         ENDIF
1508      END DO
1509      !
1510      ! South segments
1511      DO ib = 1, nbdysegs
1512         ! get mask at boundary extremities:
1513         ztestmask(1:2)=0.
1514         DO ji = 1, jpi
1515            DO jj = 1, jpj             
1516              IF (((jj + njmpp - 1) == jpjsob(ib)).AND. & 
1517               &  ((ji + nimpp - 1) == jpisdt(ib))) ztestmask(1)=tmask(ji,jj,1)
1518              IF (((jj + njmpp - 1) == jpjsob(ib)).AND. & 
1519               &  ((ji + nimpp - 1) == jpisft(ib))) ztestmask(2)=tmask(ji,jj,1) 
1520            END DO
1521         END DO
1522         CALL mpp_sum( 'bdyini', ztestmask, 2 )   ! sum over the global domain
1523
1524         IF ((ztestmask(1)==1).AND.(icorns(ib,1)==0)) THEN
1525            WRITE(ctmp1,*) ' E R R O R : Open boundary segment ', npckgs(ib)
1526            WRITE(ctmp2,*) ' ==========  does not start on land or on a corner'                                                 
1527            CALL ctl_stop( ' ', ctmp1, ctmp2, ' ' )
1528         ENDIF
1529         IF ((ztestmask(2)==1).AND.(icorns(ib,2)==0)) THEN
1530            WRITE(ctmp1,*) ' E R R O R : Open boundary segment ', npckgs(ib)
1531            WRITE(ctmp2,*) ' ==========  does not end on land or on a corner'                                                 
1532            CALL ctl_stop( ' ', ctmp1, ctmp2, ' ' )
1533         ENDIF
1534      END DO
1535      !
1536      ! North segments
1537      DO ib = 1, nbdysegn
1538         ! get mask at boundary extremities:
1539         ztestmask(1:2)=0.
1540         DO ji = 1, jpi
1541            DO jj = 1, jpj             
1542              IF (((jj + njmpp - 1) == jpjnob(ib)+1).AND. & 
1543               &  ((ji + nimpp - 1) == jpindt(ib))) ztestmask(1)=tmask(ji,jj,1)
1544              IF (((jj + njmpp - 1) == jpjnob(ib)+1).AND. & 
1545               &  ((ji + nimpp - 1) == jpinft(ib))) ztestmask(2)=tmask(ji,jj,1) 
1546            END DO
1547         END DO
1548         CALL mpp_sum( 'bdyini', ztestmask, 2 )   ! sum over the global domain
1549
1550         IF ((ztestmask(1)==1).AND.(icornn(ib,1)==0)) THEN
1551            WRITE(ctmp1,*) ' E R R O R : Open boundary segment ', npckgn(ib)
1552            WRITE(ctmp2,*) ' ==========  does not start on land'                                                 
1553            CALL ctl_stop( ' ', ctmp1, ctmp2, ' ' )
1554         ENDIF
1555         IF ((ztestmask(2)==1).AND.(icornn(ib,2)==0)) THEN
1556            WRITE(ctmp1,*) ' E R R O R : Open boundary segment ', npckgn(ib)
1557            WRITE(ctmp2,*) ' ==========  does not end on land'                                                 
1558            CALL ctl_stop( ' ', ctmp1, ctmp2, ' ' )
1559         ENDIF
1560      END DO
1561      !
1562      IF ((itest==0).AND.(lwp)) WRITE(numout,*) 'NO open boundary corner found'
1563      !
1564      ! Other tests TBD:
1565      ! segments completly on land
1566      ! optimized open boundary array length according to landmask
1567      ! Nudging layers that overlap with interior domain
1568      !
1569   END SUBROUTINE bdy_ctl_seg
1570
1571
1572   SUBROUTINE bdy_ctl_corn( ib1, ib2 )
1573      !!----------------------------------------------------------------------
1574      !!                 ***  ROUTINE bdy_ctl_corn  ***
1575      !!
1576      !! ** Purpose :   Check numerical schemes consistency between
1577      !!                segments having a common corner
1578      !!
1579      !! ** Method  :   
1580      !!----------------------------------------------------------------------
1581      INTEGER, INTENT(in)  ::   ib1, ib2
1582      INTEGER :: itest
1583      !!----------------------------------------------------------------------
1584      itest = 0
1585
1586      IF( cn_dyn2d(ib1) /= cn_dyn2d(ib2) )   itest = itest + 1
1587      IF( cn_dyn3d(ib1) /= cn_dyn3d(ib2) )   itest = itest + 1
1588      IF( cn_tra  (ib1) /= cn_tra  (ib2) )   itest = itest + 1
1589      !
1590      IF( nn_dyn2d_dta(ib1) /= nn_dyn2d_dta(ib2) )   itest = itest + 1
1591      IF( nn_dyn3d_dta(ib1) /= nn_dyn3d_dta(ib2) )   itest = itest + 1
1592      IF( nn_tra_dta  (ib1) /= nn_tra_dta  (ib2) )   itest = itest + 1
1593      !
1594      IF( nn_rimwidth(ib1) /= nn_rimwidth(ib2) )   itest = itest + 1   
1595      !
1596      IF( itest>0 ) THEN
1597         WRITE(ctmp1,*) ' E R R O R : Segments ', ib1, 'and ', ib2
1598         WRITE(ctmp2,*) ' ==========  have different open bdy schemes'                                                 
1599         CALL ctl_stop( ' ', ctmp1, ctmp2, ' ' )
1600      ENDIF
1601      !
1602   END SUBROUTINE bdy_ctl_corn
1603
1604   !!=================================================================================
1605END MODULE bdyini
Note: See TracBrowser for help on using the repository browser.