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 5989 for branches/2014/dev_r4650_UKMO10_Tidally_Meaned_Diagnostics/NEMOGCM/NEMO/OPA_SRC/DIA/diawri.F90 – NEMO

Ignore:
Timestamp:
2015-12-03T09:10:32+01:00 (8 years ago)
Author:
deazer
Message:

Merging TMB and 25h diagnostics to head of trunk
added brief documentation

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/2014/dev_r4650_UKMO10_Tidally_Meaned_Diagnostics/NEMOGCM/NEMO/OPA_SRC/DIA/diawri.F90

    r5260 r5989  
    1717   !!                 ! 2005-11  (V. Garnier) Surface pressure gradient organization 
    1818   !!            3.2  ! 2008-11  (B. Lemaire) creation from old diawri 
     19   !!            3.7  ! 2014-01  (G. Madec) remove eddy induced velocity from no-IOM output 
     20   !!                 !                     change name of output variables in dia_wri_state 
    1921   !!---------------------------------------------------------------------- 
    2022 
     
    2729   USE dynadv, ONLY: ln_dynadv_vec 
    2830   USE zdf_oce         ! ocean vertical physics 
    29    USE ldftra_oce      ! ocean active tracers: lateral physics 
    30    USE ldfdyn_oce      ! ocean dynamics: lateral physics 
    31    USE traldf_iso_grif, ONLY : psix_eiv, psiy_eiv 
    32    USE sol_oce         ! solver variables 
     31   USE ldftra          ! lateral physics: eddy diffusivity coef. 
    3332   USE sbc_oce         ! Surface boundary condition: ocean fields 
    3433   USE sbc_ice         ! Surface boundary condition: ice fields 
     
    4847   USE iom 
    4948   USE ioipsl 
     49 
    5050#if defined key_lim2 
    5151   USE limwri_2  
     
    127127      !! 
    128128      INTEGER                      ::   ji, jj, jk              ! dummy loop indices 
     129      INTEGER                      ::   jkbot                   ! 
    129130      REAL(wp)                     ::   zztmp, zztmpx, zztmpy   !  
    130131      !! 
     
    150151         CALL iom_put( "e3w" , fse3w_n(:,:,:) ) 
    151152      ENDIF 
     153 
     154      CALL iom_put( "ssh" , sshn )                 ! sea surface height 
     155      if( iom_use('ssh2') )   CALL iom_put( "ssh2", sshn(:,:) * sshn(:,:) )   ! square of sea surface height 
    152156       
    153157      CALL iom_put( "toce", tsn(:,:,:,jp_tem) )    ! 3D temperature 
     
    156160         DO jj = 1, jpj 
    157161            DO ji = 1, jpi 
    158                z2d(ji,jj) = tsn(ji,jj,MAX(mbathy(ji,jj),1),jp_tem) 
     162               jkbot = mbkt(ji,jj) 
     163               z2d(ji,jj) = tsn(ji,jj,jkbot,jp_tem) 
    159164            END DO 
    160165         END DO 
     
    167172         DO jj = 1, jpj 
    168173            DO ji = 1, jpi 
    169                z2d(ji,jj) = tsn(ji,jj,MAX(mbathy(ji,jj),1),jp_sal) 
     174               jkbot = mbkt(ji,jj) 
     175               z2d(ji,jj) = tsn(ji,jj,jkbot,jp_sal) 
    170176            END DO 
    171177         END DO 
    172178         CALL iom_put( "sbs", z2d )                ! bottom salinity 
     179      ENDIF 
     180 
     181      IF ( iom_use("taubot") ) THEN                ! bottom stress 
     182         z2d(:,:) = 0._wp 
     183         DO jj = 2, jpjm1 
     184            DO ji = fs_2, fs_jpim1   ! vector opt. 
     185               zztmpx = (  bfrua(ji  ,jj) * un(ji  ,jj,mbku(ji  ,jj))  & 
     186                      &  + bfrua(ji-1,jj) * un(ji-1,jj,mbku(ji-1,jj))  )       
     187               zztmpy = (  bfrva(ji,  jj) * vn(ji,jj  ,mbkv(ji,jj  ))  & 
     188                      &  + bfrva(ji,jj-1) * vn(ji,jj-1,mbkv(ji,jj-1))  )  
     189               z2d(ji,jj) = rau0 * SQRT( zztmpx * zztmpx + zztmpy * zztmpy ) * tmask(ji,jj,1)  
     190               ! 
     191            ENDDO 
     192         ENDDO 
     193         CALL lbc_lnk( z2d, 'T', 1. ) 
     194         CALL iom_put( "taubot", z2d )            
    173195      ENDIF 
    174196          
     
    178200         DO jj = 1, jpj 
    179201            DO ji = 1, jpi 
    180                z2d(ji,jj) = un(ji,jj,MAX(mbathy(ji,jj),1)) 
     202               jkbot = mbku(ji,jj) 
     203               z2d(ji,jj) = un(ji,jj,jkbot) 
    181204            END DO 
    182205         END DO 
    183206         CALL iom_put( "sbu", z2d )                ! bottom i-current 
     207      ENDIF 
     208 
     209      IF ( ln_dynspg_ts ) THEN 
     210         CALL iom_put(  "ubar", un_adv(:,:)      )    ! barotropic i-current 
     211      ELSE 
     212         CALL iom_put(  "ubar", un_b(:,:)        )    ! barotropic i-current 
    184213      ENDIF 
    185214       
     
    189218         DO jj = 1, jpj 
    190219            DO ji = 1, jpi 
    191                z2d(ji,jj) = vn(ji,jj,MAX(mbathy(ji,jj),1)) 
     220               jkbot = mbkv(ji,jj) 
     221               z2d(ji,jj) = vn(ji,jj,jkbot) 
    192222            END DO 
    193223         END DO 
    194224         CALL iom_put( "sbv", z2d )                ! bottom j-current 
     225      ENDIF 
     226 
     227      IF ( ln_dynspg_ts ) THEN 
     228         CALL iom_put(  "vbar", vn_adv(:,:)      )    ! barotropic j-current 
     229      ELSE 
     230         CALL iom_put(  "vbar", vn_b(:,:)        )    ! barotropic j-current 
     231      ENDIF 
     232 
     233      CALL iom_put( "woce", wn )                   ! vertical velocity 
     234      IF( iom_use('w_masstr') .OR. iom_use('w_masstr2') ) THEN   ! vertical mass transport & its square value 
     235         ! Caution: in the VVL case, it only correponds to the baroclinic mass transport. 
     236         z2d(:,:) = rau0 * e1e2t(:,:) 
     237         DO jk = 1, jpk 
     238            z3d(:,:,jk) = wn(:,:,jk) * z2d(:,:) 
     239         END DO 
     240         CALL iom_put( "w_masstr" , z3d )   
     241         IF( iom_use('w_masstr2') )   CALL iom_put( "w_masstr2", z3d(:,:,:) * z3d(:,:,:) ) 
    195242      ENDIF 
    196243 
     
    202249         DO jj = 2, jpjm1                                    ! sst gradient 
    203250            DO ji = fs_2, fs_jpim1   ! vector opt. 
    204                zztmp      = tsn(ji,jj,1,jp_tem) 
    205                zztmpx     = ( tsn(ji+1,jj  ,1,jp_tem) - zztmp ) / e1u(ji,jj) + ( zztmp - tsn(ji-1,jj  ,1,jp_tem) ) / e1u(ji-1,jj  ) 
    206                zztmpy     = ( tsn(ji  ,jj+1,1,jp_tem) - zztmp ) / e2v(ji,jj) + ( zztmp - tsn(ji  ,jj-1,1,jp_tem) ) / e2v(ji  ,jj-1) 
     251               zztmp  = tsn(ji,jj,1,jp_tem) 
     252               zztmpx = ( tsn(ji+1,jj,1,jp_tem) - zztmp ) * r1_e1u(ji,jj) + ( zztmp - tsn(ji-1,jj  ,1,jp_tem) ) * r1_e1u(ji-1,jj) 
     253               zztmpy = ( tsn(ji,jj+1,1,jp_tem) - zztmp ) * r1_e2v(ji,jj) + ( zztmp - tsn(ji  ,jj-1,1,jp_tem) ) * r1_e2v(ji,jj-1) 
    207254               z2d(ji,jj) = 0.25 * ( zztmpx * zztmpx + zztmpy * zztmpy )   & 
    208255                  &              * umask(ji,jj,1) * umask(ji-1,jj,1) * vmask(ji,jj,1) * umask(ji,jj-1,1) 
     
    365412      !!      Each nwrite time step, output the instantaneous or mean fields 
    366413      !!---------------------------------------------------------------------- 
    367       !! 
    368       INTEGER, INTENT( in ) ::   kt      ! ocean time-step index 
    369       !! 
     414      INTEGER, INTENT( in ) ::   kt   ! ocean time-step index 
     415      ! 
    370416      LOGICAL ::   ll_print = .FALSE.                        ! =T print and flush numout 
    371417      CHARACTER (len=40) ::   clhstnam, clop, clmx           ! local names 
     
    376422      INTEGER  ::   jn, ierror                               ! local integers 
    377423      REAL(wp) ::   zsto, zout, zmax, zjulian, zdt           ! local scalars 
    378       !! 
     424      ! 
    379425      REAL(wp), POINTER, DIMENSION(:,:)   :: zw2d       ! 2D workspace 
    380426      REAL(wp), POINTER, DIMENSION(:,:,:) :: zw3d       ! 3D workspace 
     
    383429      IF( nn_timing == 1 )   CALL timing_start('dia_wri') 
    384430      ! 
    385       CALL wrk_alloc( jpi , jpj      , zw2d ) 
    386       IF ( ln_traldf_gdia .OR. lk_vvl )  call wrk_alloc( jpi , jpj , jpk  , zw3d ) 
     431                     CALL wrk_alloc( jpi,jpj      , zw2d ) 
     432      IF( lk_vvl )   CALL wrk_alloc( jpi,jpj,jpk  , zw3d ) 
    387433      ! 
    388434      ! Output the initial state and forcings 
     
    402448      zdt = rdt 
    403449      IF( nacc == 1 ) zdt = rdtmin 
    404       IF( ln_mskland )   THEN   ;   clop = "only(x)"   ! put 1.e+20 on land (very expensive!!) 
    405       ELSE                      ;   clop = "x"         ! no use of the mask value (require less cpu time) 
    406       ENDIF 
     450      clop = "x"         ! no use of the mask value (require less cpu time and otherwise the model crashes) 
    407451#if defined key_diainstant 
    408452      zsto = nwrite * zdt 
     
    604648         ENDIF 
    605649 
    606          IF( .NOT. lk_cpl ) THEN 
     650         IF( .NOT. ln_cpl ) THEN 
    607651            CALL histdef( nid_T, "sohefldp", "Surface Heat Flux: Damping"         , "W/m2"   ,   &  ! qrp 
    608652               &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
     
    613657         ENDIF 
    614658 
    615          IF( lk_cpl .AND. nn_ice <= 1 ) THEN 
     659         IF( ln_cpl .AND. nn_ice <= 1 ) THEN 
    616660            CALL histdef( nid_T, "sohefldp", "Surface Heat Flux: Damping"         , "W/m2"   ,   &  ! qrp 
    617661               &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
     
    623667          
    624668         clmx ="l_max(only(x))"    ! max index on a period 
    625          CALL histdef( nid_T, "sobowlin", "Bowl Index"                         , "W-point",   &  ! bowl INDEX  
    626             &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clmx, zsto, zout ) 
     669!         CALL histdef( nid_T, "sobowlin", "Bowl Index"                         , "W-point",   &  ! bowl INDEX  
     670!            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clmx, zsto, zout ) 
    627671#if defined key_diahth 
    628672         CALL histdef( nid_T, "sothedep", "Thermocline Depth"                  , "m"      ,   & ! hth 
     
    636680#endif 
    637681 
    638          IF( lk_cpl .AND. nn_ice == 2 ) THEN 
     682         IF( ln_cpl .AND. nn_ice == 2 ) THEN 
    639683            CALL histdef( nid_T,"soicetem" , "Ice Surface Temperature"            , "K"      ,   &  ! tn_ice 
    640684               &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
     
    648692         CALL histdef( nid_U, "vozocrtx", "Zonal Current"                      , "m/s"    ,   &  ! un 
    649693            &          jpi, jpj, nh_U, ipk, 1, ipk, nz_U, 32, clop, zsto, zout ) 
    650          IF( ln_traldf_gdia ) THEN 
    651             CALL histdef( nid_U, "vozoeivu", "Zonal EIV Current"                  , "m/s"    ,   &  ! u_eiv 
    652                  &          jpi, jpj, nh_U, ipk, 1, ipk, nz_U, 32, clop, zsto, zout ) 
    653          ELSE 
    654 #if defined key_diaeiv 
    655             CALL histdef( nid_U, "vozoeivu", "Zonal EIV Current"                  , "m/s"    ,   &  ! u_eiv 
    656             &          jpi, jpj, nh_U, ipk, 1, ipk, nz_U, 32, clop, zsto, zout ) 
    657 #endif 
    658          END IF 
    659694         !                                                                                      !!! nid_U : 2D 
    660695         CALL histdef( nid_U, "sozotaux", "Wind Stress along i-axis"           , "N/m2"   ,   &  ! utau 
     
    666701         CALL histdef( nid_V, "vomecrty", "Meridional Current"                 , "m/s"    ,   &  ! vn 
    667702            &          jpi, jpj, nh_V, ipk, 1, ipk, nz_V, 32, clop, zsto, zout ) 
    668          IF( ln_traldf_gdia ) THEN 
    669             CALL histdef( nid_V, "vomeeivv", "Meridional EIV Current"             , "m/s"    ,   &  ! v_eiv 
    670                  &          jpi, jpj, nh_V, ipk, 1, ipk, nz_V, 32, clop, zsto, zout ) 
    671          ELSE  
    672 #if defined key_diaeiv 
    673             CALL histdef( nid_V, "vomeeivv", "Meridional EIV Current"             , "m/s"    ,   &  ! v_eiv 
    674             &          jpi, jpj, nh_V, ipk, 1, ipk, nz_V, 32, clop, zsto, zout ) 
    675 #endif 
    676          END IF 
    677703         !                                                                                      !!! nid_V : 2D 
    678704         CALL histdef( nid_V, "sometauy", "Wind Stress along j-axis"           , "N/m2"   ,   &  ! vtau 
     
    684710         CALL histdef( nid_W, "vovecrtz", "Vertical Velocity"                  , "m/s"    ,   &  ! wn 
    685711            &          jpi, jpj, nh_W, ipk, 1, ipk, nz_W, 32, clop, zsto, zout ) 
    686          IF( ln_traldf_gdia ) THEN 
    687             CALL histdef( nid_W, "voveeivw", "Vertical EIV Velocity"              , "m/s"    ,   &  ! w_eiv 
    688                  &          jpi, jpj, nh_W, ipk, 1, ipk, nz_W, 32, clop, zsto, zout ) 
    689          ELSE 
    690 #if defined key_diaeiv 
    691             CALL histdef( nid_W, "voveeivw", "Vertical EIV Velocity"              , "m/s"    ,   &  ! w_eiv 
    692                  &          jpi, jpj, nh_W, ipk, 1, ipk, nz_W, 32, clop, zsto, zout ) 
    693 #endif 
    694          END IF 
    695712         CALL histdef( nid_W, "votkeavt", "Vertical Eddy Diffusivity"          , "m2/s"   ,   &  ! avt 
    696713            &          jpi, jpj, nh_W, ipk, 1, ipk, nz_W, 32, clop, zsto, zout ) 
     
    703720         ENDIF 
    704721         !                                                                                      !!! nid_W : 2D 
    705 #if defined key_traldf_c2d 
    706          CALL histdef( nid_W, "soleahtw", "lateral eddy diffusivity"           , "m2/s"   ,   &  ! ahtw 
    707             &          jpi, jpj, nh_W, 1  , 1, 1  , - 99, 32, clop, zsto, zout ) 
    708 # if defined key_traldf_eiv  
    709             CALL histdef( nid_W, "soleaeiw", "eddy induced vel. coeff. at w-point", "m2/s",   &  ! aeiw 
    710                &       jpi, jpj, nh_W, 1  , 1, 1  , - 99, 32, clop, zsto, zout ) 
    711 # endif 
    712 #endif 
    713  
    714722         CALL histend( nid_W, snc4chunks=snc4set ) 
    715723 
     
    791799      ENDIF 
    792800 
    793       IF( .NOT. lk_cpl ) THEN 
     801      IF( .NOT. ln_cpl ) THEN 
    794802         CALL histwrite( nid_T, "sohefldp", it, qrp           , ndim_hT, ndex_hT )   ! heat flux damping 
    795803         CALL histwrite( nid_T, "sowafldp", it, erp           , ndim_hT, ndex_hT )   ! freshwater flux damping 
     
    797805         CALL histwrite( nid_T, "sosafldp", it, zw2d          , ndim_hT, ndex_hT )   ! salt flux damping 
    798806      ENDIF 
    799       IF( lk_cpl .AND. nn_ice <= 1 ) THEN 
     807      IF( ln_cpl .AND. nn_ice <= 1 ) THEN 
    800808         CALL histwrite( nid_T, "sohefldp", it, qrp           , ndim_hT, ndex_hT )   ! heat flux damping 
    801809         CALL histwrite( nid_T, "sowafldp", it, erp           , ndim_hT, ndex_hT )   ! freshwater flux damping 
     
    813821#endif 
    814822 
    815       IF( lk_cpl .AND. nn_ice == 2 ) THEN 
     823      IF( ln_cpl .AND. nn_ice == 2 ) THEN 
    816824         CALL histwrite( nid_T, "soicetem", it, tn_ice(:,:,1) , ndim_hT, ndex_hT )   ! surf. ice temperature 
    817825         CALL histwrite( nid_T, "soicealb", it, alb_ice(:,:,1), ndim_hT, ndex_hT )   ! ice albedo 
     
    819827 
    820828      CALL histwrite( nid_U, "vozocrtx", it, un            , ndim_U , ndex_U )    ! i-current 
    821       IF( ln_traldf_gdia ) THEN 
    822          IF (.not. ALLOCATED(psix_eiv))THEN 
    823             ALLOCATE( psix_eiv(jpi,jpj,jpk) , psiy_eiv(jpi,jpj,jpk) , STAT=ierr ) 
    824             IF( lk_mpp   )   CALL mpp_sum ( ierr ) 
    825             IF( ierr > 0 )   CALL ctl_stop('STOP', 'diawri: unable to allocate psi{x,y}_eiv') 
    826             psix_eiv(:,:,:) = 0.0_wp 
    827             psiy_eiv(:,:,:) = 0.0_wp 
    828          ENDIF 
    829          DO jk=1,jpkm1 
    830             zw3d(:,:,jk) = (psix_eiv(:,:,jk+1) - psix_eiv(:,:,jk))/fse3u(:,:,jk)  ! u_eiv = -dpsix/dz 
    831          END DO 
    832          zw3d(:,:,jpk) = 0._wp 
    833          CALL histwrite( nid_U, "vozoeivu", it, zw3d, ndim_U , ndex_U )           ! i-eiv current 
    834       ELSE 
    835 #if defined key_diaeiv 
    836          CALL histwrite( nid_U, "vozoeivu", it, u_eiv, ndim_U , ndex_U )          ! i-eiv current 
    837 #endif 
    838       ENDIF 
    839829      CALL histwrite( nid_U, "sozotaux", it, utau          , ndim_hU, ndex_hU )   ! i-wind stress 
    840830 
    841831      CALL histwrite( nid_V, "vomecrty", it, vn            , ndim_V , ndex_V  )   ! j-current 
    842       IF( ln_traldf_gdia ) THEN 
    843          DO jk=1,jpk-1 
    844             zw3d(:,:,jk) = (psiy_eiv(:,:,jk+1) - psiy_eiv(:,:,jk))/fse3v(:,:,jk)  ! v_eiv = -dpsiy/dz 
    845          END DO 
    846          zw3d(:,:,jpk) = 0._wp 
    847          CALL histwrite( nid_V, "vomeeivv", it, zw3d, ndim_V , ndex_V )           ! j-eiv current 
    848       ELSE 
    849 #if defined key_diaeiv 
    850          CALL histwrite( nid_V, "vomeeivv", it, v_eiv, ndim_V , ndex_V )          ! j-eiv current 
    851 #endif 
    852       ENDIF 
    853832      CALL histwrite( nid_V, "sometauy", it, vtau          , ndim_hV, ndex_hV )   ! j-wind stress 
    854833 
    855834      CALL histwrite( nid_W, "vovecrtz", it, wn             , ndim_T, ndex_T )    ! vert. current 
    856       IF( ln_traldf_gdia ) THEN 
    857          DO jk=1,jpk-1 
    858             DO jj = 2, jpjm1 
    859                DO ji = fs_2, fs_jpim1  ! vector opt. 
    860                   zw3d(ji,jj,jk) = (psiy_eiv(ji,jj,jk) - psiy_eiv(ji,jj-1,jk))/e2v(ji,jj) + & 
    861                        &    (psix_eiv(ji,jj,jk) - psix_eiv(ji-1,jj,jk))/e1u(ji,jj) ! w_eiv = dpsiy/dy + dpsiy/dx 
    862                END DO 
    863             END DO 
    864          END DO 
    865          zw3d(:,:,jpk) = 0._wp 
    866          CALL histwrite( nid_W, "voveeivw", it, zw3d          , ndim_T, ndex_T )    ! vert. eiv current 
    867       ELSE 
    868 #   if defined key_diaeiv 
    869          CALL histwrite( nid_W, "voveeivw", it, w_eiv          , ndim_T, ndex_T )    ! vert. eiv current 
    870 #   endif 
    871       ENDIF 
    872835      CALL histwrite( nid_W, "votkeavt", it, avt            , ndim_T, ndex_T )    ! T vert. eddy diff. coef. 
    873836      CALL histwrite( nid_W, "votkeavm", it, avmu           , ndim_T, ndex_T )    ! T vert. eddy visc. coef. 
     
    875838         CALL histwrite( nid_W, "voddmavs", it, fsavs(:,:,:), ndim_T, ndex_T )    ! S vert. eddy diff. coef. 
    876839      ENDIF 
    877 #if defined key_traldf_c2d 
    878       CALL histwrite( nid_W, "soleahtw", it, ahtw          , ndim_hT, ndex_hT )   ! lateral eddy diff. coef. 
    879 # if defined key_traldf_eiv 
    880          CALL histwrite( nid_W, "soleaeiw", it, aeiw       , ndim_hT, ndex_hT )   ! EIV coefficient at w-point 
    881 # endif 
    882 #endif 
    883840 
    884841      ! 3. Close all files 
     
    891848      ENDIF 
    892849      ! 
    893       CALL wrk_dealloc( jpi , jpj      , zw2d ) 
    894       IF ( ln_traldf_gdia .OR. lk_vvl )  call wrk_dealloc( jpi , jpj , jpk  , zw3d ) 
     850                     CALL wrk_dealloc( jpi , jpj        , zw2d ) 
     851      IF( lk_vvl )   CALL wrk_dealloc( jpi , jpj , jpk  , zw3d ) 
    895852      ! 
    896853      IF( nn_timing == 1 )   CALL timing_stop('dia_wri') 
     
    924881      !!---------------------------------------------------------------------- 
    925882      !  
    926 !     IF( nn_timing == 1 )   CALL timing_start('dia_wri_state') ! not sure this works for routines not called in first timestep 
    927  
    928883      ! 0. Initialisation 
    929884      ! ----------------- 
     
    984939         CALL histdef( id_i, "vovvldep", "T point depth"         , "m"      ,   &   ! t-point depth 
    985940            &          jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout ) 
    986       END IF 
     941      ENDIF 
    987942 
    988943#if defined key_lim2 
     
    1008963      CALL histwrite( id_i, "vomecrty", kt, vn               , jpi*jpj*jpk, idex )    ! now j-velocity 
    1009964      CALL histwrite( id_i, "vovecrtz", kt, wn               , jpi*jpj*jpk, idex )    ! now k-velocity 
    1010       CALL histwrite( id_i, "sowaflup", kt, (emp-rnf )       , jpi*jpj    , idex )    ! freshwater budget 
     965      CALL histwrite( id_i, "sowaflup", kt, emp-rnf          , jpi*jpj    , idex )    ! freshwater budget 
    1011966      CALL histwrite( id_i, "sohefldo", kt, qsr + qns        , jpi*jpj    , idex )    ! total heat flux 
    1012967      CALL histwrite( id_i, "soshfldo", kt, qsr              , jpi*jpj    , idex )    ! solar heat flux 
     
    1026981      ENDIF 
    1027982#endif 
    1028         
    1029 !     IF( nn_timing == 1 )   CALL timing_stop('dia_wri_state') ! not sure this works for routines not called in first timestep 
    1030983      !  
    1031  
    1032984   END SUBROUTINE dia_wri_state 
    1033985   !!====================================================================== 
Note: See TracChangeset for help on using the changeset viewer.