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/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/BDY – NEMO

source: branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/BDY/bdyini.F90 @ 3211

Last change on this file since 3211 was 3211, checked in by spickles2, 12 years ago

Stephen Pickles, 11 Dec 2011

Commit to bring the rest of the DCSE NEMO development branch
in line with the latest development version. This includes
array index re-ordering of all OPA_SRC/.

  • Property svn:keywords set to Id
File size: 20.7 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   !!----------------------------------------------------------------------
13#if defined key_bdy
14   !!----------------------------------------------------------------------
15   !!   'key_bdy'                     Unstructured Open Boundary Conditions
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 obc_par         ! ocean open boundary conditions
22   USE bdy_oce         ! unstructured open boundary conditions
23   USE bdydta, ONLY: bdy_dta_alloc ! open boundary data
24   USE bdytides        ! tides at open boundaries initialization (tide_init routine)
25   USE in_out_manager  ! I/O units
26   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
27   USE lib_mpp         ! for mpp_sum 
28   USE iom             ! I/O
29
30   IMPLICIT NONE
31   PRIVATE
32
33   PUBLIC   bdy_init   ! routine called by opa.F90
34
35   !! * Control permutation of array indices
36#  include "oce_ftrans.h90"
37#  include "dom_oce_ftrans.h90"
38
39   !!----------------------------------------------------------------------
40   !! NEMO/OPA 4.0 , NEMO Consortium (2011)
41   !! $Id$
42   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
43   !!----------------------------------------------------------------------
44CONTAINS
45   
46   SUBROUTINE bdy_init
47      !!----------------------------------------------------------------------
48      !!                 ***  ROUTINE bdy_init  ***
49      !!         
50      !! ** Purpose :   Initialization of the dynamics and tracer fields with
51      !!              unstructured open boundaries.
52      !!
53      !! ** Method  :   Read initialization arrays (mask, indices) to identify
54      !!              an unstructured open boundary
55      !!
56      !! ** Input   :  bdy_init.nc, input file for unstructured open boundaries
57      !!----------------------------------------------------------------------     
58      INTEGER  ::   ii, ij, ik, igrd, ib, ir   ! dummy loop indices
59      INTEGER  ::   icount, icountr, ib_len, ibr_max   ! local integers
60      INTEGER  ::   iw, ie, is, in, inum, id_dummy     !   -       -
61      INTEGER  ::   igrd_start, igrd_end               !   -       -
62      REAL(wp) ::   zefl, zwfl, znfl, zsfl              ! local scalars
63      INTEGER, DIMENSION (2)             ::   kdimsz
64      INTEGER, DIMENSION(jpbdta, jpbgrd) ::   nbidta, nbjdta   ! Index arrays: i and j indices of bdy dta
65      INTEGER, DIMENSION(jpbdta, jpbgrd) ::   nbrdta           ! Discrete distance from rim points
66      REAL(wp), DIMENSION(jpidta,jpjdta) ::   zmask            ! global domain mask
67      REAL(wp), DIMENSION(jpbdta,1)      ::   zdta             ! temporary array
68      CHARACTER(LEN=80),DIMENSION(6)     ::   clfile
69      !!
70      NAMELIST/nambdy/cn_mask, cn_dta_frs_T, cn_dta_frs_U, cn_dta_frs_V,   &
71         &            cn_dta_fla_T, cn_dta_fla_U, cn_dta_fla_V,            &
72         &            ln_tides, ln_clim, ln_vol, ln_mask,                  &
73         &            ln_dyn_fla, ln_dyn_frs, ln_tra_frs,ln_ice_frs,       &
74         &            nn_dtactl, nn_rimwidth, nn_volctl
75      !!----------------------------------------------------------------------
76
77      IF(lwp) WRITE(numout,*)
78      IF(lwp) WRITE(numout,*) 'bdy_init : initialization of unstructured open boundaries'
79      IF(lwp) WRITE(numout,*) '~~~~~~~~'
80      !
81      !                                      ! allocate bdy_oce arrays
82      IF( bdy_oce_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'bdy_init : unable to allocate oce arrays' )
83      IF( bdy_dta_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'bdy_init : unable to allocate dta arrays' )
84
85      IF( jperio /= 0 )   CALL ctl_stop( 'Cyclic or symmetric,',   &
86         &                               ' and unstructured open boundary condition are not compatible' )
87
88      IF( lk_obc      )   CALL ctl_stop( 'Straight open boundaries,',   &
89         &                               ' and unstructured open boundaries are not compatible' )
90
91      ! ---------------------------
92      REWIND( numnam )                    ! Read namelist parameters
93      READ  ( numnam, nambdy )
94
95      !                                   ! control prints
96      IF(lwp) WRITE(numout,*) '         nambdy'
97
98      !                                         ! check type of data used (nn_dtactl value)
99      IF(lwp) WRITE(numout,*) 'nn_dtactl =', nn_dtactl     
100      IF(lwp) WRITE(numout,*)
101      SELECT CASE( nn_dtactl )                   !
102      CASE( 0 )      ;   IF(lwp) WRITE(numout,*) '      initial state used for bdy data'       
103      CASE( 1 )      ;   IF(lwp) WRITE(numout,*) '      boundary data taken from file'
104      CASE DEFAULT   ;   CALL ctl_stop( 'nn_dtactl must be 0 or 1' )
105      END SELECT
106
107      IF(lwp) WRITE(numout,*)
108      IF(lwp) WRITE(numout,*) 'Boundary rim width for the FRS nn_rimwidth = ', nn_rimwidth
109
110      IF(lwp) WRITE(numout,*)
111      IF(lwp) WRITE(numout,*) '      nn_volctl = ', nn_volctl
112
113      IF( ln_vol ) THEN                     ! check volume conservation (nn_volctl value)
114         SELECT CASE ( nn_volctl )
115         CASE( 1 )      ;   IF(lwp) WRITE(numout,*) '      The total volume will be constant'
116         CASE( 0 )      ;   IF(lwp) WRITE(numout,*) '      The total volume will vary according to the surface E-P flux'
117         CASE DEFAULT   ;   CALL ctl_stop( 'nn_volctl must be 0 or 1' )
118         END SELECT
119         IF(lwp) WRITE(numout,*)
120      ELSE
121         IF(lwp) WRITE(numout,*) 'No volume correction with unstructured open boundaries'
122         IF(lwp) WRITE(numout,*)
123      ENDIF
124
125      IF( ln_tides ) THEN
126        IF(lwp) WRITE(numout,*) 'Tidal harmonic forcing at unstructured open boundaries'
127        IF(lwp) WRITE(numout,*)
128      ENDIF
129
130      IF( ln_dyn_fla ) THEN
131        IF(lwp) WRITE(numout,*) 'Flather condition on U, V at unstructured open boundaries'
132        IF(lwp) WRITE(numout,*)
133      ENDIF
134
135      IF( ln_dyn_frs ) THEN
136        IF(lwp) WRITE(numout,*) 'FRS condition on U and V at unstructured open boundaries'
137        IF(lwp) WRITE(numout,*)
138      ENDIF
139
140      IF( ln_tra_frs ) THEN
141        IF(lwp) WRITE(numout,*) 'FRS condition on T & S fields at unstructured open boundaries'
142        IF(lwp) WRITE(numout,*)
143      ENDIF
144
145      IF( ln_ice_frs ) THEN
146        IF(lwp) WRITE(numout,*) 'FRS condition on ice fields at unstructured open boundaries'
147        IF(lwp) WRITE(numout,*)
148      ENDIF
149
150      IF( ln_tides )   CALL tide_init      ! Read tides namelist
151
152
153      ! Read arrays defining unstructured open boundaries
154      ! -------------------------------------------------
155
156      ! Read global 2D mask at T-points: bdytmask
157      ! *****************************************
158      ! bdytmask = 1  on the computational domain AND on open boundaries
159      !          = 0  elsewhere   
160 
161      IF( cp_cfg == "eel" .AND. jp_cfg == 5 ) THEN          ! EEL configuration at 5km resolution
162         zmask(         :                ,:) = 0.e0
163         zmask(jpizoom+1:jpizoom+jpiglo-2,:) = 1.e0         
164      ELSE IF( ln_mask ) THEN
165         CALL iom_open( cn_mask, inum )
166         CALL iom_get ( inum, jpdom_data, 'bdy_msk', zmask(:,:) )
167         CALL iom_close( inum )
168      ELSE
169         zmask(:,:) = 1.e0
170      ENDIF
171
172      DO ij = 1, nlcj      ! Save mask over local domain     
173         DO ii = 1, nlci
174            bdytmask(ii,ij) = zmask( mig(ii), mjg(ij) )
175         END DO
176      END DO
177
178      ! Derive mask on U and V grid from mask on T grid
179      bdyumask(:,:) = 0.e0
180      bdyvmask(:,:) = 0.e0
181      DO ij=1, jpjm1
182         DO ii=1, jpim1
183            bdyumask(ii,ij)=bdytmask(ii,ij)*bdytmask(ii+1, ij )
184            bdyvmask(ii,ij)=bdytmask(ii,ij)*bdytmask(ii  ,ij+1) 
185         END DO
186      END DO
187      CALL lbc_lnk( bdyumask(:,:), 'U', 1. )   ;   CALL lbc_lnk( bdyvmask(:,:), 'V', 1. )      ! Lateral boundary cond.
188
189
190      ! Read discrete distance and mapping indices
191      ! ******************************************
192      nbidta(:,:) = 0.e0
193      nbjdta(:,:) = 0.e0
194      nbrdta(:,:) = 0.e0
195
196      IF( cp_cfg == "eel" .AND. jp_cfg == 5 ) THEN
197         icount = 0
198         DO ir = 1, nn_rimwidth                  ! Define west boundary (from ii=2 to ii=1+nn_rimwidth):
199            DO ij = 3, jpjglo-2
200               icount = icount + 1
201               nbidta(icount,:) = ir + 1 + (jpizoom-1)
202               nbjdta(icount,:) = ij     + (jpjzoom-1) 
203               nbrdta(icount,:) = ir
204            END DO
205         END DO
206         !
207         DO ir = 1, nn_rimwidth                  ! Define east boundary (from ii=jpiglo-1 to ii=jpiglo-nn_rimwidth):
208            DO ij=3,jpjglo-2
209               icount = icount + 1
210               nbidta(icount,:) = jpiglo-ir + (jpizoom-1)
211               nbidta(icount,2) = jpiglo-ir-1 + (jpizoom-1) ! special case for u points
212               nbjdta(icount,:) = ij + (jpjzoom-1)
213               nbrdta(icount,:) = ir
214            END DO
215         END DO
216         !       
217      ELSE            ! Read indices and distances in unstructured boundary data files
218         !
219         IF( ln_tides ) THEN             ! Read tides input files for preference in case there are no bdydata files
220            clfile(4) = TRIM(filtide)//TRIM(tide_cpt(1))//'_grid_T.nc'
221            clfile(5) = TRIM(filtide)//TRIM(tide_cpt(1))//'_grid_U.nc'
222            clfile(6) = TRIM(filtide)//TRIM(tide_cpt(1))//'_grid_V.nc'
223         ENDIF
224         IF( ln_dyn_fla .AND. .NOT. ln_tides ) THEN
225            clfile(4) = cn_dta_fla_T
226            clfile(5) = cn_dta_fla_U
227            clfile(6) = cn_dta_fla_V
228         ENDIF
229
230         IF( ln_tra_frs ) THEN
231            clfile(1) = cn_dta_frs_T
232            IF( .NOT. ln_dyn_frs ) THEN
233               clfile(2) = cn_dta_frs_T     ! Dummy read re read T file for sake of 6 files
234               clfile(3) = cn_dta_frs_T     !
235            ENDIF
236         ENDIF         
237         IF( ln_dyn_frs ) THEN
238            IF( .NOT. ln_tra_frs )   clfile(1) = cn_dta_frs_U      ! Dummy Read
239            clfile(2) = cn_dta_frs_U
240            clfile(3) = cn_dta_frs_V 
241         ENDIF
242
243         !                                   ! how many files are we to read in?
244         IF(ln_tides .OR. ln_dyn_fla)   igrd_start = 4
245         !
246         IF(ln_tra_frs    ) THEN   ;   igrd_start = 1
247         ELSEIF(ln_dyn_frs) THEN   ;   igrd_start = 2
248         ENDIF
249         !
250         IF( ln_tra_frs   )   igrd_end = 1
251         !
252         IF(ln_dyn_fla .OR. ln_tides) THEN   ;   igrd_end = 6
253         ELSEIF( ln_dyn_frs             ) THEN   ;   igrd_end = 3
254         ENDIF
255
256         DO igrd = igrd_start, igrd_end
257            CALL iom_open( clfile(igrd), inum )
258            id_dummy = iom_varid( inum, 'nbidta', kdimsz=kdimsz ) 
259            IF(lwp) WRITE(numout,*) 'kdimsz : ',kdimsz
260            ib_len = kdimsz(1)
261            IF( ib_len > jpbdta)   CALL ctl_stop(  'Boundary data array in file too long.',                  &
262                &                                  'File :', TRIM(clfile(igrd)),'increase parameter jpbdta.' )
263
264            CALL iom_get( inum, jpdom_unknown, 'nbidta', zdta(1:ib_len,:) )
265            DO ii = 1,ib_len
266               nbidta(ii,igrd) = INT( zdta(ii,1) )
267            END DO
268            CALL iom_get( inum, jpdom_unknown, 'nbjdta', zdta(1:ib_len,:) )
269            DO ii = 1,ib_len
270               nbjdta(ii,igrd) = INT( zdta(ii,1) )
271            END DO
272            CALL iom_get( inum, jpdom_unknown, 'nbrdta', zdta(1:ib_len,:) )
273            DO ii = 1,ib_len
274               nbrdta(ii,igrd) = INT( zdta(ii,1) )
275            END DO
276            CALL iom_close( inum )
277
278            IF( igrd < 4) THEN            ! Check that rimwidth in file is big enough for Frs case(barotropic is one):
279               ibr_max = MAXVAL( nbrdta(:,igrd) )
280               IF(lwp) WRITE(numout,*)
281               IF(lwp) WRITE(numout,*) ' Maximum rimwidth in file is ', ibr_max
282               IF(lwp) WRITE(numout,*) ' nn_rimwidth from namelist is ', nn_rimwidth
283               IF (ibr_max < nn_rimwidth)   CALL ctl_stop( 'nn_rimwidth is larger than maximum rimwidth in file' )
284            ENDIF !Check igrd < 4
285            !
286         END DO
287         !
288      ENDIF 
289
290      ! Dispatch mapping indices and discrete distances on each processor
291      ! *****************************************************************
292     
293      iw = mig(1) + 1            ! if monotasking and no zoom, iw=2
294      ie = mig(1) + nlci-1 - 1   ! if monotasking and no zoom, ie=jpim1
295      is = mjg(1) + 1            ! if monotasking and no zoom, is=2
296      in = mjg(1) + nlcj-1 - 1   ! if monotasking and no zoom, in=jpjm1
297
298      DO igrd = igrd_start, igrd_end
299         icount  = 0
300         icountr = 0
301         nblen   (igrd) = 0
302         nblenrim(igrd) = 0
303         nblendta(igrd) = 0
304         DO ir=1, nn_rimwidth
305            DO ib = 1, jpbdta
306               ! check if point is in local domain and equals ir
307               IF(  nbidta(ib,igrd) >= iw .AND. nbidta(ib,igrd) <= ie .AND.   &
308                  & nbjdta(ib,igrd) >= is .AND. nbjdta(ib,igrd) <= in .AND.   &
309                  & nbrdta(ib,igrd) == ir  ) THEN
310                  !
311                  icount = icount  + 1
312                  !
313                  IF( ir == 1 )   icountr = icountr+1
314                  IF (icount > jpbdim) THEN
315                     IF(lwp) WRITE(numout,*) 'bdy_ini: jpbdim too small'
316                     nstop = nstop + 1
317                  ELSE
318                     nbi(icount, igrd)  = nbidta(ib,igrd)- mig(1)+1
319                     nbj(icount, igrd)  = nbjdta(ib,igrd)- mjg(1)+1
320                     nbr(icount, igrd)  = nbrdta(ib,igrd)
321                     nbmap(icount,igrd) = ib
322                  ENDIF           
323               ENDIF
324            END DO
325         END DO
326         nblenrim(igrd) = icountr !: length of rim boundary data on each proc
327         nblen   (igrd) = icount  !: length of boundary data on each proc       
328      END DO 
329
330      ! Compute rim weights
331      ! -------------------
332      DO igrd = igrd_start, igrd_end
333         DO ib = 1, nblen(igrd)
334            nbw(ib,igrd) = 1.- TANH( FLOAT( nbr(ib,igrd) - 1 ) *0.5 )                     ! tanh formulation
335!           nbw(ib,igrd) = (FLOAT(nn_rimwidth+1-nbr(ib,igrd))/FLOAT(nn_rimwidth))**2      ! quadratic
336!           nbw(ib,igrd) =  FLOAT(nn_rimwidth+1-nbr(ib,igrd))/FLOAT(nn_rimwidth)          ! linear
337         END DO
338      END DO 
339   
340      ! Mask corrections
341      ! ----------------
342#if defined key_z_first
343      DO ij = 1, jpj
344         DO ii = 1, jpi
345            DO ik = 1, jpkm1
346#else
347      DO ik = 1, jpkm1
348         DO ij = 1, jpj
349            DO ii = 1, jpi
350#endif
351               tmask(ii,ij,ik) = tmask(ii,ij,ik) * bdytmask(ii,ij)
352               umask(ii,ij,ik) = umask(ii,ij,ik) * bdyumask(ii,ij)
353               vmask(ii,ij,ik) = vmask(ii,ij,ik) * bdyvmask(ii,ij)
354               bmask(ii,ij)    = bmask(ii,ij)    * bdytmask(ii,ij)
355            END DO     
356         END DO
357      END DO
358
359#if defined key_z_first
360      DO ij = 2, jpjm1
361         DO ii = 2, jpim1
362            DO ik = 1, jpkm1
363#else
364      DO ik = 1, jpkm1
365         DO ij = 2, jpjm1
366            DO ii = 2, jpim1
367#endif
368               fmask(ii,ij,ik) = fmask(ii,ij,ik) * bdytmask(ii,ij  ) * bdytmask(ii+1,ij  )   &
369                  &                              * bdytmask(ii,ij+1) * bdytmask(ii+1,ij+1)
370            END DO     
371         END DO
372      END DO
373
374#if defined key_z_first
375      bdytmask(:,:) = tmask(:,:,1)
376      tmask_i (:,:) = bdytmask(:,:) * tmask_i(:,:)             
377#else
378      tmask_i (:,:) = tmask(:,:,1) * tmask_i(:,:)             
379      bdytmask(:,:) = tmask(:,:,1)
380#endif
381
382      ! bdy masks and bmask are now set to zero on boundary points:
383      igrd = 1       ! In the free surface case, bmask is at T-points
384      DO ib = 1, nblenrim(igrd)     
385        bmask(nbi(ib,igrd), nbj(ib,igrd)) = 0.e0
386      END DO
387      !
388      igrd = 1
389      DO ib = 1, nblenrim(igrd)     
390        bdytmask(nbi(ib,igrd), nbj(ib,igrd)) = 0.e0
391      END DO
392      !
393      igrd = 2
394      DO ib = 1, nblenrim(igrd)
395        bdyumask(nbi(ib,igrd), nbj(ib,igrd)) = 0.e0
396      END DO
397      !
398      igrd = 3
399      DO ib = 1, nblenrim(igrd)
400        bdyvmask(nbi(ib,igrd), nbj(ib,igrd)) = 0.e0
401      END DO
402
403      ! Lateral boundary conditions
404      CALL lbc_lnk( fmask        , 'F', 1. )   ;   CALL lbc_lnk( bdytmask(:,:), 'T', 1. )
405      CALL lbc_lnk( bdyumask(:,:), 'U', 1. )   ;   CALL lbc_lnk( bdyvmask(:,:), 'V', 1. )
406
407      IF( ln_vol .OR. ln_dyn_fla ) THEN      ! Indices and directions of rim velocity components
408         !
409         !flagu = -1 : u component is normal to the dynamical boundary but its direction is outward
410         !flagu =  0 : u is tangential
411         !flagu =  1 : u is normal to the boundary and is direction is inward
412         icount = 0 
413         flagu(:) = 0.e0
414 
415         igrd = 2      ! u-component
416         DO ib = 1, nblenrim(igrd) 
417            zefl=bdytmask(nbi(ib,igrd)  , nbj(ib,igrd))
418            zwfl=bdytmask(nbi(ib,igrd)+1, nbj(ib,igrd))
419            IF( zefl + zwfl ==2 ) THEN
420               icount = icount +1
421            ELSE
422               flagu(ib)=-zefl+zwfl
423            ENDIF
424         END DO
425
426         !flagv = -1 : u component is normal to the dynamical boundary but its direction is outward
427         !flagv =  0 : u is tangential
428         !flagv =  1 : u is normal to the boundary and is direction is inward
429         flagv(:) = 0.e0
430
431         igrd = 3      ! v-component
432         DO ib = 1, nblenrim(igrd) 
433            znfl = bdytmask(nbi(ib,igrd), nbj(ib,igrd))
434            zsfl = bdytmask(nbi(ib,igrd), nbj(ib,igrd)+1)
435            IF( znfl + zsfl ==2 ) THEN
436               icount = icount + 1
437            ELSE
438               flagv(ib) = -znfl + zsfl
439            END IF
440         END DO
441 
442         IF( icount /= 0 ) THEN
443            IF(lwp) WRITE(numout,*)
444            IF(lwp) WRITE(numout,*) ' E R R O R : Some data velocity points,',   &
445               ' are not boundary points. Check nbi, nbj, indices.'
446            IF(lwp) WRITE(numout,*) ' ========== '
447            IF(lwp) WRITE(numout,*)
448            nstop = nstop + 1
449         ENDIF
450   
451      ENDIF
452
453      ! Compute total lateral surface for volume correction:
454      ! ----------------------------------------------------
455      bdysurftot = 0.e0 
456      IF( ln_vol ) THEN 
457         igrd = 2      ! Lateral surface at U-points
458         DO ib = 1, nblenrim(igrd)
459            bdysurftot = bdysurftot + hu     (nbi(ib,igrd)  ,nbj(ib,igrd))                      &
460               &                    * e2u    (nbi(ib,igrd)  ,nbj(ib,igrd)) * ABS( flagu(ib) )   &
461               &                    * tmask_i(nbi(ib,igrd)  ,nbj(ib,igrd))                      &
462               &                    * tmask_i(nbi(ib,igrd)+1,nbj(ib,igrd))                   
463         END DO
464
465         igrd=3 ! Add lateral surface at V-points
466         DO ib = 1, nblenrim(igrd)
467            bdysurftot = bdysurftot + hv     (nbi(ib,igrd),nbj(ib,igrd)  )                      &
468               &                    * e1v    (nbi(ib,igrd),nbj(ib,igrd)  ) * ABS( flagv(ib) )   &
469               &                    * tmask_i(nbi(ib,igrd),nbj(ib,igrd)  )                      &
470               &                    * tmask_i(nbi(ib,igrd),nbj(ib,igrd)+1)
471         END DO
472         !
473         IF( lk_mpp )   CALL mpp_sum( bdysurftot )      ! sum over the global domain
474      END IF   
475
476      ! Initialise bdy data arrays
477      ! --------------------------
478      tbdy(:,:) = 0.e0
479      sbdy(:,:) = 0.e0
480      ubdy(:,:) = 0.e0
481      vbdy(:,:) = 0.e0
482      sshbdy(:) = 0.e0
483      ubtbdy(:) = 0.e0
484      vbtbdy(:) = 0.e0
485#if defined key_lim2
486      frld_bdy(:) = 0.e0
487      hicif_bdy(:) = 0.e0
488      hsnif_bdy(:) = 0.e0
489#endif
490
491      ! Read in tidal constituents and adjust for model start time
492      ! ----------------------------------------------------------
493      IF( ln_tides )   CALL tide_data
494      !
495   END SUBROUTINE bdy_init
496
497#else
498   !!---------------------------------------------------------------------------------
499   !!   Dummy module                                   NO unstructured open boundaries
500   !!---------------------------------------------------------------------------------
501CONTAINS
502   SUBROUTINE bdy_init      ! Dummy routine
503   END SUBROUTINE bdy_init
504#endif
505
506   !!=================================================================================
507END MODULE bdyini
Note: See TracBrowser for help on using the repository browser.