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.
bdydyn.F90 in branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/BDY – NEMO

source: branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/BDY/bdydyn.F90 @ 4409

Last change on this file since 4409 was 3432, checked in by trackstand2, 12 years ago

Merge branch 'ksection_partition'

  • Property svn:keywords set to Id
File size: 10.0 KB
RevLine 
[911]1MODULE bdydyn
[1125]2   !!======================================================================
[911]3   !!                       ***  MODULE  bdydyn  ***
[1125]4   !! Unstructured Open Boundary Cond. :   Flow relaxation scheme on velocities
5   !!======================================================================
6   !! History :  1.0  !  2005-02  (J. Chanut, A. Sellar)  Original code
7   !!             -   !  2007-07  (D. Storkey) Move Flather implementation to separate routine.
8   !!            3.0  !  2008-04  (NEMO team)  add in the reference version
[1502]9   !!            3.2  !  2008-04  (R. Benshila) consider velocity instead of transport
[2528]10   !!            3.3  !  2010-09  (E.O'Dea) modifications for Shelf configurations
11   !!            3.3  !  2010-09  (D.Storkey) add ice boundary conditions
[1125]12   !!----------------------------------------------------------------------
13#if defined key_bdy 
14   !!----------------------------------------------------------------------
15   !!   'key_bdy' :                    Unstructured Open Boundary Condition
16   !!----------------------------------------------------------------------
[911]17   !!   bdy_dyn_frs    : relaxation of velocities on unstructured open boundary
18   !!   bdy_dyn_fla    : Flather condition for barotropic solution
[1125]19   !!----------------------------------------------------------------------
[911]20   USE oce             ! ocean dynamics and tracers
21   USE dom_oce         ! ocean space and time domain
22   USE bdy_oce         ! ocean open boundary conditions
23   USE dynspg_oce      ! for barotropic variables
24   USE phycst          ! physical constants
25   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
26   USE bdytides        ! for tidal harmonic forcing at boundary
[1125]27   USE in_out_manager  !
[3432]28   USE timing,  ONLY: timing_start, timing_stop
[911]29
30   IMPLICIT NONE
31   PRIVATE
32
[1125]33   PUBLIC   bdy_dyn_frs   ! routine called in dynspg_flt (free surface case ONLY)
34# if defined key_dynspg_exp || defined key_dynspg_ts
35   PUBLIC   bdy_dyn_fla   ! routine called in dynspg_flt (free surface case ONLY)
36# endif
[911]37
[3211]38   !! * Control permutation of array indices
39#  include "oce_ftrans.h90"
40#  include "dom_oce_ftrans.h90"
41
[1125]42   !!----------------------------------------------------------------------
[2528]43   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
[1146]44   !! $Id$
[2528]45   !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
[1125]46   !!----------------------------------------------------------------------
[911]47CONTAINS
48
[1125]49   SUBROUTINE bdy_dyn_frs( kt )
50      !!----------------------------------------------------------------------
51      !!                  ***  SUBROUTINE bdy_dyn_frs  ***
52      !!
[911]53      !! ** Purpose : - Apply the Flow Relaxation Scheme for dynamic in the 
54      !!                case of unstructured open boundaries.
55      !!
56      !! References :- Engedahl H., 1995: Use of the flow relaxation scheme in
57      !!               a three-dimensional baroclinic ocean model with realistic
58      !!               topography. Tellus, 365-382.
[1125]59      !!----------------------------------------------------------------------
60      INTEGER, INTENT( in ) ::   kt   ! Main time step counter
[911]61      !!
[2528]62      INTEGER  ::   jb, jk         ! dummy loop indices
63      INTEGER  ::   ii, ij, igrd   ! local integers
64      REAL(wp) ::   zwgt           ! boundary weight
[1125]65      !!----------------------------------------------------------------------
66      !
[2528]67      IF(ln_dyn_frs) THEN       ! If this is false, then this routine does nothing.
68         !
[1125]69         IF( kt == nit000 ) THEN
70            IF(lwp) WRITE(numout,*)
[2528]71            IF(lwp) WRITE(numout,*) 'bdy_dyn_frs : Flow Relaxation Scheme on momentum'
[1125]72            IF(lwp) WRITE(numout,*) '~~~~~~~'
73         ENDIF
74         !
75         igrd = 2                      ! Relaxation of zonal velocity
[2528]76         DO jb = 1, nblen(igrd)
77            DO jk = 1, jpkm1
78               ii   = nbi(jb,igrd)
79               ij   = nbj(jb,igrd)
80               zwgt = nbw(jb,igrd)
81               ua(ii,ij,jk) = ( ua(ii,ij,jk) * ( 1.- zwgt ) + ubdy(jb,jk) * zwgt ) * umask(ii,ij,jk)
[1125]82            END DO
83         END DO
84         !
85         igrd = 3                      ! Relaxation of meridional velocity
[2528]86         DO jb = 1, nblen(igrd)
87            DO jk = 1, jpkm1
88               ii   = nbi(jb,igrd)
89               ij   = nbj(jb,igrd)
90               zwgt = nbw(jb,igrd)
91               va(ii,ij,jk) = ( va(ii,ij,jk) * ( 1.- zwgt ) + vbdy(jb,jk) * zwgt ) * vmask(ii,ij,jk)
[1125]92            END DO
93         END DO
[2528]94         CALL lbc_lnk( ua, 'U', -1. )   ;   CALL lbc_lnk( va, 'V', -1. )   ! Boundary points should be updated
[1125]95         !
[2528]96      ENDIF ! ln_dyn_frs
97      !
[1125]98   END SUBROUTINE bdy_dyn_frs
[911]99
100
[2528]101# if defined   key_dynspg_exp   ||   defined key_dynspg_ts
102   !!----------------------------------------------------------------------
103   !!   'key_dynspg_exp'        OR              explicit sea surface height
104   !!   'key_dynspg_ts '                  split-explicit sea surface height
105   !!----------------------------------------------------------------------
106   
[911]107!! Option to use Flather with dynspg_flt not coded yet...
[2528]108
[1502]109   SUBROUTINE bdy_dyn_fla( pssh )
[1125]110      !!----------------------------------------------------------------------
111      !!                 ***  SUBROUTINE bdy_dyn_fla  ***
112      !!             
[911]113      !!              - Apply Flather boundary conditions on normal barotropic velocities
[2528]114      !!                (ln_dyn_fla=.true. or ln_tides=.true.)
[911]115      !!
116      !! ** WARNINGS about FLATHER implementation:
117      !!1. According to Palma and Matano, 1998 "after ssh" is used.
118      !!   In ROMS and POM implementations, it is "now ssh". In the current
119      !!   implementation (tested only in the EEL-R5 conf.), both cases were unstable.
120      !!   So I use "before ssh" in the following.
121      !!
122      !!2. We assume that the normal ssh gradient at the bdy is zero. As a matter of
123      !!   fact, the model ssh just inside the dynamical boundary is used (the outside 
124      !!   ssh in the code is not updated).
125      !!
[1125]126      !! References:  Flather, R. A., 1976: A tidal model of the northwest European
127      !!              continental shelf. Mem. Soc. R. Sci. Liege, Ser. 6,10, 141-164.     
128      !!----------------------------------------------------------------------
[3432]129      USE exchmod, Only: add_exch, bound_exch_list
[1502]130      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   pssh
131
[2528]132      INTEGER  ::   jb, igrd                         ! dummy loop indices
[1125]133      INTEGER  ::   ii, ij, iim1, iip1, ijm1, ijp1   ! 2D addresses
134      REAL(wp) ::   zcorr                            ! Flather correction
[1502]135      REAL(wp) ::   zforc                            ! temporary scalar
[1125]136      !!----------------------------------------------------------------------
[911]137
[3432]138      CALL timing_start('bdy_dyn_fla')
139
[911]140      ! ---------------------------------!
141      ! Flather boundary conditions     :!
142      ! ---------------------------------!
[3211]143
[2528]144      IF(ln_dyn_fla .OR. ln_tides) THEN ! If these are both false, then this routine does nothing.
[1125]145
[1502]146         ! Fill temporary array with ssh data (here spgu):
[2528]147         igrd = 4
[1502]148         spgu(:,:) = 0.0
[2528]149         DO jb = 1, nblenrim(igrd)
150            ii = nbi(jb,igrd)
151            ij = nbj(jb,igrd)
152            IF( ln_dyn_fla ) spgu(ii, ij) = sshbdy(jb)
153            IF( ln_tides )   spgu(ii, ij) = spgu(ii, ij) + sshtide(jb)
[1502]154         END DO
[1125]155         !
[2528]156         igrd = 5      ! Flather bc on u-velocity;
[1502]157         !             ! remember that flagu=-1 if normal velocity direction is outward
158         !             ! I think we should rather use after ssh ?
[2528]159         DO jb = 1, nblenrim(igrd)
160            ii  = nbi(jb,igrd)
161            ij  = nbj(jb,igrd) 
162            iim1 = ii + MAX( 0, INT( flagu(jb) ) )   ! T pts i-indice inside the boundary
163            iip1 = ii - MIN( 0, INT( flagu(jb) ) )   ! T pts i-indice outside the boundary
[1502]164            !
[2528]165            zcorr = - flagu(jb) * SQRT( grav * hur_e(ii, ij) ) * ( pssh(iim1, ij) - spgu(iip1,ij) )
166            zforc = ubtbdy(jb) + utide(jb)
[3432]167#if defined key_z_first
168            ua_e(ii,ij) = zforc + zcorr * umask_1(ii,ij)
169#else
170            ua_e(ii,ij) = zforc + zcorr * umask(ii,ij,1)
171#endif
[1502]172         END DO
[1125]173         !
[2528]174         igrd = 6      ! Flather bc on v-velocity
[1502]175         !             ! remember that flagv=-1 if normal velocity direction is outward
[2528]176         DO jb = 1, nblenrim(igrd)
177            ii  = nbi(jb,igrd)
178            ij  = nbj(jb,igrd) 
179            ijm1 = ij + MAX( 0, INT( flagv(jb) ) )   ! T pts j-indice inside the boundary
180            ijp1 = ij - MIN( 0, INT( flagv(jb) ) )   ! T pts j-indice outside the boundary
[1502]181            !
[2528]182            zcorr = - flagv(jb) * SQRT( grav * hvr_e(ii, ij) ) * ( pssh(ii, ijm1) - spgu(ii,ijp1) )
183            zforc = vbtbdy(jb) + vtide(jb)
[3432]184#if defined key_z_first
185            va_e(ii,ij) = zforc + zcorr * vmask_1(ii,ij)
186#else
[1502]187            va_e(ii,ij) = zforc + zcorr * vmask(ii,ij,1)
[3432]188#endif
[1502]189         END DO
[3432]190
191         !
192#if defined key_mpp_rkpart
193         ! Is quicker to pack halo-swaps for 2D fields than to do them
194         ! individually
195         CALL add_exch(jpreci,'U',Iminus,Jminus,Iplus,Jplus,ua_e,isgn=-1)
196         CALL add_exch(jpreci,'V',Iminus,Jminus,Iplus,Jplus,va_e,isgn=-1)
197         CALL bound_exch_list()
198#else
[2528]199         CALL lbc_lnk( ua_e, 'U', -1. )   ! Boundary points should be updated
200         CALL lbc_lnk( va_e, 'V', -1. )   !
[3432]201#endif
[1502]202         !
[2528]203      ENDIF ! ln_dyn_fla .or. ln_tides
[1125]204      !
[3432]205      CALL timing_stop('bdy_dyn_fla','section')
206      !
[911]207   END SUBROUTINE bdy_dyn_fla
208#endif
209
210#else
[1125]211   !!----------------------------------------------------------------------
212   !!   Dummy module                   NO Unstruct Open Boundary Conditions
213   !!----------------------------------------------------------------------
[911]214CONTAINS
[1125]215   SUBROUTINE bdy_dyn_frs( kt )      ! Empty routine
216      WRITE(*,*) 'bdy_dyn_frs: You should not have seen this print! error?', kt
[911]217   END SUBROUTINE bdy_dyn_frs
[1559]218   SUBROUTINE bdy_dyn_fla( pssh )    ! Empty routine
219      REAL :: pssh(:,:)
220      WRITE(*,*) 'bdy_dyn_fla: You should not have seen this print! error?', pssh(1,1)
[911]221   END SUBROUTINE bdy_dyn_fla
222#endif
223
[1125]224   !!======================================================================
[911]225END MODULE bdydyn
Note: See TracBrowser for help on using the repository browser.