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.
obcini.F90 on Ticket #648 – Attachment – NEMO

Ticket #648: obcini.F90

File obcini.F90, 30.3 KB (added by molines, 14 years ago)

obcini.F90 with bug fix for bmask and nbobc

Line 
1 MODULE obcini
2   !!======================================================================
3   !!                       ***  MODULE  obcini  ***
4   !! OBC initial state :  Open boundary initial state
5   !!======================================================================
6   !! History :  8.0  !  97-07  (J.M. Molines, G. Madec)  Original code
7   !!   NEMO     1.0  !  02-11  (C. Talandier, A-M. Treguier) Free surface, F90
8   !!            2.0  !  05-11  (V. Garnier) Surface pressure gradient organization
9   !!----------------------------------------------------------------------
10#if defined key_obc
11   !!----------------------------------------------------------------------
12   !!   'key_obc'                                  Open Boundary Conditions
13   !!----------------------------------------------------------------------
14   !!   obc_init       : initialization for the open boundary condition
15   !!----------------------------------------------------------------------
16   USE oce             ! ocean dynamics and tracers variables
17   USE dom_oce         ! ocean space and time domain variables
18   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
19   USE phycst          ! physical constants
20   USE obc_oce         ! ocean open boundary conditions
21   USE lib_mpp         ! for mpp_sum
22   USE in_out_manager  ! I/O units
23   USE dynspg_oce      ! flag lk_dynspg_flt
24
25   IMPLICIT NONE
26   PRIVATE
27
28   PUBLIC   obc_init   ! routine called by opa.F90
29
30   !! * Substitutions
31#  include "obc_vectopt_loop_substitute.h90"
32   !!----------------------------------------------------------------------
33   !! NEMO/OPA 9.0 , LOCEAN-IPSL (2005)
34   !! $Id: obcini.F90 1633 2009-09-22 15:46:58Z rblod $
35   !! software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt)
36   !!----------------------------------------------------------------------
37
38CONTAINS
39   
40   SUBROUTINE obc_init
41      !!----------------------------------------------------------------------
42      !!                 ***  ROUTINE obc_init  ***
43      !!         
44      !! ** Purpose :   Initialization of the dynamics and tracer fields at
45      !!              the open boundaries.
46      !!
47      !! ** Method  :   initialization of open boundary variables
48      !!      (u, v) over 3 time step and 3 rows
49      !!      (t, s) over 2 time step and 2 rows
50      !!      if ln_rstart = .FALSE. : no restart, fields set to zero
51      !!      if ln_rstart = .TRUE.  : restart, fields are read in a file
52      !!      if rdpxxx = 0 then lfbc is set true for this boundary.
53      !!
54      !! ** Input   :   restart.obc file, restart file for open boundaries
55      !!----------------------------------------------------------------------
56      USE obcrst,   ONLY :   obc_rst_read   ! Make obc_rst_read routine available
57      !!
58      INTEGER  ::   ji, jj, istop , inumfbc
59      INTEGER, DIMENSION(4) ::   icorner
60      REAL(wp), DIMENSION(2) ::   ztestmask
61      !!
62      NAMELIST/namobc/ rn_dpein, rn_dpwin, rn_dpnin, rn_dpsin,       &
63         &             rn_dpeob, rn_dpwob, rn_dpnob, rn_dpsob,       &
64         &             rn_volemp, nn_obcdta, cn_obcdta, rn_volemp,   &
65         &             ln_obc_clim, ln_vol_cst, ln_obc_fla
66      !!----------------------------------------------------------------------
67
68      REWIND( numnam )              ! Namelist namobc : open boundaries
69      READ  ( numnam, namobc )
70
71      ! convert DOCTOR namelist name into the OLD names
72      nbobc    = 0
73      nobc_dta = nn_obcdta
74      cffile   = cn_obcdta
75      rdpein   = rn_dpein
76      rdpwin   = rn_dpwin
77      rdpsin   = rn_dpsin
78      rdpnin   = rn_dpnin
79      rdpeob   = rn_dpeob
80      rdpwob   = rn_dpwob
81      rdpsob   = rn_dpsob
82      rdpnob   = rn_dpnob
83      volemp   = rn_volemp
84     
85
86
87      ! By security we set rdpxin and rdpxob respectively to 1. and 15. if the corresponding OBC is not activated
88      IF( .NOT.lp_obc_east  ) THEN   ;   rdpein = 1.   ;   rdpeob = 15.   ;   END IF
89      IF( .NOT.lp_obc_west  ) THEN   ;   rdpwin = 1.   ;   rdpwob = 15.   ;   END IF
90      IF( .NOT.lp_obc_north ) THEN   ;   rdpnin = 1.   ;   rdpnob = 15.   ;   END IF
91      IF( .NOT.lp_obc_south ) THEN   ;   rdpsin = 1.   ;   rdpsob = 15.   ;   END IF
92
93      ! number of open boudaries and open boundary indicators
94      nbobc = 0
95      IF( lp_obc_east  )   nbobc = nbobc + 1
96      IF( lp_obc_west  )   nbobc = nbobc + 1
97      IF( lp_obc_north )   nbobc = nbobc + 1
98      IF( lp_obc_south )   nbobc = nbobc + 1
99
100      IF(lwp) WRITE(numout,*)
101      IF(lwp) WRITE(numout,*) 'obc_init : initialization of open boundaries'
102      IF(lwp) WRITE(numout,*) '~~~~~~~~'
103      IF(lwp) WRITE(numout,*) '   Number of open boundaries    nbobc = ', nbobc
104      IF(lwp) WRITE(numout,*)
105
106      ! control prints
107      IF(lwp) WRITE(numout,*) '   Namelist namobc'
108      IF(lwp) WRITE(numout,*) '      data in file (=1) or initial state used (=0)   nn_obcdta   = ', nn_obcdta
109      IF(lwp) WRITE(numout,*) '      climatology (true) or not                      ln_obc_clim = ', ln_obc_clim
110      IF(lwp) WRITE(numout,*) '      vol_cst (true) or not:                         ln_vol_cst  = ', ln_vol_cst
111      IF(lwp) WRITE(numout,*) ' '
112      IF(lwp) WRITE(numout,*) '   WARNING                                                  '
113      IF(lwp) WRITE(numout,*) '      Flather"s algorithm is applied with explicit free surface scheme                 '
114      IF(lwp) WRITE(numout,*) '      or with free surface time-splitting scheme                                       '
115      IF(lwp) WRITE(numout,*) '      Nor radiation neither relaxation is allowed with explicit free surface scheme:   '
116      IF(lwp) WRITE(numout,*) '      Radiation and/or relaxation is allowed with free surface time-splitting scheme '
117      IF(lwp) WRITE(numout,*) '      depending of the choice of rdpXin = rdpXob  = 0. for open boundaries             '
118      IF(lwp) WRITE(numout,*)
119      IF(lwp) WRITE(numout,*) '      For the filtered free surface case,                                              '
120      IF(lwp) WRITE(numout,*) '      radiation, relaxation or presciption of data can be applied                      '
121      IF(lwp) WRITE(numout,*)
122
123      IF( lwp.AND.lp_obc_east ) THEN
124         WRITE(numout,*) '      East open boundary :'
125         WRITE(numout,*) '         i index                    jpieob   = ', jpieob
126         WRITE(numout,*) '         damping time scale (days)  rn_dpeob = ', rn_dpeob
127         WRITE(numout,*) '         damping time scale (days)  rn_dpein = ', rn_dpein
128      ENDIF
129
130      IF( lwp.AND.lp_obc_west ) THEN
131         WRITE(numout,*) '      West open boundary :'
132         WRITE(numout,*) '         i index                    jpiwob   = ', jpiwob
133         WRITE(numout,*) '         damping time scale (days)  rn_dpwob = ', rn_dpwob
134         WRITE(numout,*) '         damping time scale (days)  rn_dpwin = ', rn_dpwin
135      ENDIF
136
137      IF( lwp.AND.lp_obc_north ) THEN
138         WRITE(numout,*) '      North open boundary :'
139         WRITE(numout,*) '         j index                    jpjnob   = ', jpjnob
140         WRITE(numout,*) '         damping time scale (days)  rn_dpnob = ', rn_dpnob
141         WRITE(numout,*) '         damping time scale (days)  rn_dpnin = ', rn_dpnin
142      ENDIF
143
144      IF( lwp.AND.lp_obc_south ) THEN
145         WRITE(numout,*) '      South open boundary :'
146         WRITE(numout,*) '         j index                    jpjsob   = ', jpjsob
147         WRITE(numout,*) '         damping time scale (days)  rn_dpsob = ', rn_dpsob
148         WRITE(numout,*) '         damping time scale (days)  rn_dpsin = ', rn_dpsin
149         WRITE(numout,*)
150      ENDIF
151
152      IF( nbobc >= 2 .AND. jperio /= 0 )   &
153         &   CALL ctl_stop( ' Cyclic or symmetric, and open boundary condition are not compatible' )
154
155      ! 1. Initialisation of constants
156      ! ------------------------------
157      ! ...                          convert rdp$ob in seconds
158      ! Fixed Bdy flag              inbound                outbound
159      lfbceast  = .FALSE.   ;   rdpein = rdpein * rday    ;   rdpeob = rdpeob * rday
160      lfbcwest  = .FALSE.   ;   rdpwin = rdpwin * rday    ;   rdpwob = rdpwob * rday
161      lfbcnorth = .FALSE.   ;   rdpnin = rdpnin * rday    ;   rdpnob = rdpnob * rday
162      lfbcsouth = .FALSE.   ;   rdpsin = rdpsin * rday    ;   rdpsob = rdpsob * rday
163      inumfbc = 0
164      ! ... look for Fixed Boundaries (rdp = 0 )
165      ! ... When specified, lbcxxx flags are set to TRUE and rdpxxx are set to
166      ! ...  a small arbitrary value, (to avoid division by zero further on).
167      ! ...  rdpxxx is not used anymore.
168      IF( lp_obc_east )  THEN
169         IF( (rdpein+rdpeob) == 0 )  THEN
170            lfbceast = .TRUE.   ;   rdpein = 1e-3   ;   rdpeob = 1e-3
171            inumfbc = inumfbc+1
172         ELSEIF ( (rdpein*rdpeob) == 0 )  THEN
173            CALL ctl_stop( 'obc_init : rn_dpein & rn_dpeob must be both zero or non zero' )
174         END IF
175      END IF
176
177      IF( lp_obc_west )  THEN
178         IF( (rdpwin + rdpwob) == 0 )  THEN
179            lfbcwest = .TRUE.     ;     rdpwin = 1e-3     ;     rdpwob = 1e-3
180            inumfbc = inumfbc+1
181         ELSEIF ( (rdpwin*rdpwob) == 0 )  THEN
182            CALL ctl_stop( 'obc_init : rn_dpwin & rn_dpwob must be both zero or non zero' )
183         END IF
184      END IF
185      IF( lp_obc_north )  THEN
186         IF( (rdpnin + rdpnob) == 0 )  THEN
187            lfbcnorth = .TRUE.     ;     rdpnin = 1e-3     ;     rdpnob = 1e-3
188            inumfbc = inumfbc+1
189         ELSEIF ( (rdpnin*rdpnob) == 0 )  THEN
190            CALL ctl_stop( 'obc_init : rn_dpnin & rn_dpnob must be both zero or non zero' )
191         END IF
192      END IF
193      IF( lp_obc_south )  THEN
194         IF( (rdpsin + rdpsob) == 0 )  THEN
195            lfbcsouth = .TRUE.   ;   rdpsin = 1e-3   ;   rdpsob = 1e-3
196            inumfbc = inumfbc+1
197         ELSEIF ( (rdpsin*rdpsob) == 0 )  THEN
198            CALL ctl_stop( 'obc_init : rn_dpsin & rn_dpsob must be both zero or non zero' )
199         END IF
200      END IF
201
202      ! 2.  Clever mpp indices for loops on the open boundaries.
203      !     The loops will be performed only on the processors
204      !     that contain a given open boundary.
205      ! --------------------------------------------------------
206
207      IF( lp_obc_east ) THEN
208         ! ...   mpp initialization
209         nie0   = max( 1, min(jpieob   - nimpp+1, jpi     ) )
210         nie1   = max( 0, min(jpieob   - nimpp+1, jpi - 1 ) )
211         nie0p1 = max( 1, min(jpieob+1 - nimpp+1, jpi     ) )
212         nie1p1 = max( 0, min(jpieob+1 - nimpp+1, jpi - 1 ) )
213         nie0m1 = max( 1, min(jpieob-1 - nimpp+1, jpi     ) )
214         nie1m1 = max( 0, min(jpieob-1 - nimpp+1, jpi - 1 ) )
215         nje0   = max( 2, min(jpjed    - njmpp+1, jpj     ) )
216         nje1   = max( 0, min(jpjef    - njmpp+1, jpj - 1 ) )
217         nje0p1 = max( 1, min(jpjedp1  - njmpp+1, jpj     ) )
218         nje0m1 = max( 1, min(jpjed    - njmpp+1, jpj     ) )
219         nje1m1 = max( 0, min(jpjefm1  - njmpp+1, jpj - 1 ) )
220         nje1m2 = max( 0, min(jpjefm1-1- njmpp+1, jpj - 1 ) )
221         IF(lwp) THEN
222            IF( lfbceast ) THEN
223               WRITE(numout,*)'     '
224               WRITE(numout,*)'         Specified East Open Boundary'
225            ELSE
226               WRITE(numout,*)'     '
227               WRITE(numout,*)'         Radiative East Open Boundary'
228            END IF
229         END IF
230      END IF
231
232      IF( lp_obc_west ) THEN
233         ! ...   mpp initialization
234         niw0   = max( 1, min(jpiwob   - nimpp+1, jpi     ) )
235         niw1   = max( 0, min(jpiwob   - nimpp+1, jpi - 1 ) )
236         niw0p1 = max( 1, min(jpiwob+1 - nimpp+1, jpi     ) )
237         niw1p1 = max( 0, min(jpiwob+1 - nimpp+1, jpi - 1 ) )
238         njw0   = max( 2, min(jpjwd    - njmpp+1, jpj     ) )
239         njw1   = max( 0, min(jpjwf    - njmpp+1, jpj - 1 ) )
240         njw0p1 = max( 1, min(jpjwdp1  - njmpp+1, jpj     ) )
241         njw0m1 = max( 1, min(jpjwd    - njmpp+1, jpj     ) )
242         njw1m1 = max( 0, min(jpjwfm1  - njmpp+1, jpj - 1 ) )
243         njw1m2 = max( 0, min(jpjwfm1-1- njmpp+1, jpj - 1 ) )
244         IF(lwp) THEN
245            IF( lfbcwest ) THEN
246               WRITE(numout,*)'     '
247               WRITE(numout,*)'         Specified West Open Boundary'
248            ELSE
249               WRITE(numout,*)'     '
250               WRITE(numout,*)'         Radiative West Open Boundary'
251            END IF
252         END IF
253      END IF
254 
255      IF( lp_obc_north ) THEN
256         ! ...   mpp initialization
257         nin0   = max( 2, min(jpind    - nimpp+1, jpi     ) )
258         nin1   = max( 0, min(jpinf    - nimpp+1, jpi - 1 ) )
259         nin0p1 = max( 1, min(jpindp1  - nimpp+1, jpi     ) )
260         nin0m1 = max( 1, min(jpind    - nimpp+1, jpi     ) )
261         nin1m1 = max( 0, min(jpinfm1  - nimpp+1, jpi - 1 ) )
262         nin1m2 = max( 0, min(jpinfm1-1- nimpp+1, jpi - 1 ) )
263         njn0   = max( 1, min(jpjnob   - njmpp+1, jpj     ) )
264         njn1   = max( 0, min(jpjnob   - njmpp+1, jpj - 1 ) )
265         njn0p1 = max( 1, min(jpjnob+1 - njmpp+1, jpj     ) )
266         njn1p1 = max( 0, min(jpjnob+1 - njmpp+1, jpj - 1 ) )
267         njn0m1 = max( 1, min(jpjnob-1 - njmpp+1, jpj     ) )
268         njn1m1 = max( 0, min(jpjnob-1 - njmpp+1, jpj - 1 ) )
269         IF(lwp) THEN
270            IF( lfbcnorth ) THEN
271               WRITE(numout,*)'     '
272               WRITE(numout,*)'         Specified North Open Boundary'
273            ELSE
274               WRITE(numout,*)'     '
275               WRITE(numout,*)'         Radiative North Open Boundary'
276            END IF
277         END IF
278      END IF
279
280      IF( lp_obc_south ) THEN
281         ! ...   mpp initialization
282         nis0   = max( 2, min(jpisd    - nimpp+1, jpi     ) )
283         nis1   = max( 0, min(jpisf    - nimpp+1, jpi - 1 ) )
284         nis0p1 = max( 1, min(jpisdp1  - nimpp+1, jpi     ) )
285         nis0m1 = max( 1, min(jpisd    - nimpp+1, jpi     ) )
286         nis1m1 = max( 0, min(jpisfm1  - nimpp+1, jpi - 1 ) )
287         nis1m2 = max( 0, min(jpisfm1-1- nimpp+1, jpi - 1 ) )
288         njs0   = max( 1, min(jpjsob   - njmpp+1, jpj     ) )
289         njs1   = max( 0, min(jpjsob   - njmpp+1, jpj - 1 ) )
290         njs0p1 = max( 1, min(jpjsob+1 - njmpp+1, jpj     ) )
291         njs1p1 = max( 0, min(jpjsob+1 - njmpp+1, jpj - 1 ) )
292         IF(lwp) THEN
293            IF( lfbcsouth ) THEN
294               WRITE(numout,*)'     '
295               WRITE(numout,*)'         Specified South Open Boundary'
296            ELSE
297               WRITE(numout,*)'     '
298               WRITE(numout,*)'         Radiative South Open Boundary'
299            END IF
300         END IF
301      END IF
302
303      ! 3. mask correction for OBCs
304      ! ---------------------------
305
306      IF( lp_obc_east ) THEN
307         !... (jpjed,jpjefm1),jpieob
308         bmask(nie0p1:nie1p1,nje0:nje1m1) = 0.e0
309
310         ! ... initilization to zero
311         uemsk(:,:) = 0.e0   ;   vemsk(:,:) = 0.e0   ;   temsk(:,:) = 0.e0
312
313         ! ... set 2D mask on East OBC,  Vopt
314         DO ji = fs_nie0, fs_nie1
315            DO jj = nje0, nje1
316               uemsk(jj,:) = umask(ji,  jj,:) * tmask_i(ji,jj)   * tmask_i(ji+1,jj)
317               vemsk(jj,:) = vmask(ji+1,jj,:) * tmask_i(ji+1,jj) 
318               temsk(jj,:) = tmask(ji+1,jj,:) * tmask_i(ji+1,jj) 
319            END DO
320         END DO
321
322      END IF
323
324      IF( lp_obc_west ) THEN
325         ! ... (jpjwd,jpjwfm1),jpiwob
326         bmask(niw0:niw1,njw0:njw1m1) = 0.e0
327
328         ! ... initilization to zero
329         uwmsk(:,:) = 0.e0   ;   vwmsk(:,:) = 0.e0   ;   twmsk(:,:) = 0.e0 
330
331         ! ... set 2D mask on West OBC,  Vopt
332         DO ji = fs_niw0, fs_niw1
333            DO jj = njw0, njw1
334               uwmsk(jj,:) = umask(ji,jj,:) * tmask_i(ji,jj)   * tmask_i(ji+1,jj)
335               vwmsk(jj,:) = vmask(ji,jj,:) * tmask_i(ji,jj) 
336               twmsk(jj,:) = tmask(ji,jj,:) * tmask_i(ji,jj)
337            END DO
338         END DO
339
340      END IF
341
342      IF( lp_obc_north ) THEN
343         ! ... jpjnob,(jpind,jpisfm1)
344         bmask(nin0:nin1m1,njn0p1:njn1p1) = 0.e0
345
346         ! ... initilization to zero
347         unmsk(:,:) = 0.e0   ;   vnmsk(:,:) = 0.e0   ;   tnmsk(:,:) = 0.e0
348
349         ! ... set 2D mask on North OBC,  Vopt
350         DO jj = fs_njn0, fs_njn1
351            DO ji = nin0, nin1
352               unmsk(ji,:) = umask(ji,jj+1,:) * tmask_i(ji,jj+1) 
353               vnmsk(ji,:) = vmask(ji,jj  ,:) * tmask_i(ji,jj)   * tmask_i(ji,jj+1)
354               tnmsk(ji,:) = tmask(ji,jj+1,:) * tmask_i(ji,jj+1)
355            END DO
356         END DO
357
358      END IF
359
360      IF( lp_obc_south ) THEN 
361         ! ... jpjsob,(jpisd,jpisfm1)
362         bmask(nis0:nis1m1,njs0:njs1) = 0.e0
363
364         ! ... initilization to zero
365         usmsk(:,:) = 0.e0   ;   vsmsk(:,:) = 0.e0   ;   tsmsk(:,:) = 0.e0
366
367         ! ... set 2D mask on South OBC,  Vopt
368         DO jj = fs_njs0, fs_njs1 
369            DO ji = nis0, nis1
370               usmsk(ji,:) = umask(ji,jj,:) * tmask_i(ji,jj) 
371               vsmsk(ji,:) = vmask(ji,jj,:) * tmask_i(ji,jj) * tmask_i(ji,jj+1)
372               tsmsk(ji,:) = tmask(ji,jj,:) * tmask_i(ji,jj)
373            END DO
374         END DO
375
376      END IF
377
378      IF ( ln_vol_cst .OR. lk_dynspg_flt ) THEN
379        ! ... Initialize obcumask and obcvmask for the Force filtering
380        !     boundary condition in dynspg_flt
381        obcumask(:,:) = umask(:,:,1)
382        obcvmask(:,:) = vmask(:,:,1)
383
384        ! ... Initialize obctmsk on overlap region and obcs. This mask
385        !     is used in obcvol.F90 to calculate cumulate flux E-P.
386        !     obc Tracer point are outside the domain ( U/V obc points) ==> masked by obctmsk
387        !     - no flux E-P on obcs and overlap region (jpreci = jprecj = 1)
388        obctmsk(:,:) = tmask_i(:,:)     
389
390        IF( lp_obc_east ) THEN
391           ! ... East obc Force filtering mask for the grad D
392           obcumask(nie0  :nie1  ,nje0p1:nje1m1) = 0.e0
393           obcvmask(nie0p1:nie1p1,nje0p1:nje1m1) = 0.e0
394           ! ... set to 0 on East OBC
395           obctmsk(nie0p1:nie1p1,nje0:nje1) = 0.e0
396        END IF
397 
398        IF( lp_obc_west ) THEN
399           ! ... West obc Force filtering mask for the grad D
400           obcumask(niw0:niw1,njw0:njw1) = 0.e0
401           obcvmask(niw0:niw1,njw0:njw1) = 0.e0
402           ! ... set to 0 on West OBC
403           obctmsk(niw0:niw1,njw0:njw1) = 0.e0
404        END IF
405 
406        IF( lp_obc_north ) THEN
407           ! ... North obc Force filtering mask for the grad D
408           obcumask(nin0p1:nin1m1,njn0p1:njn1p1) = 0.e0
409           obcvmask(nin0p1:nin1m1,njn0  :njn1  ) = 0.e0
410           ! ... set to 0 on North OBC
411           obctmsk(nin0:nin1,njn0p1:njn1p1) = 0.e0
412        END IF
413 
414        IF( lp_obc_south ) THEN
415           ! ... South obc Force filtering mask for the grad D
416           obcumask(nis0p1:nis1m1,njs0:njs1) = 0.e0
417           obcvmask(nis0p1:nis1m1,njs0:njs1) = 0.e0
418           ! ... set to 0 on South OBC
419           obctmsk(nis0:nis1,njs0:njs1) = 0.e0
420        END IF
421      ENDIF
422
423      IF ( ln_vol_cst .OR. lk_dynspg_flt  ) THEN
424
425         ! 3.1 Total lateral surface
426         ! -------------------------
427         obcsurftot = 0.e0
428 
429         IF( lp_obc_east ) THEN ! ... East open boundary lateral surface
430            DO ji = nie0, nie1
431               DO jj = 1, jpj 
432                  obcsurftot = obcsurftot+hu(ji,jj)*e2u(ji,jj)*uemsk(jj,1) * MAX(obctmsk(ji,jj),obctmsk(ji+1,jj) )
433               END DO
434            END DO
435         END IF
436 
437         IF( lp_obc_west ) THEN ! ... West open boundary lateral surface
438            DO ji = niw0, niw1
439               DO jj = 1, jpj 
440                  obcsurftot = obcsurftot+hu(ji,jj)*e2u(ji,jj)*uwmsk(jj,1) * MAX(obctmsk(ji,jj),obctmsk(ji+1,jj) )
441               END DO
442            END DO
443         END IF
444         IF ( cp_cfg == 'natl' .AND. jp_cfg == 025 ) THEN
445            ! we do not apply the volume correction on the shallow north boundary for NATL025
446         ELSE
447         IF( lp_obc_north ) THEN ! ... North open boundary lateral surface
448            DO jj = njn0, njn1
449               DO ji = 1, jpi
450                  obcsurftot = obcsurftot+hv(ji,jj)*e1v(ji,jj)*vnmsk(ji,1) * MAX(obctmsk(ji,jj),obctmsk(ji,jj+1) )
451               END DO
452            END DO
453         END IF
454         END IF
455 
456         IF( lp_obc_south ) THEN ! ... South open boundary lateral surface
457            DO jj = njs0, njs1
458               DO ji = 1, jpi
459                  obcsurftot = obcsurftot+hv(ji,jj)*e1v(ji,jj)*vsmsk(ji,1) * MAX(obctmsk(ji,jj),obctmsk(ji,jj+1) )
460               END DO
461            END DO
462         END IF
463 
464         IF( lk_mpp )   CALL mpp_sum( obcsurftot )   ! sum over the global domain
465      ENDIF
466
467      ! 5. Control print on mask
468      !    The extremities of the open boundaries must be in land
469      !    or else correspond to an "ocean corner" between two open boundaries.
470      !    corner 1 is southwest, 2 is south east, 3 is northeast, 4 is northwest.
471      ! --------------------------------------------------------------------------
472
473      icorner(:)=0
474
475      ! ... control of the west boundary
476      IF( lp_obc_west ) THEN
477         IF( jpiwob < 2 .OR.  jpiwob >= jpiglo-2 ) THEN
478            WRITE(ctmp1,*) ' jpiwob exceed ', jpiglo-2, 'or less than 2'
479            CALL ctl_stop( ctmp1 )
480         END IF
481         ztestmask(:)=0.
482         DO ji=niw0,niw1
483            IF( (njw0 + njmpp - 1) == jpjwd ) ztestmask(1)=ztestmask(1)+ tmask(ji,njw0,1)
484            IF( (njw1 + njmpp - 1) == jpjwf ) ztestmask(2)=ztestmask(2)+ tmask(ji,njw1,1)
485         END DO
486         IF( lk_mpp )   CALL mpp_sum( ztestmask, 2 )   ! sum over the global domain
487
488         IF( ztestmask(1) /= 0. ) icorner(1)=icorner(1)+1
489         IF( ztestmask(2) /= 0. ) icorner(4)=icorner(4)+1
490      END IF
491
492      ! ... control of the east boundary
493      IF( lp_obc_east ) THEN
494         IF( jpieob < 4 .OR.  jpieob >= jpiglo ) THEN
495            WRITE(ctmp1,*) ' jpieob exceed ', jpiglo, ' or less than 4'
496            CALL ctl_stop( ctmp1 )
497         END IF
498         ztestmask(:)=0.
499         DO ji=nie0p1,nie1p1
500            IF( (nje0 + njmpp - 1) == jpjed ) ztestmask(1)=ztestmask(1)+ tmask(ji,nje0,1)
501            IF( (nje1 + njmpp - 1) == jpjef ) ztestmask(2)=ztestmask(2)+ tmask(ji,nje1,1)
502         END DO
503         IF( lk_mpp )   CALL mpp_sum( ztestmask, 2 )   ! sum over the global domain
504
505        IF( ztestmask(1) /= 0. ) icorner(2)=icorner(2)+1
506        IF( ztestmask(2) /= 0. ) icorner(3)=icorner(3)+1
507      END IF
508
509      ! ... control of the north boundary
510      IF( lp_obc_north ) THEN
511         IF( jpjnob < 4 .OR.  jpjnob >= jpjglo ) THEN
512            WRITE(ctmp1,*) 'jpjnob exceed ', jpjglo, ' or less than 4'
513            CALL ctl_stop( ctmp1 )
514         END IF
515         ztestmask(:)=0.
516         DO jj=njn0p1,njn1p1
517            IF( (nin0 + nimpp - 1) == jpind ) ztestmask(1)=ztestmask(1)+ tmask(nin0,jj,1)
518            IF( (nin1 + nimpp - 1) == jpinf ) ztestmask(2)=ztestmask(2)+ tmask(nin1,jj,1)
519         END DO
520         IF( lk_mpp )   CALL mpp_sum( ztestmask, 2 )   ! sum over the global domain
521
522         IF( ztestmask(1) /= 0. ) icorner(4)=icorner(4)+1
523         IF( ztestmask(2) /= 0. ) icorner(3)=icorner(3)+1
524      END IF
525
526      ! ... control of the south boundary
527      IF( lp_obc_south ) THEN
528         IF( jpjsob < 2 .OR.  jpjsob >= jpjglo-2 ) THEN
529            WRITE(ctmp1,*) ' jpjsob exceed ', jpjglo-2, ' or less than 2'
530            CALL ctl_stop( ctmp1 )
531         END IF
532         ztestmask(:)=0.
533         DO jj=njs0,njs1
534            IF( (nis0 + nimpp - 1) == jpisd ) ztestmask(1)=ztestmask(1)+ tmask(nis0,jj,1)
535            IF( (nis1 + nimpp - 1) == jpisf ) ztestmask(2)=ztestmask(2)+ tmask(nis1,jj,1)
536         END DO
537         IF( lk_mpp )   CALL mpp_sum( ztestmask, 2 )   ! sum over the global domain
538
539         IF( ztestmask(1) /= 0. ) icorner(1)=icorner(1)+1
540         IF( ztestmask(2) /= 0. ) icorner(2)=icorner(2)+1
541      END IF
542
543      IF( icorner(1) == 2 ) THEN
544         IF(lwp) WRITE(numout,*)
545         IF(lwp) WRITE(numout,*) ' South West ocean corner, two open boudaries'
546         IF(lwp) WRITE(numout,*) ' ========== '
547         IF(lwp) WRITE(numout,*)
548         IF( jpisd /= jpiwob.OR.jpjsob /= jpjwd ) &
549              &   CALL ctl_stop( ' Open boundaries do not fit, we stop' )
550
551      ELSE IF( icorner(1) == 1 ) THEN
552         CALL ctl_stop( ' Open boundaries do not fit at SW corner, we stop' )
553      END IF
554
555      IF( icorner(2) == 2 ) THEN
556          IF(lwp) WRITE(numout,*)
557          IF(lwp) WRITE(numout,*) ' South East ocean corner, two open boudaries'
558          IF(lwp) WRITE(numout,*) ' ========== '
559          IF(lwp) WRITE(numout,*)
560          IF( jpisf /= jpieob+1.OR.jpjsob /= jpjed ) &
561               &   CALL ctl_stop( ' Open boundaries do not fit, we stop' )
562      ELSE IF( icorner(2) == 1 ) THEN
563         CALL ctl_stop( ' Open boundaries do not fit at SE corner, we stop' )
564      END IF
565
566      IF( icorner(3) == 2 ) THEN
567         IF(lwp) WRITE(numout,*)
568         IF(lwp) WRITE(numout,*) ' North East ocean corner, two open boudaries'
569         IF(lwp) WRITE(numout,*) ' ========== '
570         IF(lwp) WRITE(numout,*)
571         IF( jpinf /= jpieob+1 .OR. jpjnob+1 /= jpjef ) &
572              &   CALL ctl_stop( ' Open boundaries do not fit, we stop' )
573       ELSE IF( icorner(3) == 1 ) THEN
574          CALL ctl_stop( ' Open boundaries do not fit at NE corner, we stop' )
575       END IF
576
577      IF( icorner(4) == 2 ) THEN
578         IF(lwp) WRITE(numout,*)
579         IF(lwp) WRITE(numout,*) ' North West ocean corner, two open boudaries'
580         IF(lwp) WRITE(numout,*) ' ========== '
581         IF(lwp) WRITE(numout,*)
582         IF( jpind /= jpiwob.OR.jpjnob+1 /= jpjwf ) &
583              &   CALL ctl_stop( ' Open boundaries do not fit, we stop' )
584       ELSE IF( icorner(4) == 1 ) THEN
585          CALL ctl_stop( ' Open boundaries do not fit at NW corner, we stop' )
586       END IF 
587
588      ! 6. Initialization of open boundary variables (u, v, t, s)
589      ! --------------------------------------------------------------
590      !   only if at least one boundary is  radiative
591      IF ( inumfbc < nbobc .AND.  ln_rstart ) THEN
592         !  Restart from restart.obc
593         CALL obc_rst_read
594      ELSE
595
596!         ! ... Initialization to zero of radiation arrays.
597!         !     Those have dimensions of local subdomains
598
599          uebnd(:,:,:,:) = 0.e0   ;   unbnd(:,:,:,:) = 0.e0
600          vebnd(:,:,:,:) = 0.e0   ;   vnbnd(:,:,:,:) = 0.e0
601          tebnd(:,:,:,:) = 0.e0   ;   tnbnd(:,:,:,:) = 0.e0
602          sebnd(:,:,:,:) = 0.e0   ;   snbnd(:,:,:,:) = 0.e0
603
604          uwbnd(:,:,:,:) = 0.e0   ;   usbnd(:,:,:,:) = 0.e0
605          vwbnd(:,:,:,:) = 0.e0   ;   vsbnd(:,:,:,:) = 0.e0
606          twbnd(:,:,:,:) = 0.e0   ;   tsbnd(:,:,:,:) = 0.e0
607          swbnd(:,:,:,:) = 0.e0   ;   ssbnd(:,:,:,:) = 0.e0
608
609      END IF
610
611      ! 7. Control print
612      ! -----------------------------------------------------------------
613
614      ! ... control of the east boundary
615      IF( lp_obc_east ) THEN
616         istop = 0
617         IF( jpieob < 4 .OR.  jpieob >= jpiglo ) THEN
618            IF(lwp) WRITE(numout,cform_err)
619            IF(lwp) WRITE(numout,*) '            jpieob exceed ', jpim1, ' or less than 4'
620            istop = istop + 1
621         END IF
622
623         IF( lk_mpp ) THEN
624            ! ...
625            IF( nimpp > jpieob-5) THEN
626               IF(lwp) WRITE(numout,cform_err)
627               IF(lwp) WRITE(numout,*) '        A sub-domain is too close to the East OBC'
628               IF(lwp) WRITE(numout,*) '        nimpp must be < jpieob-5'
629               istop = istop + 1
630            ENDIF
631         ELSE
632
633            ! ... stop if  e r r o r (s)   detected
634            IF( istop /= 0 ) THEN
635               WRITE(ctmp1,*) istop,' obcini : E R R O R (S) detected : stop'
636               CALL ctl_stop( ctmp1 )
637            ENDIF
638         ENDIF
639      ENDIF
640
641      ! ... control of the west boundary
642      IF( lp_obc_west ) THEN
643         istop = 0
644         IF( jpiwob < 2 .OR.  jpiwob >= jpiglo ) THEN
645            IF(lwp) WRITE(numout,cform_err)
646            IF(lwp) WRITE(numout,*) '            jpiwob exceed ', jpim1, ' or less than 2'
647            istop = istop + 1
648         END IF
649
650         IF( lk_mpp ) THEN
651            IF( (nimpp < jpiwob+5) .AND. (nimpp > 1) ) THEN
652               IF(lwp) WRITE(numout,cform_err)
653               IF(lwp) WRITE(numout,*) '        A sub-domain is too close to the West OBC'
654               IF(lwp) WRITE(numout,*) '        nimpp must be > jpiwob-5 or =1'
655               istop = istop + 1
656            ENDIF
657         ELSE
658   
659            ! ... stop if  e r r o r (s)   detected
660            IF( istop /= 0 ) THEN
661               WRITE(ctmp1,*) istop,' obcini : E R R O R (S) detected : stop'
662               CALL ctl_stop( ctmp1 )
663            ENDIF
664         ENDIF
665      ENDIF
666
667      ! control of the north boundary
668      IF( lp_obc_north ) THEN
669         istop = 0
670         IF( jpjnob < 4 .OR.  jpjnob >= jpjglo ) THEN
671            IF(lwp) WRITE(numout,cform_err)
672            IF(lwp) WRITE(numout,*) '          jpjnob exceed ', jpjm1,' or less than 4'
673            istop = istop + 1
674         END IF
675
676         IF( lk_mpp ) THEN
677            IF( njmpp > jpjnob-5) THEN
678               IF(lwp) WRITE(numout,cform_err)
679               IF(lwp) WRITE(numout,*) '        A sub-domain is too close to the North OBC'
680               IF(lwp) WRITE(numout,*) '        njmpp must be < jpjnob-5'
681               istop = istop + 1
682            ENDIF
683         ELSE
684   
685            ! ... stop if  e r r o r (s)   detected
686            IF( istop /= 0 ) THEN
687                WRITE(ctmp1,*) istop,' obcini : E R R O R (S) detected : stop'
688               CALL ctl_stop( ctmp1 )
689           ENDIF
690         ENDIF
691      ENDIF
692
693      ! control of the south boundary
694      IF( lp_obc_south ) THEN
695         istop = 0
696         IF( jpjsob < 2 .OR. jpjsob >= jpjglo ) THEN
697            IF(lwp) WRITE(numout,cform_err)
698            IF(lwp) WRITE(numout,*) '          jpjsob exceed ', jpjm1,' or less than 2'
699            istop = istop + 1
700         END IF
701
702         IF( lk_mpp ) THEN
703            IF( (njmpp < jpjsob+5) .AND. (njmpp > 1) ) THEN
704               IF(lwp) WRITE(numout,cform_err)
705               IF(lwp) WRITE(numout,*) '        A sub-domain is too close to the South OBC'
706               IF(lwp) WRITE(numout,*) '        njmpp must be > jpjsob+5 or =1'
707               istop = istop + 1
708            ENDIF
709         ELSE
710   
711            ! ... stop if  e r r o r (s)   detected
712            IF( istop /= 0 ) THEN
713               WRITE(ctmp1,*) istop,' obcini : E R R O R (S) detected : stop'
714               CALL ctl_stop( ctmp1 )
715            ENDIF
716         ENDIF
717      ENDIF
718
719   END SUBROUTINE obc_init
720
721#else
722   !!---------------------------------------------------------------------------------
723   !!   Dummy module                                                NO open boundaries
724   !!---------------------------------------------------------------------------------
725CONTAINS
726   SUBROUTINE obc_init      ! Dummy routine
727   END SUBROUTINE obc_init
728#endif
729
730   !!=================================================================================
731END MODULE obcini