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 1125 for trunk/NEMO/OPA_SRC/BDY/bdydyn.F90 – NEMO

Ignore:
Timestamp:
2008-06-23T11:05:02+02:00 (16 years ago)
Author:
ctlod
Message:

trunk: BDY package code review (coding rules), see ticket: #214

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/NEMO/OPA_SRC/BDY/bdydyn.F90

    r911 r1125  
    11MODULE bdydyn 
    2    !!================================================================================= 
     2   !!====================================================================== 
    33   !!                       ***  MODULE  bdydyn  *** 
    4    !! Ocean dynamics:   Flow relaxation scheme of velocities on unstruc. open boundary 
    5    !!================================================================================= 
    6  
    7    !!--------------------------------------------------------------------------------- 
     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 
     9   !!---------------------------------------------------------------------- 
     10#if defined key_bdy  
     11   !!---------------------------------------------------------------------- 
     12   !!   'key_bdy' :                    Unstructured Open Boundary Condition 
     13   !!---------------------------------------------------------------------- 
    814   !!   bdy_dyn_frs    : relaxation of velocities on unstructured open boundary 
    915   !!   bdy_dyn_fla    : Flather condition for barotropic solution 
    10    !!--------------------------------------------------------------------------------- 
    11 #if defined key_bdy || defined key_bdy_tides 
    12    !!--------------------------------------------------------------------------------- 
    13    !! * Modules used 
     16   !!---------------------------------------------------------------------- 
    1417   USE oce             ! ocean dynamics and tracers  
    1518   USE dom_oce         ! ocean space and time domain 
     
    1922   USE lbclnk          ! ocean lateral boundary conditions (or mpp link) 
    2023   USE bdytides        ! for tidal harmonic forcing at boundary 
    21    USE in_out_manager 
     24   USE in_out_manager  ! 
    2225 
    2326   IMPLICIT NONE 
    2427   PRIVATE 
    2528 
    26    !! * Accessibility 
    27    PUBLIC bdy_dyn_frs  ! routine called in dynspg_flt (free surface case ONLY) 
    28 #if defined key_dynspg_exp || defined key_dynspg_ts 
    29    PUBLIC bdy_dyn_fla  ! routine called in dynspg_flt (free surface case ONLY) 
    30 #endif 
     29   PUBLIC   bdy_dyn_frs   ! routine called in dynspg_flt (free surface case ONLY) 
     30# if defined key_dynspg_exp || defined key_dynspg_ts 
     31   PUBLIC   bdy_dyn_fla   ! routine called in dynspg_flt (free surface case ONLY) 
     32# endif 
    3133 
    32    !!--------------------------------------------------------------------------------- 
     34   !!---------------------------------------------------------------------- 
     35   !! NEMO/OPA 3.0 , LOCEAN-IPSL (2008)  
     36   !! $Id: $  
     37   !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt) 
     38   !!---------------------------------------------------------------------- 
    3339 
    3440CONTAINS 
    3541 
    36    SUBROUTINE bdy_dyn_frs ( kt ) 
    37       !!------------------------------------------------------------------------------ 
    38       !!                      SUBROUTINE bdy_dyn_frs 
    39       !!                     *********************** 
     42   SUBROUTINE bdy_dyn_frs( kt ) 
     43      !!---------------------------------------------------------------------- 
     44      !!                  ***  SUBROUTINE bdy_dyn_frs  *** 
     45      !! 
    4046      !! ** Purpose : - Apply the Flow Relaxation Scheme for dynamic in the   
    4147      !!                case of unstructured open boundaries. 
     
    4450      !!               a three-dimensional baroclinic ocean model with realistic 
    4551      !!               topography. Tellus, 365-382. 
     52      !!---------------------------------------------------------------------- 
     53      INTEGER, INTENT( in ) ::   kt   ! Main time step counter 
    4654      !! 
    47       !! History : 
    48       !!   9.0  !  05-02 (J. Chanut, A. Sellar) Original 
    49       !!        !  07-07 (D. Storkey) Move Flather implementation to separate routine. 
    50       !!------------------------------------------------------------------------------ 
    51       !! * Arguments 
    52       INTEGER, INTENT( in ) ::   kt    ! Main time step counter 
     55      INTEGER  ::   ib, ik, igrd      ! dummy loop indices 
     56      INTEGER  ::   ii, ij            ! 2D addresses 
     57      REAL(wp) ::   zwgt              ! boundary weight 
     58      !!---------------------------------------------------------------------- 
     59      ! 
     60      IF(ln_bdy_dyn_frs) THEN ! If this is false, then this routine does nothing.  
    5361 
    54       !! * Local declarations 
    55       INTEGER ::   jb, jk, jgrd        ! dummy loop indices 
    56       INTEGER ::   ii, ij              ! 2D addresses 
    57       REAL(wp) :: zwgt                 ! boundary weight 
    58  
    59       !!------------------------------------------------------------------------------ 
    60       !!  OPA 9.0, LODYC-IPSL (2003) 
    61       !!------------------------------------------------------------------------------ 
    62  
    63       ! ---------------------------! 
    64       ! FRS on the total velocity :! 
    65       ! ---------------------------!   
    66    
    67       jgrd=2 !: Relaxation of zonal velocity 
    68       DO jb = 1, nblen(jgrd) 
    69         DO jk = 1, jpkm1 
    70           ii = nbi(jb,jgrd) 
    71           ij = nbj(jb,jgrd) 
    72           zwgt = nbw(jb,jgrd) 
    73  
    74           ua(ii,ij,jk) = ( ua(ii,ij,jk)*(1.-zwgt) + ubdy(jb,jk)*zwgt ) & 
    75                                                         * umask(ii,ij,jk) 
    76         END DO 
    77       END DO 
    78  
    79       jgrd=3 !: Relaxation of meridional velocity 
    80       DO jb = 1, nblen(jgrd) 
    81         DO jk = 1, jpkm1 
    82           ii = nbi(jb,jgrd) 
    83           ij = nbj(jb,jgrd) 
    84           zwgt = nbw(jb,jgrd) 
    85  
    86           va(ii,ij,jk) = ( va(ii,ij,jk)*(1.-zwgt) + vbdy(jb,jk)*zwgt ) & 
    87                                                         * vmask(ii,ij,jk) 
    88         END DO 
    89       END DO  
    90  
    91       CALL lbc_lnk( ua, 'U', 1. ) ! Boundary points should be updated 
    92       CALL lbc_lnk( va, 'V', 1. ) ! 
     62         IF( kt == nit000 ) THEN 
     63            IF(lwp) WRITE(numout,*) 
     64            IF(lwp) WRITE(numout,*) 'bdy_dyn : Flow Relaxation Scheme on momentum' 
     65            IF(lwp) WRITE(numout,*) '~~~~~~~' 
     66         ENDIF 
     67         ! 
     68         igrd = 2                      ! Relaxation of zonal velocity 
     69         DO ib = 1, nblen(igrd) 
     70            DO ik = 1, jpkm1 
     71               ii = nbi(ib,igrd) 
     72               ij = nbj(ib,igrd) 
     73               zwgt = nbw(ib,igrd) 
     74               ua(ii,ij,ik) = ( ua(ii,ij,ik) * ( 1.- zwgt ) + ubdy(ib,ik) * zwgt ) * umask(ii,ij,ik) 
     75            END DO 
     76         END DO 
     77         ! 
     78         igrd = 3                      ! Relaxation of meridional velocity 
     79         DO ib = 1, nblen(igrd) 
     80            DO ik = 1, jpkm1 
     81               ii = nbi(ib,igrd) 
     82               ij = nbj(ib,igrd) 
     83               zwgt = nbw(ib,igrd) 
     84               va(ii,ij,ik) = ( va(ii,ij,ik) * ( 1.- zwgt ) + vbdy(ib,ik) * zwgt ) * vmask(ii,ij,ik) 
     85            END DO 
     86         END DO  
     87         ! 
     88         CALL lbc_lnk( ua, 'U', 1. )   ! Boundary points should be updated 
     89         CALL lbc_lnk( va, 'V', 1. )   ! 
     90         ! 
     91      ENDIF ! ln_bdy_dyn_frs 
    9392 
    9493   END SUBROUTINE bdy_dyn_frs 
     94 
    9595 
    9696#if defined key_dynspg_exp || defined key_dynspg_ts 
    9797!! Option to use Flather with dynspg_flt not coded yet... 
    9898   SUBROUTINE bdy_dyn_fla 
    99       !!------------------------------------------------------------------------------ 
    100       !!                      SUBROUTINE bdy_dyn_fla 
    101       !!                     *********************** 
     99      !!---------------------------------------------------------------------- 
     100      !!                 ***  SUBROUTINE bdy_dyn_fla  *** 
     101      !!              
    102102      !!              - Apply Flather boundary conditions on normal barotropic velocities  
    103       !!                (ln_bdy_fla=.true.) 
     103      !!                (ln_bdy_dyn_fla=.true. or ln_bdy_tides=.true.) 
    104104      !! 
    105105      !! ** WARNINGS about FLATHER implementation: 
     
    113113      !!   ssh in the code is not updated). 
    114114      !! 
    115       !!             - Flather, R. A., 1976: A tidal model of the northwest European 
    116       !!               continental shelf. Mem. Soc. R. Sci. Liege, Ser. 6,10, 141-164.      
    117       !! History : 
    118       !!   9.0  !  05-02 (J. Chanut, A. Sellar) Original 
    119       !!        !  07-07 (D. Storkey) Flather algorithm in separate routine. 
    120       !!------------------------------------------------------------------------------ 
    121       !! * Local declarations 
    122       INTEGER ::   jb, jk, jgrd, ji, jj, jim1, jip1, jjm1, jjp1   ! dummy loop indices 
    123       INTEGER ::   ii, ij                    ! 2D addresses 
    124       REAL(wp) ::  corr ! Flather correction 
    125       REAL(wp) :: zwgt                       ! boundary weight 
    126       REAL(wp) :: elapsed 
    127  
    128       !!------------------------------------------------------------------------------ 
    129       !!  OPA 9.0, LODYC-IPSL (2003) 
    130       !!------------------------------------------------------------------------------ 
     115      !! References:  Flather, R. A., 1976: A tidal model of the northwest European 
     116      !!              continental shelf. Mem. Soc. R. Sci. Liege, Ser. 6,10, 141-164.      
     117      !!---------------------------------------------------------------------- 
     118      INTEGER  ::   ib, igrd                         ! dummy loop indices 
     119      INTEGER  ::   ii, ij, iim1, iip1, ijm1, ijp1   ! 2D addresses 
     120      REAL(wp) ::   zcorr                            ! Flather correction 
     121      !!---------------------------------------------------------------------- 
    131122 
    132123      ! ---------------------------------! 
     
    134125      ! ---------------------------------!  
    135126      
     127      IF(ln_bdy_dyn_fla .OR. ln_bdy_tides) THEN ! If these are both false, then this routine does nothing.  
     128 
    136129      ! Fill temporary array with ssh data (here spgu): 
    137       jgrd = 1 
    138       DO jb = 1, nblenrim(jgrd) 
    139         ji = nbi(jb,jgrd) 
    140         jj = nbj(jb,jgrd) 
    141         spgu(ji, jj) = sshbdy(jb) 
    142 #if defined key_bdy_tides 
    143         spgu(ji, jj) = spgu(ji, jj) + sshtide(jb) 
    144 #endif 
     130      igrd = 1 
     131      spgu(:,:) = 0.0 
     132      DO ib = 1, nblenrim(igrd) 
     133         ii = nbi(ib,igrd) 
     134         ij = nbj(ib,igrd) 
     135         IF( ln_bdy_dyn_fla ) spgu(ii, ij) = sshbdy(ib) 
     136         IF( ln_bdy_tides )   spgu(ii, ij) = spgu(ii, ij) + sshtide(ib) 
    145137      END DO 
    146  
    147       jgrd = 2  !: Flather bc on u-velocity;  
    148                 ! remember that flagu=-1 if normal velocity direction is outward 
    149                 ! I think we should rather use after ssh ? 
    150  
    151       DO jb = 1, nblenrim(jgrd) 
    152         ji  = nbi(jb,jgrd) 
    153         jj  = nbj(jb,jgrd)  
    154         jim1 = ji+MAX(0, INT(flagu(jb))) ! T pts i-indice inside the boundary 
    155         jip1 = ji-MIN(0, INT(flagu(jb))) ! T pts i-indice outside the boundary  
    156       
    157         corr = - flagu(jb) * sqrt (grav / (hu_e(ji, jj) + 1.e-20) )  & 
    158                                       * ( sshn_e(jim1, jj) - spgu(jip1,jj) ) 
    159         ua_e(ji, jj) = ( ubtbdy(jb) + utide(jb) ) * hu_e(ji,jj) 
    160         if ( ln_bdy_fla ) then 
    161           ua_e(ji,jj) = ua_e(ji,jj) + corr * umask(ji,jj,1) * hu_e(ji,jj) 
    162         endif 
    163  
     138      ! 
     139      igrd = 2      ! Flather bc on u-velocity;  
     140      !             ! remember that flagu=-1 if normal velocity direction is outward 
     141      !             ! I think we should rather use after ssh ? 
     142      DO ib = 1, nblenrim(igrd) 
     143         ii  = nbi(ib,igrd) 
     144         ij  = nbj(ib,igrd)  
     145         iim1 = ii + MAX( 0, INT( flagu(ib) ) )   ! T pts i-indice inside the boundary 
     146         iip1 = ii - MIN( 0, INT( flagu(ib) ) )   ! T pts i-indice outside the boundary  
     147         ! 
     148         zcorr = - flagu(ib) * SQRT( grav / (hu_e(ii, ij) + 1.e-20) ) * ( sshn_e(iim1, ij) - spgu(iip1,ij) ) 
     149         ua_e(ii, ij) = ( ubtbdy(ib) + utide(ib) ) * hu_e(ii,ij) 
     150         ua_e(ii,ij) = ua_e(ii,ij) + zcorr * umask(ii,ij,1) * hu_e(ii,ij) 
    164151      END DO 
    165         
    166       jgrd = 3  !: Flather bc on v-velocity 
    167                 ! remember that flagv=-1 if normal velocity direction is outward 
    168                   
    169       DO jb = 1, nblenrim(jgrd) 
    170         ji  = nbi(jb,jgrd) 
    171         jj  = nbj(jb,jgrd)  
    172         jjm1 = jj+MAX(0, INT(flagv(jb))) ! T pts j-indice inside the boundary 
    173         jjp1 = jj-MIN(0, INT(flagv(jb))) ! T pts j-indice outside the boundary  
    174       
    175         corr = - flagv(jb) * sqrt (grav / (hv_e(ji, jj) + 1.e-20) )  & 
    176                                       * ( sshn_e(ji, jjm1) - spgu(ji,jjp1) ) 
    177         va_e(ji, jj) = ( vbtbdy(jb) + vtide(jb) ) * hv_e(ji,jj) 
    178         if ( ln_bdy_fla ) then 
    179           va_e(ji,jj) = va_e(ji,jj) + corr * vmask(ji,jj,1) * hv_e(ji,jj) 
    180         endif 
    181  
     152      ! 
     153      igrd = 3      ! Flather bc on v-velocity 
     154      !             ! remember that flagv=-1 if normal velocity direction is outward 
     155      DO ib = 1, nblenrim(igrd) 
     156         ii  = nbi(ib,igrd) 
     157         ij  = nbj(ib,igrd)  
     158         ijm1 = ij + MAX( 0, INT( flagv(ib) ) )   ! T pts j-indice inside the boundary 
     159         ijp1 = ij - MIN( 0, INT( flagv(ib) ) )   ! T pts j-indice outside the boundary  
     160         ! 
     161         zcorr = - flagv(ib) * SQRT( grav / (hv_e(ii, ij) + 1.e-20) ) * ( sshn_e(ii, ijm1) - spgu(ii,ijp1) ) 
     162         va_e(ii, ij) = ( vbtbdy(ib) + vtide(ib) ) * hv_e(ii,ij) 
     163         va_e(ii,ij) = va_e(ii,ij) + zcorr * vmask(ii,ij,1) * hv_e(ii,ij) 
    182164      END DO 
    183  
     165      ! 
    184166      CALL lbc_lnk( ua_e, 'U', 1. ) ! Boundary points should be updated 
    185167      CALL lbc_lnk( va_e, 'V', 1. ) ! 
     168      ! 
     169      ENDIF ! ln_bdy_dyn_fla .or. ln_bdy_tides 
    186170 
    187171   END SUBROUTINE bdy_dyn_fla 
     
    189173 
    190174#else 
    191    !!================================================================================= 
    192    !!                       ***  MODULE  bdydyn  *** 
    193    !! Ocean dynamics:   Flow relaxation scheme of velocities on unstruc. open boundary 
    194    !!================================================================================= 
     175   !!---------------------------------------------------------------------- 
     176   !!   Dummy module                   NO Unstruct Open Boundary Conditions 
     177   !!---------------------------------------------------------------------- 
    195178CONTAINS 
    196  
    197    SUBROUTINE bdy_dyn_frs 
    198                               ! No Unstructured open boundaries ==> empty routine 
     179   SUBROUTINE bdy_dyn_frs( kt )      ! Empty routine 
     180      WRITE(*,*) 'bdy_dyn_frs: You should not have seen this print! error?', kt 
    199181   END SUBROUTINE bdy_dyn_frs 
    200  
    201    SUBROUTINE bdy_dyn_fla 
    202                               ! No Unstructured open boundaries ==> empty routine 
     182   SUBROUTINE bdy_dyn_fla            ! Empty routine 
     183      WRITE(*,*) 'bdy_dyn_fla: You should not have seen this print! error?' 
    203184   END SUBROUTINE bdy_dyn_fla 
    204  
    205185#endif 
    206186 
     187   !!====================================================================== 
    207188END MODULE bdydyn 
Note: See TracChangeset for help on using the changeset viewer.