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

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

Initial revision

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