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

Changeset 1334 for trunk/NEMO/OPA_SRC


Ignore:
Timestamp:
2009-03-03T15:07:48+01:00 (15 years ago)
Author:
smasson
Message:

complete work on time origin in outputs (ticket:335) + downward vertical axis (ticket:357)

Location:
trunk/NEMO/OPA_SRC
Files:
6 edited

Legend:

Unmodified
Added
Removed
  • trunk/NEMO/OPA_SRC/DIA/diagap.F90

    r1317 r1334  
    8282      !! * local declarations 
    8383      INTEGER ::   ji, jj, jk   ! dummy loop indices 
    84       INTEGER ::   it 
     84      INTEGER ::   it, itmod    ! time step indices 
    8585      CHARACTER (len=40) ::   clhstnam 
    8686      CHARACTER (len=40) ::   clop 
     
    170170         ! Horizontal grid : zphi() 
    171171         CALL histbeg( clhstnam, 1, zfoo, 1, zloo,   & 
    172                        1, 1, 1, 1, 0, zjulian, zdt, nhoridg, numgap, domain_id=nidom ) 
     172                       1, 1, 1, 1, nit000-1, zjulian, zdt, nhoridg, numgap, domain_id=nidom ) 
    173173         ! Vertical grids : gdept, gdepw 
    174174         CALL histvert( numgap, "deptht", "Vertical T levels",   & 
    175                          "m", jpk, gdept_0, ndepidg ) 
     175                         "m", jpk, gdept_0, ndepidg, "down" ) 
    176176 
    177177         ! define fields to be stored 
     
    198198      ! ---------------------- 
    199199 
    200       it = kt - nit000 + 1       ! define time axis 
    201       IF( MOD( it, ngap ) == 0 ) THEN 
     200      itmod = kt - nit000 + 1       ! define time axis 
     201      it = kt  
     202      IF( MOD( itmod, ngap ) == 0 ) THEN 
    202203 
    203204         ! initialization 
     
    240241          ! ----------------------------====== 
    241242 
    242           IF( MOD( it, nprg ) == 0 ) THEN 
     243          IF( MOD( itmod, nprg ) == 0 ) THEN 
    243244              IF(lwp) THEN 
    244245                  WRITE(numout,*) 'dia_gap: time step = ', kt, 'model - data' 
  • trunk/NEMO/OPA_SRC/DIA/diaptr.F90

    r1317 r1334  
    418418 
    419419      CHARACTER (len=40)       ::   clhstnam, clop                   ! temporary names 
    420       INTEGER                  ::   iline, it, ji                    ! 
     420      INTEGER                  ::   iline, it, ji, itmod             ! 
    421421      REAL(wp)                 ::   zsto, zout, zdt, zjulian   ! temporary scalars 
    422422      REAL(wp), DIMENSION(jpj) ::   zphi, zfoo 
     
    424424       
    425425      ! define time axis 
    426       it = kt - nit000 + 1 
     426      it = kt 
     427      itmod = kt - nit000 + 1 
    427428 
    428429      ! Initialization 
     
    481482         ! Horizontal grid : zphi() 
    482483         CALL histbeg(clhstnam, 1, zfoo, jpj, zphi,   & 
    483             1, 1, 1, jpj, 0, zjulian, zdt, nhoridz, numptr, domain_id=nidom ) 
     484            1, 1, 1, jpj, nit000-1, zjulian, zdt, nhoridz, numptr, domain_id=nidom ) 
    484485         ! Vertical grids : gdept_0, gdepw_0 
    485486         CALL histvert( numptr, "deptht", "Vertical T levels",   & 
    486             "m", jpk, gdept_0, ndepidzt ) 
     487            "m", jpk, gdept_0, ndepidzt, "down" ) 
    487488         CALL histvert( numptr, "depthw", "Vertical W levels",   & 
    488             "m", jpk, gdepw_0, ndepidzw ) 
     489            "m", jpk, gdepw_0, ndepidzw, "down" ) 
    489490          
    490491         !  Zonal mean T and S 
     
    555556      ENDIF 
    556557 
    557       IF( MOD( it, nf_ptr ) == 0 ) THEN 
     558      IF( MOD( itmod, nf_ptr ) == 0 ) THEN 
    558559 
    559560         IF(lwp) THEN 
  • trunk/NEMO/OPA_SRC/DIA/diawri.F90

    r1318 r1334  
    109109      INTEGER ::   inum = 11             ! temporary logical unit 
    110110      INTEGER ::   & 
    111          iimi, iima, ipk, it,         &  ! temporary integers 
     111         iimi, iima, ipk, it, itmod,  &  ! temporary integers 
    112112         ijmi, ijma                      !    "          " 
    113113      REAL(wp) ::   & 
     
    148148 
    149149      ! define time axis 
    150       it = kt - nit000 + 1 
     150      it = kt 
     151      itmod = kt - nit000 + 1 
    151152 
    152153 
     
    184185         CALL histbeg( clhstnam, jpi, glamt, jpj, gphit,           &  ! Horizontal grid: glamt and gphit 
    185186            &          iimi, iima-iimi+1, ijmi, ijma-ijmi+1,       & 
    186             &          0, zjulian, zdt, nh_T, nid_T, domain_id=nidom ) 
     187            &          nit000-1, zjulian, zdt, nh_T, nid_T, domain_id=nidom ) 
    187188         CALL histvert( nid_T, "deptht", "Vertical T levels",      &  ! Vertical grid: gdept 
    188             &           "m", ipk, gdept_0, nz_T ) 
     189            &           "m", ipk, gdept_0, nz_T, "down" ) 
    189190         !                                                            ! Index of ocean points 
    190191         CALL wheneq( jpi*jpj*ipk, tmask, 1, 1., ndex_T , ndim_T  )      ! volume 
     
    197198         CALL histbeg( clhstnam, jpi, glamu, jpj, gphiu,           &  ! Horizontal grid: glamu and gphiu 
    198199            &          iimi, iima-iimi+1, ijmi, ijma-ijmi+1,       & 
    199             &          0, zjulian, zdt, nh_U, nid_U, domain_id=nidom ) 
     200            &          nit000-1, zjulian, zdt, nh_U, nid_U, domain_id=nidom ) 
    200201         CALL histvert( nid_U, "depthu", "Vertical U levels",      &  ! Vertical grid: gdept 
    201             &           "m", ipk, gdept_0, nz_U ) 
     202            &           "m", ipk, gdept_0, nz_U, "down" ) 
    202203         !                                                            ! Index of ocean points 
    203204         CALL wheneq( jpi*jpj*ipk, umask, 1, 1., ndex_U , ndim_U  )      ! volume 
     
    210211         CALL histbeg( clhstnam, jpi, glamv, jpj, gphiv,           &  ! Horizontal grid: glamv and gphiv 
    211212            &          iimi, iima-iimi+1, ijmi, ijma-ijmi+1,       & 
    212             &          0, zjulian, zdt, nh_V, nid_V, domain_id=nidom ) 
     213            &          nit000-1, zjulian, zdt, nh_V, nid_V, domain_id=nidom ) 
    213214         CALL histvert( nid_V, "depthv", "Vertical V levels",      &  ! Vertical grid : gdept 
    214             &          "m", ipk, gdept_0, nz_V ) 
     215            &          "m", ipk, gdept_0, nz_V, "down" ) 
    215216         !                                                            ! Index of ocean points 
    216217         CALL wheneq( jpi*jpj*ipk, vmask, 1, 1., ndex_V , ndim_V  )      ! volume 
     
    223224         CALL histbeg( clhstnam, jpi, glamt, jpj, gphit,           &  ! Horizontal grid: glamt and gphit 
    224225            &          iimi, iima-iimi+1, ijmi, ijma-ijmi+1,       & 
    225             &          0, zjulian, zdt, nh_W, nid_W, domain_id=nidom ) 
     226            &          nit000-1, zjulian, zdt, nh_W, nid_W, domain_id=nidom ) 
    226227         CALL histvert( nid_W, "depthw", "Vertical W levels",      &  ! Vertical grid: gdepw 
    227             &          "m", ipk, gdepw_0, nz_W ) 
     228            &          "m", ipk, gdepw_0, nz_W, "down" ) 
    228229 
    229230 
     
    406407      ! donne le nombre d'elements, et ndex la liste des indices a sortir 
    407408 
    408       IF( lwp .AND. MOD( it, nwrite ) == 0 ) THEN  
     409      IF( lwp .AND. MOD( itmod, nwrite ) == 0 ) THEN  
    409410         WRITE(numout,*) 'dia_wri : write model outputs in NetCDF files at ', kt, 'time-step' 
    410411         WRITE(numout,*) '~~~~~~ ' 
     
    512513 
    513514      !  Create an output files (output.abort.nc) if S < 0 or u > 20 m/s 
    514       IF( kindic < 0 )   CALL dia_wri_state( 'output.abort' ) 
     515      IF( kindic < 0 )   CALL dia_wri_state( 'output.abort', kt ) 
    515516 
    516517      ! 3. Close all files 
     
    526527 
    527528 
    528    SUBROUTINE dia_wri_state( cdfile_name ) 
     529   SUBROUTINE dia_wri_state( cdfile_name, kt ) 
    529530      !!--------------------------------------------------------------------- 
    530531      !!                 ***  ROUTINE dia_wri_state  *** 
     
    546547      !!---------------------------------------------------------------------- 
    547548      !! * Arguments 
    548       CHARACTER (len=* ), INTENT( in ) ::   & 
    549          cdfile_name      ! name of the file created 
     549      CHARACTER (len=* ), INTENT( in ) ::   cdfile_name      ! name of the file created 
     550      INTEGER           , INTENT( in ) ::   kt               ! ocean time-step index 
    550551 
    551552      !! * Local declarations 
     
    588589      zjulian = zjulian - adatrj   !   set calendar origin to the beginning of the experiment 
    589590      CALL histbeg( clname, jpi, glamt, jpj, gphit,   & 
    590           1, jpi, 1, jpj, 0, zjulian, zdt, nh_i, id_i, domain_id=nidom )          ! Horizontal grid : glamt and gphit 
     591          1, jpi, 1, jpj, nit000-1, zjulian, zdt, nh_i, id_i, domain_id=nidom )          ! Horizontal grid : glamt and gphit 
    591592      CALL histvert( id_i, "deptht", "Vertical T levels",   &    ! Vertical grid : gdept 
    592           "m", jpk, gdept_0, nz_i) 
     593          "m", jpk, gdept_0, nz_i, "down") 
    593594 
    594595      ! Declare all the output fields as NetCDF variables 
     
    634635 
    635636      ! Write all fields on T grid 
    636       CALL histwrite( id_i, "votemper", 1, tn      , jpi*jpj*jpk, idex )    ! now temperature 
    637       CALL histwrite( id_i, "vosaline", 1, sn      , jpi*jpj*jpk, idex )    ! now salinity 
    638 #if defined key_dynspg_rl 
    639       CALL histwrite( id_i, "sobarstf", 1, bsfn     , jpi*jpj    , idex )    ! barotropic streamfunction 
     637      CALL histwrite( id_i, "votemper", kt, tn      , jpi*jpj*jpk, idex )    ! now temperature 
     638      CALL histwrite( id_i, "vosaline", kt, sn      , jpi*jpj*jpk, idex )    ! now salinity 
     639#if defined key_dynspg_rl 
     640      CALL histwrite( id_i, "sobarstf", kt, bsfn     , jpi*jpj    , idex )    ! barotropic streamfunction 
    640641#else 
    641       CALL histwrite( id_i, "sossheig", 1, sshn     , jpi*jpj    , idex )    ! sea surface height 
    642 #endif 
    643       CALL histwrite( id_i, "vozocrtx", 1, un       , jpi*jpj*jpk, idex )    ! now i-velocity 
    644       CALL histwrite( id_i, "vomecrty", 1, vn       , jpi*jpj*jpk, idex )    ! now j-velocity 
    645       CALL histwrite( id_i, "vovecrtz", 1, wn       , jpi*jpj*jpk, idex )    ! now k-velocity 
    646       CALL histwrite( id_i, "sowaflup", 1, emp      , jpi*jpj    , idex )    ! freshwater budget 
    647       CALL histwrite( id_i, "sohefldo", 1, qsr + qns, jpi*jpj    , idex )    ! total heat flux 
    648       CALL histwrite( id_i, "soshfldo", 1, qsr      , jpi*jpj    , idex )    ! solar heat flux 
    649       CALL histwrite( id_i, "soicecov", 1, fr_i     , jpi*jpj    , idex )    ! ice fraction 
    650       CALL histwrite( id_i, "sozotaux", 1, utau     , jpi*jpj    , idex )    ! i-wind stress 
    651       CALL histwrite( id_i, "sometauy", 1, vtau     , jpi*jpj    , idex )    ! j-wind stress 
     642      CALL histwrite( id_i, "sossheig", kt, sshn     , jpi*jpj    , idex )    ! sea surface height 
     643#endif 
     644      CALL histwrite( id_i, "vozocrtx", kt, un       , jpi*jpj*jpk, idex )    ! now i-velocity 
     645      CALL histwrite( id_i, "vomecrty", kt, vn       , jpi*jpj*jpk, idex )    ! now j-velocity 
     646      CALL histwrite( id_i, "vovecrtz", kt, wn       , jpi*jpj*jpk, idex )    ! now k-velocity 
     647      CALL histwrite( id_i, "sowaflup", kt, emp      , jpi*jpj    , idex )    ! freshwater budget 
     648      CALL histwrite( id_i, "sohefldo", kt, qsr + qns, jpi*jpj    , idex )    ! total heat flux 
     649      CALL histwrite( id_i, "soshfldo", kt, qsr      , jpi*jpj    , idex )    ! solar heat flux 
     650      CALL histwrite( id_i, "soicecov", kt, fr_i     , jpi*jpj    , idex )    ! ice fraction 
     651      CALL histwrite( id_i, "sozotaux", kt, utau     , jpi*jpj    , idex )    ! i-wind stress 
     652      CALL histwrite( id_i, "sometauy", kt, vtau     , jpi*jpj    , idex )    ! j-wind stress 
    652653 
    653654      ! 3. Close the file 
  • trunk/NEMO/OPA_SRC/TRD/trdmld.F90

    r1317 r1334  
    230230      INTEGER, INTENT( in ) ::   kt   ! ocean time-step index 
    231231      !! 
    232       INTEGER :: ji, jj, jk, jl, ik, it 
     232      INTEGER :: ji, jj, jk, jl, ik, it, itmod 
    233233      LOGICAL :: lldebug = .TRUE. 
    234234      REAL(wp) :: zavt, zfn, zfn2 
     
    389389      smltrd(:,:,:) = smltrd(:,:,:) * ucf   !  is no longer used, and is reset to 0. at next time step) 
    390390       
    391       it = kt - nit000 + 1 
    392  
    393       MODULO_NTRD : IF( MOD( it, ntrd ) == 0 ) THEN        ! nitend MUST be multiple of ntrd 
     391      ! define time axis 
     392      it = kt 
     393      itmod = kt - nit000 + 1 
     394 
     395      MODULO_NTRD : IF( MOD( itmod, ntrd ) == 0 ) THEN        ! nitend MUST be multiple of ntrd 
    394396         ! 
    395397         ztmltot (:,:) = 0.e0   ;   zsmltot (:,:) = 0.e0   ! reset arrays to zero 
     
    576578#if defined key_dimgout 
    577579 
    578       IF( MOD( it, ntrd ) == 0 ) THEN 
     580      IF( MOD( itmod, ntrd ) == 0 ) THEN 
    579581         iyear =  ndastp/10000 
    580582         imon  = (ndastp-iyear*10000)/100 
     
    593595      ! ---------------------------------- 
    594596 
    595       IF( lwp .AND. MOD( it , ntrd ) == 0 ) THEN 
     597      IF( lwp .AND. MOD( itmod , ntrd ) == 0 ) THEN 
    596598         WRITE(numout,*) ' ' 
    597599         WRITE(numout,*) 'trd_mld : write trends in the NetCDF file :' 
     
    683685#endif 
    684686 
    685       IF( MOD( it, ntrd ) == 0 ) THEN 
     687      IF( MOD( itmod, ntrd ) == 0 ) THEN 
    686688         ! 
    687689         ! III.5 Reset cumulative arrays to zero 
     
    876878      IF(lwp) WRITE(numout,*) ' Name of NETCDF file ', clhstnam 
    877879      CALL histbeg( clhstnam, jpi, glamt, jpj, gphit,                                            & 
    878       &             1, jpi, 1, jpj, 0, zjulian, rdt, nh_t, nidtrd, domain_id=nidom ) 
     880      &             1, jpi, 1, jpj, nit000-1, zjulian, rdt, nh_t, nidtrd, domain_id=nidom ) 
    879881 
    880882      !-- Define the ML depth variable 
  • trunk/NEMO/OPA_SRC/TRD/trdvor.F90

    r1317 r1334  
    312312      INTEGER, INTENT( in ) ::   kt   ! ocean time-step index 
    313313      !! 
    314       INTEGER  ::   ji, jj, jk, jl, it 
     314      INTEGER  ::   ji, jj, jk, jl, it, itmod 
    315315      REAL(wp) ::   zmean 
    316316      REAL(wp), DIMENSION(jpi,jpj) ::   zun, zvn 
     
    406406 
    407407      ! define time axis 
    408       it= kt - nit000 + 1 
     408      it = kt 
     409      itmod = kt - nit000 + 1 
    409410 
    410411      IF( MOD( it, ntrd ) == 0 ) THEN 
     
    455456      IF( kt >=  nit000+1 ) THEN 
    456457 
    457          IF( lwp .AND. MOD( it, ntrd ) == 0 ) THEN 
     458         IF( lwp .AND. MOD( itmod, ntrd ) == 0 ) THEN 
    458459            WRITE(numout,*) '' 
    459460            WRITE(numout,*) 'trd_vor : write trends in the NetCDF file at kt = ', kt 
     
    568569      IF(lwp) WRITE(numout,*) ' Name of NETCDF file ', clhstnam 
    569570      CALL histbeg( clhstnam, jpi, glamf, jpj, gphif,1, jpi,   &  ! Horizontal grid : glamt and gphit 
    570          &          1, jpj, 0, zjulian, rdt, nh_t, nidvor, domain_id=nidom ) 
     571         &          1, jpj, nit000-1, zjulian, rdt, nh_t, nidvor, domain_id=nidom ) 
    571572      CALL wheneq( jpi*jpj, fmask, 1, 1., ndexvor1, ndimvor1 )    ! surface 
    572573 
  • trunk/NEMO/OPA_SRC/step.F90

    r1246 r1334  
    195195 
    196196      IF( ninist == 1 ) THEN                       ! Output the initial state and forcings 
    197                         CALL dia_wri_state( 'output.init' ) 
     197                        CALL dia_wri_state( 'output.init', kstp ) 
    198198                        ninist = 0 
    199199      ENDIF 
Note: See TracChangeset for help on using the changeset viewer.