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 trunk/NEMO/OPA_SRC/OBC – NEMO

source: trunk/NEMO/OPA_SRC/OBC/obcini.F90 @ 78

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

CT : UPDATE052 : change logical lpXXXobc to lp_obc_XXX for Open Boundaries case

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