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 12724 for NEMO/branches/2020/dev_r12377_KERNEL-06_techene_e3/src/OCE/IOM – NEMO

Ignore:
Timestamp:
2020-04-08T21:37:59+02:00 (4 years ago)
Author:
techene
Message:

branch KERNEL-06 : merge with trunk@12698 #2385 - in duplcated files : changes to comply to the new trunk variables and some loop bug fixes

Location:
NEMO/branches/2020/dev_r12377_KERNEL-06_techene_e3
Files:
5 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/2020/dev_r12377_KERNEL-06_techene_e3

    • Property svn:externals
      •  

        old new  
        33^/utils/build/mk@HEAD         mk 
        44^/utils/tools@HEAD            tools 
        5 ^/vendors/AGRIF/dev_r11615_ENHANCE-04_namelists_as_internalfiles_agrif@HEAD      ext/AGRIF 
         5^/vendors/AGRIF/dev@HEAD      ext/AGRIF 
        66^/vendors/FCM@HEAD            ext/FCM 
        77^/vendors/IOIPSL@HEAD         ext/IOIPSL 
         8 
         9# SETTE 
         10^/utils/CI/sette@HEAD         sette 
  • NEMO/branches/2020/dev_r12377_KERNEL-06_techene_e3/src/OCE/IOM/iom.F90

    r12377 r12724  
    111111      CHARACTER(len=lc) :: clname 
    112112      INTEGER             :: irefyear, irefmonth, irefday 
    113       INTEGER           :: ji, jkmin 
     113      INTEGER           :: ji 
    114114      LOGICAL :: llrst_context              ! is context related to restart 
    115115      ! 
     
    220220           
    221221          ! Add vertical grid bounds 
    222           jkmin = MIN(2,jpk)  ! in case jpk=1 (i.e. sas2D) 
    223           zt_bnds(2,:        ) = gdept_1d(:) 
    224           zt_bnds(1,jkmin:jpk) = gdept_1d(1:jpkm1) 
    225           zt_bnds(1,1        ) = gdept_1d(1) - e3w_1d(1) 
    226           zw_bnds(1,:        ) = gdepw_1d(:) 
    227           zw_bnds(2,1:jpkm1  ) = gdepw_1d(jkmin:jpk) 
    228           zw_bnds(2,jpk:     ) = gdepw_1d(jpk) + e3t_1d(jpk) 
     222          zt_bnds(2,:      ) = gdept_1d(:) 
     223          zt_bnds(1,2:jpk  ) = gdept_1d(1:jpkm1) 
     224          zt_bnds(1,1      ) = gdept_1d(1) - e3w_1d(1) 
     225          zw_bnds(1,:      ) = gdepw_1d(:) 
     226          zw_bnds(2,1:jpkm1) = gdepw_1d(2:jpk) 
     227          zw_bnds(2,jpk:   ) = gdepw_1d(jpk) + e3t_1d(jpk) 
    229228          CALL iom_set_axis_attr(  "deptht", bounds=zw_bnds ) 
    230229          CALL iom_set_axis_attr(  "depthu", bounds=zw_bnds ) 
     
    274273      ! 
    275274      ! set time step length 
    276       dtime%second = rdt 
     275      dtime%second = rn_Dt 
    277276      CALL xios_set_timestep( dtime ) 
    278277      ! 
     
    410409   IF(cdmdl == "OPA") THEN 
    411410!from restart.F90 
    412    CALL iom_set_rstw_var_active("rdt") 
     411   CALL iom_set_rstw_var_active("rn_Dt") 
    413412   IF ( .NOT. ln_diurnal_only ) THEN 
    414413        CALL iom_set_rstw_var_active('ub'  ) 
     
    448447 
    449448        i = 0 
    450         i = i + 1; fields(i)%vname="rdt";            fields(i)%grid="grid_scalar" 
     449        i = i + 1; fields(i)%vname="rn_Dt";            fields(i)%grid="grid_scalar" 
    451450        i = i + 1; fields(i)%vname="un";             fields(i)%grid="grid_N_3D" 
    452451        i = i + 1; fields(i)%vname="ub";             fields(i)%grid="grid_N_3D" 
     
    665664 
    666665 
    667    SUBROUTINE iom_open( cdname, kiomid, ldwrt, kdom, ldstop, ldiof, kdlev ) 
     666   SUBROUTINE iom_open( cdname, kiomid, ldwrt, kdom, ldstop, ldiof, kdlev, cdcomp ) 
    668667      !!--------------------------------------------------------------------- 
    669668      !!                   ***  SUBROUTINE  iom_open  *** 
     
    678677      LOGICAL         , INTENT(in   ), OPTIONAL ::   ldiof    ! Interp On the Fly, needed for AGRIF (default = .FALSE.) 
    679678      INTEGER         , INTENT(in   ), OPTIONAL ::   kdlev    ! number of vertical levels 
     679      CHARACTER(len=3), INTENT(in   ), OPTIONAL ::   cdcomp   ! name of component calling iom_nf90_open 
    680680      ! 
    681681      CHARACTER(LEN=256)    ::   clname    ! the name of the file based on cdname [[+clcpu]+clcpu] 
     
    823823      ENDIF 
    824824      IF( istop == nstop ) THEN   ! no error within this routine 
    825          CALL iom_nf90_open( clname, kiomid, llwrt, llok, idompar, kdlev = kdlev ) 
     825         CALL iom_nf90_open( clname, kiomid, llwrt, llok, idompar, kdlev = kdlev, cdcomp = cdcomp ) 
    826826      ENDIF 
    827827      ! 
     
    23582358            idx = INDEX(clname,'@startdate@') + INDEX(clname,'@STARTDATE@') 
    23592359            DO WHILE ( idx /= 0 )  
    2360                cldate = iom_sdate( fjulday - rdt / rday ) 
     2360               cldate = iom_sdate( fjulday - rn_Dt / rday ) 
    23612361               clname = clname(1:idx-1)//TRIM(cldate)//clname(idx+11:LEN_TRIM(clname)) 
    23622362               idx = INDEX(clname,'@startdate@') + INDEX(clname,'@STARTDATE@') 
     
    23652365            idx = INDEX(clname,'@startdatefull@') + INDEX(clname,'@STARTDATEFULL@') 
    23662366            DO WHILE ( idx /= 0 )  
    2367                cldate = iom_sdate( fjulday - rdt / rday, ldfull = .TRUE. ) 
     2367               cldate = iom_sdate( fjulday - rn_Dt / rday, ldfull = .TRUE. ) 
    23682368               clname = clname(1:idx-1)//TRIM(cldate)//clname(idx+15:LEN_TRIM(clname)) 
    23692369               idx = INDEX(clname,'@startdatefull@') + INDEX(clname,'@STARTDATEFULL@') 
     
    23722372            idx = INDEX(clname,'@enddate@') + INDEX(clname,'@ENDDATE@') 
    23732373            DO WHILE ( idx /= 0 )  
    2374                cldate = iom_sdate( fjulday + rdt / rday * REAL( nitend - nit000, wp ), ld24 = .TRUE. ) 
     2374               cldate = iom_sdate( fjulday + rn_Dt / rday * REAL( nitend - nit000, wp ), ld24 = .TRUE. ) 
    23752375               clname = clname(1:idx-1)//TRIM(cldate)//clname(idx+9:LEN_TRIM(clname)) 
    23762376               idx = INDEX(clname,'@enddate@') + INDEX(clname,'@ENDDATE@') 
     
    23792379            idx = INDEX(clname,'@enddatefull@') + INDEX(clname,'@ENDDATEFULL@') 
    23802380            DO WHILE ( idx /= 0 )  
    2381                cldate = iom_sdate( fjulday + rdt / rday * REAL( nitend - nit000, wp ), ld24 = .TRUE., ldfull = .TRUE. ) 
     2381               cldate = iom_sdate( fjulday + rn_Dt / rday * REAL( nitend - nit000, wp ), ld24 = .TRUE., ldfull = .TRUE. ) 
    23822382               clname = clname(1:idx-1)//TRIM(cldate)//clname(idx+13:LEN_TRIM(clname)) 
    23832383               idx = INDEX(clname,'@enddatefull@') + INDEX(clname,'@ENDDATEFULL@') 
  • NEMO/branches/2020/dev_r12377_KERNEL-06_techene_e3/src/OCE/IOM/iom_def.F90

    r12377 r12724  
    5050   TYPE, PUBLIC ::   file_descriptor 
    5151      CHARACTER(LEN=240)                        ::   name     !: name of the file 
     52      CHARACTER(LEN=3  )                        ::   comp     !: name of component opening the file ('OCE', 'ICE'...) 
    5253      INTEGER                                   ::   nfid     !: identifier of the file (0 if closed) 
    5354                                                              !: jpioipsl option has been removed) 
     
    6465      REAL(kind=wp), DIMENSION(jpmax_vars)      ::   scf      !: scale_factor of the variables 
    6566      REAL(kind=wp), DIMENSION(jpmax_vars)      ::   ofs      !: add_offset of the variables 
    66       INTEGER                                   ::   nlev     ! number of vertical levels 
    6767   END TYPE file_descriptor 
    6868   TYPE(file_descriptor), DIMENSION(jpmax_files), PUBLIC ::   iom_file !: array containing the info for all opened files 
  • NEMO/branches/2020/dev_r12377_KERNEL-06_techene_e3/src/OCE/IOM/iom_nf90.F90

    r12377 r12724  
    1919   !!---------------------------------------------------------------------- 
    2020   USE dom_oce         ! ocean space and time domain 
    21    USE sbc_oce, ONLY: jpka, ght_abl ! abl vertical level number and height 
     21   USE sbc_oce, ONLY: ght_abl ! abl vertical level number and height 
    2222   USE lbclnk          ! lateal boundary condition / mpp exchanges 
    2323   USE iom_def         ! iom variables definitions 
     
    4646CONTAINS 
    4747 
    48    SUBROUTINE iom_nf90_open( cdname, kiomid, ldwrt, ldok, kdompar, kdlev ) 
     48   SUBROUTINE iom_nf90_open( cdname, kiomid, ldwrt, ldok, kdompar, kdlev, cdcomp ) 
    4949      !!--------------------------------------------------------------------- 
    5050      !!                   ***  SUBROUTINE  iom_open  *** 
     
    5858      INTEGER, DIMENSION(2,5), INTENT(in   ), OPTIONAL ::   kdompar     ! domain parameters:  
    5959      INTEGER                , INTENT(in   ), OPTIONAL ::   kdlev       ! size of the ice/abl third dimension 
     60      CHARACTER(len=3)       , INTENT(in   ), OPTIONAL ::   cdcomp      ! name of component calling iom_nf90_open 
    6061 
    6162      CHARACTER(LEN=256) ::   clinfo           ! info character 
    6263      CHARACTER(LEN=256) ::   cltmp            ! temporary character 
     64      CHARACTER(LEN=3  ) ::   clcomp           ! name of component calling iom_nf90_open 
    6365      INTEGER            ::   iln              ! lengths of character 
    6466      INTEGER            ::   istop            ! temporary storage of nstop 
     
    7072      INTEGER            ::   ihdf5            ! local variable for retrieval of value for NF90_HDF5 
    7173      LOGICAL            ::   llclobber        ! local definition of ln_clobber 
    72       INTEGER            ::   ilevels          ! vertical levels 
    7374      !--------------------------------------------------------------------- 
    7475      ! 
     
    7778      ! 
    7879      !                 !number of vertical levels 
    79       IF( PRESENT(kdlev) )   THEN   ;   ilevels = kdlev    ! use input value (useful for sea-ice and abl) 
    80       ELSE                          ;   ilevels = jpk      ! by default jpk 
     80      IF( PRESENT(cdcomp) )   THEN 
     81         IF( .NOT. PRESENT(kdlev) ) CALL ctl_stop( 'iom_nf90_open: cdcomp and kdlev must both be present' ) 
     82         clcomp = cdcomp    ! use input value 
     83      ELSE 
     84         clcomp = 'OCE'     ! by default  
    8185      ENDIF 
    8286      ! 
     
    125129            CALL iom_nf90_check(NF90_SET_FILL( if90id, NF90_NOFILL,                   idmy ), clinfo) 
    126130            ! define dimensions 
    127             CALL iom_nf90_check(NF90_DEF_DIM( if90id,            'x',   kdompar(1,1), idmy ), clinfo) 
    128             CALL iom_nf90_check(NF90_DEF_DIM( if90id,            'y',   kdompar(2,1), idmy ), clinfo) 
    129             IF( PRESENT(kdlev) ) THEN 
    130               IF( kdlev == jpka ) THEN 
    131                  CALL iom_nf90_check(NF90_DEF_DIM( if90id, 'nav_lev',          kdlev, idmy ), clinfo) 
    132                  CALL iom_nf90_check(NF90_DEF_DIM( if90id, 'time_counter', NF90_UNLIMITED, idmy ), clinfo) 
    133               ELSE 
    134                  CALL iom_nf90_check(NF90_DEF_DIM( if90id, 'nav_lev',            jpk, idmy ), clinfo) 
    135                  CALL iom_nf90_check(NF90_DEF_DIM( if90id, 'time_counter', NF90_UNLIMITED, idmy ), clinfo) 
    136                  CALL iom_nf90_check(NF90_DEF_DIM( if90id,  'numcat',          kdlev, idmy ), clinfo) 
    137               ENDIF 
    138             ELSE 
    139                CALL iom_nf90_check(NF90_DEF_DIM( if90id, 'nav_lev',            jpk, idmy ), clinfo) 
    140                CALL iom_nf90_check(NF90_DEF_DIM( if90id, 'time_counter', NF90_UNLIMITED, idmy ), clinfo) 
    141             ENDIF 
     131                               CALL iom_nf90_check(NF90_DEF_DIM( if90id,            'x',   kdompar(1,1), idmy ), clinfo) 
     132                               CALL iom_nf90_check(NF90_DEF_DIM( if90id,            'y',   kdompar(2,1), idmy ), clinfo) 
     133            SELECT CASE (clcomp) 
     134            CASE ('OCE')   ;   CALL iom_nf90_check(NF90_DEF_DIM( if90id,      'nav_lev',            jpk, idmy ), clinfo) 
     135            CASE ('ICE')   ;   CALL iom_nf90_check(NF90_DEF_DIM( if90id,       'numcat',          kdlev, idmy ), clinfo) 
     136            CASE ('ABL')   ;   CALL iom_nf90_check(NF90_DEF_DIM( if90id,      'nav_lev',          kdlev, idmy ), clinfo) 
     137            CASE ('SED')   ;   CALL iom_nf90_check(NF90_DEF_DIM( if90id,       'numsed',          kdlev, idmy ), clinfo) 
     138            CASE DEFAULT   ;   CALL ctl_stop( 'iom_nf90_open unknown component type' ) 
     139            END SELECT 
     140                               CALL iom_nf90_check(NF90_DEF_DIM( if90id, 'time_counter', NF90_UNLIMITED, idmy ), clinfo) 
    142141            ! global attributes 
    143142            CALL iom_nf90_check(NF90_PUT_ATT( if90id, NF90_GLOBAL, 'DOMAIN_number_total'   , jpnij              ), clinfo) 
     
    165164         ENDDO 
    166165         iom_file(kiomid)%name   = TRIM(cdname) 
     166         iom_file(kiomid)%comp   = clcomp 
    167167         iom_file(kiomid)%nfid   = if90id 
    168168         iom_file(kiomid)%nvars  = 0 
    169169         iom_file(kiomid)%irec   = -1   ! useless for NetCDF files, used to know if the file is in define mode  
    170          iom_file(kiomid)%nlev   = ilevels 
    171170         CALL iom_nf90_check(NF90_Inquire(if90id, unlimitedDimId = iom_file(kiomid)%iduld), clinfo) 
    172171         IF( iom_file(kiomid)%iduld .GE. 0 ) THEN 
     
    529528      INTEGER, DIMENSION(4) :: idimid               ! dimensions id 
    530529      CHARACTER(LEN=256)    :: clinfo               ! info character 
    531       CHARACTER(LEN= 12), DIMENSION(5) :: cltmp     ! temporary character 
    532530      INTEGER               :: if90id               ! nf90 file identifier 
    533       INTEGER               :: idmy                 ! dummy variable 
    534531      INTEGER               :: itype                ! variable type 
    535532      INTEGER, DIMENSION(4) :: ichunksz             ! NetCDF4 chunk sizes. Will be computed using 
     
    540537      !                                             ! when appropriate (currently chunking is applied to 4d fields only) 
    541538      INTEGER               :: idlv                 ! local variable 
    542       INTEGER               :: idim3                ! id of the third dimension 
    543539      !--------------------------------------------------------------------- 
    544540      ! 
     
    554550         ENDIF 
    555551         ! define the dimension variables if it is not already done 
    556          ! Warning: we must use the same character length in an array constructor (at least for gcc compiler) 
    557          cltmp = (/ 'nav_lon     ', 'nav_lat     ', 'nav_lev     ', 'time_counter', 'numcat      ' /)    
    558          CALL iom_nf90_check(NF90_DEF_VAR( if90id, TRIM(cltmp(1)), NF90_FLOAT , (/ 1, 2 /), iom_file(kiomid)%nvid(1) ), clinfo) 
    559          CALL iom_nf90_check(NF90_DEF_VAR( if90id, TRIM(cltmp(2)), NF90_FLOAT , (/ 1, 2 /), iom_file(kiomid)%nvid(2) ), clinfo) 
    560          CALL iom_nf90_check(NF90_DEF_VAR( if90id, TRIM(cltmp(3)), NF90_FLOAT , (/ 3    /), iom_file(kiomid)%nvid(3) ), clinfo) 
    561          CALL iom_nf90_check(NF90_DEF_VAR( if90id, TRIM(cltmp(4)), NF90_DOUBLE, (/ 4    /), iom_file(kiomid)%nvid(4) ), clinfo) 
     552         DO jd = 1, 2 
     553            CALL iom_nf90_check(NF90_INQUIRE_DIMENSION(if90id,jd,iom_file(kiomid)%cn_var(jd),iom_file(kiomid)%dimsz(jd,jd)),clinfo) 
     554            CALL iom_nf90_check(NF90_DEF_VAR( if90id, TRIM(iom_file(kiomid)%cn_var(jd)), NF90_FLOAT , (/ 1, 2 /),   & 
     555               &                              iom_file(kiomid)%nvid(jd) ), clinfo) 
     556         END DO 
     557         iom_file(kiomid)%dimsz(2,1) = iom_file(kiomid)%dimsz(2,2)   ! second dim of first  variable 
     558         iom_file(kiomid)%dimsz(1,2) = iom_file(kiomid)%dimsz(1,1)   ! first  dim of second variable 
     559         DO jd = 3, 4 
     560            CALL iom_nf90_check(NF90_INQUIRE_DIMENSION(if90id,jd,iom_file(kiomid)%cn_var(jd),iom_file(kiomid)%dimsz(1,jd)), clinfo) 
     561            CALL iom_nf90_check(NF90_DEF_VAR( if90id, TRIM(iom_file(kiomid)%cn_var(jd)), NF90_FLOAT , (/ jd   /),   & 
     562               &                              iom_file(kiomid)%nvid(jd) ), clinfo) 
     563         END DO 
    562564         ! update informations structure related the dimension variable we just added... 
    563565         iom_file(kiomid)%nvars       = 4 
    564566         iom_file(kiomid)%luld(1:4)   = (/ .FALSE., .FALSE., .FALSE., .TRUE. /) 
    565          iom_file(kiomid)%cn_var(1:4) = cltmp(1:4) 
    566567         iom_file(kiomid)%ndims(1:4)  = (/ 2, 2, 1, 1 /) 
    567          IF( NF90_INQ_DIMID( if90id, 'numcat', idmy ) == nf90_noerr ) THEN   ! add a 5th variable corresponding to the 5th dimension 
    568             CALL iom_nf90_check(NF90_DEF_VAR( if90id, TRIM(cltmp(5)), NF90_FLOAT , (/ 5 /), iom_file(kiomid)%nvid(5) ), clinfo) 
    569             iom_file(kiomid)%nvars     = 5 
    570             iom_file(kiomid)%luld(5)   = .FALSE. 
    571             iom_file(kiomid)%cn_var(5) = cltmp(5) 
    572             iom_file(kiomid)%ndims(5)  = 1 
    573          ENDIF 
    574          ! trick: defined to 0 to say that dimension variables are defined but not yet written 
    575          iom_file(kiomid)%dimsz(1, 1)  = 0    
    576568         IF(lwp) WRITE(numout,*) TRIM(clinfo)//' define dimension variables done' 
    577569      ENDIF 
     
    594586         IF(     PRESENT(pv_r0d) ) THEN   ;   idims = 0 
    595587         ELSEIF( PRESENT(pv_r1d) ) THEN 
    596             IF(( SIZE(pv_r1d,1) == jpk ).OR.( SIZE(pv_r1d,1) == jpka )) THEN   ;   idim3 = 3 
    597             ELSE                                                               ;   idim3 = 5 
    598             ENDIF 
    599                                               idims = 2   ;   idimid(1:idims) = (/idim3,4/) 
    600          ELSEIF( PRESENT(pv_r2d) ) THEN   ;   idims = 3   ;   idimid(1:idims) = (/1,2  ,4/) 
     588                                              idims = 2   ;   idimid(1:idims) = (/3,4/) 
     589         ELSEIF( PRESENT(pv_r2d) ) THEN   ;   idims = 3   ;   idimid(1:idims) = (/1,2,4/) 
    601590         ELSEIF( PRESENT(pv_r3d) ) THEN 
    602             IF(( SIZE(pv_r3d,3) == jpk ).OR.( SIZE(pv_r3d,3) == jpka )) THEN   ;   idim3 = 3 
    603             ELSE                                                               ;   idim3 = 5 
    604             ENDIF 
    605                                               idims = 4   ;   idimid(1:idims) = (/1,2,idim3,4/) 
     591                                              idims = 4   ;   idimid(1:idims) = (/1,2,3,4/) 
    606592         ENDIF 
    607593         IF( PRESENT(ktype) ) THEN   ! variable external type 
     
    678664            ! ============= 
    679665            ! trick: is defined to 0 => dimension variable are defined but not yet written 
    680             IF( iom_file(kiomid)%dimsz(1, 1) == 0 ) THEN 
    681                CALL iom_nf90_check( NF90_INQ_VARID( if90id, 'nav_lon'     , idmy )         , clinfo ) 
    682                CALL iom_nf90_check( NF90_PUT_VAR  ( if90id, idmy, glamt(ix1:ix2, iy1:iy2) ), clinfo ) 
    683                CALL iom_nf90_check( NF90_INQ_VARID( if90id, 'nav_lat'     , idmy )         , clinfo ) 
    684                CALL iom_nf90_check( NF90_PUT_VAR  ( if90id, idmy, gphit(ix1:ix2, iy1:iy2) ), clinfo ) 
    685                CALL iom_nf90_check( NF90_INQ_VARID( if90id, 'nav_lev'     , idmy ), clinfo ) 
    686                IF (iom_file(kiomid)%nlev == jpka) THEN   ;   CALL iom_nf90_check( NF90_PUT_VAR  ( if90id, idmy,  ght_abl), clinfo ) 
    687                ELSE                                      ;   CALL iom_nf90_check( NF90_PUT_VAR  ( if90id, idmy, gdept_1d), clinfo ) 
    688                ENDIF 
    689                IF( NF90_INQ_VARID( if90id, 'numcat', idmy ) == nf90_noerr ) THEN 
    690                   CALL iom_nf90_check( NF90_PUT_VAR  ( if90id, idmy, (/ (idlv, idlv = 1,iom_file(kiomid)%nlev) /)), clinfo ) 
    691                ENDIF 
    692                ! +++ WRONG VALUE: to be improved but not really useful... 
    693                CALL iom_nf90_check( NF90_INQ_VARID( if90id, 'time_counter', idmy ), clinfo ) 
    694                CALL iom_nf90_check( NF90_PUT_VAR( if90id, idmy, kt                      ), clinfo )    
    695                ! update the values of the variables dimensions size 
    696                CALL iom_nf90_check( NF90_INQUIRE_DIMENSION( if90id, 1, len = iom_file(kiomid)%dimsz(1,1) ), clinfo ) 
    697                CALL iom_nf90_check( NF90_INQUIRE_DIMENSION( if90id, 2, len = iom_file(kiomid)%dimsz(2,1) ), clinfo ) 
    698                iom_file(kiomid)%dimsz(1:2, 2) = iom_file(kiomid)%dimsz(1:2, 1) 
    699                CALL iom_nf90_check( NF90_INQUIRE_DIMENSION( if90id, 3, len = iom_file(kiomid)%dimsz(1,3) ), clinfo ) 
    700                iom_file(kiomid)%dimsz(1  , 4) = 1   ! unlimited dimension 
     666            IF( iom_file(kiomid)%dimsz(1, 4) == 0 ) THEN   ! time_counter = 0 
     667               CALL iom_nf90_check(    NF90_PUT_VAR( if90id, 1,                            glamt(ix1:ix2, iy1:iy2) ), clinfo ) 
     668               CALL iom_nf90_check(    NF90_PUT_VAR( if90id, 2,                            gphit(ix1:ix2, iy1:iy2) ), clinfo ) 
     669               SELECT CASE (iom_file(kiomid)%comp) 
     670               CASE ('OCE')   
     671                  CALL iom_nf90_check( NF90_PUT_VAR( if90id, 3,                                           gdept_1d ), clinfo ) 
     672               CASE ('ABL') 
     673                  CALL iom_nf90_check( NF90_PUT_VAR( if90id, 3,                                            ght_abl ), clinfo ) 
     674               CASE DEFAULT 
     675                  CALL iom_nf90_check( NF90_PUT_VAR( if90id, 3, (/ (idlv, idlv = 1,iom_file(kiomid)%dimsz(1,3)) /) ), clinfo ) 
     676               END SELECT 
     677               ! "wrong" value: to be improved but not really useful... 
     678               CALL iom_nf90_check(   NF90_PUT_VAR( if90id, 4,                                                  kt ), clinfo )    
     679               ! update the size of the variable corresponding to the unlimited dimension 
     680               iom_file(kiomid)%dimsz(1, 4) = 1   ! so we don't enter this IF case any more... 
    701681               IF(lwp) WRITE(numout,*) TRIM(clinfo)//' write dimension variables done' 
    702682            ENDIF 
  • NEMO/branches/2020/dev_r12377_KERNEL-06_techene_e3/src/OCE/IOM/restart.F90

    r12680 r12724  
    144144      !!---------------------------------------------------------------------- 
    145145                     IF(lwxios) CALL iom_swap(      cwxios_context          ) 
    146                      CALL iom_rstput( kt, nitrst, numrow, 'rdt'    , rdt       , ldxios = lwxios)   ! dynamics time step 
     146                     CALL iom_rstput( kt, nitrst, numrow, 'rdt'    , rn_Dt       , ldxios = lwxios)   ! dynamics time step 
    147147                     CALL iom_delay_rst( 'WRITE', 'OCE', numrow )   ! save only ocean delayed global communication variables 
    148148 
     
    247247      IF( iom_varid( numror, 'rdt', ldstop = .FALSE. ) > 0 )   THEN 
    248248         CALL iom_get( numror, 'rdt', zrdt, ldxios = lrxios ) 
    249          IF( zrdt /= rdt )   neuler = 0 
     249         IF( zrdt /= rn_Dt ) THEN 
     250            IF(lwp) WRITE( numout,*) 
     251            IF(lwp) WRITE( numout,*) 'rst_read:  rdt not equal to the read one' 
     252            IF(lwp) WRITE( numout,*) 
     253            IF(lwp) WRITE( numout,*) '      ==>>>   forced euler first time-step' 
     254            l_1st_euler =  .TRUE. 
     255         ENDIF 
    250256      ENDIF 
    251257 
     
    256262      IF ( ln_diurnal_only ) THEN  
    257263         IF(lwp) WRITE( numout, * ) & 
    258          &   "rst_read:- ln_diurnal_only set, setting rhop=rau0"  
    259          rhop = rau0 
     264         &   "rst_read:- ln_diurnal_only set, setting rhop=rho0"  
     265         rhop = rho0 
    260266         CALL iom_get( numror, jpdom_autoglo, 'tn'     , w3d, ldxios = lrxios )  
    261267         ts(:,:,1,jp_tem,Kmm) = w3d(:,:,1) 
     
    270276         CALL iom_get( numror, jpdom_autoglo, 'sshb'   ,ssh(:,:         ,Kbb), ldxios = lrxios ) 
    271277      ELSE 
    272          neuler = 0 
     278         l_1st_euler =  .TRUE.      ! before field not found, forced euler 1st time-step 
    273279      ENDIF 
    274280      ! 
     
    284290      ENDIF 
    285291      ! 
    286       IF( neuler == 0 ) THEN                                  ! Euler restart (neuler=0) 
     292      IF( l_1st_euler ) THEN                                  ! Euler restart  
    287293         ts (:,:,:,:,Kbb) = ts (:,:,:,:,Kmm)                  ! all before fields set to now values 
    288294         uu (:,:,:  ,Kbb) = uu (:,:,:  ,Kmm) 
Note: See TracChangeset for help on using the changeset viewer.