source: trunk/NEMOGCM/NEMO/OPA_SRC/OBC/obcini.F90 @ 23

Last change on this file since 23 was 23, checked in by rblod, 12 years ago

Add possibility to choose where to apply the barotropic correction

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