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.
obcspg.F90 in trunk/NEMO/OPA_SRC/OBC – NEMO

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

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

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

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 25.5 KB
Line 
1MODULE obcspg 
2   !!======================================================================
3   !!                       ***  MODULE  obcspg  ***
4   !! Open Boundaries  :   Radiation of barotropic stream function on each
5   !!                      open boundary
6   !!======================================================================
7#if   defined key_obc   &&   defined key_dynspg_rl
8   !!----------------------------------------------------------------------
9   !!   'key_obc'    and                            Open Boundary Condition
10   !!   'key_dynspg_rl'                                 Rigid-Lid
11   !!----------------------------------------------------------------------
12   !!   obc_spg       : call the subroutine for each open boundary
13   !!   obc_spg_east  : radiation of the east open boundary streamfunction
14   !!   obc_spg_west  : radiation of the west open boundary streamfunction
15   !!   obc_spg_north : radiation of the north open boundary streamfunction
16   !!   obc_spg_south : radiation of the south open boundary streamfunction
17   !!----------------------------------------------------------------------
18   !! * Modules used
19   USE oce             ! ocean dynamics and tracers variables
20   USE dom_oce         ! ocean space and time domain variables
21   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
22   USE phycst          ! physical constants
23   USE obc_oce         ! ocean open boundary conditions
24   USE lib_mpp         ! for mppobc
25   USE in_out_manager  ! I/O manager
26
27   IMPLICIT NONE
28   PRIVATE
29
30   !! * Accessibility
31   PUBLIC obc_spg     ! routine called in step.F90 (rigid lid case)
32
33   !! * Module variables
34   INTEGER ::   ji, jj, jk, jnic   ! dummy loop indices
35
36   INTEGER ::      & ! ... boundary space indices
37      nib   = 1,   & ! nib   = boundary point
38      nibm  = 2,   & ! nibm  = 1st interior point
39      nibm2 = 3,   & ! nibm2 = 2nd interior point
40      !              ! ... boundary time indices
41      nit   = 1,   & ! nit    = now
42      nitm  = 2,   & ! nitm   = before
43      nitm2 = 3      ! nitm2  = before-before
44
45   REAL(wp) ::   rtaue  , rtauw  , rtaun  , rtaus  ,   &  !
46                 rtauein, rtauwin, rtaunin, rtausin       !
47
48   !! * Substitutions
49#  include "obc_vectopt_loop_substitute.h90"
50   !!----------------------------------------------------------------------
51   !!   OPA 9.0 , LODYC-IPSL  (2003)
52   !!----------------------------------------------------------------------
53
54CONTAINS
55
56   SUBROUTINE obc_spg ( kt )
57      !!----------------------------------------------------------------------
58      !!                    ***  ROUTINE obc_spg  ***
59      !!
60      !! **  Purpose :
61      !!       Compute now barotropic stream function at the open boundaries.
62      !!       (lp_obc_east, and/or lp_obc_west, and/or lp_obc_north, and/or lp_obc_south).
63      !!       Deduce the correct bsf trend on the open boundaries and isolated
64      !!       coastlines previous to the call to the barotropic solver.
65      !!
66      !! ** Method :
67      !!      In case of open boundaries, there can be a net barotropic flow
68      !!      through the boundaries, hence the potential on the coastlines
69      !!      on each side of the OBC is different.
70      !!      This routine:
71      !!           1. compute the contribution of the isolated coastlines to the
72      !!              rhs of the barotropic equation
73      !!           2. compute the contribution of the OBC to the rhs of the
74      !!              barotropic equation using a radiation equation as explained
75      !!              in the OBC routine.
76      !!
77      !! Reference :
78      !!   Marchesiello P., 1995, these de l'universite J. Fourier, Grenoble, France.
79      !! History :
80      !!        ! 95-03 (J.-M. Molines) Original, SPEM
81      !!        ! 97-07 (G. Madec, J.-M. Molines) additions
82      !!   8.5  ! 02-10 (C. Talandier, A-M. Treguier) F90
83      !!----------------------------------------------------------------------
84      !! * Arguments
85      INTEGER, INTENT( in ) ::   kt
86      !!----------------------------------------------------------------------
87
88      IF( kt == nit000 .OR. ln_rstart ) THEN      ! Initialization
89         ! ... Boundary restoring coefficient
90         rtaue = 2. * rdt / rdpeob
91         rtauw = 2. * rdt / rdpwob
92         rtaun = 2. * rdt / rdpnob
93         rtaus = 2. * rdt / rdpsob
94         ! ... Boundary restoring coefficient for inflow ( all boundaries)
95         rtauein = 2. * rdt / rdpein 
96         rtauwin = 2. * rdt / rdpwin
97         rtaunin = 2. * rdt / rdpnin
98         rtausin = 2. * rdt / rdpsin 
99      ENDIF
100
101      ! right hand side of the barotropic elliptic equation
102      ! ---------------------------------------------------
103
104      ! Isolated coastline contribution to the RHS of the barotropic Eq.
105      gcbob(:,:) = 0.e0
106      DO jnic = 1, nbobc-1
107         gcbob(:,:) = gcbob(:,:) + gcfobc(:,:,jnic) * gcbic(jnic)
108      END DO
109
110      IF( lp_obc_east  )   CALL obc_spg_east ( kt )    ! East open boundary
111
112      IF( lp_obc_west  )   CALL obc_spg_west ( kt )    ! West open boundary
113
114      IF( lp_obc_north )   CALL obc_spg_north( kt )    ! North open boundary
115
116      IF( lp_obc_south )   CALL obc_spg_south( kt )    ! South open boundary
117
118      IF( lk_mpp )   CALL lbc_lnk( gcbob, 'G', 1. )
119 
120   END SUBROUTINE obc_spg
121
122
123   SUBROUTINE obc_spg_east ( kt )
124      !!------------------------------------------------------------------------------
125      !!                ***  SUBROUTINE obc_spg_east  ***
126      !!                 
127      !! ** Purpose :   Apply the radiation algorithm on east OBC stream function.
128      !!      If lfbceast=T , there is no radiation but only fixed OBC
129      !!
130      !!  History :
131      !!         ! 95-03 (J.-M. Molines) Original from SPEM
132      !!         ! 97-07 (G. Madec, J.-M. Molines) additions
133      !!         ! 97-12 (M. Imbard) Mpp adaptation
134      !!         ! 00-06 (J.-M. Molines)
135      !!    8.5  ! 02-10 (C. Talandier, A-M Treguier) F90
136      !!------------------------------------------------------------------------------
137      !! * Arguments
138      INTEGER, INTENT( in ) ::   kt
139
140      !! * Local declarations
141      INTEGER ::   ij
142      REAL(wp) ::   z2dtr, ztau, zin
143      REAL(wp) ::   z05cx, zdt, z4nor2, z2dx, z2dy
144      !!------------------------------------------------------------------------------
145
146      ! 1. First three time steps and more if lfbceast is .TRUE.
147      !    In that case open boundary conditions are FIXED.
148      ! --------------------------------------------------------
149
150      IF( ( kt < nit000 + 3 .AND. .NOT.ln_rstart) .OR. lfbceast ) THEN
151
152         ! 1.1 Fixed barotropic stream function
153         ! ------------------------------------
154         DO jj = nje0m1, nje1 
155            ij = jj -1 + njmpp
156            bsfeob(jj)=bfoe(ij)
157         END DO
158
159      ELSE
160
161      ! 2. Beyond the fourth time step if lfbceast is .FALSE.
162      ! -----------------------------------------------------
163
164         ! 2.1. Barotropic stream function radiation
165         ! ----------------------------------------
166         !
167         !          nibm2      nibm      nib
168         !            |   nibm  |   nib   |///
169         !            |    |    |    |    |///
170         !  jj-line --f----v----f----v----f---
171         !            |         |         |///
172         !            |    |    |    |    |///
173         !         jpieob-2   jpieob-1    jpieob
174         !                 |         |       
175         !              jpieob-1    jpieob     
176         !
177         ! ... radiative conditions plus restoring term toward climatology
178         ! ... Initialize bsfeob to clim in any case, at the first step
179         !     to ensure proper values at the ends of the open line.
180         ! ... Take care that the j indices starts at nje0 (jpjed) and finish
181         !     at nje1 (jpjef) to be sure that jpjefm1 and jpjef are set OK.
182         DO ji = fs_nie0-1, fs_nie1-1 ! Vector opt.
183            DO jj = nje0p1, nje1m2 
184               ij = jj -1 + njmpp
185         ! ... 2* gradi(bsf) (v-point i=nibm, time mean)
186               z2dx = ( bebnd(jj,nibm ,nit) + bebnd(jj,nibm ,nitm2) - 2.*bebnd(jj,nibm2,nitm) ) &
187                      / e1v(ji,jj)
188         ! ... 2* gradj(bsf) (f-point i=nibm, time nitm)
189               z2dy = ( bebnd(jj+1,nibm,nitm) - bebnd(jj-1,nibm,nitm) ) / e2f(ji,jj)
190         ! ... square of the norm of grad(bsf)
191               z4nor2 = z2dx * z2dx + z2dy * z2dy
192         ! ... minus time derivative (leap-frog) at nibm, without / 2 dt
193               zdt = bebnd(jj,nibm,nitm2) - bebnd(jj,nibm,nit)
194         ! ... i-phase speed ratio (bounded by 1) and MASKED!
195               IF( z4nor2 == 0 ) THEN
196                  IF(lwp) WRITE(numout,*)' PB dans obc_spg_east au pt ',jj,' : z4nor=0'
197                  z4nor2 = 0.001
198               ENDIF
199               z05cx = zdt * z2dx / z4nor2 * bmask(ji,jj)
200               z05cx = z05cx / e1v(ji+1,jj)
201               z05cx = min( z05cx, 1. )
202         ! ... z05cx < 0, inflow  zin=0, ztau=1 
203         !          => 0, outflow zin=1, ztau=rtaue
204               zin = sign( 1., z05cx )
205               zin = 0.5*( zin + abs(zin) )
206         ! ... Modification JM:  We maintain a restoring term toward
207         !                   bsfb even in case of inflow
208         ! But restoring is stronger when in flow (10 days) (ztau in set in obcspg.F)
209               ztau = (1.-zin ) * rtauein + zin * rtaue
210               z05cx = z05cx * zin
211         ! ... update bsfn with radiative or climatological bsf (not mask!)
212               bsfeob(jj) = ( ( 1. - z05cx - ztau ) * bebnd(jj,nib ,nitm) + 2.*z05cx  &
213                               * bebnd(jj,nibm,nit) + ztau * bfoe (ij) )              &
214                            / (1. + z05cx)
215            END DO
216         END DO
217
218      ENDIF
219      IF( lk_mpp )   CALL mppobc(bsfeob,jpjed,jpjef,jpieob-1,1,2,jpj)
220
221
222      ! 3. right hand side of the barotropic elliptic equation
223      ! ------------------------------------------------------
224 
225      IF( ( neuler == 0 ) .AND. ( kt == nit000 ) ) THEN
226         z2dtr = 1.0 / rdt
227      ELSE
228         z2dtr = 0.5 / rdt
229      ENDIF
230      DO ji = fs_nie0-1, fs_nie1-1 ! Vector opt.
231         DO jj = nje0m1, nje1 
232            gcbob(ji,jj) = gcbob(ji,jj) - hvr(ji+1,jj) * e2v(ji+1,jj) / e1v(ji+1,jj) &
233                           * ( bsfeob(jj) - bsfb(ji+1,jj) ) * z2dtr * bmask(ji,jj) 
234         END DO
235      END DO
236
237   END SUBROUTINE obc_spg_east
238
239   SUBROUTINE obc_spg_west ( kt )
240      !!------------------------------------------------------------------------------
241      !!                  ***  SUBROUTINE obc_spg_west  ***
242      !!                   
243      !! ** Purpose :
244      !!      Apply the radiation algorithm on west OBC stream function.
245      !!      If the logical lfbcwest is .TRUE., there is no radiation but only fixed OBC
246      !!
247      !!  History :
248      !!         ! 95-03 (J.-M. Molines) Original from SPEM
249      !!         ! 97-07 (G. Madec, J.-M. Molines) additions
250      !!         ! 97-12 (M. Imbard) Mpp adaptation
251      !!         ! 00-06 (J.-M. Molines)
252      !!    8.5  ! 02-10 (C. Talandier, A-M Treguier) F90
253      !!------------------------------------------------------------------------------
254      !! * Arguments
255      INTEGER, INTENT( in ) ::   kt
256
257      !! * Local declarations
258      INTEGER ::   ij
259
260      REAL(wp) ::   z2dtr, ztau, zin
261      REAL(wp) ::   z05cx, zdt, z4nor2, z2dx, z2dy
262
263      !!------------------------------------------------------------------------------
264      !!  OPA 8.5, LODYC-IPSL (2002)
265      !!------------------------------------------------------------------------------
266
267      ! 1. First three time steps and more if lfbcwest is .TRUE.
268      !    In that case open boundary conditions are FIXED.
269      ! --------------------------------------------------------
270
271      IF( ( kt < nit000 + 3 .AND. .NOT.ln_rstart ) .OR. lfbcwest ) THEN
272
273         ! 1.1 Fixed barotropic stream function
274         ! ------------------------------------
275         DO jj = njw0m1, njw1
276            ij = jj -1 + njmpp
277            bsfwob(jj)=bfow(ij)
278         END DO
279
280      ELSE
281
282      ! 2. Beyond the fourth time step if lfbcwest is .FALSE.
283      ! -----------------------------------------------------
284
285         ! 2.1. Barotropic stream function radiation
286         ! ----------------------------------------
287         !
288         !         nib       nibm     nibm2
289         !       ///|   nib   |   nibm  |
290         !       ///|    |    |    |    |
291         !       ---f----v----f----v----f-- jj-line
292         !       ///|         |         |
293         !       ///|    |    |    |    |
294         !        jpiwob    jpiwob+1    jpiwob+2
295         !               |         |       
296         !             jpiwob+1    jpiwob+2     
297         !
298         ! ... radiative conditions plus restoring term toward climatology
299         ! ... Initialize bsfwob to clim in any case, at the first step
300         !     to ensure proper values at the ends of the open line.
301         ! ... Take care that the j indices starts at njw0 (jpjwd) and finish
302         ! ... at njw1 (jpjwf) to be sure that jpjwfm1 and jpjwf are set OK.
303         DO ji = fs_niw0+1, fs_niw1+1 ! Vector opt.
304            DO jj = njw0p1, njw1m2
305               ij = jj -1 + njmpp
306         ! ... 2* gradi(bsf) (v-point i=nibm, time mean)
307               z2dx = ( -  bwbnd(jj,nibm ,nit ) - bwbnd(jj,nibm ,nitm2) + 2.*bwbnd(jj,nibm2,nitm ) ) &
308                      / e1v(ji+1,jj)
309         ! ... 2* gradj(bsf) (f-point i=nibm, time nitm)
310               z2dy = ( bwbnd(jj+1,nibm,nitm) - bwbnd(jj-1,nibm,nitm) ) / e2f(ji,jj)
311         ! ... square of the norm of grad(bsf)
312               z4nor2 = z2dx * z2dx + z2dy * z2dy
313         ! ... minus time derivative (leap-frog) at nibm, without / 2 dt
314               zdt = bwbnd(jj,nibm,nitm2) - bwbnd(jj,nibm,nit)
315         ! ... i-phase speed ratio (bounded by 1) and MASKED!
316               IF( z4nor2 == 0 ) THEN
317                  IF(lwp) WRITE(numout,*)' PB dans obc_spg_west au pt ',jj,' : z4nor =0'
318                  z4nor2=0.0001
319               ENDIF
320               z05cx = zdt * z2dx / z4nor2 * bmask(ji,jj)
321               z05cx = z05cx / e1v(ji,jj)
322               z05cx = max( z05cx, -1. )
323         ! ... z05cx => 0, inflow  zin=0, ztau=1 
324         !           <  0, outflow zin=1, ztau=rtauw
325               zin = sign( 1., -1. * z05cx )
326               zin = 0.5*( zin + abs(zin) )
327               ztau = (1.-zin )*rtauwin + zin * rtauw
328               z05cx = z05cx * zin
329         !  ... update bsfn with radiative or climatological bsf (not mask!)
330               bsfwob(jj) = ( ( 1. + z05cx - ztau ) * bwbnd(jj,nib ,nitm) - 2.*z05cx &
331                               * bwbnd(jj,nibm,nit) + ztau * bfow (ij) )             &
332                            / (1. - z05cx)
333            END DO
334         END DO
335
336      ENDIF
337      IF( lk_mpp )   CALL mppobc(bsfwob,jpjwd,jpjwf,jpiwob+1,1,2,jpj) 
338
339
340      ! 3. right hand side of the barotropic elliptic equation
341      ! -------------------------------------------------------
342
343      IF( ( neuler == 0 ) .AND. ( kt == nit000 ) ) THEN
344         z2dtr = 1.0 / rdt
345      ELSE
346         z2dtr = 0.5 / rdt
347      ENDIF
348      DO ji = fs_niw0+1, fs_niw1+1 ! Vector opt.
349         DO jj = njw0m1, njw1
350            gcbob(ji,jj) = gcbob(ji,jj) - hvr(ji,jj) * e2v(ji,jj) / e1v(ji,jj)     &
351                           * ( bsfwob(jj) - bsfb(ji-1,jj) ) * z2dtr * bmask(ji,jj)
352         END DO
353      END DO
354
355   END SUBROUTINE obc_spg_west
356
357   SUBROUTINE obc_spg_north ( kt )
358      !!------------------------------------------------------------------------------
359      !!                 ***  SUBROUTINE obc_spg_north  ***
360      !!
361      !! ** Purpose :   Apply the radiation algorithm on north OBC stream function.
362      !!      If lfbcnorth=T, there is no radiation but only fixed OBC
363      !!
364      !!  History :
365      !!         ! 95-03 (J.-M. Molines) Original from SPEM
366      !!         ! 97-07 (G. Madec, J.-M. Molines) additions
367      !!         ! 97-12 (M. Imbard) Mpp adaptation
368      !!         ! 00-06 (J.-M. Molines)
369      !!    8.5  ! 02-10 (C. Talandier, A-M Treguier) F90
370      !!------------------------------------------------------------------------------
371      !! * Arguments
372      INTEGER, INTENT( in ) ::   kt
373
374      !! * Local declarations
375      INTEGER ::   ii
376      REAL(wp) ::   z2dtr, ztau, zin
377      REAL(wp) ::   z05cx, zdt, z4nor2, z2dx, z2dy
378      !!------------------------------------------------------------------------------
379
380      ! 1. First three time steps and more if lfbcnorth is .TRUE.
381      !    In that case open boundary conditions are FIXED.
382      ! --------------------------------------------------------
383
384      IF( ( kt < nit000 + 3 .AND. .NOT.ln_rstart ) .OR. lfbcnorth ) THEN
385
386         ! 1.1 Fixed barotropic stream function
387         ! ------------------------------------
388         DO ji = nin0m1, nin1
389            ii = ji -1 + nimpp
390            bsfnob(ji)=bfon(ii)
391         END DO
392
393      ELSE     
394
395      ! 2. Beyond the fourth time step if lfbcnorth is .FALSE.
396      ! -----------------------------------------------------
397
398         ! 2.1. Barotropic stream function radiation
399         ! -----------------------------------------
400         !
401         !           ji-row
402         !             |
403         !        ////////////
404         !        ////////////
405         !   nib  -----f------  jpjnob
406         !             |   
407         !      nib--  u   ---- jpjnob
408         !             |       
409         !  nibm  -----f-----   jpjnob-1
410         !             |       
411         !     nibm--  u   ---- jpjnob-1
412         !             |       
413         !  nibm2 -----f-----   jpjnob-2
414         !             |
415         !
416         ! ... radiative conditions plus restoring term toward climatology
417         ! ... z05cx is always the cross boundary phase velocity
418         ! ... Initialize bsfnob to clim in any case, at the first step
419         !     to ensure proper values at the ends of the open line.
420         ! ... Take care that the i indices starts at nin0 (jpind) and finish
421         ! ... at nin1 (jpinf) to be sure that jpinfm1 and jpinf are set OK.
422         DO jj = fs_njn0-1, fs_njn1-1 ! Vector opt.
423            DO ji = nin0p1, nin1m2
424               ii = ji -1 + nimpp
425         ! ... 2* gradj(bsf) (u-point i=nibm, time mean)
426               z2dx = ( bnbnd(ji,nibm ,nit) + bnbnd(ji,nibm ,nitm2) - 2.*bnbnd(ji,nibm2,nitm) ) &
427                      / e2u(ji,jj)
428         ! ... 2* gradi(bsf) (f-point i=nibm, time nitm)
429               z2dy = ( bnbnd(ji+1,nibm,nitm) - bnbnd(ji-1,nibm,nitm) ) / e1f(ji,jj)
430         ! ... square of the norm of grad(bsf)
431               z4nor2 = z2dx * z2dx + z2dy * z2dy
432         ! ... minus time derivative (leap-frog) at nibm, without / 2 dt
433               zdt = bnbnd(ji,nibm,nitm2) - bnbnd(ji,nibm,nit)
434         ! ... j-phase speed ratio (bounded by 1) and MASKED!
435               IF( z4nor2 == 0 ) THEN
436                  IF(lwp) WRITE(numout,*)' PB dans obc_spg_north au pt',ji,' : z4nor =0'
437               ENDIF
438               z05cx = zdt * z2dx / z4nor2 * bmask(ji,jj)
439               z05cx = z05cx / e2u(ji,jj+1)
440               z05cx = min( z05cx, 1. )
441         ! ... z05cx < 0, inflow  zin=0, ztau=1
442         !          => 0, outflow zin=1, ztau=rtaun
443               zin = sign( 1., z05cx )
444               zin = 0.5*( zin + abs(zin) )
445               ztau = (1.-zin ) * rtaunin + zin * rtaun
446               z05cx = z05cx * zin
447         ! ... update bsfn with radiative or climatological bsf (not mask!)
448               bsfnob(ji) = ( ( 1. - z05cx - ztau ) * bnbnd(ji,nib ,nitm) + 2.*z05cx  &
449                               * bnbnd(ji,nibm,nit) + ztau * bfon (ii) )              &
450                            / (1. + z05cx)
451            END DO
452         END DO
453
454      ENDIF
455      IF( lk_mpp )   CALL mppobc(bsfnob,jpind,jpinf,jpjnob-1,1,1,jpi)
456
457
458      ! 3. right hand side of the barotropic elliptic equation
459      !-------------------------------------------------------
460
461      IF( ( neuler == 0 ) .AND. ( kt == nit000 ) ) THEN
462         z2dtr = 1.0 / rdt
463      ELSE
464         z2dtr = 0.5 / rdt
465      ENDIF
466      DO jj = fs_njn0-1, fs_njn1-1 ! Vector opt.
467         DO ji = nin0m1, nin1
468            gcbob(ji,jj) = gcbob(ji,jj) - hur(ji,jj+1) *  e1u(ji,jj+1) / e2u(ji,jj+1) &
469                           * ( bsfnob(ji) - bsfb(ji,jj+1) ) * z2dtr * bmask(ji,jj)
470         END DO
471      END DO
472
473   END SUBROUTINE obc_spg_north
474
475
476   SUBROUTINE obc_spg_south ( kt )
477      !!------------------------------------------------------------------------------
478      !!                  ***  SUBROUTINE obc_spg_south  ***
479      !!               
480      !! ** Purpose :   Apply the radiation algorithm on south OBC stream function.
481      !!      If lfbcsouth=T, there is no radiation but only fixed OBC
482      !!
483      !!  History :
484      !!         ! 95-03 (J.-M. Molines) Original from SPEM
485      !!         ! 97-07 (G. Madec, J.-M. Molines) additions
486      !!         ! 97-12 (M. Imbard) Mpp adaptation
487      !!         ! 00-06 (J.-M. Molines)
488      !!    8.5  ! 02-10 (C. Talandier, A-M Treguier) F90
489      !!------------------------------------------------------------------------------
490      !! * Arguments
491      INTEGER, INTENT(in) ::   kt  ! ocean time-step index
492
493      !! * Local declarations
494      INTEGER  ::   ii             ! temporary integers
495      REAL(wp) ::   &
496         z2dtr, ztau, zin   ,   &  ! temporary scalars
497         z05cx, zdt , z4nor2,   &  !    "         "
498         z2dx , z2dy               !    "         "
499      !!------------------------------------------------------------------------------
500
501      ! 1. First three time steps and more if lfbcsouth is .TRUE.
502      !    In that case open boundary conditions are FIXED.
503      ! --------------------------------------------------------
504     
505      IF( ( kt < nit000 + 3 .AND. .NOT.ln_rstart ) .OR. lfbcsouth ) THEN
506
507         ! 1.1 Fixed barotropic stream function
508         ! ------------------------------------
509         DO ji = nis0m1, nis1 
510            ii = ji -1 + nimpp 
511            bsfsob(ji)=bfos(ii)
512         END DO
513
514      ELSE
515
516      ! 2. Beyond the fourth time step if lfbcsouth is .FALSE.
517      ! -----------------------------------------------------
518
519         ! 2.1. Barotropic stream function radiation
520         ! -----------------------------------------
521         !
522         !           ji-row
523         !             |
524         ! nibm2  -----f------  jpjsob + 2
525         !             |   
526         !   nibm  --  u  ----- jpjsob + 2
527         !             |       
528         !  nibm  -----f-----   jpjsob + 1
529         !             |       
530         !    nib  --  u  ----- jpjsob + 1
531         !             |       
532         !    nib -----f-----   jpjsob
533         !        ///////////     
534         !        ///////////
535         !
536         ! ... radiative conditions plus restoring term toward climatology
537         ! ... z05cx is always the cross boundary phase velocity
538         ! ... Initialize bsfsob to clim in any case, at the first step
539         !     to ensure proper values at the ends of the open line.
540         ! ... Take care that the i indices starts at nis0 (jpisd) and finish
541         ! ... at nis1 (jpisf) to be sure that jpisfm1 and jpisf are set OK.
542         DO jj = fs_njs0+1, fs_njs1+1 ! Vector opt.
543            DO ji = nis0p1, nis1m2
544               ii = ji -1 + nimpp
545         ! ... 2* gradj(bsf) (u-point i=nibm, time mean)
546               z2dx = ( - bsbnd(ji,nibm ,nit) - bsbnd(ji,nibm ,nitm2) + 2.*bsbnd(ji,nibm2,nitm) ) &
547                      / e2u(ji,jj+1)
548         ! ... 2* gradi(bsf) (f-point i=nibm, time nitm)
549               z2dy = ( bsbnd(ji+1,nibm,nitm) - bsbnd(ji-1,nibm,nitm) ) / e1f(ji,jj)
550         ! ... square of the norm of grad(bsf)
551               z4nor2 = z2dx * z2dx + z2dy * z2dy
552         ! ... minus time derivative (leap-frog) at nibm, without / 2 dt
553               zdt = bsbnd(ji,nibm,nitm2) - bsbnd(ji,nibm,nit)
554         ! ... j-phase speed ratio (bounded by -1) and MASKED!
555               IF( z4nor2 == 0 ) THEN
556                  IF(lwp) WRITE(numout,*)' PB dans obc_spg_south au pt ',ji,' : z4nor =0'
557               ENDIF
558               z05cx = zdt * z2dx / z4nor2 * bmask(ji,jj)
559               z05cx = z05cx / e2u(ji,jj)
560               z05cx = max( z05cx, -1. )
561         ! ... z05cx => 0, inflow  zin=0, ztau=1
562         !           <  0, outflow zin=1, ztau=rtaus
563               zin = sign( 1., -1. * z05cx )
564               zin = 0.5*( zin + abs(zin) )
565               ztau = (1.-zin ) *rtausin  + zin * rtaus
566               z05cx = z05cx * zin
567         ! ... update bsfn with radiative or climatological bsf (not mask!)
568               bsfsob(ji) = ( ( 1. + z05cx - ztau ) * bsbnd(ji,nib ,nitm) - 2.*z05cx  & 
569                               * bsbnd(ji,nibm,nit) + ztau * bfos (ii) )              &
570                            / (1. - z05cx)
571            END DO
572         END DO
573
574      ENDIF
575      IF( lk_mpp )   CALL mppobc(bsfsob,jpisd,jpisf,jpjsob+1,1,1,jpi)
576
577 
578      ! 3. right hand side of the barotropic elliptic equation
579      ! -------------------------------------------------------
580
581      IF( ( neuler == 0 ) .AND. ( kt == nit000 ) ) THEN
582         z2dtr = 1.0 / rdt
583      ELSE
584         z2dtr = 0.5 / rdt
585      ENDIF
586      DO jj = fs_njs0+1, fs_njs1+1 ! Vector opt.
587         DO ji = nis0m1, nis1 
588            gcbob(ji,jj) = gcbob(ji,jj) - hur(ji,jj) * e1u(ji,jj) / e2u(ji,jj) &
589                           * ( bsfsob(ji) - bsfb(ji,jj-1) ) * z2dtr * bmask(ji,jj)
590         END DO
591      END DO
592
593   END SUBROUTINE obc_spg_south
594
595#else
596   !!----------------------------------------------------------------------
597   !!   Default case :                                         Empty module
598   !!----------------------------------------------------------------------
599CONTAINS
600   SUBROUTINE obc_spg( kt )        ! Empty routine
601      INTEGER, INTENT( in ) :: kt
602      WRITE(*,*) 'obc_spg: You should not have seen this print! error?', kt
603   END SUBROUTINE obc_spg
604#endif
605
606   !!======================================================================
607END MODULE obcspg
Note: See TracBrowser for help on using the repository browser.