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 911 for trunk/NEMO/OPA_SRC – NEMO

Changeset 911 for trunk/NEMO/OPA_SRC


Ignore:
Timestamp:
2008-04-28T11:31:32+02:00 (16 years ago)
Author:
ctlod
Message:

Implementation of the BDY package, see ticket: #126

Location:
trunk/NEMO/OPA_SRC
Files:
9 added
15 edited

Legend:

Unmodified
Added
Removed
  • trunk/NEMO/OPA_SRC/DOM/domvvl.F90

    r888 r911  
    55   !!====================================================================== 
    66   !! History :   9.0  !  06-06  (B. Levier, L. Marie)  original code 
     7   !!              "   !  07-07  (D. Storkey) Bug fixes and code for BDY option. 
    78   !!---------------------------------------------------------------------- 
    89 
     
    2425   USE lib_mpp         ! distributed memory computing library 
    2526   USE lbclnk          ! ocean lateral boundary conditions (or mpp link) 
     27   USE bdy_oce         ! unstructured open boundary conditions 
    2628 
    2729   IMPLICIT NONE 
     
    116118            zpk = mbathy(ji,jj) - 1 
    117119            DO jk = 1, zpk 
     120#if defined key_sigma_vvl 
     121               mut(ji,jj,jk) = 1./SUM( fsve3t(ji,jj,1:zpk) )  
     122               muu(ji,jj,jk) = 1./SUM( fsve3u(ji,jj,1:zpk) )  
     123               muv(ji,jj,jk) = 1./SUM( fsve3v(ji,jj,1:zpk) )  
     124               muf(ji,jj,jk) = 1./SUM( fsve3f(ji,jj,1:zpk) )  
     125#else 
    118126               mut(ji,jj,jk) = SUM( fsve3t(ji,jj,jk:zpk) ) / zmut(ji,jj) 
    119127               muu(ji,jj,jk) = SUM( fsve3u(ji,jj,jk:zpk) ) / zmuu(ji,jj) 
    120128               muv(ji,jj,jk) = SUM( fsve3v(ji,jj,jk:zpk) ) / zmuv(ji,jj) 
    121129               muf(ji,jj,jk) = SUM( fsve3f(ji,jj,jk:zpk) ) / zmuf(ji,jj) 
     130#endif 
    122131            END DO 
    123132         END DO 
     
    191200      END DO 
    192201 
     202      IF( MINVAL(fse3t(:,:,:)) < 0.0 ) THEN 
     203         CALL ctl_stop('E R R O R : dom_vvl','Level thickness fse3t less than zero.') 
     204         nstop = nstop+1 
     205      ENDIF 
     206 
    193207      ! Local depth or Inverse of the local depth of the water column at u- and v-points 
    194208      ! ------------------------------ 
     
    288302#endif 
    289303 
     304#if defined key_bdy || key_bdy_tides 
     305      DO jj = 1, jpj 
     306         DO ji = 1, jpi 
     307            zhdiv(ji,jj) = zhdiv(ji,jj)*bdytmask(ji,jj) 
     308         END DO 
     309      END DO 
     310#endif 
     311 
    290312      CALL lbc_lnk( zhdiv, 'T', 1. ) 
    291313 
  • trunk/NEMO/OPA_SRC/DYN/divcur.F90

    r780 r911  
    1414   USE in_out_manager  ! I/O manager 
    1515   USE obc_oce         ! ocean lateral open boundary condition 
     16   USE bdy_oce        ! Unstructured open boundaries variables 
    1617   USE lbclnk          ! ocean lateral boundary conditions (or mpp link) 
    1718 
     
    7778      !!   8.2  !  00-03  (G. Madec)  no slip accurate 
    7879      !!   9.0  !  03-08  (G. Madec)  merged of cur and div, free form, F90 
     80      !!        !  05-01  (J. Chanut, A. Sellar) unstructured open boundaries 
    7981      !!---------------------------------------------------------------------- 
    8082      !! * Arguments 
     
    134136         ENDIF 
    135137#endif 
     138#endif          
     139#if defined key_bdy || key_bdy_tides 
     140         ! unstructured open boundaries (div must be zero behind the open boundary) 
     141         DO jj = 1, jpj 
     142           DO ji = 1, jpi 
     143             hdivn(ji,jj,jk)=hdivn(ji,jj,jk)*bdytmask(ji,jj) 
     144           END DO 
     145         END DO 
    136146#endif          
    137147#if defined key_agrif 
     
    291301      !!   8.1  !  97-08  (J.M. Molines)  Open boundaries 
    292302      !!   9.0  !  02-09  (G. Madec, E. Durand)  Free form, F90 
     303      !!        !  05-01  (J. Chanut) Unstructured open boundaries 
    293304      !!---------------------------------------------------------------------- 
    294305      !! * Arguments 
     
    344355#endif 
    345356#endif          
     357#if defined key_bdy || key_bdy_tides 
     358         ! unstructured open boundaries (div must be zero behind the open boundary) 
     359         DO jj = 1, jpj 
     360           DO ji = 1, jpi 
     361             hdivn(ji,jj,jk)=hdivn(ji,jj,jk)*bdytmask(ji,jj) 
     362           END DO 
     363         END DO 
     364#endif         
    346365#if defined key_agrif 
    347366         if ( .NOT. AGRIF_Root() ) then 
  • trunk/NEMO/OPA_SRC/DYN/dynnxt.F90

    r782 r911  
    1616   USE obcdyn_bt       ! 2D open boundary condition for momentum (obc_dyn_bt routine) 
    1717   USE obcvol          ! ocean open boundary condition (obc_vol routines) 
     18   USE bdy_oce         ! unstructured open boundary conditions 
     19   USE bdydta          ! unstructured open boundary conditions 
     20   USE bdydyn          ! unstructured open boundary conditions 
    1821   USE dynspg_oce      ! type of surface pressure gradient 
    1922   USE lbclnk          ! lateral boundary condition (or mpp link) 
     
    6669      !!        !  02-10  (C. Talandier, A-M. Treguier) Open boundary cond. 
    6770      !!   9.0  !  05-11  (V. Garnier) Surface pressure gradient organization 
     71      !!    "   !  07-07  (D. Storkey) Calls to BDY routines.  
    6872      !!---------------------------------------------------------------------- 
    6973      !! * Arguments 
     
    186190            END DO 
    187191         END DO 
     192 
     193         IF( lk_vvl ) THEN 
     194            ! Unweight velocities prior to updating open boundaries. 
     195 
     196            DO jj = 1, jpj                      ! caution: don't use (:,:) for this loop 
     197               DO ji = 1, jpi                   ! it causes optimization problems on NEC in auto-tasking 
     198                  ua(ji,jj,jk) = ua(ji,jj,jk) / fse3u(ji,jj,jk) 
     199                  va(ji,jj,jk) = va(ji,jj,jk) / fse3v(ji,jj,jk) 
     200 
     201                  un(ji,jj,jk) = un(ji,jj,jk) / fse3u(ji,jj,jk) 
     202                  vn(ji,jj,jk) = vn(ji,jj,jk) / fse3v(ji,jj,jk) 
     203 
     204                  ub(ji,jj,jk) = ub(ji,jj,jk) / zfse3ub(ji,jj,jk) 
     205                  vb(ji,jj,jk) = vb(ji,jj,jk) / zfse3vb(ji,jj,jk) 
     206               END DO 
     207            END DO 
     208 
     209         ENDIF 
     210 
    188211# if defined key_obc 
    189212         !                                             ! =============== 
     
    209232 
    210233         IF ( ln_vol_cst ) CALL obc_vol( kt ) 
     234 
     235      ENDIF 
     236 
     237      !                                                ! =============== 
     238      DO jk = 1, jpkm1                                 ! Horizontal slab 
     239         !                                             ! =============== 
     240# elif defined key_bdy || key_bdy_tides 
     241         !                                             ! =============== 
     242      END DO                                           !   End of slab 
     243      !                                                ! =============== 
     244      ! Update (ua,va) along open boundaries (for exp or ts options). 
     245      IF ( lk_dynspg_exp .or. lk_dynspg_ts ) THEN 
     246   
     247         CALL bdy_dyn_frs( kt ) 
     248 
     249         IF ( ln_bdy_fla ) THEN 
     250 
     251            ua_e(:,:)=0.0 
     252            va_e(:,:)=0.0 
     253 
     254            ! Set these variables for use in bdy_dyn_fla 
     255            hu_e(:,:) = hu(:,:) 
     256            hv_e(:,:) = hv(:,:) 
     257 
     258            DO jk = 1, jpkm1 
     259               !! Vertically integrated momentum trends 
     260               ua_e(:,:) = ua_e(:,:) + fse3u(:,:,jk) * umask(:,:,jk) * ua(:,:,jk) 
     261               va_e(:,:) = va_e(:,:) + fse3v(:,:,jk) * vmask(:,:,jk) * va(:,:,jk) 
     262            END DO 
     263 
     264            DO jk = 1 , jpkm1 
     265               ua(:,:,jk) = ua(:,:,jk) - ua_e(:,:) * hur(:,:) 
     266               va(:,:,jk) = va(:,:,jk) - va_e(:,:) * hvr(:,:) 
     267            END DO 
     268 
     269            CALL bdy_dta_bt( kt+1, 0) 
     270            CALL bdy_dyn_fla 
     271 
     272         ENDIF 
     273 
     274         DO jk = 1 , jpkm1 
     275            ua(:,:,jk) = ua(:,:,jk) + ua_e(:,:) * hur(:,:) 
     276            va(:,:,jk) = va(:,:,jk) + va_e(:,:) * hvr(:,:) 
     277         END DO 
    211278 
    212279      ENDIF 
  • trunk/NEMO/OPA_SRC/DYN/dynspg_flt.F90

    r888 r911  
    1111   !!            " "  !  05-11  (V. Garnier) Surface pressure gradient organization 
    1212   !!            " "  !  06-07  (S. Masson)  distributed restart using iom 
     13   !!            " "  !  05-01  (J.Chanut, A.Sellar) Calls to BDY routines.  
    1314   !!---------------------------------------------------------------------- 
    1415#if defined key_dynspg_flt   ||   defined key_esopa   
     
    3637   USE obcdyn          ! ocean open boundary condition (obc_dyn routines) 
    3738   USE obcvol          ! ocean open boundary condition (obc_vol routines) 
     39   USE bdy_oce         ! Unstructured open boundaries condition 
     40   USE bdydyn          ! Unstructured open boundaries condition (bdy_dyn routine)  
     41   USE bdyvol          ! Unstructured open boundaries condition (bdy_vol routine) 
    3842   USE cla_dynspg      ! cross land advection 
    3943   USE in_out_manager  ! I/O manager 
     
    249253      CALL obc_vol( kt )      ! Correction of the barotropic componant velocity to control the volume of the system 
    250254#endif 
     255#if defined key_bdy 
     256      ! Update velocities on unstructured boundary using the Flow Relaxation Scheme 
     257      CALL bdy_dyn_frs( kt ) 
     258 
     259      IF (ln_bdy_vol) THEN 
     260      ! Correction of the barotropic component velocity to control the volume of the system 
     261        CALL bdy_vol( kt ) 
     262      ENDIF 
     263#endif 
    251264#if defined key_agrif 
    252265      CALL Agrif_dyn( kt )    ! Update velocities on each coarse/fine interfaces  
     
    377390            spgu(ji,jj) = z2dt * ztdgu * obcumask(ji,jj) 
    378391            spgv(ji,jj) = z2dt * ztdgv * obcvmask(ji,jj) 
     392#elif defined key_bdy 
     393            ! caution : grad D = 0 along open boundaries 
     394            ! Remark: The filtering force could be reduced here in the FRS zone 
     395            !         by multiplying spgu/spgv by (1-alpha) ??   
     396            spgu(ji,jj) = z2dt * ztdgu * bdyumask(ji,jj) 
     397            spgv(ji,jj) = z2dt * ztdgv * bdyvmask(ji,jj)            
    379398#else 
    380399            spgu(ji,jj) = z2dt * ztdgu 
  • trunk/NEMO/OPA_SRC/DYN/dynspg_oce.F90

    r719 r911  
    4646         un_b  , vn_b                     ! vertically integrated horizontal velocities (now) 
    4747      REAL(wp), PUBLIC, DIMENSION(jpi,jpj) :: & ! variables of the explicit barotropic loop 
    48          sshn_e, ssha_e,               &  ! sea surface heigth (now,after) 
    49          ua_e  , va_e                     ! vertically integrated horizontal velocities (after) 
     48         sshb_e, sshn_e, ssha_e,       &  ! sea surface heigth (before,now,after) 
     49         ua_e  , va_e,                 &  ! vertically integrated horizontal velocities (after) 
     50         hu_e  , hv_e                     ! depth arrays for the barotropic solution 
    5051#endif 
    5152 
  • trunk/NEMO/OPA_SRC/DYN/dynspg_ts.F90

    r888 r911  
    77   !!             " "  !  05-11  (V. Garnier, G. Madec)  optimization 
    88   !!             9.0  !  06-08  (S. Masson)  distributed restart using iom 
    9    !!             9.0  !  08-01  (R. Benshila)  change averaging method 
     9   !!              "   !  08-01  (R. Benshila)  change averaging method 
     10   !!              "   !  07-07  (D. Storkey) calls to BDY routines 
    1011   !!--------------------------------------------------------------------- 
    1112#if defined key_dynspg_ts   ||   defined key_esopa 
     
    3031   USE obc_oce         ! Lateral open boundary condition 
    3132   USE obc_par         ! open boundary condition parameters 
     33   USE bdy_oce         ! unstructured open boundaries 
     34   USE bdy_par         ! unstructured open boundaries 
     35   USE bdydta          ! unstructured open boundaries 
     36   USE bdydyn          ! unstructured open boundaries 
     37   USE bdytides        ! tidal forcing at unstructured open boundaries. 
    3238   USE lib_mpp         ! distributed memory computing library 
    3339   USE lbclnk          ! ocean lateral boundary conditions (or mpp link) 
     
    102108         zua, zva, zub, zvb,                &  !     "        " 
    103109         zssha_b, zua_b, zva_b,             &  !     "        " 
    104          zsshb_e, zub_e, zvb_e,             &  !     "        " 
     110         zub_e, zvb_e,                      &  !     "        " 
    105111         zun_e, zvn_e                          !     "        " 
    106112      !! Variable volume 
    107113      REAL(wp), DIMENSION(jpi,jpj)     ::   & 
    108          zspgu_1, zspgv_1, zsshun_e, zsshvn_e, hu_e, hv_e         ! 2D workspace 
     114         zspgu_1, zspgv_1, zsshun_e, zsshvn_e                     ! 2D workspace 
    109115      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zfse3un_e, zfse3vn_e  ! 3D workspace 
    110116      !!---------------------------------------------------------------------- 
     
    130136         ua_e(:,:)   = un_b(:,:) 
    131137         va_e(:,:)   = vn_b(:,:) 
     138         hu_e(:,:)   = hu(:,:) 
     139         hv_e(:,:)   = hv(:,:) 
    132140 
    133141         IF( ln_dynvor_een ) THEN 
     
    300308 
    301309      ! variables for the barotropic equations 
    302       zsshb_e(:,:) = sshn_b(:,:)       ! (barotropic) sea surface height (before and now) 
     310      sshb_e (:,:) = sshn_b(:,:)       ! (barotropic) sea surface height (before and now) 
    303311      sshn_e (:,:) = sshn_b(:,:) 
    304312      zub_e  (:,:) = un_b  (:,:)       ! barotropic transports issued from the barotropic equations (before and now) 
     
    309317      zua_b  (:,:) = 0.e0 
    310318      zva_b  (:,:) = 0.e0 
     319      hu_e   (:,:) = hu   (:,:)     ! (barotropic) ocean depth at u-point 
     320      hv_e   (:,:) = hv   (:,:)     ! (barotropic) ocean depth at v-point 
    311321      IF( lk_vvl ) THEN 
    312322         zsshun_e(:,:)    = sshu (:,:)     ! (barotropic) sea surface height (now) at u-point 
    313323         zsshvn_e(:,:)    = sshv (:,:)     ! (barotropic) sea surface height (now) at v-point 
    314          hu_e(:,:)        = hu   (:,:)     ! (barotropic) ocean depth at u-point 
    315          hv_e(:,:)        = hv   (:,:)     ! (barotropic) ocean depth at v-point 
    316324         zfse3un_e(:,:,:) = fse3u(:,:,:)   ! (barotropic) scale factors  at u-point 
    317325         zfse3un_e(:,:,:) = fse3v(:,:,:)   ! (barotropic) scale factors  at v-point 
     
    338346         ! Time interpolation of open boundary condition data 
    339347         IF( lk_obc )   CALL obc_dta_bt( kt, jit ) 
     348         IF( lk_bdy .or. lk_bdy_tides)   CALL bdy_dta_bt( kt, jit+1 ) 
    340349 
    341350         ! Horizontal divergence of barotropic transports 
    342351         !-------------------------------------------------- 
     352         zhdiv(:,:) = 0.e0 
    343353         DO jj = 2, jpjm1 
    344354            DO ji = fs_2, fs_jpim1   ! vector opt. 
     
    360370#endif 
    361371 
     372         IF( lk_bdy .or. lk_bdy_tides ) THEN 
     373            DO jj = 1, jpj 
     374               DO ji = 1, jpi 
     375                  zhdiv(ji,jj) = zhdiv(ji,jj)*bdytmask(ji,jj) 
     376               END DO 
     377            END DO 
     378         ENDIF 
     379 
    362380         ! Sea surface height from the barotropic system 
    363381         !---------------------------------------------- 
    364382         DO jj = 2, jpjm1 
    365383            DO ji = fs_2, fs_jpim1   ! vector opt. 
    366                ssha_e(ji,jj) = ( zsshb_e(ji,jj) - z2dt_e *  ( zraur * emp(ji,jj)  & 
     384               ssha_e(ji,jj) = ( sshb_e(ji,jj) - z2dt_e *  ( zraur * emp(ji,jj)  & 
    367385            &            +  zhdiv(ji,jj) ) ) * tmask(ji,jj,1) 
    368386            END DO 
     
    457475         !         - Correct the barotropic transports 
    458476         IF( lk_obc )   CALL obc_fla_ts 
    459  
     477         IF( lk_bdy .or. lk_bdy_tides )   CALL bdy_dyn_fla 
    460478 
    461479         ! ... Boundary conditions on ua_e, va_e, ssha_e 
     
    475493         ! --------------------------------------- 
    476494         IF( neuler == 0 .AND. jit == 1 ) THEN   ! Euler (forward) time stepping 
    477             zsshb_e(:,:) = sshn_e(:,:) 
     495            sshb_e (:,:) = sshn_e(:,:) 
    478496            zub_e  (:,:) = zun_e (:,:) 
    479497            zvb_e  (:,:) = zvn_e (:,:) 
     
    482500            zvn_e  (:,:) = va_e  (:,:) 
    483501         ELSE                                        ! Asselin filtering 
    484             zsshb_e(:,:) = atfp * ( zsshb_e(:,:) + ssha_e(:,:) ) + atfp1 * sshn_e(:,:) 
     502            sshb_e (:,:) = atfp * ( sshb_e(:,:) + ssha_e(:,:) ) + atfp1 * sshn_e(:,:) 
    485503            zub_e  (:,:) = atfp * ( zub_e  (:,:) + ua_e  (:,:) ) + atfp1 * zun_e  (:,:) 
    486504            zvb_e  (:,:) = atfp * ( zvb_e  (:,:) + va_e  (:,:) ) + atfp1 * zvn_e  (:,:) 
     
    571589      IF( lp_obc_south )   zhdiv(nis0  :nis1  ,njs0  :njs1)   = 0.e0    ! south 
    572590#endif 
     591 
     592      IF ( lk_bdy .or. lk_bdy_tides ) THEN 
     593         DO jj = 1, jpj 
     594           DO ji = 1, jpi 
     595             zhdiv(ji,jj) = zhdiv(ji,jj)*bdytmask(ji,jj) 
     596           END DO 
     597         END DO 
     598      ENDIF 
    573599 
    574600      ! sea surface height 
  • trunk/NEMO/OPA_SRC/DYN/wzvmod.F90

    r888 r911  
    77   !!             7.0  !  96-01  (G. Madec)  Statement function for e3 
    88   !!             8.5  !  02-07  (G. Madec)  Free form, F90 
     9   !!              "   !  07-07  (D. Storkey) Zero zhdiv at open boundary (BDY)  
    910   !!---------------------------------------------------------------------- 
    1011   !!   wzv        : Compute the vertical velocity 
     
    1819   USE prtctl          ! Print control 
    1920   USE phycst 
     21   USE bdy_oce         ! unstructured open boundaries 
    2022   USE lbclnk          ! ocean lateral boundary condition (or mpp link) 
    2123 
     
    5456 
    5557      !! * Local declarations 
     58      INTEGER  ::     jgrd, jb           ! temporary scalars 
    5659      INTEGER  ::           jk           ! dummy loop indices 
    5760      !! Variable volume 
     
    8891         ! Horizontal divergence of barotropic transports 
    8992         !-------------------------------------------------- 
     93         zhdiv(:,:) = 0.e0 
    9094         DO jj = 2, jpjm1 
    9195            DO ji = 2, jpim1   ! vector opt. 
     
    105109         IF( lp_obc_north )   zhdiv(nin0  :nin1  ,njn0p1:njn1p1) = 0.e0    ! north 
    106110         IF( lp_obc_south )   zhdiv(nis0  :nis1  ,njs0  :njs1)   = 0.e0    ! south 
     111#endif 
     112 
     113#if defined key_bdy || defined key_bdy_tides 
     114         jgrd=1 !: tracer grid. 
     115         DO jb = 1, nblenrim(jgrd) 
     116           ji = nbi(jb,jgrd) 
     117           jj = nbj(jb,jgrd) 
     118           zhdiv(ji, jj) = 0.e0 
     119           zhdiv(ji, jj) = 0.e0 
     120           zhdiv(ji, jj) = 0.e0 
     121         END DO 
    107122#endif 
    108123 
  • trunk/NEMO/OPA_SRC/IOM/iom.F90

    r752 r911  
    66   !! History :  9.0  ! 05 12  (J. Belier) Original code 
    77   !!            9.0  ! 06 02  (S. Masson) Adaptation to NEMO 
     8   !!             "   ! 07 07  (D. Storkey) Changes to iom_gettime 
    89   !!-------------------------------------------------------------------- 
    910   !!gm  caution add !DIR nec: improved performance to be checked as well as no result changes 
     
    646647 
    647648 
    648    SUBROUTINE iom_gettime( kiomid, cdvar, ptime ) 
     649   SUBROUTINE iom_gettime( kiomid, ptime, cdvar, kntime, cdunits, cdcalendar ) 
    649650      !!-------------------------------------------------------------------- 
    650651      !!                   ***  SUBROUTINE iom_gettime  *** 
     
    652653      !! ** Purpose : read the time axis cdvar in the file  
    653654      !!-------------------------------------------------------------------- 
    654       INTEGER               , INTENT(in   ) ::   kiomid   ! file Identifier 
    655       CHARACTER(len=*)      , INTENT(in   ) ::   cdvar    ! time axis name 
    656       REAL(wp), DIMENSION(:), INTENT(  out) ::   ptime    ! the time axis 
    657       ! 
     655      INTEGER                    , INTENT(in   ) ::   kiomid     ! file Identifier 
     656      REAL(wp), DIMENSION(:)     , INTENT(  out) ::   ptime      ! the time axis 
     657      CHARACTER(len=*), OPTIONAL , INTENT(in   ) ::   cdvar      ! time axis name 
     658      INTEGER         , OPTIONAL , INTENT(  out) ::   kntime     ! number of times in file 
     659      CHARACTER(len=*), OPTIONAL , INTENT(  out) ::   cdunits    ! units attribute of time coordinate 
     660      CHARACTER(len=*), OPTIONAL , INTENT(  out) ::   cdcalendar ! calendar attribute of  
     661      ! 
     662      INTEGER, DIMENSION(1) :: kdimsz 
    658663      INTEGER            ::   idvar    ! id of the variable 
     664      CHARACTER(LEN=32)  ::   tname    ! local name of time coordinate 
    659665      CHARACTER(LEN=100) ::   clinfo   ! info character 
    660666      !--------------------------------------------------------------------- 
    661667      ! 
    662       IF( kiomid > 0 ) THEN 
    663          clinfo = 'iom_gettime, file: '//trim(iom_file(kiomid)%name)//', var: '//trim(cdvar) 
    664          idvar = iom_varid( kiomid, cdvar ) 
     668      IF ( PRESENT(cdvar) ) THEN 
     669         tname = cdvar 
     670      ELSE 
     671         tname = iom_file(kiomid)%uldname 
     672      ENDIF 
     673      IF( kiomid > 0 ) THEN 
     674         clinfo = 'iom_gettime, file: '//trim(iom_file(kiomid)%name)//', var: '//trim(tname) 
     675         IF ( PRESENT(kntime) ) THEN 
     676            idvar  = iom_varid( kiomid, tname, kdimsz = kdimsz ) 
     677            kntime = kdimsz(1) 
     678         ELSE 
     679            idvar = iom_varid( kiomid, tname ) 
     680         ENDIF 
    665681         ! 
    666682         ptime(:) = 0. ! default definition 
     
    670686                  IF( iom_file(kiomid)%dimsz(1,idvar) == size(ptime) ) THEN 
    671687                     SELECT CASE (iom_file(kiomid)%iolib) 
    672                      CASE (jpioipsl )   ;   CALL iom_ioipsl_gettime( kiomid, idvar, ptime ) 
    673                      CASE (jpnf90   )   ;   CALL iom_nf90_gettime(   kiomid, idvar, ptime ) 
     688                     CASE (jpioipsl )   ;   CALL iom_ioipsl_gettime( kiomid, idvar, ptime, cdunits, cdcalendar ) 
     689                     CASE (jpnf90   )   ;   CALL iom_nf90_gettime(   kiomid, idvar, ptime, cdunits, cdcalendar ) 
    674690                     CASE (jprstdimg)   ;   CALL ctl_stop( TRIM(clinfo)//' case IO library == jprstdimg not coded...' ) 
    675691                     CASE DEFAULT     
  • trunk/NEMO/OPA_SRC/IOM/iom_def.F90

    r839 r911  
    55   !!==================================================================== 
    66   !! History :  9.0  ! 06 09  (S. Masson) Original code 
     7   !!             "   ! 07 07  (D. Storkey) Add uldname 
    78   !!-------------------------------------------------------------------- 
    89   !!--------------------------------------------------------------------------------- 
     
    5758      INTEGER                                   ::   iduld    !: id of the unlimited dimension 
    5859      INTEGER                                   ::   irec     !: writing record position   
     60      CHARACTER(LEN=32)                         ::   uldname  !: name of the unlimited dimension 
    5961      CHARACTER(LEN=32), DIMENSION(jpmax_vars)  ::   cn_var   !: names of the variables 
    6062      INTEGER, DIMENSION(jpmax_vars)            ::   nvid     !: id of the variables 
  • trunk/NEMO/OPA_SRC/IOM/iom_ioipsl.F90

    r745 r911  
    66   !! History :  9.0  ! 05 12  (J. Belier) Original code 
    77   !!            9.0  ! 06 02  (S. Masson) Adaptation to NEMO 
     8   !!             "   ! 07 07  (D. Storkey) Changes to iom_ioipsl_gettime 
    89   !!-------------------------------------------------------------------- 
    910   !!gm  caution add !DIR nec: improved performance to be checked as well as no result changes 
     
    254255 
    255256 
    256    SUBROUTINE iom_ioipsl_gettime( kiomid, kvid, ptime ) 
     257   SUBROUTINE iom_ioipsl_gettime( kiomid, kvid, ptime, cdunits, cdcalendar ) 
    257258      !!-------------------------------------------------------------------- 
    258259      !!                   ***  SUBROUTINE iom_gettime  *** 
     
    260261      !! ** Purpose : read the time axis kvid in the file with IOIPSL (only fliocom module) 
    261262      !!-------------------------------------------------------------------- 
    262       INTEGER               , INTENT(in   ) ::   kiomid   ! file Identifier 
    263       INTEGER               , INTENT(in   ) ::   kvid     ! variable id 
    264       REAL(wp), DIMENSION(:), INTENT(  out) ::   ptime    ! the time axis 
     263      INTEGER                   , INTENT(in   ) ::   kiomid     ! file Identifier 
     264      INTEGER                   , INTENT(in   ) ::   kvid       ! variable id 
     265      REAL(wp), DIMENSION(:)    , INTENT(  out) ::   ptime      ! the time axis 
     266      CHARACTER(len=*), OPTIONAL, INTENT(  out) ::   cdunits    ! units attribute 
     267      CHARACTER(len=*), OPTIONAL, INTENT(  out) ::   cdcalendar ! calendar attribute 
    265268      !--------------------------------------------------------------------- 
    266269      ! 
    267270      CALL fliogetv( iom_file(kiomid)%nfid, TRIM(iom_file(kiomid)%cn_var(kvid)), ptime(:),   & 
    268271            &         start=(/ 1 /), count=(/ iom_file(kiomid)%dimsz(1, kvid) /) ) 
     272      IF ( PRESENT(cdunits) ) THEN  
     273         CALL fliogeta( iom_file(kiomid)%nfid, TRIM(iom_file(kiomid)%cn_var(kvid)), "units", cdunits ) 
     274      ENDIF 
     275      IF ( PRESENT(cdcalendar) ) THEN  
     276         CALL fliogeta( iom_file(kiomid)%nfid, TRIM(iom_file(kiomid)%cn_var(kvid)), "calendar", cdcalendar ) 
     277      ENDIF 
    269278      ! 
    270279   END SUBROUTINE iom_ioipsl_gettime 
  • trunk/NEMO/OPA_SRC/IOM/iom_nf90.F90

    r745 r911  
    66   !! History :  9.0  ! 05 12  (J. Belier) Original code 
    77   !!            9.0  ! 06 02  (S. Masson) Adaptation to NEMO 
     8   !!             "   ! 07 07  (D. Storkey) Changes to iom_nf90_gettime 
    89   !!-------------------------------------------------------------------- 
    910   !!gm  caution add !DIR nec: improved performance to be checked as well as no result changes 
     
    121122         iom_file(kiomid)%irec   = -1   ! useless for NetCDF files, used to know if the file is in define mode  
    122123         CALL iom_nf90_check(NF90_Inquire(if90id, unlimitedDimId = iom_file(kiomid)%iduld), clinfo) 
     124         IF ( iom_file(kiomid)%iduld .GE. 0 ) THEN 
     125           CALL iom_nf90_check(NF90_Inquire_Dimension(if90id, iom_file(kiomid)%iduld,   & 
     126        &                                               name = iom_file(kiomid)%uldname), clinfo) 
     127         ENDIF 
    123128         IF(lwp) WRITE(numout,*) '                   ---> '//TRIM(cdname)//' OK' 
    124129      ELSE 
     
    272277 
    273278 
    274    SUBROUTINE iom_nf90_gettime( kiomid, kvid, ptime ) 
     279   SUBROUTINE iom_nf90_gettime( kiomid, kvid, ptime, cdunits, cdcalendar ) 
    275280      !!-------------------------------------------------------------------- 
    276281      !!                   ***  SUBROUTINE iom_gettime  *** 
     
    278283      !! ** Purpose : read the time axis kvid in the file with NF90 
    279284      !!-------------------------------------------------------------------- 
    280       INTEGER               , INTENT(in   ) ::   kiomid   ! file Identifier 
    281       INTEGER               , INTENT(in   ) ::   kvid     ! variable id 
    282       REAL(wp), DIMENSION(:), INTENT(  out) ::   ptime    ! the time axis 
     285      INTEGER                   , INTENT(in   ) ::   kiomid     ! file Identifier 
     286      INTEGER                   , INTENT(in   ) ::   kvid       ! variable id 
     287      REAL(wp), DIMENSION(:)    , INTENT(  out) ::   ptime      ! the time axis 
     288      CHARACTER(len=*), OPTIONAL, INTENT(  out) ::   cdunits    ! units attribute 
     289      CHARACTER(len=*), OPTIONAL, INTENT(  out) ::   cdcalendar ! calendar attribute 
    283290      ! 
    284291      CHARACTER(LEN=100) ::   clinfo     ! info character 
     
    287294      CALL iom_nf90_check(NF90_GET_VAR(iom_file(kiomid)%nfid, iom_file(kiomid)%nvid(kvid), ptime(:),   & 
    288295            &                           start=(/ 1 /), count=(/ iom_file(kiomid)%dimsz(1, kvid) /)), clinfo) 
     296      IF ( PRESENT(cdunits) ) THEN  
     297         CALL iom_nf90_check(NF90_GET_ATT(iom_file(kiomid)%nfid, iom_file(kiomid)%nvid(kvid), "units", & 
     298            &                           values=cdunits), clinfo) 
     299      ENDIF 
     300      IF ( PRESENT(cdcalendar) ) THEN  
     301         CALL iom_nf90_check(NF90_GET_ATT(iom_file(kiomid)%nfid, iom_file(kiomid)%nvid(kvid), "calendar", & 
     302            &                           values=cdcalendar), clinfo) 
     303      ENDIF 
    289304      ! 
    290305   END SUBROUTINE iom_nf90_gettime 
  • trunk/NEMO/OPA_SRC/OBC/obcdta.F90

    r719 r911  
    88   !!   8.5  !  02-10 (C. Talandier, A-M. Treguier) Free surface, F90 
    99   !!   9.0  !  04-06 (F. Durand, A-M. Treguier) Netcdf BC files on input 
     10   !!    "   ! 07-07  (D. Storkey) Change to arguments for iom_gettime 
    1011   !!------------------------------------------------------------------------------ 
    1112#if defined key_obc 
     
    159160               IF ( lp_obc_east ) THEN 
    160161                  IF ( llnot_done ) THEN 
    161                      CALL iom_gettime ( id_e, TRIM(cl_vname), tcobc ) 
     162                     CALL iom_gettime ( id_e, tcobc, TRIM(cl_vname) ) 
    162163                     llnot_done = .FALSE. 
    163164                  ENDIF 
     
    166167               IF ( lp_obc_west ) THEN 
    167168                  IF ( llnot_done ) THEN 
    168                      CALL iom_gettime ( id_w, TRIM(cl_vname), tcobc ) 
     169                     CALL iom_gettime ( id_w, tcobc, TRIM(cl_vname) ) 
    169170                     llnot_done = .FALSE. 
    170171                 ENDIF 
     
    173174               IF ( lp_obc_north ) THEN 
    174175                  IF ( llnot_done ) THEN 
    175                      CALL iom_gettime ( id_n, TRIM(cl_vname), tcobc ) 
     176                     CALL iom_gettime ( id_n, tcobc, TRIM(cl_vname) ) 
    176177                     llnot_done = .FALSE. 
    177178                  ENDIF 
     
    180181               IF ( lp_obc_south ) THEN 
    181182                  IF ( llnot_done ) THEN 
    182                      CALL iom_gettime ( id_s, TRIM(cl_vname), tcobc ) 
     183                     CALL iom_gettime ( id_s, tcobc, TRIM(cl_vname) ) 
    183184                     llnot_done = .FALSE. 
    184185                  ENDIF 
  • trunk/NEMO/OPA_SRC/TRA/tranxt.F90

    r888 r911  
    2727   USE domvvl          ! variable volume 
    2828   USE obctra          ! open boundary condition (obc_tra routine) 
     29   USE bdytra          ! Unstructured open boundary condition (bdy_tra routine) 
    2930   USE in_out_manager  ! I/O manager 
    3031   USE lbclnk          ! ocean lateral boundary conditions (or mpp link) 
     
    159160      ! Update tracers on open boundaries. 
    160161      CALL obc_tra( kt ) 
     162 
     163      !                                                ! =============== 
     164      DO jk = 1, jpkm1                                 ! Horizontal slab 
     165         !                                             ! =============== 
     166#elif defined key_bdy 
     167         !                                             ! =============== 
     168      END DO                                           !   End of slab 
     169      !                                                ! =============== 
     170 
     171      ! Update tracers on open boundaries. 
     172      CALL bdy_tra( kt ) 
     173 
    161174      !                                                ! =============== 
    162175      DO jk = 1, jpkm1                                 ! Horizontal slab 
  • trunk/NEMO/OPA_SRC/opa.F90

    r900 r911  
    3434   !!    "   !  06-03  (L. Debreu, C. Mazauric)  Agrif implementation 
    3535   !!    "   !  06-04  (G. Madec, R. Benshila)  Step reorganization 
     36   !!    "   !  07-07  (J. Chanut, A. Sellar) Unstructured open boundaries (BDY) 
    3637   !!---------------------------------------------------------------------- 
    3738   !! * Modules used 
     
    5051   USE obc_par         ! open boundary cond. parameters 
    5152   USE obcini          ! open boundary cond. initialization (obc_ini routine) 
     53   USE bdy_par         ! unstructured open boundary cond. parameters 
     54   USE bdyini          ! unstructured open boundary cond. initialization (bdy_init routine) 
     55   USE bdytides        ! tides at open boundaries initialization (lk_bdy_tides) 
    5256   USE istate          ! initial state setting          (istate_init routine) 
    5357   USE eosbn2          ! equation of state            (eos bn2 routine) 
     
    269273      IF( lk_obc    )   CALL obc_init       ! Open boundaries  
    270274 
     275      IF( lk_bdy .OR. lk_bdy_tides )   & 
     276         &              CALL bdy_init       ! Unstructured open boundaries 
     277 
    271278      CALL istate_init                      ! ocean initial state (Dynamics and tracers) 
    272279 
  • trunk/NEMO/OPA_SRC/step.F90

    r888 r911  
    2020   !!             " "  !  06-07  (S. Masson)  restart using iom 
    2121   !!             " "  !  06-08  (G. Madec)  surface module  
     22   !!             " "  !  07-07  (J. Chanut, A. Sellar) Unstructured open boundaries (BDY) 
    2223   !!---------------------------------------------------------------------- 
    2324 
     
    7475   USE obcrad          ! open boundary cond. radiation    (obc_rad routine) 
    7576   USE obcspg          ! open boundary cond  spg          (obc_spg routine) 
     77 
     78   USE bdy_par         ! unstructured open boundary data variables 
     79   USE bdydta          ! unstructured open boundary data  (bdy_dta routine) 
    7680 
    7781   USE divcur          ! hor. divergence and curl      (div & cur routines) 
     
    187191      IF( lk_obc     )   CALL obc_rad( kstp )         ! compute phase velocities at open boundaries 
    188192 
     193      IF( lk_bdy     )   CALL bdy_dta( kstp )         ! update dynamic and tracer data at unstructured open boundary 
     194 
    189195      IF( ninist == 1 ) THEN                       ! Output the initial state and forcings 
    190196                        CALL dia_wri_state( 'output.init' ) 
Note: See TracChangeset for help on using the changeset viewer.