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 12489 for NEMO/trunk/src/OCE/DIA – NEMO

Ignore:
Timestamp:
2020-02-28T16:55:11+01:00 (4 years ago)
Author:
davestorkey
Message:

Preparation for new timestepping scheme #2390.
Main changes:

  1. Initial euler timestep now handled in stp and not in TRA/DYN routines.
  2. Renaming of all timestep parameters. In summary, the namelist parameter is now rn_Dt and the current timestep is rDt (and rDt_ice, rDt_trc etc).
  3. Renaming of a few miscellaneous parameters, eg. atfp -> rn_atfp (namelist parameter used everywhere) and rau0 -> rho0.

This version gives bit-comparable results to the previous version of the trunk.

Location:
NEMO/trunk/src/OCE/DIA
Files:
10 edited

Legend:

Unmodified
Added
Removed
  • NEMO/trunk/src/OCE/DIA/dia25h.F90

    r12377 r12489  
    140140      ! ----------------- 
    141141      ! Define frequency of summing to create 25 h mean 
    142       IF( MOD( 3600,NINT(rdt) ) == 0 ) THEN 
    143          i_steps = 3600/NINT(rdt) 
     142      IF( MOD( 3600,NINT(rn_Dt) ) == 0 ) THEN 
     143         i_steps = 3600/NINT(rn_Dt) 
    144144      ELSE 
    145          CALL ctl_stop('STOP', 'dia_wri_tide: timestep must give MOD(3600,rdt) = 0 otherwise no hourly values are possible') 
     145         CALL ctl_stop('STOP', 'dia_wri_tide: timestep must give MOD(3600,rn_Dt) = 0 otherwise no hourly values are possible') 
    146146      ENDIF 
    147147 
  • NEMO/trunk/src/OCE/DIA/diaar5.F90

    r12377 r12489  
    103103         END DO 
    104104         CALL iom_put( 'volcello'  , zrhd(:,:,:)  )  ! WARNING not consistent with CMIP DR where volcello is at ca. 2000 
    105          CALL iom_put( 'masscello' , rau0 * e3t(:,:,:,Kmm) * tmask(:,:,:) )  ! ocean mass 
     105         CALL iom_put( 'masscello' , rho0 * e3t(:,:,:,Kmm) * tmask(:,:,:) )  ! ocean mass 
    106106      ENDIF  
    107107      ! 
     
    181181         CALL iom_put( 'sshsteric', zssh_steric ) 
    182182         !                                         ! ocean bottom pressure 
    183          zztmp = rau0 * grav * 1.e-4_wp               ! recover pressure from pressure anomaly and cover to dbar = 1.e4 Pa 
     183         zztmp = rho0 * grav * 1.e-4_wp               ! recover pressure from pressure anomaly and cover to dbar = 1.e4 Pa 
    184184         zbotpres(:,:) = zztmp * ( zbotpres(:,:) + ssh(:,:,Kmm) + thick0(:,:) ) 
    185185         CALL iom_put( 'botpres', zbotpres ) 
     
    213213         ztemp = glob_sum( 'diaar5', ztsn(:,:,1,jp_tem) ) 
    214214         zsal  = glob_sum( 'diaar5', ztsn(:,:,1,jp_sal) ) 
    215          zmass = rau0 * ( zarho + zvol )       
     215         zmass = rho0 * ( zarho + zvol )       
    216216         ! 
    217217         CALL iom_put( 'masstot', zmass ) 
     
    251251               z2d(:,:) = 0._wp 
    252252               DO_3D_11_11( 1, jpkm1 ) 
    253                   z2d(ji,jj) = z2d(ji,jj) + rau0 * e3t(ji,jj,jk,Kmm) *  ztpot(ji,jj,jk) 
     253                  z2d(ji,jj) = z2d(ji,jj) + rho0 * e3t(ji,jj,jk,Kmm) *  ztpot(ji,jj,jk) 
    254254               END_3D 
    255255               CALL iom_put( 'tosmint_pot', z2d )  
     
    285285          ELSE 
    286286            DO_3D_11_11( 1, jpk ) 
    287                zpe(ji,jj) = zpe(ji,jj) + avt(ji,jj,jk) * MIN(0._wp,rn2(ji,jj,jk)) * rau0 * e3w(ji,jj,jk,Kmm) 
     287               zpe(ji,jj) = zpe(ji,jj) + avt(ji,jj,jk) * MIN(0._wp,rn2(ji,jj,jk)) * rho0 * e3w(ji,jj,jk,Kmm) 
    288288            END_3D 
    289289         ENDIF 
     
    325325       CALL lbc_lnk( 'diaar5', z2d, 'U', -1. ) 
    326326       IF( cptr == 'adv' ) THEN 
    327           IF( ktra == jp_tem ) CALL iom_put( 'uadv_heattr' , rau0_rcp * z2d )  ! advective heat transport in i-direction 
    328           IF( ktra == jp_sal ) CALL iom_put( 'uadv_salttr' , rau0     * z2d )  ! advective salt transport in i-direction 
     327          IF( ktra == jp_tem ) CALL iom_put( 'uadv_heattr' , rho0_rcp * z2d )  ! advective heat transport in i-direction 
     328          IF( ktra == jp_sal ) CALL iom_put( 'uadv_salttr' , rho0     * z2d )  ! advective salt transport in i-direction 
    329329       ENDIF 
    330330       IF( cptr == 'ldf' ) THEN 
    331           IF( ktra == jp_tem ) CALL iom_put( 'udiff_heattr' , rau0_rcp * z2d ) ! diffusive heat transport in i-direction 
    332           IF( ktra == jp_sal ) CALL iom_put( 'udiff_salttr' , rau0     * z2d ) ! diffusive salt transport in i-direction 
     331          IF( ktra == jp_tem ) CALL iom_put( 'udiff_heattr' , rho0_rcp * z2d ) ! diffusive heat transport in i-direction 
     332          IF( ktra == jp_sal ) CALL iom_put( 'udiff_salttr' , rho0     * z2d ) ! diffusive salt transport in i-direction 
    333333       ENDIF 
    334334       ! 
     
    339339       CALL lbc_lnk( 'diaar5', z2d, 'V', -1. ) 
    340340       IF( cptr == 'adv' ) THEN 
    341           IF( ktra == jp_tem ) CALL iom_put( 'vadv_heattr' , rau0_rcp * z2d )  ! advective heat transport in j-direction 
    342           IF( ktra == jp_sal ) CALL iom_put( 'vadv_salttr' , rau0     * z2d )  ! advective salt transport in j-direction 
     341          IF( ktra == jp_tem ) CALL iom_put( 'vadv_heattr' , rho0_rcp * z2d )  ! advective heat transport in j-direction 
     342          IF( ktra == jp_sal ) CALL iom_put( 'vadv_salttr' , rho0     * z2d )  ! advective salt transport in j-direction 
    343343       ENDIF 
    344344       IF( cptr == 'ldf' ) THEN 
    345           IF( ktra == jp_tem ) CALL iom_put( 'vdiff_heattr' , rau0_rcp * z2d ) ! diffusive heat transport in j-direction 
    346           IF( ktra == jp_sal ) CALL iom_put( 'vdiff_salttr' , rau0     * z2d ) ! diffusive salt transport in j-direction 
     345          IF( ktra == jp_tem ) CALL iom_put( 'vdiff_heattr' , rho0_rcp * z2d ) ! diffusive heat transport in j-direction 
     346          IF( ktra == jp_sal ) CALL iom_put( 'vdiff_salttr' , rho0     * z2d ) ! diffusive salt transport in j-direction 
    347347       ENDIF 
    348348           
  • NEMO/trunk/src/OCE/DIA/diacfl.F90

    r12377 r12489  
    5252      ! 
    5353      INTEGER                          ::   ji, jj, jk                       ! dummy loop indices 
    54       REAL(wp)                         ::   z2dt, zCu_max, zCv_max, zCw_max  ! local scalars 
     54      REAL(wp)                         ::   zCu_max, zCv_max, zCw_max        ! local scalars 
    5555      INTEGER , DIMENSION(3)           ::   iloc_u , iloc_v , iloc_w , iloc  ! workspace 
    5656      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zCu_cfl, zCv_cfl, zCw_cfl        ! workspace 
     
    5959      IF( ln_timing )   CALL timing_start('dia_cfl') 
    6060      ! 
    61       !                       ! setup timestep multiplier to account for initial Eulerian timestep 
    62       IF( neuler == 0 .AND. kt == nit000 ) THEN   ;    z2dt = rdt 
    63       ELSE                                        ;    z2dt = rdt * 2._wp 
    64       ENDIF 
    65       ! 
    66       !                 
    6761      DO_3D_11_11( 1, jpk ) 
    68          zCu_cfl(ji,jj,jk) = ABS( uu(ji,jj,jk,Kmm) ) * z2dt / e1u  (ji,jj)      ! for i-direction 
    69          zCv_cfl(ji,jj,jk) = ABS( vv(ji,jj,jk,Kmm) ) * z2dt / e2v  (ji,jj)      ! for j-direction 
    70          zCw_cfl(ji,jj,jk) = ABS( ww(ji,jj,jk) ) * z2dt / e3w(ji,jj,jk,Kmm)   ! for k-direction 
     62         zCu_cfl(ji,jj,jk) = ABS( uu(ji,jj,jk,Kmm) ) * rDt / e1u  (ji,jj)      ! for i-direction 
     63         zCv_cfl(ji,jj,jk) = ABS( vv(ji,jj,jk,Kmm) ) * rDt / e2v  (ji,jj)      ! for j-direction 
     64         zCw_cfl(ji,jj,jk) = ABS( ww(ji,jj,jk) ) * rDt / e3w(ji,jj,jk,Kmm)   ! for k-direction 
    7165      END_3D 
    7266      ! 
     
    118112         WRITE(numcfl,*) '******************************************' 
    119113         WRITE(numcfl,FMT='(3x,a12,6x,f7.4,1x,i4,1x,i4,1x,i4)') 'Run Max Cu', rCu_max, nCu_loc(1), nCu_loc(2), nCu_loc(3) 
    120          WRITE(numcfl,FMT='(3x,a8,11x,f15.1)') ' => dt/C', z2dt/rCu_max 
     114         WRITE(numcfl,FMT='(3x,a8,11x,f15.1)') ' => dt/C', rDt/rCu_max 
    121115         WRITE(numcfl,*) '******************************************' 
    122116         WRITE(numcfl,FMT='(3x,a12,6x,f7.4,1x,i4,1x,i4,1x,i4)') 'Run Max Cv', rCv_max, nCv_loc(1), nCv_loc(2), nCv_loc(3) 
    123          WRITE(numcfl,FMT='(3x,a8,11x,f15.1)') ' => dt/C', z2dt/rCv_max 
     117         WRITE(numcfl,FMT='(3x,a8,11x,f15.1)') ' => dt/C', rDt/rCv_max 
    124118         WRITE(numcfl,*) '******************************************' 
    125119         WRITE(numcfl,FMT='(3x,a12,6x,f7.4,1x,i4,1x,i4,1x,i4)') 'Run Max Cw', rCw_max, nCw_loc(1), nCw_loc(2), nCw_loc(3) 
    126          WRITE(numcfl,FMT='(3x,a8,11x,f15.1)') ' => dt/C', z2dt/rCw_max 
     120         WRITE(numcfl,FMT='(3x,a8,11x,f15.1)') ' => dt/C', rDt/rCw_max 
    127121         CLOSE( numcfl )  
    128122         ! 
     
    131125         WRITE(numout,*) 'dia_cfl : Maximum Courant number information for the run ' 
    132126         WRITE(numout,*) '~~~~~~~' 
    133          WRITE(numout,*) '   Max Cu = ', rCu_max, ' at (i,j,k) = (',nCu_loc(1),nCu_loc(2),nCu_loc(3),') => dt/C = ', z2dt/rCu_max 
    134          WRITE(numout,*) '   Max Cv = ', rCv_max, ' at (i,j,k) = (',nCv_loc(1),nCv_loc(2),nCv_loc(3),') => dt/C = ', z2dt/rCv_max 
    135          WRITE(numout,*) '   Max Cw = ', rCw_max, ' at (i,j,k) = (',nCw_loc(1),nCw_loc(2),nCw_loc(3),') => dt/C = ', z2dt/rCw_max 
     127         WRITE(numout,*) '   Max Cu = ', rCu_max, ' at (i,j,k) = (',nCu_loc(1),nCu_loc(2),nCu_loc(3),') => dt/C = ', rDt/rCu_max 
     128         WRITE(numout,*) '   Max Cv = ', rCv_max, ' at (i,j,k) = (',nCv_loc(1),nCv_loc(2),nCv_loc(3),') => dt/C = ', rDt/rCv_max 
     129         WRITE(numout,*) '   Max Cw = ', rCw_max, ' at (i,j,k) = (',nCw_loc(1),nCw_loc(2),nCw_loc(3),') => dt/C = ', rDt/rCw_max 
    136130      ENDIF 
    137131      ! 
  • NEMO/trunk/src/OCE/DIA/diadct.F90

    r12377 r12489  
    676676                  zsn   = interp(Kmm,k%I,k%J,jk,'V',ts(:,:,:,jp_sal,Kmm) )  
    677677                  zrhop = interp(Kmm,k%I,k%J,jk,'V',rhop)  
    678                   zrhoi = interp(Kmm,k%I,k%J,jk,'V',rhd*rau0+rau0)  
     678                  zrhoi = interp(Kmm,k%I,k%J,jk,'V',rhd*rho0+rho0)  
    679679                  zsshn =  0.5*( ssh(k%I,k%J,Kmm) + ssh(k%I,k%J+1,Kmm)    ) * vmask(k%I,k%J,1)  
    680680               CASE(2,3)  
     
    682682                  zsn   = interp(Kmm,k%I,k%J,jk,'U',ts(:,:,:,jp_sal,Kmm) )  
    683683                  zrhop = interp(Kmm,k%I,k%J,jk,'U',rhop)  
    684                   zrhoi = interp(Kmm,k%I,k%J,jk,'U',rhd*rau0+rau0)  
     684                  zrhoi = interp(Kmm,k%I,k%J,jk,'U',rhd*rho0+rho0)  
    685685                  zsshn =  0.5*( ssh(k%I,k%J,Kmm) + ssh(k%I+1,k%J,Kmm)    ) * umask(k%I,k%J,1)   
    686686               END SELECT  
     
    849849                 zsn   = interp(Kmm,k%I,k%J,jk,'V',ts(:,:,:,jp_sal,Kmm) )  
    850850                 zrhop = interp(Kmm,k%I,k%J,jk,'V',rhop)  
    851                  zrhoi = interp(Kmm,k%I,k%J,jk,'V',rhd*rau0+rau0)  
     851                 zrhoi = interp(Kmm,k%I,k%J,jk,'V',rhd*rho0+rho0)  
    852852 
    853853              CASE(2,3)  
     
    855855                 zsn   = interp(Kmm,k%I,k%J,jk,'U',ts(:,:,:,jp_sal,Kmm) )  
    856856                 zrhop = interp(Kmm,k%I,k%J,jk,'U',rhop)  
    857                  zrhoi = interp(Kmm,k%I,k%J,jk,'U',rhd*rau0+rau0)  
     857                 zrhoi = interp(Kmm,k%I,k%J,jk,'U',rhd*rho0+rho0)  
    858858                 zsshn =  0.5*( ssh(k%I,k%J,Kmm)    + ssh(k%I+1,k%J,Kmm)    ) * umask(k%I,k%J,1)   
    859859              END SELECT  
  • NEMO/trunk/src/OCE/DIA/diadetide.F90

    r12377 r12489  
    99   USE in_out_manager , ONLY :   lwp, numout 
    1010   USE iom            , ONLY :   iom_put 
    11    USE dom_oce        , ONLY :   rdt, nsec_day 
     11   USE dom_oce        , ONLY :   rn_Dt, nsec_day 
    1212   USE phycst         , ONLY :   rpi 
    1313   USE tide_mod 
     
    100100      zwght = 0.0_wp 
    101101      DO jn = 1, ndiadetide 
    102          ztmp = ( tdiadetide(jn) - REAL( nsec_day, KIND=wp ) ) / rdt 
     102         ztmp = ( tdiadetide(jn) - REAL( nsec_day, KIND=wp ) ) / rn_Dt 
    103103         IF ( ( ztmp < 0.5_wp ).AND.( ztmp >= -0.5_wp ) ) THEN 
    104104            zwght = zwght + 1.0_wp / REAL( ndiadetide, KIND=wp ) 
  • NEMO/trunk/src/OCE/DIA/diahsb.F90

    r12377 r12489  
    9191      ! 1 - Trends due to forcing ! 
    9292      ! ------------------------- ! 
    93       z_frc_trd_v = r1_rau0 * glob_sum( 'diahsb', - ( emp(:,:) - rnf(:,:) + fwfisf_cav(:,:) + fwfisf_par(:,:) ) * surf(:,:) )   ! volume fluxes 
     93      z_frc_trd_v = r1_rho0 * glob_sum( 'diahsb', - ( emp(:,:) - rnf(:,:) + fwfisf_cav(:,:) + fwfisf_par(:,:) ) * surf(:,:) )   ! volume fluxes 
    9494      z_frc_trd_t =           glob_sum( 'diahsb', sbc_tsc(:,:,jp_tem) * surf(:,:) )                       ! heat fluxes 
    9595      z_frc_trd_s =           glob_sum( 'diahsb', sbc_tsc(:,:,jp_sal) * surf(:,:) )                       ! salt fluxes 
     
    101101         &                          + glob_sum( 'diahsb', ( risf_cav_tsc(:,:,jp_tem) + risf_par_tsc(:,:,jp_tem) ) * surf(:,:) ) 
    102102      !                    ! Add penetrative solar radiation 
    103       IF( ln_traqsr )   z_frc_trd_t = z_frc_trd_t + r1_rau0_rcp * glob_sum( 'diahsb', qsr     (:,:) * surf(:,:) ) 
     103      IF( ln_traqsr )   z_frc_trd_t = z_frc_trd_t + r1_rho0_rcp * glob_sum( 'diahsb', qsr     (:,:) * surf(:,:) ) 
    104104      !                    ! Add geothermal heat flux 
    105105      IF( ln_trabbc )   z_frc_trd_t = z_frc_trd_t +               glob_sum( 'diahsb', qgh_trd0(:,:) * surf(:,:) ) 
     
    121121      ENDIF 
    122122 
    123       frc_v = frc_v + z_frc_trd_v * rdt 
    124       frc_t = frc_t + z_frc_trd_t * rdt 
    125       frc_s = frc_s + z_frc_trd_s * rdt 
     123      frc_v = frc_v + z_frc_trd_v * rn_Dt 
     124      frc_t = frc_t + z_frc_trd_t * rn_Dt 
     125      frc_s = frc_s + z_frc_trd_s * rn_Dt 
    126126      !                                          ! Advection flux through fixed surface (z=0) 
    127127      IF( ln_linssh ) THEN 
    128          frc_wn_t = frc_wn_t + z_wn_trd_t * rdt 
    129          frc_wn_s = frc_wn_s + z_wn_trd_s * rdt 
     128         frc_wn_t = frc_wn_t + z_wn_trd_t * rn_Dt 
     129         frc_wn_s = frc_wn_s + z_wn_trd_s * rn_Dt 
    130130      ENDIF 
    131131 
     
    197197 
    198198      CALL iom_put(   'bgfrcvol' , frc_v    * 1.e-9    )              ! vol - surface forcing (km3)  
    199       CALL iom_put(   'bgfrctem' , frc_t    * rau0 * rcp * 1.e-20 )   ! hc  - surface forcing (1.e20 J)  
    200       CALL iom_put(   'bgfrchfx' , frc_t    * rau0 * rcp /  &         ! hc  - surface forcing (W/m2)  
    201          &                       ( surf_tot * kt * rdt )        ) 
     199      CALL iom_put(   'bgfrctem' , frc_t    * rho0 * rcp * 1.e-20 )   ! hc  - surface forcing (1.e20 J)  
     200      CALL iom_put(   'bgfrchfx' , frc_t    * rho0 * rcp /  &         ! hc  - surface forcing (W/m2)  
     201         &                       ( surf_tot * kt * rn_Dt )        ) 
    202202      CALL iom_put(   'bgfrcsal' , frc_s    * 1.e-9    )              ! sc  - surface forcing (psu*km3)  
    203203 
     
    205205         CALL iom_put( 'bgtemper' , zdiff_hc / zvol_tot )              ! Temperature drift     (C)  
    206206         CALL iom_put( 'bgsaline' , zdiff_sc / zvol_tot )              ! Salinity    drift     (PSU) 
    207          CALL iom_put( 'bgheatco' , zdiff_hc * 1.e-20 * rau0 * rcp )   ! Heat content drift    (1.e20 J)  
    208          CALL iom_put( 'bgheatfx' , zdiff_hc * rau0 * rcp /  &         ! Heat flux drift       (W/m2)  
    209             &                       ( surf_tot * kt * rdt )        ) 
     207         CALL iom_put( 'bgheatco' , zdiff_hc * 1.e-20 * rho0 * rcp )   ! Heat content drift    (1.e20 J)  
     208         CALL iom_put( 'bgheatfx' , zdiff_hc * rho0 * rcp /  &         ! Heat flux drift       (W/m2)  
     209            &                       ( surf_tot * kt * rn_Dt )        ) 
    210210         CALL iom_put( 'bgsaltco' , zdiff_sc * 1.e-9    )              ! Salt content drift    (psu*km3) 
    211211         CALL iom_put( 'bgvolssh' , zdiff_v1 * 1.e-9    )              ! volume ssh drift      (km3)   
     
    225225         CALL iom_put( 'bgtemper' , zdiff_hc1 / zvol_tot)              ! Heat content drift    (C)  
    226226         CALL iom_put( 'bgsaline' , zdiff_sc1 / zvol_tot)              ! Salt content drift    (PSU) 
    227          CALL iom_put( 'bgheatco' , zdiff_hc1 * 1.e-20 * rau0 * rcp )  ! Heat content drift    (1.e20 J)  
    228          CALL iom_put( 'bgheatfx' , zdiff_hc1 * rau0 * rcp /  &        ! Heat flux drift       (W/m2)  
    229             &                       ( surf_tot * kt * rdt )         ) 
     227         CALL iom_put( 'bgheatco' , zdiff_hc1 * 1.e-20 * rho0 * rcp )  ! Heat content drift    (1.e20 J)  
     228         CALL iom_put( 'bgheatfx' , zdiff_hc1 * rho0 * rcp /  &        ! Heat flux drift       (W/m2)  
     229            &                       ( surf_tot * kt * rn_Dt )         ) 
    230230         CALL iom_put( 'bgsaltco' , zdiff_sc1 * 1.e-9    )             ! Salt content drift    (psu*km3) 
    231231         CALL iom_put( 'bgvolssh' , zdiff_v1 * 1.e-9    )              ! volume ssh drift      (km3)   
  • NEMO/trunk/src/OCE/DIA/diahth.F90

    r12377 r12489  
    261261            zzdep = 300. 
    262262            CALL  dia_hth_htc( Kmm, zzdep, ts(:,:,:,jp_tem,Kmm), htc3 ) 
    263             CALL iom_put( 'hc300', rau0_rcp * htc3 )  ! vertically integrated heat content (J/m2) 
     263            CALL iom_put( 'hc300', rho0_rcp * htc3 )  ! vertically integrated heat content (J/m2) 
    264264         ENDIF 
    265265         ! 
     
    270270            zzdep = 700. 
    271271            CALL  dia_hth_htc( Kmm, zzdep, ts(:,:,:,jp_tem,Kmm), htc7 ) 
    272             CALL iom_put( 'hc700', rau0_rcp * htc7 )  ! vertically integrated heat content (J/m2) 
     272            CALL iom_put( 'hc700', rho0_rcp * htc7 )  ! vertically integrated heat content (J/m2) 
    273273   
    274274         ENDIF 
     
    280280            zzdep = 2000. 
    281281            CALL  dia_hth_htc( Kmm, zzdep, ts(:,:,:,jp_tem,Kmm), htc20 ) 
    282             CALL iom_put( 'hc2000', rau0_rcp * htc20 )  ! vertically integrated heat content (J/m2)   
     282            CALL iom_put( 'hc2000', rho0_rcp * htc20 )  ! vertically integrated heat content (J/m2)   
    283283         ENDIF 
    284284         ! 
  • NEMO/trunk/src/OCE/DIA/dianam.F90

    r10068 r12489  
    7272 
    7373      IF( llfsec .OR. kfreq < 0 ) THEN   ;   inbsec = kfreq                       ! output frequency already in seconds 
    74       ELSE                               ;   inbsec = kfreq * NINT( rdt )   ! from time-step to seconds 
     74      ELSE                               ;   inbsec = kfreq * NINT( rn_Dt )   ! from time-step to seconds 
    7575      ENDIF 
    7676      iddss = NINT( rday          )                                         ! number of seconds in 1 day 
     
    116116      ! date of the beginning and the end of the run 
    117117 
    118       zdrun = rdt / rday * REAL( nitend - nit000, wp )                ! length of the run in days 
    119       zjul  = fjulday - rdt / rday 
     118      zdrun = rn_Dt / rday * REAL( nitend - nit000, wp )                ! length of the run in days 
     119      zjul  = fjulday - rn_Dt / rday 
    120120      CALL ju2ymds( zjul        , iyear1, imonth1, iday1, zsec1 )           ! year/month/day of the beginning of run 
    121121      CALL ju2ymds( zjul + zdrun, iyear2, imonth2, iday2, zsec2 )           ! year/month/day of the end       of run 
  • NEMO/trunk/src/OCE/DIA/diaptr.F90

    r12377 r12489  
    5050 
    5151   REAL(wp) ::   rc_sv    = 1.e-6_wp   ! conversion from m3/s to Sverdrup 
    52    REAL(wp) ::   rc_pwatt = 1.e-15_wp  ! conversion from W    to PW (further x rau0 x Cp) 
    53    REAL(wp) ::   rc_ggram = 1.e-9_wp   ! conversion from g    to Gg  (further x rau0) 
     52   REAL(wp) ::   rc_pwatt = 1.e-15_wp  ! conversion from W    to PW (further x rho0 x Cp) 
     53   REAL(wp) ::   rc_ggram = 1.e-9_wp   ! conversion from g    to Gg  (further x rho0) 
    5454 
    5555   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: btmsk   ! T-point basin interior masks 
     
    346346         IF( dia_ptr_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'dia_ptr_init : unable to allocate arrays' ) 
    347347 
    348          rc_pwatt = rc_pwatt * rau0_rcp          ! conversion from K.s-1 to PetaWatt 
    349          rc_ggram = rc_ggram * rau0              ! conversion from m3/s to Gg/s 
     348         rc_pwatt = rc_pwatt * rho0_rcp          ! conversion from K.s-1 to PetaWatt 
     349         rc_ggram = rc_ggram * rho0              ! conversion from m3/s to Gg/s 
    350350 
    351351         IF( lk_mpp )   CALL mpp_ini_znl( numout )     ! Define MPI communicator for zonal sum 
  • NEMO/trunk/src/OCE/DIA/diawri.F90

    r12377 r12489  
    173173 
    174174      IF ( iom_use("taubot") ) THEN                ! bottom stress 
    175          zztmp = rau0 * 0.25 
     175         zztmp = rho0 * 0.25 
    176176         z2d(:,:) = 0._wp 
    177177         DO_2D_00_00 
     
    212212      IF( iom_use('w_masstr') .OR. iom_use('w_masstr2') ) THEN   ! vertical mass transport & its square value 
    213213         ! Caution: in the VVL case, it only correponds to the baroclinic mass transport. 
    214          z2d(:,:) = rau0 * e1e2t(:,:) 
     214         z2d(:,:) = rho0 * e1e2t(:,:) 
    215215         DO jk = 1, jpk 
    216216            z3d(:,:,jk) = ww(:,:,jk) * z2d(:,:) 
     
    249249            z2d(ji,jj) = z2d(ji,jj) + e3t(ji,jj,jk,Kmm) * ts(ji,jj,jk,jp_tem,Kmm) * tmask(ji,jj,jk) 
    250250         END_3D 
    251          CALL iom_put( "heatc", rau0_rcp * z2d )   ! vertically integrated heat content (J/m2) 
     251         CALL iom_put( "heatc", rho0_rcp * z2d )   ! vertically integrated heat content (J/m2) 
    252252      ENDIF 
    253253 
     
    257257            z2d(ji,jj) = z2d(ji,jj) + e3t(ji,jj,jk,Kmm) * ts(ji,jj,jk,jp_sal,Kmm) * tmask(ji,jj,jk) 
    258258         END_3D 
    259          CALL iom_put( "saltc", rau0 * z2d )          ! vertically integrated salt content (PSU*kg/m2) 
     259         CALL iom_put( "saltc", rho0 * z2d )          ! vertically integrated salt content (PSU*kg/m2) 
    260260      ENDIF 
    261261      ! 
     
    279279         z2d(:,:) = 0.e0 
    280280         DO jk = 1, jpkm1 
    281             z3d(:,:,jk) = rau0 * uu(:,:,jk,Kmm) * e2u(:,:) * e3u(:,:,jk,Kmm) * umask(:,:,jk) 
     281            z3d(:,:,jk) = rho0 * uu(:,:,jk,Kmm) * e2u(:,:) * e3u(:,:,jk,Kmm) * umask(:,:,jk) 
    282282            z2d(:,:) = z2d(:,:) + z3d(:,:,jk) 
    283283         END DO 
     
    308308         z3d(:,:,jpk) = 0.e0 
    309309         DO jk = 1, jpkm1 
    310             z3d(:,:,jk) = rau0 * vv(:,:,jk,Kmm) * e1v(:,:) * e3v(:,:,jk,Kmm) * vmask(:,:,jk) 
     310            z3d(:,:,jk) = rho0 * vv(:,:,jk,Kmm) * e1v(:,:) * e3v(:,:,jk,Kmm) * vmask(:,:,jk) 
    311311         END DO 
    312312         CALL iom_put( "v_masstr", z3d )              ! mass transport in j-direction 
     
    337337         END_3D 
    338338         CALL lbc_lnk( 'diawri', z2d, 'T', -1. ) 
    339          CALL iom_put( "tosmint", rau0 * z2d )        ! Vertical integral of temperature 
     339         CALL iom_put( "tosmint", rho0 * z2d )        ! Vertical integral of temperature 
    340340      ENDIF 
    341341      IF( iom_use("somint") ) THEN 
     
    345345         END_3D 
    346346         CALL lbc_lnk( 'diawri', z2d, 'T', -1. ) 
    347          CALL iom_put( "somint", rau0 * z2d )         ! Vertical integral of salinity 
     347         CALL iom_put( "somint", rho0 * z2d )         ! Vertical integral of salinity 
    348348      ENDIF 
    349349 
     
    432432      clop = "x"         ! no use of the mask value (require less cpu time and otherwise the model crashes) 
    433433#if defined key_diainstant 
    434       zsto = nn_write * rdt 
     434      zsto = nn_write * rn_Dt 
    435435      clop = "inst("//TRIM(clop)//")" 
    436436#else 
    437       zsto=rdt 
     437      zsto=rn_Dt 
    438438      clop = "ave("//TRIM(clop)//")" 
    439439#endif 
    440       zout = nn_write * rdt 
    441       zmax = ( nitend - nit000 + 1 ) * rdt 
     440      zout = nn_write * rn_Dt 
     441      zmax = ( nitend - nit000 + 1 ) * rn_Dt 
    442442 
    443443      ! Define indices of the horizontal output zoom and vertical limit storage 
     
    460460 
    461461         ! Compute julian date from starting date of the run 
    462          CALL ymds2ju( nyear, nmonth, nday, rdt, zjulian ) 
     462         CALL ymds2ju( nyear, nmonth, nday, rn_Dt, zjulian ) 
    463463         zjulian = zjulian - adatrj   !   set calendar origin to the beginning of the experiment 
    464464         IF(lwp)WRITE(numout,*) 
     
    482482         CALL histbeg( clhstnam, jpi, glamt, jpj, gphit,           &  ! Horizontal grid: glamt and gphit 
    483483            &          iimi, iima-iimi+1, ijmi, ijma-ijmi+1,       & 
    484             &          nit000-1, zjulian, rdt, nh_T, nid_T, domain_id=nidom, snc4chunks=snc4set ) 
     484            &          nit000-1, zjulian, rn_Dt, nh_T, nid_T, domain_id=nidom, snc4chunks=snc4set ) 
    485485         CALL histvert( nid_T, "deptht", "Vertical T levels",      &  ! Vertical grid: gdept 
    486486            &           "m", ipk, gdept_1d, nz_T, "down" ) 
     
    518518         CALL histbeg( clhstnam, jpi, glamu, jpj, gphiu,           &  ! Horizontal grid: glamu and gphiu 
    519519            &          iimi, iima-iimi+1, ijmi, ijma-ijmi+1,       & 
    520             &          nit000-1, zjulian, rdt, nh_U, nid_U, domain_id=nidom, snc4chunks=snc4set ) 
     520            &          nit000-1, zjulian, rn_Dt, nh_U, nid_U, domain_id=nidom, snc4chunks=snc4set ) 
    521521         CALL histvert( nid_U, "depthu", "Vertical U levels",      &  ! Vertical grid: gdept 
    522522            &           "m", ipk, gdept_1d, nz_U, "down" ) 
     
    531531         CALL histbeg( clhstnam, jpi, glamv, jpj, gphiv,           &  ! Horizontal grid: glamv and gphiv 
    532532            &          iimi, iima-iimi+1, ijmi, ijma-ijmi+1,       & 
    533             &          nit000-1, zjulian, rdt, nh_V, nid_V, domain_id=nidom, snc4chunks=snc4set ) 
     533            &          nit000-1, zjulian, rn_Dt, nh_V, nid_V, domain_id=nidom, snc4chunks=snc4set ) 
    534534         CALL histvert( nid_V, "depthv", "Vertical V levels",      &  ! Vertical grid : gdept 
    535535            &          "m", ipk, gdept_1d, nz_V, "down" ) 
     
    544544         CALL histbeg( clhstnam, jpi, glamt, jpj, gphit,           &  ! Horizontal grid: glamt and gphit 
    545545            &          iimi, iima-iimi+1, ijmi, ijma-ijmi+1,       & 
    546             &          nit000-1, zjulian, rdt, nh_W, nid_W, domain_id=nidom, snc4chunks=snc4set ) 
     546            &          nit000-1, zjulian, rn_Dt, nh_W, nid_W, domain_id=nidom, snc4chunks=snc4set ) 
    547547         CALL histvert( nid_W, "depthw", "Vertical W levels",      &  ! Vertical grid: gdepw 
    548548            &          "m", ipk, gdepw_1d, nz_W, "down" ) 
     
    554554            CALL histbeg( clhstnam, jpi, glamt, jpj, gphit,           &  ! Horizontal grid: glamt and gphit 
    555555               &          iimi, iima-iimi+1, ijmi, ijma-ijmi+1,       & 
    556                &          nit000-1, zjulian, rdt, nh_A, nid_A, domain_id=nidom, snc4chunks=snc4set ) 
     556               &          nit000-1, zjulian, rn_Dt, nh_A, nid_A, domain_id=nidom, snc4chunks=snc4set ) 
    557557            CALL histvert( nid_A, "ght_abl", "Vertical T levels",      &  ! Vertical grid: gdept 
    558558               &           "m", ipka, ght_abl(2:jpka), nz_A, "up" ) 
Note: See TracChangeset for help on using the changeset viewer.