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.
Changeset 2800 for branches/2011/UKMO_MERCATOR_obc_bdy_merge/NEMOGCM/NEMO/OPA_SRC/OBC/obcdyn2d.F90 – NEMO

Ignore:
Timestamp:
2011-07-13T15:31:05+02:00 (13 years ago)
Author:
davestorkey
Message:
  1. Application of boundary conditions to barotropic and baroclinic velocities clearly separated.
  2. Option to input full velocities in boundary data (default expects barotropic and baroclinic velocities separately).
  3. Option to use initial conditions as boundary conditions coded.
File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/2011/UKMO_MERCATOR_obc_bdy_merge/NEMOGCM/NEMO/OPA_SRC/OBC/obcdyn2d.F90

    r2797 r2800  
    3939CONTAINS 
    4040 
    41    SUBROUTINE obc_dyn2d( kt, pssh ) 
     41   SUBROUTINE obc_dyn2d( kt ) 
    4242      !!---------------------------------------------------------------------- 
    4343      !!                  ***  SUBROUTINE obc_dyn2d  *** 
     
    4747      !!---------------------------------------------------------------------- 
    4848      INTEGER,                      INTENT(in) ::   kt   ! Main time step counter 
    49       REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   pssh 
    5049      !! 
    5150      INTEGER                                  ::   ib_obc ! Loop counter 
     
    5655         CASE(jp_none) 
    5756            CYCLE 
     57         CASE(jp_frs) 
     58            CALL obc_dyn2d_frs( idx_obc(ib_obc), dta_obc(ib_obc), kt ) 
    5859         CASE(jp_flather) 
    59             CALL obc_dyn2d_fla( idx_obc(ib_obc), dta_obc(ib_obc), pssh ) 
     60            CALL obc_dyn2d_fla( idx_obc(ib_obc), dta_obc(ib_obc) ) 
    6061         CASE DEFAULT 
    6162            CALL ctl_stop( 'obc_dyn3d : unrecognised option for open boundaries for barotropic variables' ) 
     
    6566   END SUBROUTINE obc_dyn2d 
    6667 
    67    SUBROUTINE obc_dyn2d_fla( idx, dta, pssh ) 
    68       !!---------------------------------------------------------------------- 
    69       !!                 ***  SUBROUTINE obc_dyn_fla  *** 
     68   SUBROUTINE obc_dyn2d_frs( idx, dta, kt ) 
     69      !!---------------------------------------------------------------------- 
     70      !!                  ***  SUBROUTINE obc_dyn2d_frs  *** 
     71      !! 
     72      !! ** Purpose : - Apply the Flow Relaxation Scheme for barotropic velocities 
     73      !!                at open boundaries. 
     74      !! 
     75      !! References :- Engedahl H., 1995: Use of the flow relaxation scheme in  
     76      !!               a three-dimensional baroclinic ocean model with realistic 
     77      !!               topography. Tellus, 365-382. 
     78      !!---------------------------------------------------------------------- 
     79      INTEGER,         INTENT(in) :: kt 
     80      TYPE(OBC_INDEX), INTENT(in) ::   idx  ! OBC indices 
     81      TYPE(OBC_DATA),  INTENT(in) ::   dta  ! OBC external data 
     82      !! 
     83      INTEGER  ::   jb, jk         ! dummy loop indices 
     84      INTEGER  ::   ii, ij, igrd   ! local integers 
     85      REAL(wp) ::   zwgt           ! boundary weight 
     86      !!---------------------------------------------------------------------- 
     87      ! 
     88      ! 
     89      igrd = 2                      ! Relaxation of zonal velocity 
     90      DO jb = 1, idx%nblen(igrd) 
     91         ii   = idx%nbi(jb,igrd) 
     92         ij   = idx%nbj(jb,igrd) 
     93         zwgt = idx%nbw(jb,igrd) 
     94         pu2d(ii,ij) = ( pu2d(ii,ij) + zwgt * ( dta%u2d(jb) - pu2d(ii,ij) ) ) * umask(ii,ij,1) 
     95      END DO 
     96      ! 
     97      igrd = 3                      ! Relaxation of meridional velocity 
     98      DO jb = 1, idx%nblen(igrd) 
     99         ii   = idx%nbi(jb,igrd) 
     100         ij   = idx%nbj(jb,igrd) 
     101         zwgt = idx%nbw(jb,igrd) 
     102         pv2d(ii,ij) = ( pv2d(ii,ij) + zwgt * ( dta%v2d(jb) - pv2d(ii,ij) ) ) * vmask(ii,ij,1) 
     103      END DO  
     104      CALL lbc_lnk( pu2d, 'U', -1. )  
     105      CALL lbc_lnk( pv2d, 'V', -1. )   ! Boundary points should be updated 
     106      ! 
     107 
     108   END SUBROUTINE obc_dyn2d_frs 
     109 
     110 
     111   SUBROUTINE obc_dyn2d_fla( idx, dta ) 
     112      !!---------------------------------------------------------------------- 
     113      !!                 ***  SUBROUTINE obc_dyn2d_fla  *** 
    70114      !!              
    71115      !!              - Apply Flather boundary conditions on normal barotropic velocities  
     
    100144      
    101145!!!!!!!!!!!! SOME WORK TO DO ON THE TIDES HERE !!!!!!!!!!!!! 
     146 
     147!!! REPLACE spgu with nemo_wrk work space 
    102148 
    103149      ! Fill temporary array with ssh data (here spgu): 
     
    120166         iip1 = ii - MIN( 0, INT( idx%flagu(jb) ) )   ! T pts i-indice outside the boundary  
    121167         ! 
    122          zcorr = - idx%flagu(jb) * SQRT( grav * hur_e(ii, ij) ) * ( pssh(iim1, ij) - spgu(iip1,ij) ) 
     168         zcorr = - idx%flagu(jb) * SQRT( grav * phur(ii, ij) ) * ( pssh(iim1, ij) - spgu(iip1,ij) ) 
    123169!!$         zforc = dta%u2d(jb) + utide(jb) 
    124170         zforc = dta%u2d(jb) 
    125          ua_e(ii,ij) = zforc + zcorr * umask(ii,ij,1)  
     171         pu2d(ii,ij) = zforc + zcorr * umask(ii,ij,1)  
    126172      END DO 
    127173      ! 
     
    134180         ijp1 = ij - MIN( 0, INT( idx%flagv(jb) ) )   ! T pts j-indice outside the boundary  
    135181         ! 
    136          zcorr = - idx%flagv(jb) * SQRT( grav * hvr_e(ii, ij) ) * ( pssh(ii, ijm1) - spgu(ii,ijp1) ) 
     182         zcorr = - idx%flagv(jb) * SQRT( grav * phvr(ii, ij) ) * ( pssh(ii, ijm1) - spgu(ii,ijp1) ) 
    137183!!$         zforc = dta%v2d(jb) + vtide(jb) 
    138184         zforc = dta%v2d(jb) 
    139          va_e(ii,ij) = zforc + zcorr * vmask(ii,ij,1) 
    140       END DO 
    141       CALL lbc_lnk( ua_e, 'U', -1. )   ! Boundary points should be updated 
    142       CALL lbc_lnk( va_e, 'V', -1. )   ! 
     185         pv2d(ii,ij) = zforc + zcorr * vmask(ii,ij,1) 
     186      END DO 
     187      CALL lbc_lnk( pu2d, 'U', -1. )   ! Boundary points should be updated 
     188      CALL lbc_lnk( pv2d, 'V', -1. )   ! 
    143189      ! 
    144190      ! 
Note: See TracChangeset for help on using the changeset viewer.