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 in tags/nemo_dev_x1/NEMO/OPA_SRC/OBC – NEMO

source: tags/nemo_dev_x1/NEMO/OPA_SRC/OBC/obcini.F90 @ 6474

Last change on this file since 6474 was 25, checked in by opalod, 20 years ago

CT : BUGFIX008 : Running problem for EEL5 configuration is solved

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 33.3 KB
Line 
1 MODULE obcini
2   !!=================================================================================
3   !!                       ***  MODULE  obcini  ***
4   !! OBC initial state :  Open boundary initial state
5   !!=================================================================================
6#if defined key_obc
7   !!---------------------------------------------------------------------------------
8   !!   'key_obc'                                             Open Boundary Conditions
9   !!---------------------------------------------------------------------------------
10   !! * Modules used
11   USE oce             ! ocean dynamics and tracers variables
12   USE dom_oce         ! ocean space and time domain variables
13   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
14   USE phycst          ! physical constants
15   USE obc_oce         ! ocean open boundary conditions
16   USE lib_mpp         ! for mpp_sum
17   USE in_out_manager  ! I/O units
18
19   IMPLICIT NONE
20   PRIVATE
21
22   !! * Routine accessibility
23   PUBLIC obc_init        ! routine called by opa.F90
24
25   !! * Substitutions
26#  include "obc_vectopt_loop_substitute.h90"
27   !!---------------------------------------------------------------------------------
28   !!   OPA 9.0 , LODYC-IPSL  (2003)
29   !!---------------------------------------------------------------------------------
30
31CONTAINS
32   
33   SUBROUTINE obc_init
34      !!----------------------------------------------------------------------
35      !!                 ***  ROUTINE obc_init  ***
36      !!         
37      !! ** Purpose :   Initialization of the dynamics and tracer fields at
38      !!      the open boundaries.
39      !!
40      !! ** Method  :   initialization of open boundary variables
41      !!      (u, v, bsf) over 3 time step and 3 rows
42      !!      (t, s)      over 2 time step and 2 rows
43      !!      if ln_rstart = .FALSE. : no restart, fields set to zero
44      !!      if ln_rstart = .TRUE.  : restart, fields are read in a file
45      !!      if rdpxxx = 0 then lfbc is set true for this boundary.
46      !!
47      !! ** Input   :   restart.obc file, restart file for open boundaries
48      !!
49      !! History :
50      !!   8.0  !  97-07  (G. Madec)  Original code
51      !!        !  97-11  (J.M. Molines)
52      !!   8.5  !  02-11  (C. Talandier, A-M. Treguier) Free surface, F90
53      !!----------------------------------------------------------------------
54      !! * Modules used
55      USE obcrst,   ONLY :   obc_rst_lec   ! Make obc_rst_lec routine available
56      USE obcdom,   ONLY :   obc_dom       ! Make obc_dom routine available
57
58      !! * Local declarations
59      INTEGER  ::   ji, jj, istop
60      INTEGER, DIMENSION(4) ::   icorner
61      REAL(wp) ::   zbsic1, zbsic2, zbsic3
62      REAL(wp), DIMENSION(2) ::   ztestmask
63
64      NAMELIST/namobc/ rdpein, rdpwin, rdpnin, rdpsin,   &
65         &             rdpeob, rdpwob, rdpnob, rdpsob,   &
66         &             zbsic1, zbsic2, zbsic3,           &
67         &             nbic, volemp, nobc_dta
68      !!----------------------------------------------------------------------
69
70      IF(lwp) WRITE(numout,*)
71      IF(lwp) WRITE(numout,*) 'obc_init : initialization of open boundaries'
72      IF(lwp) WRITE(numout,*) '~~~~~~~~'
73
74
75      ! 0. read namelist parameters
76      ! ---------------------------
77      ! default values already set except:
78      zbsic1 = 0.e0
79      zbsic2 = 0.e0
80      zbsic3 = 0.e0
81
82      ! Namelist namobc : open boundaries
83      REWIND( numnam )
84      READ  ( numnam, namobc )
85
86      bsfic0(1) = zbsic1
87      bsfic (2) = zbsic2
88      bsfic (3) = zbsic3
89
90      ! By security we set rdpxin and rdpxob respectively
91      ! to 1. and 15. if the corresponding OBC is not activated
92      IF( .NOT. lpeastobc ) THEN
93         rdpein = 1.
94         rdpeob = 15.
95      END IF
96      IF( .NOT. lpwestobc ) THEN
97         rdpwin = 1.
98         rdpwob = 15.
99      END IF
100      IF( .NOT. lpnorthobc ) THEN
101         rdpnin = 1.
102         rdpnob = 15.
103      END IF
104      IF( .NOT. lpsouthobc ) THEN
105         rdpsin = 1.
106         rdpsob = 15.
107      END IF
108
109      ! number of open boudaries and open boundary indicators
110      nbobc = 0
111      IF( lpeastobc ) THEN
112         nbobc=nbobc+1
113      END IF
114      IF( lpwestobc ) THEN
115         nbobc=nbobc+1
116      END IF
117      IF( lpnorthobc ) THEN
118         nbobc=nbobc+1
119      END IF
120      IF( lpsouthobc ) THEN
121         nbobc=nbobc+1
122      END IF
123
124      IF(lwp) WRITE(numout,*) '         Number of open boundaries    nbobc = ',nbobc
125      IF(lwp) WRITE(numout,*)
126      IF( nbobc /= 0 .AND. jperio /= 0 ) THEN
127         IF(lwp) WRITE(numout,*)
128         IF(lwp) WRITE(numout,*) ' E R R O R : Cyclic or symmetric,',   &
129            ' and open boundary condition are not compatible'
130         IF(lwp) WRITE(numout,*) ' ========== '
131         IF(lwp) WRITE(numout,*)
132         nstop = nstop + 1
133      END IF
134
135      ! control prints
136      IF(lwp) WRITE(numout,*) '         namobc'
137      IF(lwp) WRITE(numout,*) ' '
138      IF(lwp) WRITE(numout,*) '         data in file (=1) or     nobc_dta = ', nobc_dta
139      IF(lwp) WRITE(numout,*) '         initial state used (=0)             '
140      IF( lwp.AND.lpeastobc ) THEN
141         WRITE(numout,*) '         East open boundary :'
142         WRITE(numout,*) '              i index                    jpieob = ', jpieob
143         WRITE(numout,*) '              damping time scale (days)  rdpeob = ', rdpeob
144         WRITE(numout,*) '              damping time scale (days)  rdpein = ', rdpein
145      END IF
146
147      IF( lwp.AND.lpwestobc ) THEN
148         WRITE(numout,*) '         West open boundary :'
149         WRITE(numout,*) '              i index                    jpiwob = ', jpiwob
150         WRITE(numout,*) '              damping time scale (days)  rdpwob = ', rdpwob
151         WRITE(numout,*) '              damping time scale (days)  rdpwin = ', rdpwin
152      END IF
153
154      IF( lwp.AND.lpnorthobc ) THEN
155         WRITE(numout,*) '         North open boundary :'
156         WRITE(numout,*) '               j index                    jpjnob = ', jpjnob
157         WRITE(numout,*) '               damping time scale (days)  rdpnob = ', rdpnob
158         WRITE(numout,*) '               damping time scale (days)  rdpnin = ', rdpnin
159      END IF
160
161      IF( lwp.AND.lpsouthobc ) THEN
162         WRITE(numout,*) '         South open boundary :'
163         WRITE(numout,*) '               j index                    jpjsob = ', jpjsob
164         WRITE(numout,*) '               damping time scale (days)  rdpsob = ', rdpsob
165         WRITE(numout,*) '               damping time scale (days)  rdpsin = ', rdpsin
166         WRITE(numout,*) ' '
167      END IF
168
169      ! 1. Initialisation of constants
170      ! ------------------------------
171
172      ! ... convert rdp$ob in seconds
173      rdpein = rdpein * rday
174      rdpwin = rdpwin * rday
175      rdpnin = rdpnin * rday
176      rdpsin = rdpsin * rday
177      rdpeob = rdpeob * rday
178      rdpwob = rdpwob * rday
179      rdpnob = rdpnob * rday
180      rdpsob = rdpsob * rday
181      lfbceast  = .FALSE.
182      lfbcwest  = .FALSE.
183      lfbcnorth = .FALSE.
184      lfbcsouth = .FALSE.
185      ! ... look for Fixed Boundaries (rdp = 0 )
186      ! ... When specified, lbcxxx flags are set to TRUE and rdpxxx are set to
187      ! ...  a small arbitrary value, (to avoid division by zero further on).
188      ! ...  rdpxxx is not used anymore.
189      IF( lpeastobc .AND. ( rdpein == 0 .AND. rdpeob == 0 ) )  THEN
190         lfbceast = .TRUE.
191         rdpein = 1e-3
192         rdpeob = 1e-3
193      END IF
194      IF( lpwestobc .AND. ( rdpwin == 0 .AND. rdpwob == 0 ) )  THEN
195         lfbcwest = .TRUE.
196         rdpwin = 1e-3
197         rdpwob = 1e-3
198      END IF
199      IF( lpnorthobc .AND. ( rdpnin == 0 .AND. rdpnob == 0 ) )  THEN
200         lfbcnorth = .TRUE.
201         rdpnin = 1e-3
202         rdpnob = 1e-3
203      END IF
204      IF( lpsouthobc .AND. ( rdpsin == 0 .AND. rdpsob == 0 ) )  THEN
205         lfbcsouth = .TRUE.
206         rdpsin = 1e-3
207         rdpsob = 1e-3
208      END IF
209
210      ! 2.  Clever mpp indices for loops on the open boundaries.
211      !     The loops will be performed only on the processors
212      !     that contain a given open boundary.
213      ! --------------------------------------------------------
214
215      IF( lpeastobc ) THEN
216         ! ...   mpp initialization
217         nie0   = max( 1, min(jpieob   - nimpp+1, jpi     ) )
218         nie1   = max( 0, min(jpieob   - nimpp+1, jpi - 1 ) )
219         nie0p1 = max( 1, min(jpieob+1 - nimpp+1, jpi     ) )
220         nie1p1 = max( 0, min(jpieob+1 - nimpp+1, jpi - 1 ) )
221         nie0m1 = max( 1, min(jpieob-1 - nimpp+1, jpi     ) )
222         nie1m1 = max( 0, min(jpieob-1 - nimpp+1, jpi - 1 ) )
223         nje0   = max( 2, min(jpjed    - njmpp+1, jpj     ) )
224         nje1   = max( 0, min(jpjef    - njmpp+1, jpj - 1 ) )
225         nje0p1 = max( 1, min(jpjedp1  - njmpp+1, jpj     ) )
226         nje0m1 = max( 1, min(jpjed    - njmpp+1, jpj     ) )
227         nje1m1 = max( 0, min(jpjefm1  - njmpp+1, jpj - 1 ) )
228         nje1m2 = max( 0, min(jpjefm1-1- njmpp+1, jpj - 1 ) )
229         IF(lwp) THEN
230            IF( lfbceast ) THEN
231               WRITE(numout,*)'     '
232               WRITE(numout,*)'         Specified East Open Boundary'
233            ELSE
234               WRITE(numout,*)'     '
235               WRITE(numout,*)'         Radiative East Open Boundary'
236            END IF
237         END IF
238      END IF
239
240      IF( lpwestobc ) THEN
241         ! ...   mpp initialization
242         niw0   = max( 1, min(jpiwob   - nimpp+1, jpi     ) )
243         niw1   = max( 0, min(jpiwob   - nimpp+1, jpi - 1 ) )
244         niw0p1 = max( 1, min(jpiwob+1 - nimpp+1, jpi     ) )
245         niw1p1 = max( 0, min(jpiwob+1 - nimpp+1, jpi - 1 ) )
246         njw0   = max( 2, min(jpjwd    - njmpp+1, jpj     ) )
247         njw1   = max( 0, min(jpjwf    - njmpp+1, jpj - 1 ) )
248         njw0p1 = max( 1, min(jpjwdp1  - njmpp+1, jpj     ) )
249         njw0m1 = max( 1, min(jpjwd    - njmpp+1, jpj     ) )
250         njw1m1 = max( 0, min(jpjwfm1  - njmpp+1, jpj - 1 ) )
251         njw1m2 = max( 0, min(jpjwfm1-1- njmpp+1, jpj - 1 ) )
252         IF(lwp) THEN
253            IF( lfbcwest ) THEN
254               WRITE(numout,*)'     '
255               WRITE(numout,*)'         Specified West Open Boundary'
256            ELSE
257               WRITE(numout,*)'     '
258               WRITE(numout,*)'         Radiative West Open Boundary'
259            END IF
260         END IF
261      END IF
262 
263      IF( lpnorthobc ) THEN
264         ! ...   mpp initialization
265         nin0   = max( 2, min(jpind    - nimpp+1, jpi     ) )
266         nin1   = max( 0, min(jpinf    - nimpp+1, jpi - 1 ) )
267         nin0p1 = max( 1, min(jpindp1  - nimpp+1, jpi     ) )
268         nin0m1 = max( 1, min(jpind    - nimpp+1, jpi     ) )
269         nin1m1 = max( 0, min(jpinfm1  - nimpp+1, jpi - 1 ) )
270         nin1m2 = max( 0, min(jpinfm1-1- nimpp+1, jpi - 1 ) )
271         njn0   = max( 1, min(jpjnob   - njmpp+1, jpj     ) )
272         njn1   = max( 0, min(jpjnob   - njmpp+1, jpj - 1 ) )
273         njn0p1 = max( 1, min(jpjnob+1 - njmpp+1, jpj     ) )
274         njn1p1 = max( 0, min(jpjnob+1 - njmpp+1, jpj - 1 ) )
275         njn0m1 = max( 1, min(jpjnob-1 - njmpp+1, jpj     ) )
276         njn1m1 = max( 0, min(jpjnob-1 - njmpp+1, jpj - 1 ) )
277         IF(lwp) THEN
278            IF( lfbcnorth ) THEN
279               WRITE(numout,*)'     '
280               WRITE(numout,*)'         Specified North Open Boundary'
281            ELSE
282               WRITE(numout,*)'     '
283               WRITE(numout,*)'         Radiative North Open Boundary'
284            END IF
285         END IF
286      END IF
287
288      IF( lpsouthobc ) THEN
289         ! ...   mpp initialization
290         nis0   = max( 2, min(jpisd    - nimpp+1, jpi     ) )
291         nis1   = max( 0, min(jpisf    - nimpp+1, jpi - 1 ) )
292         nis0p1 = max( 1, min(jpisdp1  - nimpp+1, jpi     ) )
293         nis0m1 = max( 1, min(jpisd    - nimpp+1, jpi     ) )
294         nis1m1 = max( 0, min(jpisfm1  - nimpp+1, jpi - 1 ) )
295         nis1m2 = max( 0, min(jpisfm1-1- nimpp+1, jpi - 1 ) )
296         njs0   = max( 1, min(jpjsob   - njmpp+1, jpj     ) )
297         njs1   = max( 0, min(jpjsob   - njmpp+1, jpj - 1 ) )
298         njs0p1 = max( 1, min(jpjsob+1 - njmpp+1, jpj     ) )
299         njs1p1 = max( 0, min(jpjsob+1 - njmpp+1, jpj - 1 ) )
300         IF(lwp) THEN
301            IF( lfbcsouth ) THEN
302               WRITE(numout,*)'     '
303               WRITE(numout,*)'         Specified South Open Boundary'
304            ELSE
305               WRITE(numout,*)'     '
306               WRITE(numout,*)'         Radiative South Open Boundary'
307            END IF
308         END IF
309      END IF
310
311      ! 3. mask correction for OBCs
312      ! ---------------------------
313
314      IF( lpeastobc ) THEN
315         !... (jpjed,jpjefm1),jpieob
316         DO jj = nje0, nje1m1
317# if defined key_dynspg_fsc
318            DO ji = nie0p1, nie1p1
319# else
320            DO ji = nie0, nie1
321# endif
322               bmask(ji,jj) = 0.e0
323            END DO
324         END DO
325
326         ! ... initilization to zero
327         uemsk(:,:) = 0.e0
328         vemsk(:,:) = 0.e0
329         temsk(:,:) = 0.e0
330
331         ! ... set 2D mask on East OBC,  Vopt
332         DO ji = fs_nie0, fs_nie1
333            DO jj = nje0, nje1
334               uemsk(jj,:) = umask(ji,  jj,:)
335               vemsk(jj,:) = vmask(ji+1,jj,:)
336               temsk(jj,:) = tmask(ji+1,jj,:)
337            END DO
338         END DO
339
340      END IF
341
342      IF( lpwestobc ) THEN
343         ! ... (jpjwd,jpjwfm1),jpiwob
344         DO jj = njw0, njw1m1
345            DO ji = niw0, niw1
346               bmask(ji,jj) = 0.e0
347            END DO
348         END DO
349
350         ! ... initilization to zero
351         uwmsk(:,:) = 0.e0
352         vwmsk(:,:) = 0.e0
353         twmsk(:,:) = 0.e0
354
355         ! ... set 2D mask on West OBC,  Vopt
356         DO ji = fs_niw0, fs_niw1
357            DO jj = njw0, njw1
358               uwmsk(jj,:) = umask(ji,jj,:)
359               vwmsk(jj,:) = vmask(ji,jj,:)
360               twmsk(jj,:) = tmask(ji,jj,:)
361            END DO
362         END DO
363
364      END IF
365
366      IF( lpnorthobc ) THEN
367         ! ... jpjnob,(jpind,jpisfm1)
368# if defined key_dynspg_fsc
369         DO jj = njn0p1, njn1p1
370# else
371         DO jj = njn0, njn1
372# endif
373            DO ji = nin0, nin1m1
374               bmask(ji,jj) = 0.e0
375            END DO
376         END DO
377
378         ! ... initilization to zero
379         unmsk(:,:) = 0.e0
380         vnmsk(:,:) = 0.e0
381         tnmsk(:,:) = 0.e0
382
383         ! ... set 2D mask on North OBC,  Vopt
384         DO jj = fs_njn0, fs_njn1
385            DO ji = nin0, nin1
386               unmsk(ji,:) = umask(ji,jj+1,:)
387               vnmsk(ji,:) = vmask(ji,jj  ,:)
388               tnmsk(ji,:) = tmask(ji,jj+1,:)
389            END DO
390         END DO
391
392      END IF
393
394      IF( lpsouthobc ) THEN 
395         ! ... jpjsob,(jpisd,jpisfm1)
396         DO jj = njs0, njs1
397            DO ji = nis0, nis1m1
398               bmask(ji,jj) = 0.e0
399            END DO
400         END DO
401
402         ! ... initilization to zero
403         usmsk(:,:) = 0.e0
404         vsmsk(:,:) = 0.e0
405         tsmsk(:,:) = 0.e0
406
407         ! ... set 2D mask on South OBC,  Vopt
408         DO jj = njs0, njs1
409            DO ji = nis0, nis1
410               usmsk(ji,:) = umask(ji,jj,:)
411               vsmsk(ji,:) = vmask(ji,jj,:)
412               tsmsk(ji,:) = tmask(ji,jj,:)
413            END DO
414         END DO
415
416      END IF
417
418# if defined key_dynspg_fsc
419
420      ! ... Initialize obcumask and obcvmask for the Force filtering
421      !     boundary condition in dynspg_fsc
422      obcumask(:,:) = umask(:,:,1)
423      obcvmask(:,:) = vmask(:,:,1)
424
425      ! ... Initialize obctmsk on overlap region and obcs. This mask
426      !     is used in obcvol.F90 to calculate cumulate flux E-P.
427      !     - no flux E-P on obcs and overlap region (jpereci = jprecj = 1)
428      obctmsk(:,:) = tmask(:,:,1)     
429      obctmsk(1  ,:) = 0.e0           
430      obctmsk(jpi,:) = 0.e0           
431      obctmsk(:  ,1) = 0.e0           
432      obctmsk(:,jpj) = 0.e0           
433
434      IF( lpeastobc ) THEN
435         ! ... East obc Force filtering mask for the grad D
436         DO ji = nie0, nie1
437            DO jj = nje0p1, nje1m1
438               obcumask(ji  ,jj)=0.e0
439               obcvmask(ji+1,jj)=0.e0
440            END DO
441         END DO
442
443         ! ... set to 0 on East OBC
444         DO jj = nje0p1, nje1m1
445            DO ji = nie0p1, nie1p1
446               obctmsk(ji,jj) = 0.e0
447            END DO
448         END DO
449      END IF
450
451      IF( lpwestobc ) THEN
452         ! ... West obc Force filtering mask for the grad D
453         DO ji = niw0, niw1
454            DO jj = njw0p1, njw1m1
455               obcumask(ji,jj)=0.e0
456               obcvmask(ji,jj)=0.e0
457            END DO
458         END DO
459
460         ! ... set to 0 on West OBC
461         DO jj = njw0p1, njw1m1
462            DO ji = niw0, niw1
463               obctmsk(ji,jj) = 0.e0
464            END DO
465         END DO
466      END IF
467
468      IF( lpnorthobc ) THEN
469         ! ... North obc Force filtering mask for the grad D
470         DO jj = njn0, njn1
471            DO ji = nin0p1, nin1m1
472               obcvmask(ji,jj  )=0.e0
473               obcumask(ji,jj+1)=0.e0
474            END DO
475         END DO
476
477         ! ... set to 0 on North OBC
478         DO jj = njn0p1, njn1p1
479            DO ji = nin0p1, nin1m1
480               obctmsk(ji,jj) = 0.e0
481            END DO
482         END DO
483      END IF
484
485      IF( lpsouthobc ) THEN
486         ! ... South obc Force filtering mask for the grad D
487         DO jj = njs0, njs1 
488            DO ji = nis0p1, nis1m1
489               obcumask(ji,jj)=0.e0
490               obcvmask(ji,jj)=0.e0
491            END DO
492         END DO
493
494         ! ... set to 0 on South OBC
495         DO jj = njs0, njs1
496            DO ji = nis0p1, nis1m1
497               obctmsk(ji,jj) = 0.e0
498            END DO
499         END DO
500      END IF
501
502      ! 3.1 Total lateral surface for each open boundary
503      ! ------------------------------------------------
504
505      ! ... West open boundary surface
506      IF( lpwestobc ) THEN
507         DO ji = niw0, niw1
508            DO jj = 1, jpj 
509               obcsurftot = obcsurftot+hu(ji,jj)*e2u(ji,jj)*uwmsk(jj,1)
510            END DO
511         END DO
512       END IF
513
514      ! ... East open boundary surface
515      IF( lpeastobc ) THEN
516         DO ji = nie0, nie1
517            DO jj = 1, jpj 
518               obcsurftot = obcsurftot+hu(ji,jj)*e2u(ji,jj)*uemsk(jj,1)
519            END DO
520         END DO
521      END IF
522
523      ! ... North open boundary vertical surface
524      IF( lpnorthobc ) THEN
525         DO jj = njn0, njn1
526            DO ji = 1, jpi
527               obcsurftot = obcsurftot+hv(ji,jj)*e1v(ji,jj)*vnmsk(ji,1)
528            END DO
529         END DO
530      END IF
531
532      ! ... South open boundary vertical surface
533      IF( lpsouthobc ) THEN
534         DO jj = njs0, njs1
535            DO ji = 1, jpi
536               obcsurftot = obcsurftot+hv(ji,jj)*e1v(ji,jj)*vsmsk(ji,1)
537            END DO
538         END DO
539      END IF
540      IF( lk_mpp )   CALL mpp_sum( obcsurftot )   ! sum over the global domain
541
542# endif
543
544      ! 5. Control print on mask
545      !    The extremities of the open boundaries must be in land
546      !    or else correspond to an "ocean corner" between two open boundaries.
547      !    corner 1 is southwest, 2 is south east, 3 is northeast, 4 is northwest.
548      ! --------------------------------------------------------------------------
549
550      icorner(:)=0
551
552      ! ... control of the west boundary
553      IF( lpwestobc ) THEN
554         IF( jpiwob < 2 .OR.  jpiwob >= jpiglo-2 ) THEN
555            IF(lwp) WRITE(numout,*)
556            IF(lwp) WRITE(numout,*) ' E R R O R : jpiwob exceed ', jpiglo-2, 'or less than 2'
557            IF(lwp) WRITE(numout,*) ' ========== '
558            IF(lwp) WRITE(numout,*)
559            nstop = nstop + 1
560         END IF
561         ztestmask(:)=0.
562         DO ji=niw0,niw1
563            IF( (njw0 + njmpp - 1) == jpjwd ) ztestmask(1)=ztestmask(1)+ tmask(ji,njw0,1)
564            IF( (njw1 + njmpp - 1) == jpjwf ) ztestmask(2)=ztestmask(2)+ tmask(ji,njw1,1)
565         END DO
566         IF( lk_mpp )   CALL mpp_sum( ztestmask, 2 )   ! sum over the global domain
567
568         IF( ztestmask(1) /= 0. ) icorner(1)=icorner(1)+1
569         IF( ztestmask(2) /= 0. ) icorner(4)=icorner(4)+1
570      END IF
571
572      ! ... control of the east boundary
573      IF( lpeastobc ) THEN
574         IF( jpieob < 4 .OR.  jpieob >= jpiglo ) THEN
575            IF(lwp) WRITE(numout,*)
576            IF(lwp) WRITE(numout,*) ' E R R O R : jpieob exceed ', jpiglo, ' or less than 4'
577            IF(lwp) WRITE(numout,*) ' ========== '
578            IF(lwp) WRITE(numout,*)
579            nstop = nstop + 1
580         END IF
581         ztestmask(:)=0.
582         DO ji=nie0p1,nie1p1
583            IF( (nje0 + njmpp - 1) == jpjed ) ztestmask(1)=ztestmask(1)+ tmask(ji,nje0,1)
584            IF( (nje1 + njmpp - 1) == jpjef ) ztestmask(2)=ztestmask(2)+ tmask(ji,nje1,1)
585         END DO
586         IF( lk_mpp )   CALL mpp_sum( ztestmask, 2 )   ! sum over the global domain
587
588        IF( ztestmask(1) /= 0. ) icorner(2)=icorner(2)+1
589        IF( ztestmask(2) /= 0. ) icorner(3)=icorner(3)+1
590      END IF
591
592      ! ... control of the north boundary
593      IF( lpnorthobc ) THEN
594         IF( jpjnob < 4 .OR.  jpjnob >= jpjglo ) THEN
595            IF(lwp) WRITE(numout,*)
596            IF(lwp) WRITE(numout,*) ' E R R O R : jpjnob exceed ', jpjglo, ' or less than 4'
597            IF(lwp) WRITE(numout,*) ' ========== '
598            IF(lwp) WRITE(numout,*)
599            nstop = nstop + 1
600         END IF
601         ztestmask(:)=0.
602         DO jj=njn0p1,njn1p1
603            IF( (nin0 + nimpp - 1) == jpind ) ztestmask(1)=ztestmask(1)+ tmask(nin0,jj,1)
604            IF( (nin1 + nimpp - 1) == jpinf ) ztestmask(2)=ztestmask(2)+ tmask(nin1,jj,1)
605         END DO
606         IF( lk_mpp )   CALL mpp_sum( ztestmask, 2 )   ! sum over the global domain
607
608         IF( ztestmask(1) /= 0. ) icorner(4)=icorner(4)+1
609         IF( ztestmask(2) /= 0. ) icorner(3)=icorner(3)+1
610      END IF
611
612      ! ... control of the south boundary
613      IF( lpsouthobc ) THEN
614         IF( jpjsob < 2 .OR.  jpjsob >= jpjglo-2 ) THEN
615            IF(lwp) WRITE(numout,*)
616            IF(lwp) WRITE(numout,*) ' E R R O R : jpjsob exceed ', jpjglo-2, ' or less than 2'
617            IF(lwp) WRITE(numout,*) ' ========== '
618            IF(lwp) WRITE(numout,*)
619            nstop = nstop + 1
620         END IF
621         ztestmask(:)=0.
622         DO jj=njs0,njs1
623            IF( (nis0 + nimpp - 1) == jpisd ) ztestmask(1)=ztestmask(1)+ tmask(nis0,jj,1)
624            IF( (nis1 + nimpp - 1) == jpisf ) ztestmask(2)=ztestmask(2)+ tmask(nis1,jj,1)
625         END DO
626         IF( lk_mpp )   CALL mpp_sum( ztestmask, 2 )   ! sum over the global domain
627
628         IF( ztestmask(1) /= 0. ) icorner(1)=icorner(1)+1
629         IF( ztestmask(2) /= 0. ) icorner(2)=icorner(2)+1
630      END IF
631
632      IF( icorner(1) == 2 ) THEN
633         IF(lwp) WRITE(numout,*)
634         IF(lwp) WRITE(numout,*) ' South West ocean corner, two open boudaries'
635         IF(lwp) WRITE(numout,*) ' ========== '
636         IF(lwp) WRITE(numout,*)
637         IF( jpisd /= jpiwob.OR.jpjsob /= jpjwd ) THEN
638            IF(lwp) WRITE(numout,*) ' Open boundaries do not fit, we stop'
639            nstop = nstop + 1
640         END IF
641      ELSE IF( icorner(1) == 1 ) THEN
642              IF(lwp) WRITE(numout,*) ' Open boundaries do not fit at SW corner, we stop'
643              nstop = nstop + 1
644      END IF
645
646      IF( icorner(2) == 2 ) THEN
647          IF(lwp) WRITE(numout,*)
648          IF(lwp) WRITE(numout,*) ' South East ocean corner, two open boudaries'
649          IF(lwp) WRITE(numout,*) ' ========== '
650          IF(lwp) WRITE(numout,*)
651          IF( jpisf /= jpieob+1.OR.jpjsob /= jpjed ) THEN
652             IF(lwp) WRITE(numout,*) ' Open boundaries do not fit, we stop'
653             nstop = nstop + 1
654          END IF
655      ELSE IF( icorner(2) == 1 ) THEN
656              IF(lwp) WRITE(numout,*) ' Open boundaries do not fit at SE corner, we stop'
657              nstop = nstop + 1
658      END IF
659
660      IF( icorner(3) == 2 ) THEN
661         IF(lwp) WRITE(numout,*)
662         IF(lwp) WRITE(numout,*) ' North East ocean corner, two open boudaries'
663         IF(lwp) WRITE(numout,*) ' ========== '
664         IF(lwp) WRITE(numout,*)
665         IF( jpinf /= jpieob+1 .OR. jpjnob+1 /= jpjef ) THEN
666            IF(lwp) WRITE(numout,*) ' Open boundaries do not fit, we stop'
667            nstop = nstop + 1
668         END IF
669       ELSE IF( icorner(3) == 1 ) THEN
670               IF(lwp) WRITE(numout,*) ' Open boundaries do not fit at NE corner, we stop'
671               nstop = nstop + 1
672       END IF
673
674      IF( icorner(4) == 2 ) THEN
675         IF(lwp) WRITE(numout,*)
676         IF(lwp) WRITE(numout,*) ' North West ocean corner, two open boudaries'
677         IF(lwp) WRITE(numout,*) ' ========== '
678         IF(lwp) WRITE(numout,*)
679         IF( jpind /= jpiwob.OR.jpjnob+1 /= jpjwf ) THEN
680            IF(lwp) WRITE(numout,*) ' Open boundaries do not fit, we stop'
681            nstop = nstop + 1
682         END IF
683       ELSE IF( icorner(4) == 1 ) THEN
684               IF(lwp) WRITE(numout,*) ' Open boundaries do not fit at NW corner, we stop'
685               nstop = nstop + 1
686       END IF 
687
688      ! 6. Initialization of open boundary variables (u, v, bsf, t, s)
689      ! --------------------------------------------------------------
690
691      ! ... Restart from restart.obc
692      IF( ln_rstart ) THEN
693         CALL obc_rst_lec
694      ELSE
695
696          ! ... Initialization to zero of radiation arrays.
697          !     Those have dimensions of local subdomains
698
699          bebnd(:,:,:)   = 0.e0   ;   bnbnd(:,:,:)   = 0.e0
700          uebnd(:,:,:,:) = 0.e0   ;   unbnd(:,:,:,:) = 0.e0
701          vebnd(:,:,:,:) = 0.e0   ;   vnbnd(:,:,:,:) = 0.e0
702          tebnd(:,:,:,:) = 0.e0   ;   tnbnd(:,:,:,:) = 0.e0 
703          sebnd(:,:,:,:) = 0.e0   ;   snbnd(:,:,:,:) = 0.e0
704
705          bwbnd(:,:,:)   = 0.e0   ;   bsbnd(:,:,:)   = 0.e0
706          uwbnd(:,:,:,:) = 0.e0   ;   usbnd(:,:,:,:) = 0.e0
707          vwbnd(:,:,:,:) = 0.e0   ;   vsbnd(:,:,:,:) = 0.e0
708          twbnd(:,:,:,:) = 0.e0   ;   tsbnd(:,:,:,:) = 0.e0 
709          swbnd(:,:,:,:) = 0.e0   ;   ssbnd(:,:,:,:) = 0.e0
710
711      END IF
712
713# if ! defined key_dynspg_fsc
714      ! 7. Isolated coastline arrays initialization (rigid lid case only)
715      ! -----------------------------------------------------------------
716      CALL obc_dom
717# endif
718
719      ! 8. Control print
720      ! ... control of the east boundary
721      IF( lpeastobc ) THEN
722         istop = 0
723         IF( jpieob < 4 .OR.  jpieob >= jpiglo ) THEN
724            IF(lwp) WRITE(numout,cform_err)
725            IF(lwp) WRITE(numout,*) '            jpieob exceed ', jpim1, ' or less than 4'
726            istop = istop + 1
727         END IF
728
729         IF( lk_mpp ) THEN
730            ! ...
731            IF( nimpp > jpieob-5) THEN
732               IF(lwp) WRITE(numout,cform_err)
733               IF(lwp) WRITE(numout,*) '        A sub-domain is too close to the East OBC'
734               IF(lwp) WRITE(numout,*) '        nimpp must be < jpieob-5'
735               istop = istop + 1
736            ENDIF
737         ELSE
738            IF( tmask(jpieob+1,jpjed  ,1) /= 0. .OR.    &
739                tmask(jpieob+1,jpjed+1,1) /= 1.         ) THEN
740               IF(lwp) WRITE(numout,cform_err)
741               IF(lwp) WRITE(numout,*) '           starting point is not a land T-point.'
742               IF(lwp) WRITE(numout,*) '      or   starting point + 1 is not a ocean T-point.'
743               istop = istop + 1
744            END IF
745            IF( tmask(jpieob+1,jpjef  ,1) /= 0. .OR.  &
746                tmask(jpieob+1,jpjef-1,1) /= 1.       ) THEN
747               IF(lwp) WRITE(numout,cform_err)
748               IF(lwp) WRITE(numout,*) '           ending point is not a land T-point.'
749               IF(lwp) WRITE(numout,*) '        or ending point - 1 is not a ocean T-point.'
750               istop = istop + 1
751            END IF
752
753            ! ... stop if  e r r o r (s)   detected
754            IF( istop /= 0 ) THEN
755               IF(lwp)WRITE(numout,*)
756               IF(lwp)WRITE(numout,*) istop,' E R R O R (S) detected : stop'
757               IF(lwp)WRITE(numout,*) ' =============== '
758               IF(lwp)WRITE(numout,*)
759               nstop = nstop + 1
760            ENDIF
761         ENDIF
762      ENDIF
763
764      ! ... control of the west boundary
765      IF( lpwestobc ) THEN
766         istop = 0
767         IF( jpiwob < 2 .OR.  jpiwob >= jpiglo ) THEN
768            IF(lwp) WRITE(numout,cform_err)
769            IF(lwp) WRITE(numout,*) '            jpiwob exceed ', jpim1, ' or less than 2'
770            istop = istop + 1
771         END IF
772
773         IF( lk_mpp ) THEN
774            IF( (nimpp < jpiwob+5) .AND. (nimpp > 1) ) THEN
775               IF(lwp) WRITE(numout,cform_err)
776               IF(lwp) WRITE(numout,*) '        A sub-domain is too close to the West OBC'
777               IF(lwp) WRITE(numout,*) '        nimpp must be > jpiwob-5 or =1'
778               istop = istop + 1
779            ENDIF
780         ELSE
781            IF( tmask(jpiwob,jpjwd  ,1) /= 0. .OR.    &
782                tmask(jpiwob,jpjwd+1,1) /= 1.         ) THEN
783               IF(lwp) WRITE(numout,cform_err)
784               IF(lwp) WRITE(numout,*) '           starting point is not a land T-point.'
785               IF(lwp) WRITE(numout,*) '      or   starting point + 1 is not a ocean T-point.'
786               istop = istop + 1
787            END IF
788            IF ( tmask(jpieob+1,jpjef  ,1) /= 0. .OR.  &
789                 tmask(jpieob+1,jpjef-1,1) /= 1.       ) THEN
790                IF(lwp) WRITE(numout,cform_err)
791                IF(lwp) WRITE(numout,*) '           ending point is not a land T-point.'
792                IF(lwp) WRITE(numout,*) '        or ending point - 1 is not a ocean T-point.'
793                istop = istop + 1
794            END IF
795   
796            ! ... stop if  e r r o r (s)   detected
797            IF( istop /= 0 ) THEN
798               IF(lwp)WRITE(numout,*)
799               IF(lwp)WRITE(numout,*) istop,' E R R O R (S) detected : stop'
800               IF(lwp)WRITE(numout,*) ' =============== '
801               IF(lwp)WRITE(numout,*)
802               nstop = nstop + 1
803            ENDIF
804         ENDIF
805      ENDIF
806
807      ! control of the north boundary
808      IF( lpnorthobc ) THEN
809         istop = 0
810         IF( jpjnob < 4 .OR.  jpjnob >= jpjglo ) THEN
811            IF(lwp) WRITE(numout,cform_err)
812            IF(lwp) WRITE(numout,*) '          jpjnob exceed ', jpjm1,' or less than 4'
813            istop = istop + 1
814         END IF
815
816         IF( lk_mpp ) THEN
817            IF( njmpp > jpjnob-5) THEN
818               IF(lwp) WRITE(numout,cform_err)
819               IF(lwp) WRITE(numout,*) '        A sub-domain is too close to the North OBC'
820               IF(lwp) WRITE(numout,*) '        njmpp must be < jpjnob-5'
821               istop = istop + 1
822            ENDIF
823         ELSE
824            IF( tmask(jpind  , jpjnob+1,1) /= 0. .OR.    &
825                tmask(jpind+1, jpjnob+1,1) /= 1.         ) THEN
826               IF(lwp) WRITE(numout,cform_err)
827               IF(lwp) WRITE(numout,*) '           starting point is not a land T-point.'
828               IF(lwp) WRITE(numout,*) '      or   starting point + 1 is not a ocean T-point.'
829               istop = istop + 1
830            END IF
831            IF( tmask(jpinf  ,jpjnob+1,1) /= 0. .OR.  &
832                tmask(jpinf-1,jpjnob+1,1) /= 1.       ) THEN
833               IF(lwp) WRITE(numout,cform_err)
834               IF(lwp) WRITE(numout,*) '           ending point is not a land T-point.'
835               IF(lwp) WRITE(numout,*) '        or ending point - 1 is not a ocean T-point.'
836               istop = istop + 1
837            END IF
838   
839            ! ... stop if  e r r o r (s)   detected
840            IF( istop /= 0 ) THEN
841               IF(lwp)WRITE(numout,*)
842               IF(lwp)WRITE(numout,*) istop,' E R R O R (S) detected : stop'
843               IF(lwp)WRITE(numout,*) ' =============== '
844               IF(lwp)WRITE(numout,*)
845               nstop = nstop + 1
846            ENDIF
847         ENDIF
848      ENDIF
849
850      ! control of the south boundary
851      IF( lpsouthobc ) THEN
852         istop = 0
853         IF( jpjsob < 2 .OR. jpjsob >= jpjglo ) THEN
854            IF(lwp) WRITE(numout,cform_err)
855            IF(lwp) WRITE(numout,*) '          jpjsob exceed ', jpjm1,' or less than 2'
856            istop = istop + 1
857         END IF
858
859         IF( lk_mpp ) THEN
860            IF( (njmpp < jpjsob+5) .AND. (njmpp > 1) ) THEN
861               IF(lwp) WRITE(numout,cform_err)
862               IF(lwp) WRITE(numout,*) '        A sub-domain is too close to the South OBC'
863               IF(lwp) WRITE(numout,*) '        njmpp must be > jpjsob+5 or =1'
864               istop = istop + 1
865            ENDIF
866         ELSE
867            IF( tmask(jpisd  , jpjsob,1) /= 0. .OR.    &
868                tmask(jpisd+1, jpjsob,1) /= 1.         ) THEN
869               IF(lwp) WRITE(numout,cform_err)
870               IF(lwp) WRITE(numout,*) '           starting point is not a land T-point.'
871               IF(lwp) WRITE(numout,*) '      or   starting point + 1 is not a ocean T-point.'
872               istop = istop + 1
873            END IF
874            IF( tmask(jpisf  ,jpjsob,1) /= 0. .OR.  &
875                tmask(jpisf-1,jpjsob,1) /= 1.       ) THEN
876               IF(lwp) WRITE(numout,cform_err)
877               IF(lwp) WRITE(numout,*) '           ending point is not a land T-point.'
878               IF(lwp) WRITE(numout,*) '        or ending point - 1 is not a ocean T-point.'
879               istop = istop + 1
880            END IF
881   
882            ! ... stop if  e r r o r (s)   detected
883            IF( istop /= 0 ) THEN
884               IF(lwp)WRITE(numout,*)
885               IF(lwp)WRITE(numout,*) istop,' E R R O R (S) detected : stop'
886               IF(lwp)WRITE(numout,*) ' =============== '
887               IF(lwp)WRITE(numout,*)
888               nstop = nstop + 1
889            ENDIF
890         ENDIF
891      ENDIF
892
893   END SUBROUTINE obc_init
894
895#else
896   !!---------------------------------------------------------------------------------
897   !!   Dummy module                                                NO open boundaries
898   !!---------------------------------------------------------------------------------
899CONTAINS
900   SUBROUTINE obc_init      ! Dummy routine
901   END SUBROUTINE obc_init
902#endif
903
904   !!=================================================================================
905END MODULE obcini
Note: See TracBrowser for help on using the repository browser.