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 15080 – NEMO

Changeset 15080


Ignore:
Timestamp:
2021-07-05T16:19:12+02:00 (3 years ago)
Author:
ayoung
Message:

GO8 package merge, essential changes only. Ticket #2648.

Location:
NEMO/branches/UKMO/NEMO_r4.2RC_GO8_package
Files:
19 added
9 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/UKMO/NEMO_r4.2RC_GO8_package/src/ICE/icerst.F90

    r14239 r15080  
    2727   USE in_out_manager ! I/O manager 
    2828   USE iom            ! I/O manager library 
     29   USE ioipsl  , ONLY : ju2ymds                     ! for calendar 
    2930   USE lib_mpp        ! MPP library 
    3031   USE lib_fortran    ! fortran utilities (glob_sum + no signed zero) 
     
    5253      INTEGER, INTENT(in) ::   kt       ! number of iteration 
    5354      ! 
     55      INTEGER             ::   iyear, imonth, iday 
     56      REAL (wp)           ::   zsec 
     57      REAL (wp)           ::   zfjulday 
    5458      CHARACTER(len=20)   ::   clkt     ! ocean time-step define as a character 
    5559      CHARACTER(len=50)   ::   clname   ! ice output restart file name 
     
    6872         IF( nitrst <= nitend .AND. nitrst > 0 ) THEN 
    6973            ! beware of the format used to write kt (default is i8.8, that should be large enough...) 
    70             IF( nitrst > 99999999 ) THEN   ;   WRITE(clkt, *       ) nitrst 
    71             ELSE                           ;   WRITE(clkt, '(i8.8)') nitrst 
     74            IF ( ln_rstdate ) THEN 
     75               zfjulday = fjulday + (2*nn_fsbc+1)*rdt / rday 
     76               IF( ABS(zfjulday - REAL(NINT(zfjulday),wp)) < 0.1 / rday )   zfjulday = REAL(NINT(zfjulday),wp)   ! avoid truncation error 
     77               CALL ju2ymds( zfjulday, iyear, imonth, iday, zsec )            
     78               WRITE(clkt, '(i4.4,2i2.2)') iyear, imonth, iday 
     79            ELSE 
     80               IF( nitrst > 99999999 ) THEN   ;   WRITE(clkt, *       ) nitrst 
     81               ELSE                           ;   WRITE(clkt, '(i8.8)') nitrst 
     82               ENDIF 
    7283            ENDIF 
    7384            ! create the file 
  • NEMO/branches/UKMO/NEMO_r4.2RC_GO8_package/src/ICE/icethd_pnd.F90

    r14997 r15080  
    545545      ! a_ip_frac -> apnd 
    546546 
    547       CALL ctl_stop( 'STOP', 'icethd_pnd : topographic melt ponds are still an ongoing work' ) 
     547!      CALL ctl_stop( 'STOP', 'icethd_pnd : topographic melt ponds are still an ongoing work' ) 
    548548 
    549549      !--------------------------------------------------------------- 
  • NEMO/branches/UKMO/NEMO_r4.2RC_GO8_package/src/OCE/DIA/diawri.F90

    r15017 r15080  
    5454   USE iom            !  
    5555   USE ioipsl         !  
     56   USE eosbn2 
    5657 
    5758#if defined key_si3 
     
    125126      REAL(wp), DIMENSION(A2D(     0))     ::   z2d   ! 2D workspace 
    126127      REAL(wp), DIMENSION(A2D(nn_hls),jpk) ::   z3d   ! 3D workspace 
     128      CHARACTER(len=4),SAVE :: ttype , stype           ! temperature and salinity type 
    127129      !!---------------------------------------------------------------------- 
    128130      !  
    129131      IF( ln_timing )   CALL timing_start('dia_wri') 
     132      ! 
     133      IF( kt == nit000 ) THEN 
     134         IF( ln_TEOS10 ) THEN 
     135            IF ( iom_use("toce_pot") .OR. iom_use("soce_pra") .OR. iom_use("sst_pot") .OR. iom_use("sss_pra") & 
     136                  & .OR. iom_use("sbt_pot") .OR. iom_use("sbs_pra") .OR. iom_use("sstgrad_pot") .OR. iom_use("sstgrad2_pot") & 
     137                  & .OR. iom_use("tosmint_pot") .OR. iom_use("somint_pra"))  THEN  
     138               CALL ctl_stop( 'diawri: potential temperature and practical salinity not available with ln_TEOS10' ) 
     139            ELSE 
     140               ttype='con' ; stype='abs'   ! teos-10 using conservative temperature and absolute salinity 
     141            ENDIF  
     142         ELSE IF( ln_EOS80  ) THEN 
     143            IF ( iom_use("toce_con") .OR. iom_use("soce_abs") .OR. iom_use("sst_con") .OR. iom_use("sss_abs") & 
     144                  & .OR. iom_use("sbt_con") .OR. iom_use("sbs_abs") .OR. iom_use("sstgrad_con") .OR. iom_use("sstgrad2_con") & 
     145                  & .OR. iom_use("tosmint_con") .OR. iom_use("somint_abs"))  THEN  
     146               CALL ctl_stop( 'diawri: conservative temperature and absolute salinity not available with ln_EOS80' ) 
     147            ELSE 
     148               ttype='pot' ; stype='pra'   ! eos-80 using potential temperature and practical salinity 
     149            ENDIF 
     150         ELSE IF ( ln_SEOS) THEN 
     151            ttype='seos' ; stype='seos' ! seos using Simplified Equation of state 
     152         ENDIF 
     153      ENDIF 
     154 
    130155      !  
    131156      ! Output the initial state and forcings 
     
    200225#endif 
    201226       
    202       CALL iom_put( "toce", ts(:,:,:,jp_tem,Kmm) )    ! 3D temperature 
    203       CALL iom_put(  "sst", ts(:,:,1,jp_tem,Kmm) )    ! surface temperature 
    204       IF ( iom_use("sbt") ) THEN 
     227      CALL iom_put( "toce_"//ttype, ts(:,:,:,jp_tem,Kmm) )    ! 3D temperature 
     228      CALL iom_put(  "sst_"//ttype, ts(:,:,1,jp_tem,Kmm) )    ! surface temperature 
     229      IF ( iom_use("sbt_"//ttype) ) THEN 
    205230         DO_2D( 0, 0, 0, 0 ) 
    206231            ikbot = mbkt(ji,jj) 
    207232            z2d(ji,jj) = ts(ji,jj,ikbot,jp_tem,Kmm) 
    208233         END_2D 
    209          CALL iom_put( "sbt", z2d )                ! bottom temperature 
     234         CALL iom_put( "sbt_"//ttype, z2d )                ! bottom temperature 
    210235      ENDIF 
    211236       
    212       CALL iom_put( "soce", ts(:,:,:,jp_sal,Kmm) )    ! 3D salinity 
    213       CALL iom_put(  "sss", ts(:,:,1,jp_sal,Kmm) )    ! surface salinity 
    214       IF ( iom_use("sbs") ) THEN 
     237      CALL iom_put( "soce_"//stype, ts(:,:,:,jp_sal,Kmm) )    ! 3D salinity 
     238      CALL iom_put(  "sss_"//stype, ts(:,:,1,jp_sal,Kmm) )    ! surface salinity 
     239      IF ( iom_use("sbs_"//stype) ) THEN 
    215240         DO_2D( 0, 0, 0, 0 ) 
    216241            ikbot = mbkt(ji,jj) 
    217242            z2d(ji,jj) = ts(ji,jj,ikbot,jp_sal,Kmm) 
    218243         END_2D 
    219          CALL iom_put( "sbs", z2d )                ! bottom salinity 
     244         CALL iom_put( "sbs_"//stype, z2d )                ! bottom salinity 
    220245      ENDIF 
    221246 
     
    298323      ENDIF 
    299324          
    300       IF ( iom_use("sstgrad") .OR. iom_use("sstgrad2") ) THEN 
     325      IF ( iom_use("sstgrad_"//ttype) .OR. iom_use("sstgrad2_"//ttype) ) THEN 
    301326         DO_2D( 0, 0, 0, 0 )                       ! sst gradient 
    302327            zztmp  = ts(ji,jj,1,jp_tem,Kmm) 
     
    306331               &                 * umask(ji,jj,1) * umask(ji-1,jj,1) * vmask(ji,jj,1) * vmask(ji,jj-1,1) 
    307332         END_2D 
    308          CALL iom_put( "sstgrad2",  z2d )          ! square of module of sst gradient 
    309          IF ( iom_use("sstgrad") ) THEN 
     333         CALL iom_put( "sstgrad2_"//ttype,  z2d )          ! square of module of sst gradient 
     334         IF ( iom_use("sstgrad_"//ttype) ) THEN 
    310335            DO_2D( 0, 0, 0, 0 ) 
    311336               z2d(ji,jj) = SQRT( z2d(ji,jj) ) 
    312337            END_2D 
    313             CALL iom_put( "sstgrad",  z2d )        ! module of sst gradient 
     338            CALL iom_put( "sstgrad_"//ttype,  z2d )        ! module of sst gradient 
    314339         ENDIF 
    315340      ENDIF 
     
    441466      ENDIF 
    442467 
    443       IF( iom_use("tosmint") ) THEN 
     468      IF( (.NOT.l_ldfeiv_time) .AND. ( iom_use('RossRad')  .OR. iom_use('RossRadlim') & 
     469            &                     .OR. iom_use('Tclinic_recip') .OR. iom_use('RR_GS')      & 
     470            &                     .OR. iom_use('aeiu_2d')  .OR. iom_use('aeiv_2d') ) ) THEN 
     471         CALL ldf_eiv(kt, 75.0, z2d, z3d(:,:,1), Kmm) 
     472         CALL iom_put('aeiu_2d', z2d) 
     473         CALL iom_put('aeiv_2d', z3d(:,:,1)) 
     474      ENDIF 
     475 
     476 
     477      IF( iom_use("tosmint_"//ttype) ) THEN 
    444478         z2d(:,:) = 0._wp 
    445479         DO_3D( 0, 0, 0, 0, 1, jpkm1 ) 
    446480            z2d(ji,jj) = z2d(ji,jj) + rho0 * e3t(ji,jj,jk,Kmm) * ts(ji,jj,jk,jp_tem,Kmm) 
    447481         END_3D 
    448          CALL iom_put( "tosmint", z2d )            ! Vertical integral of temperature 
    449       ENDIF 
    450       IF( iom_use("somint") ) THEN 
     482         CALL iom_put( "tosmint_"//ttype, z2d )            ! Vertical integral of temperature 
     483      ENDIF 
     484      IF( iom_use("somint_"//stype) ) THEN 
    451485         z2d(:,:) = 0._wp 
    452486         DO_3D( 0, 0, 0, 0, 1, jpkm1 ) 
    453487            z2d(ji,jj) = z2d(ji,jj) + rho0 * e3t(ji,jj,jk,Kmm) * ts(ji,jj,jk,jp_sal,Kmm) 
    454488         END_3D 
    455          CALL iom_put( "somint", z2d )             ! Vertical integral of salinity 
     489         CALL iom_put( "somint_"//stype, z2d )             ! Vertical integral of salinity 
    456490      ENDIF 
    457491 
  • NEMO/branches/UKMO/NEMO_r4.2RC_GO8_package/src/OCE/DOM/domain.F90

    r15023 r15080  
    316316         &             nn_it000, nn_itend , nn_date0    , nn_time0     , nn_leapy  , nn_istate ,     & 
    317317         &             nn_stock, nn_write , ln_mskland  , ln_clobber   , nn_chunksz, ln_1st_euler  , & 
    318          &             ln_cfmeta, ln_xios_read, nn_wxios 
     318         &             ln_cfmeta, ln_xios_read, nn_wxios, ln_rstdate, ln_rst_eos 
    319319      NAMELIST/namdom/ ln_linssh, rn_Dt, rn_atfp, ln_crs, ln_c1d, ln_meshmask 
    320320      NAMELIST/namtile/ ln_tile, nn_ltile_i, nn_ltile_j 
     
    417417#endif 
    418418         WRITE(numout,*) '      mask land points                ln_mskland      = ', ln_mskland 
     419         WRITE(numout,*) '      date-stamp restart files        ln_rstdate      = ', ln_rstdate 
    419420         WRITE(numout,*) '      additional CF standard metadata ln_cfmeta       = ', ln_cfmeta 
    420421         WRITE(numout,*) '      overwrite an existing file      ln_clobber      = ', ln_clobber 
  • NEMO/branches/UKMO/NEMO_r4.2RC_GO8_package/src/OCE/DOM/dommsk.F90

    r15014 r15080  
    3434   USE lbclnk         ! ocean lateral boundary conditions (or mpp link) 
    3535   USE lib_mpp        ! Massively Parallel Processing library 
     36   USE iom             ! For shlat2d 
     37   USE fldread         ! for sn_shlat2d 
     38 
    3639 
    3740   IMPLICIT NONE 
     
    9396      INTEGER  ::   ios, inum 
    9497      !! 
    95       NAMELIST/namlbc/ rn_shlat, ln_vorlat 
     98      REAL(wp) :: zshlat           !: locally modified shlat for some strait 
     99      REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: zshlat2d 
     100      LOGICAL                         :: ln_shlat2d 
     101      CHARACTER(len = 256)            :: cn_shlat2d_file, cn_shlat2d_var   
     102      NAMELIST/namlbc/ rn_shlat, ln_vorlat, ln_shlat2d, cn_shlat2d_file, cn_shlat2d_var 
    96103      NAMELIST/nambdy/ ln_bdy ,nb_bdy, ln_coords_file, cn_coords_file,         & 
    97104         &             ln_mask_file, cn_mask_file, cn_dyn2d, nn_dyn2d_dta,     & 
     
    118125      ! 
    119126      IF(lwp) WRITE(numout,*) 
    120       IF     (      rn_shlat == 0.               ) THEN   ;   IF(lwp) WRITE(numout,*) '   ==>>>   ocean lateral  free-slip' 
    121       ELSEIF (      rn_shlat == 2.               ) THEN   ;   IF(lwp) WRITE(numout,*) '   ==>>>   ocean lateral  no-slip' 
    122       ELSEIF ( 0. < rn_shlat .AND. rn_shlat < 2. ) THEN   ;   IF(lwp) WRITE(numout,*) '   ==>>>   ocean lateral  partial-slip' 
    123       ELSEIF ( 2. < rn_shlat                     ) THEN   ;   IF(lwp) WRITE(numout,*) '   ==>>>   ocean lateral  strong-slip' 
    124       ELSE 
    125          CALL ctl_stop( 'dom_msk: wrong value for rn_shlat (i.e. a negalive value). We stop.' ) 
     127      IF ( ln_shlat2d ) THEN 
     128         IF(lwp) WRITE(numout,*) '         READ shlat as a 2D coefficient in a file ' 
     129         ALLOCATE( zshlat2d(jpi,jpj) ) 
     130         CALL iom_open(TRIM(cn_shlat2d_file), inum) 
     131         CALL iom_get (inum, jpdom_global, TRIM(cn_shlat2d_var), zshlat2d, 1) ! 
     132         CALL iom_close(inum) 
     133       ELSE 
     134         IF     (      rn_shlat == 0.               ) THEN   ;   IF(lwp) WRITE(numout,*) '   ==>>>   ocean lateral  free-slip' 
     135         ELSEIF (      rn_shlat == 2.               ) THEN   ;   IF(lwp) WRITE(numout,*) '   ==>>>   ocean lateral  no-slip' 
     136         ELSEIF ( 0. < rn_shlat .AND. rn_shlat < 2. ) THEN   ;   IF(lwp) WRITE(numout,*) '   ==>>>   ocean lateral  partial-slip' 
     137         ELSEIF ( 2. < rn_shlat                     ) THEN   ;   IF(lwp) WRITE(numout,*) '   ==>>>   ocean lateral  strong-slip' 
     138         ELSE 
     139            CALL ctl_stop( 'dom_msk: wrong value for rn_shlat (i.e. a negalive value). We stop.' ) 
     140         ENDIF 
    126141      ENDIF 
    127142 
     
    202217      ! Lateral boundary conditions on velocity (modify fmask) 
    203218      ! ---------------------------------------   
    204       IF( rn_shlat /= 0._wp ) THEN      ! Not free-slip lateral boundary condition 
    205          ! 
    206          DO_3D( 0, 0, 0, 0, 1, jpk ) 
    207             IF( fmask(ji,jj,jk) == 0._wp ) THEN 
    208                fmask(ji,jj,jk) = rn_shlat * MIN( 1._wp , MAX( umask(ji,jj,jk), umask(ji,jj+1,jk), & 
    209                   &                                           vmask(ji,jj,jk), vmask(ji+1,jj,jk) ) ) 
    210             ENDIF 
    211          END_3D 
     219      IF( rn_shlat /= 0._wp .or. ln_shlat2d ) THEN      ! Not free-slip lateral boundary condition 
     220         ! 
     221         IF (  ln_shlat2d ) THEN 
     222            DO_3D( 0, 0, 0, 0, 1, jpk ) 
     223               IF( fmask(ji,jj,jk) == 0._wp ) THEN 
     224                  fmask(ji,jj,jk) = zshlat2d(ji,jj) * MIN( 1._wp , MAX( umask(ji,jj,jk), umask(ji,jj+1,jk), & 
     225                     &                                           vmask(ji,jj,jk), vmask(ji+1,jj,jk) ) ) 
     226               ENDIF 
     227            END_3D 
     228         ELSE 
     229            DO_3D( 0, 0, 0, 0, 1, jpk ) 
     230               IF( fmask(ji,jj,jk) == 0._wp ) THEN 
     231                  fmask(ji,jj,jk) = rn_shlat * MIN( 1._wp , MAX( umask(ji,jj,jk), umask(ji,jj+1,jk), & 
     232                     &                                           vmask(ji,jj,jk), vmask(ji+1,jj,jk) ) ) 
     233               ENDIF 
     234            END_3D 
     235         END IF 
     236         ! 
     237         IF( ln_shlat2d ) DEALLOCATE( zshlat2d ) 
     238         ! 
    212239         CALL lbc_lnk( 'dommsk', fmask, 'F', 1._wp )      ! Lateral boundary conditions on fmask 
    213240         ! 
     
    219246      ! --------------------------------  
    220247      ! 
    221       CALL usr_def_fmask( cn_cfg, nn_cfg, fmask ) 
     248      IF ( .not. ln_shlat2d ) THEN       
     249         CALL usr_def_fmask( cn_cfg, nn_cfg, fmask ) 
     250      ENDIF 
    222251      ! 
    223252   END SUBROUTINE dom_msk 
  • NEMO/branches/UKMO/NEMO_r4.2RC_GO8_package/src/OCE/ICB/icbrst.F90

    r13286 r15080  
    2525   USE netcdf         ! netcdf routines for IO 
    2626   USE iom 
     27   USE ioipsl, ONLY : ju2ymds    ! for calendar 
    2728   USE icb_oce        ! define iceberg arrays 
    2829   USE icbutl         ! iceberg utility routines 
     
    191192      INTEGER ::   idg  ! number of digits 
    192193      INTEGER ::   ix_dim, iy_dim, ik_dim, in_dim 
     194      INTEGER                ::   iyear, imonth, iday 
     195      REAL (wp)              ::   zsec 
     196      REAL (wp)              ::   zfjulday 
    193197      CHARACTER(len=256)     :: cl_path 
    194198      CHARACTER(len=256)     :: cl_filename 
    195       CHARACTER(len=8  )     :: cl_kt 
     199      CHARACTER(LEN=20)      ::   cl_kt     ! ocean time-step deine as a character 
    196200      CHARACTER(LEN=12 )     :: clfmt            ! writing format 
    197201      TYPE(iceberg), POINTER :: this 
     
    213217         ! 
    214218         ! file name 
    215          WRITE(cl_kt, '(i8.8)') kt 
    216          cl_filename = TRIM(cexper)//"_"//cl_kt//"_"//TRIM(cn_icbrst_out) 
     219         IF ( ln_rstdate ) THEN 
     220            zfjulday = fjulday + rdt / rday 
     221            IF( ABS(zfjulday - REAL(NINT(zfjulday),wp)) < 0.1 / rday )   zfjulday = REAL(NINT(zfjulday),wp)   ! avoid truncation error 
     222            CALL ju2ymds( zfjulday, iyear, imonth, iday, zsec )            
     223            WRITE(cl_kt, '(i4.4,2i2.2)') iyear, imonth, iday 
     224         ELSE 
     225            IF( kt > 999999999 ) THEN   ;   WRITE(cl_kt, *       ) kt 
     226            ELSE                        ;   WRITE(cl_kt, '(i8.8)') kt 
     227            ENDIF 
     228         ENDIF 
     229         cl_filename = TRIM(cexper)//"_icebergs_"//TRIM(ADJUSTL(cl_kt))//"_restart" 
    217230         IF( lk_mpp ) THEN 
    218231            idg = MAX( INT(LOG10(REAL(MAX(1,jpnij-1),wp))) + 1, 4 )          ! how many digits to we need to write? min=4, max=9 
  • NEMO/branches/UKMO/NEMO_r4.2RC_GO8_package/src/OCE/IOM/in_out_manager.F90

    r14553 r15080  
    2828   LOGICAL       ::   ln_rstart        !: start from (F) rest or (T) a restart file 
    2929   LOGICAL       ::   ln_rst_list      !: output restarts at list of times (T) or by frequency (F) 
     30   LOGICAL       ::   ln_rst_eos       !: check equation of state used for the restart is consistent with model 
    3031   INTEGER       ::   nn_rstctl        !: control of the time step (0, 1 or 2) 
    3132   INTEGER       ::   nn_rstssh   = 0  !: hand made initilization of ssh or not (1/0) 
     
    4041   INTEGER, DIMENSION(10) :: nn_stocklist  !: restart dump times 
    4142   LOGICAL       ::   ln_mskland       !: mask land points in NetCDF outputs (costly: + ~15%) 
     43   LOGICAL       ::   ln_rstdate       !: T=> stamp output restart files with date instead of timestep 
    4244   LOGICAL       ::   ln_cfmeta        !: output additional data to netCDF files required for compliance with the CF metadata standard 
    4345   LOGICAL       ::   ln_clobber       !: clobber (overwrite) an existing file 
  • NEMO/branches/UKMO/NEMO_r4.2RC_GO8_package/src/OCE/IOM/restart.F90

    r14834 r15080  
    3434   USE in_out_manager ! I/O manager 
    3535   USE iom            ! I/O module 
     36   USE ioipsl, ONLY : ju2ymds    ! for calendar 
    3637   USE lib_mpp        ! distribued memory computing library 
    3738 
     
    6768      INTEGER, INTENT(in) ::   kt     ! ocean time-step 
    6869      !! 
     70      INTEGER             ::   iyear, imonth, iday 
     71      REAL (wp)           ::   zsec 
     72      REAL (wp)           ::   zfjulday 
    6973      CHARACTER(LEN=20)   ::   clkt     ! ocean time-step deine as a character 
    7074      CHARACTER(LEN=50)   ::   clname   ! ocean output restart file name 
     
    98102         IF( nitrst <= nitend .AND. nitrst > 0 ) THEN 
    99103            ! beware of the format used to write kt (default is i8.8, that should be large enough...) 
    100             IF( nitrst > 999999999 ) THEN   ;   WRITE(clkt, *       ) nitrst 
    101             ELSE                            ;   WRITE(clkt, '(i8.8)') nitrst 
     104             
     105            IF ( ln_rstdate ) THEN 
     106               zfjulday = fjulday + rdt / rday 
     107               IF( ABS(zfjulday - REAL(NINT(zfjulday),wp)) < 0.1 / rday )   zfjulday = REAL(NINT(zfjulday),wp)   ! avoid truncation error 
     108               CALL ju2ymds( zfjulday, iyear, imonth, iday, zsec )            
     109               WRITE(clkt, '(i4.4,2i2.2)') iyear, imonth, iday 
     110            ELSE 
     111               IF( nitrst > 999999999 ) THEN   ;   WRITE(clkt, *       ) nitrst 
     112               ELSE                            ;   WRITE(clkt, '(i8.8)') nitrst 
     113               ENDIF 
    102114            ENDIF 
    103115            ! create the file 
     
    182194      ENDIF 
    183195 
     196      CALL iom_rstput( kt, nitrst, numrow, 'neos'    , REAL(neos))   ! equation of state 
     197      !CALL iom_rstput( kt, nitrst, numrow, 'neos'    , neos      , ktype = jp_i1, ldxios = lwxios)   ! equation of state 
     198 
    184199      IF( ln_diurnal )   CALL iom_rstput( kt, nitrst, numrow, 'Dsst', x_dsst ) 
    185200      IF( kt == nitrst ) THEN 
     
    266281      INTEGER          , INTENT(in) ::   Kbb, Kmm   ! ocean time level indices 
    267282      INTEGER  ::   jk 
     283      REAL(wp) ::   zeos 
    268284      REAL(wp), DIMENSION(jpi, jpj, jpk) :: w3d 
    269285      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:)   :: zgdept       ! 3D workspace for QCO 
     
    282298         RETURN 
    283299      ENDIF 
     300!      IF ( ln_rst_eos ) THEN 
     301!         ! Check equation of state used is consistent with the restart 
     302!         IF( iom_varid( numror, 'neos') == -1) THEN 
     303!            CALL ctl_stop( 'restart, rst_read: variable neos not found. STOP check that the equations of state in the restart file and in the namelist nameos are consistent and use ln_rst_eos=F') 
     304!         ELSE 
     305!            CALL iom_get( numror, 'neos', zeos, ldxios = lrxios ) 
     306!            IF ( INT(zeos) /= neos ) CALL ctl_stop( 'restart, rst_read: equation of state used in restart file differs from namelist nameos') 
     307!         ENDIF 
     308!      ENDIF 
     309 
    284310      ! 
    285311#if defined key_RK3 
  • NEMO/branches/UKMO/NEMO_r4.2RC_GO8_package/src/OCE/module_example.F90

    r14842 r15080  
    5252#  include "do_loop_substitute.h90" 
    5353   !for other substitutions 
    54 #  include "exampl_substitute.h90" 
    5554   !!---------------------------------------------------------------------- 
    5655   !! NEMO/OCE 4.0 , NEMO Consortium (2018) 
Note: See TracChangeset for help on using the changeset viewer.