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 1353 for trunk/NEMO/TOP_SRC/trcdia.F90 – NEMO

Ignore:
Timestamp:
2009-03-31T10:15:02+02:00 (15 years ago)
Author:
cetlod
Message:

correction of calendar in TOP outputs, see ticket:368

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/NEMO/TOP_SRC/trcdia.F90

    r1334 r1353  
    141141 
    142142      ! define time axis 
    143       it = kt 
    144143      itmod = kt - nittrc000 + 1 
     144      it    = kt 
    145145 
    146146      ! Define NETCDF files and fields at beginning of first time step 
     
    157157         IF(lwp)WRITE(numout,*)' Date 0 used :', nittrc000                         & 
    158158            &                 ,' YEAR ', nyear, ' MONTH ', nmonth, ' DAY ', nday   & 
    159             &                 ,'Julian day : ', xjulian     
     159            &                 ,'Julian day : ', xjulian   
     160   
    160161         IF(lwp) WRITE(numout,*) ' indexes of zoom = ', iimi, iima, ijmi, ijma,  & 
    161162            &                    ' limit storage in depth = ', ipk 
    162163 
    163164 
    164 ! Define the NETCDF files for passive tracer concentration 
    165  
     165         ! Define the NETCDF files for passive tracer concentration 
    166166         CALL dia_nam( clhstnam, nwritetrc, 'ptrc_T' ) 
    167167         IF(lwp)WRITE(numout,*)" Name of NETCDF file ", clhstnam 
    168 ! Horizontal grid : glamt and gphit 
     168 
     169         ! Horizontal grid : glamt and gphit 
    169170         CALL histbeg( clhstnam, jpi, glamt, jpj, gphit,     & 
    170171            &          iimi, iima-iimi+1, ijmi, ijma-ijmi+1,         &  
    171             &          nittrc000-1, xjulian, zdt, nhorit5, nit5 , domain_id=nidom) 
    172 ! Vertical grid for tracer : gdept 
    173          CALL histvert( nit5, 'deptht', 'Vertical T levels', & 
    174             &            'm', ipk, gdept_0, ndepit5, 'down') 
    175  
    176 ! Index of ocean points in 3D and 2D (surface) 
    177          CALL wheneq( jpi*jpj*ipk,tmask,1,1.,ndext50,ndimt50 ) 
    178          CALL wheneq( jpi*jpj,tmask,1,1.,ndext51,ndimt51 ) 
    179  
    180 ! Declare all the output fields as NETCDF variables 
    181  
    182 ! tracer concentrations 
     172            &          nittrc000-ndttrc, xjulian, zdt, nhorit5, nit5 , domain_id=nidom) 
     173 
     174         ! Vertical grid for tracer : gdept 
     175         CALL histvert( nit5, 'deptht', 'Vertical T levels', 'm', ipk, gdept_0, ndepit5) 
     176 
     177         ! Index of ocean points in 3D and 2D (surface) 
     178         CALL wheneq( jpi*jpj*ipk, tmask, 1, 1., ndext50, ndimt50 ) 
     179         CALL wheneq( jpi*jpj    , tmask, 1, 1., ndext51, ndimt51 ) 
     180 
     181         ! Declare all the output fields as NETCDF variables 
    183182         DO jn = 1, jptra 
    184183            IF( lutsav(jn) ) THEN 
     
    187186               cltrau = ctrcun(jn)   ! UNIT for tracer 
    188187               CALL histdef( nit5, cltra, cltral, cltrau, jpi, jpj, nhorit5,  & 
    189                   &               ipk, 1, ipk,  ndepit5, 32, clop, zsto, zout)  
     188                  &          ipk, 1, ipk,  ndepit5, 32, clop, zsto, zout )  
    190189            ENDIF 
    191190         END DO 
     
    278277 
    279278      ! define time axis 
    280       it = kt 
    281279      itmod = kt - nittrc000 + 1 
     280      it    = kt 
    282281 
    283282      ! Define the NETCDF files (one per tracer) 
     
    300299               CALL histbeg( clhstnam, jpi, glamt, jpj, gphit,       & 
    301300                  &          iimi, iima-iimi+1, ijmi, ijma-ijmi+1,   & 
    302                   &          nittrc000-1, xjulian, rdt, nhorit6(jn),           & 
     301                  &          nittrc000-ndttrc, xjulian, zdt, nhorit6(jn),  & 
    303302                  &          nit6(jn) , domain_id=nidom ) 
    304303 
    305304               ! Vertical grid for tracer trend - one per each tracer IF needed 
    306                CALL histvert( nit6(jn), 'deptht', 'Vertical T levels',   & 
    307                   &           'm', ipk, gdept_0, ndepit6(jn), 'down' )  
     305               CALL histvert( nit6(jn), 'deptht', 'Vertical T levels', 'm', ipk, gdept_0, ndepit6(jn) )  
    308306             END IF 
    309307          END DO 
     
    504502      ENDIF 
    505503#  if defined key_diainstant 
    506       zsto=nwritedia*zdt 
     504      zsto = nwritedia * zdt 
    507505      clop = "inst("//TRIM(clop)//")" 
    508506#  else 
    509       zsto=zdt 
     507      zsto = zdt 
    510508      clop = "ave("//TRIM(clop)//")" 
    511509#  endif 
    512       zout=nwritedia*zdt 
     510      zout = nwritedia * zdt 
    513511 
    514512      ! Define indices of the horizontal output zoom and vertical limit storage 
     
    518516 
    519517      ! define time axis 
    520       it = kt 
    521518      itmod = kt - nittrc000 + 1 
     519      it    = kt 
    522520 
    523521      ! 1. Define NETCDF files and fields at beginning of first time step 
     
    539537         CALL histbeg( clhstnam, jpi, glamt, jpj, gphit,             & 
    540538            &          iimi, iima-iimi+1, ijmi, ijma-ijmi+1,         & 
    541             &          nittrc000-1, xjulian, zdt, nhoritd, nitd , domain_id=nidom ) 
     539            &          nittrc000-ndttrc, xjulian, zdt, nhoritd, nitd , domain_id=nidom ) 
    542540 
    543541         ! Vertical grid for 2d and 3d arrays 
    544542 
    545          CALL histvert( nitd, 'deptht', 'Vertical T levels',   & 
    546             &           'm', ipk, gdept_0, ndepitd, 'down') 
     543         CALL histvert( nitd, 'deptht', 'Vertical T levels','m', ipk, gdept_0, ndepitd) 
    547544 
    548545         ! Declare all the output fields as NETCDF variables 
     
    655652      ENDIF 
    656653#        if defined key_diainstant 
    657       zsto=nwritebio*zdt 
     654      zsto = nwritebio * zdt 
    658655      clop = "inst("//TRIM(clop)//")" 
    659656#        else 
    660       zsto=zdt 
     657      zsto = zdt 
    661658      clop = "ave("//TRIM(clop)//")" 
    662659#        endif 
    663       zout=nwritebio*zdt 
    664  
    665       ! Define indices of the horizontal output zoom and vertical limit storage      iimi = 1      ;      iima = jpi 
     660      zout = nwritebio * zdt 
     661 
     662      ! Define indices of the horizontal output zoom and vertical limit storage 
    666663      iimi = 1      ;      iima = jpi 
    667664      ijmi = 1      ;      ijma = jpj 
     
    669666 
    670667      ! define time axis 
    671       it = kt 
    672668      itmod = kt - nittrc000 + 1 
     669      it    = kt 
    673670 
    674671      ! Define NETCDF files and fields at beginning of first time step 
     
    684681         IF(lwp)WRITE(numout,*) " Name of NETCDF file for biological trends ", clhstnam 
    685682         ! Horizontal grid : glamt and gphit 
    686     CALL histbeg(clhstnam, jpi, glamt, jpj, gphit,      & 
     683         CALL histbeg( clhstnam, jpi, glamt, jpj, gphit,      & 
    687684            &    iimi, iima-iimi+1, ijmi, ijma-ijmi+1,          & 
    688             &    nittrc000-1, xjulian, rdt, nhoritb, nitb , domain_id=nidom) 
     685            &    nittrc000-ndttrc, xjulian, zdt, nhoritb, nitb , domain_id=nidom ) 
    689686         ! Vertical grid for biological trends 
    690          CALL histvert(nitb, 'deptht', 'Vertical T levels',  & 
    691             &    'm', ipk, gdept_0, ndepitb, 'down') 
     687         CALL histvert(nitb, 'deptht', 'Vertical T levels', 'm', ipk, gdept_0, ndepitb) 
    692688 
    693689         ! Declare all the output fields as NETCDF variables 
     
    697693            cltral = ctrbil(jn)   ! long title for biological diagnostic 
    698694            cltrau = ctrbiu(jn)   ! UNIT for biological diagnostic 
    699             CALL histdef(nitb, cltra, cltral, cltrau, jpi, jpj, nhoritb,  & 
     695            CALL histdef( nitb, cltra, cltral, cltrau, jpi, jpj, nhoritb,  & 
    700696               &         ipk, 1, ipk,  ndepitb, 32, clop, zsto, zout) 
    701697         END DO 
    702698 
    703699         ! CLOSE netcdf Files 
    704           CALL histend(nitb) 
     700          CALL histend( nitb ) 
    705701 
    706702         IF(lwp) WRITE(numout,*) 
     
    720716 
    721717      DO jn = 1, jpdiabio 
    722          cltra=ctrbio(jn)  ! short title for biological diagnostic 
     718         cltra = ctrbio(jn)  ! short title for biological diagnostic 
    723719         CALL histwrite(nitb, cltra, it, trbio(:,:,:,jn), ndimt50,ndext50) 
    724720      END DO 
Note: See TracChangeset for help on using the changeset viewer.