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

Last change on this file since 247 was 247, checked in by opalod, 19 years ago

CL : Add CVS Header and CeCILL licence information

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