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 13248 for NEMO/branches/2020/dev_r12558_HPC-08_epico_Extra_Halo/src/OCE/DIA/diawri.F90 – NEMO

Ignore:
Timestamp:
2020-07-03T20:46:53+02:00 (4 years ago)
Author:
francesca
Message:

dev_r12558_HPC-08_epico_Extra_Halo: merge with trunk@13237, see #2366

File:
1 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/2020/dev_r12558_HPC-08_epico_Extra_Halo/src/OCE/DIA/diawri.F90

    r13247 r13248  
    8585   !! * Substitutions 
    8686#  include "do_loop_substitute.h90" 
     87#  include "domzgr_substitute.h90" 
    8788   !!---------------------------------------------------------------------- 
    8889   !! NEMO/OCE 4.0 , NEMO Consortium (2018) 
     
    136137      CALL iom_put("e3v_0", e3v_0(:,:,:) ) 
    137138      ! 
    138       CALL iom_put( "e3t" , e3t(:,:,:,Kmm) ) 
    139       CALL iom_put( "e3u" , e3u(:,:,:,Kmm) ) 
    140       CALL iom_put( "e3v" , e3v(:,:,:,Kmm) ) 
    141       CALL iom_put( "e3w" , e3w(:,:,:,Kmm) ) 
    142       IF( iom_use("e3tdef") )   & 
    143          CALL iom_put( "e3tdef"  , ( ( e3t(:,:,:,Kmm) - e3t_0(:,:,:) ) / e3t_0(:,:,:) * 100 * tmask(:,:,:) ) ** 2 ) 
    144  
    145       IF( ll_wd ) THEN 
    146          CALL iom_put( "ssh" , (ssh(:,:,Kmm)+ssh_ref)*tmask(:,:,1) )   ! sea surface height (brought back to the reference used for wetting and drying) 
     139      IF ( iom_use("e3t") .OR. iom_use("e3tdef") ) THEN  ! time-varying e3t 
     140         DO jk = 1, jpk 
     141            z3d(:,:,jk) =  e3t(:,:,jk,Kmm) 
     142         END DO 
     143         CALL iom_put( "e3t"     ,     z3d(:,:,:) ) 
     144         CALL iom_put( "e3tdef"  , ( ( z3d(:,:,:) - e3t_0(:,:,:) ) / e3t_0(:,:,:) * 100 * tmask(:,:,:) ) ** 2 )  
     145      ENDIF  
     146      IF ( iom_use("e3u") ) THEN                         ! time-varying e3u 
     147         DO jk = 1, jpk 
     148            z3d(:,:,jk) =  e3u(:,:,jk,Kmm) 
     149         END DO  
     150         CALL iom_put( "e3u" , z3d(:,:,:) ) 
     151      ENDIF 
     152      IF ( iom_use("e3v") ) THEN                         ! time-varying e3v 
     153         DO jk = 1, jpk 
     154            z3d(:,:,jk) =  e3v(:,:,jk,Kmm) 
     155         END DO  
     156         CALL iom_put( "e3v" , z3d(:,:,:) ) 
     157      ENDIF 
     158      IF ( iom_use("e3w") ) THEN                         ! time-varying e3w 
     159         DO jk = 1, jpk 
     160            z3d(:,:,jk) =  e3w(:,:,jk,Kmm) 
     161         END DO  
     162         CALL iom_put( "e3w" , z3d(:,:,:) ) 
     163      ENDIF 
     164 
     165      IF( ll_wd ) THEN                                   ! sea surface height (brought back to the reference used for wetting and drying) 
     166         CALL iom_put( "ssh" , (ssh(:,:,Kmm)+ssh_ref)*tmask(:,:,1) ) 
    147167      ELSE 
    148168         CALL iom_put( "ssh" , ssh(:,:,Kmm) )              ! sea surface height 
     
    172192      ENDIF 
    173193 
     194#if ! defined key_qco 
    174195      CALL iom_put( "rhop", rhop(:,:,:) )          ! 3D potential density (sigma0) 
     196#endif 
    175197 
    176198      IF ( iom_use("taubot") ) THEN                ! bottom stress 
     
    210232 
    211233      IF( ln_zad_Aimp ) ww = ww + wi               ! Recombine explicit and implicit parts of vertical velocity for diagnostic output 
    212       ! 
    213234      CALL iom_put( "woce", ww )                   ! vertical velocity 
     235 
    214236      IF( iom_use('w_masstr') .OR. iom_use('w_masstr2') ) THEN   ! vertical mass transport & its square value 
    215237         ! Caution: in the VVL case, it only correponds to the baroclinic mass transport. 
     
    417439      ! 
    418440      REAL(wp), DIMENSION(jpi,jpj)   :: zw2d       ! 2D workspace 
    419       REAL(wp), DIMENSION(jpi,jpj,jpk) :: zw3d       ! 3D workspace 
     441      REAL(wp), DIMENSION(jpi,jpj,jpk) :: zw3d, ze3t, zgdept       ! 3D workspace 
    420442      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE :: zw3d_abl   ! ABL 3D workspace 
    421443      !!---------------------------------------------------------------------- 
     
    457479      it = kt 
    458480      itmod = kt - nit000 + 1 
     481 
     482      ! store e3t for subsitute 
     483      DO jk = 1, jpk 
     484         ze3t  (:,:,jk) =  e3t  (:,:,jk,Kmm) 
     485         zgdept(:,:,jk) =  gdept(:,:,jk,Kmm) 
     486      END DO 
     487 
    459488 
    460489      ! 1. Define NETCDF files and fields at beginning of first time step 
     
    570599         DEALLOCATE(zw3d_abl) 
    571600         ENDIF 
     601         ! 
    572602 
    573603         ! Declare all the output fields as NETCDF variables 
     
    579609            &          jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout ) 
    580610         IF(  .NOT.ln_linssh  ) THEN 
    581             CALL histdef( nid_T, "vovvle3t", "Level thickness"                    , "m"      ,&  ! e3t(:,:,:,Kmm) 
     611            CALL histdef( nid_T, "vovvle3t", "Level thickness"                    , "m"      ,&  ! e3t n 
    582612            &             jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout ) 
    583             CALL histdef( nid_T, "vovvldep", "T point depth"                      , "m"      ,&  ! e3t(:,:,:,Kmm) 
     613            CALL histdef( nid_T, "vovvldep", "T point depth"                      , "m"      ,&  ! e3t n 
    584614            &             jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout ) 
    585             CALL histdef( nid_T, "vovvldef", "Squared level deformation"          , "%^2"    ,&  ! e3t(:,:,:,Kmm) 
     615            CALL histdef( nid_T, "vovvldef", "Squared level deformation"          , "%^2"    ,&  ! e3t n 
    586616            &             jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout ) 
    587617         ENDIF 
     
    767797 
    768798      IF( .NOT.ln_linssh ) THEN 
    769          CALL histwrite( nid_T, "votemper", it, ts(:,:,:,jp_tem,Kmm) * e3t(:,:,:,Kmm) , ndim_T , ndex_T  )   ! heat content 
    770          CALL histwrite( nid_T, "vosaline", it, ts(:,:,:,jp_sal,Kmm) * e3t(:,:,:,Kmm) , ndim_T , ndex_T  )   ! salt content 
    771          CALL histwrite( nid_T, "sosstsst", it, ts(:,:,1,jp_tem,Kmm) * e3t(:,:,1,Kmm) , ndim_hT, ndex_hT )   ! sea surface heat content 
    772          CALL histwrite( nid_T, "sosaline", it, ts(:,:,1,jp_sal,Kmm) * e3t(:,:,1,Kmm) , ndim_hT, ndex_hT )   ! sea surface salinity content 
     799         CALL histwrite( nid_T, "votemper", it, ts(:,:,:,jp_tem,Kmm) * ze3t(:,:,:) , ndim_T , ndex_T  )   ! heat content 
     800         CALL histwrite( nid_T, "vosaline", it, ts(:,:,:,jp_sal,Kmm) * ze3t(:,:,:) , ndim_T , ndex_T  )   ! salt content 
     801         CALL histwrite( nid_T, "sosstsst", it, ts(:,:,1,jp_tem,Kmm) * ze3t(:,:,1) , ndim_hT, ndex_hT )   ! sea surface heat content 
     802         CALL histwrite( nid_T, "sosaline", it, ts(:,:,1,jp_sal,Kmm) * ze3t(:,:,1) , ndim_hT, ndex_hT )   ! sea surface salinity content 
    773803      ELSE 
    774804         CALL histwrite( nid_T, "votemper", it, ts(:,:,:,jp_tem,Kmm) , ndim_T , ndex_T  )   ! temperature 
     
    778808      ENDIF 
    779809      IF( .NOT.ln_linssh ) THEN 
    780          zw3d(:,:,:) = ( ( e3t(:,:,:,Kmm) - e3t_0(:,:,:) ) / e3t_0(:,:,:) * 100 * tmask(:,:,:) ) ** 2 
    781          CALL histwrite( nid_T, "vovvle3t", it, e3t (:,:,:,Kmm) , ndim_T , ndex_T  )   ! level thickness 
    782          CALL histwrite( nid_T, "vovvldep", it, gdept(:,:,:,Kmm) , ndim_T , ndex_T  )   ! t-point depth 
     810         zw3d(:,:,:) = ( ( ze3t(:,:,:) - e3t_0(:,:,:) ) / e3t_0(:,:,:) * 100 * tmask(:,:,:) ) ** 2 
     811         CALL histwrite( nid_T, "vovvle3t", it, ze3t (:,:,:)    , ndim_T , ndex_T  )   ! level thickness 
     812         CALL histwrite( nid_T, "vovvldep", it, zgdept , ndim_T , ndex_T  )   ! t-point depth  
    783813         CALL histwrite( nid_T, "vovvldef", it, zw3d             , ndim_T , ndex_T  )   ! level thickness deformation 
    784814      ENDIF 
     
    919949      !! 
    920950      INTEGER :: inum, jk 
     951      REAL(wp), DIMENSION(jpi,jpj,jpk) :: ze3t, zgdept      ! 3D workspace !!st patch to use substitution 
    921952      !!---------------------------------------------------------------------- 
    922953      !  
    923       IF(lwp) WRITE(numout,*) 
    924       IF(lwp) WRITE(numout,*) 'dia_wri_state : single instantaneous ocean state' 
    925       IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~   and forcing fields file created ' 
    926       IF(lwp) WRITE(numout,*) '                and named :', cdfile_name, '...nc' 
     954      IF(lwp) THEN 
     955         WRITE(numout,*) 
     956         WRITE(numout,*) 'dia_wri_state : single instantaneous ocean state' 
     957         WRITE(numout,*) '~~~~~~~~~~~~~   and forcing fields file created ' 
     958         WRITE(numout,*) '                and named :', cdfile_name, '...nc' 
     959      ENDIF  
     960      ! 
     961      DO jk = 1, jpk 
     962         ze3t(:,:,jk) =  e3t(:,:,jk,Kmm) 
     963         zgdept(:,:,jk) =  gdept(:,:,jk,Kmm) 
     964      END DO 
    927965      ! 
    928966      CALL iom_open( TRIM(cdfile_name), inum, ldwrt = .TRUE. ) 
     
    939977      ENDIF 
    940978      CALL iom_rstput( 0, 0, inum, 'risfdep', risfdep            )    ! now k-velocity 
    941       CALL iom_rstput( 0, 0, inum, 'ht'     , ht                 )    ! now water column height 
     979      CALL iom_rstput( 0, 0, inum, 'ht'     , ht(:,:)            )    ! now water column height 
    942980      ! 
    943981      IF ( ln_isf ) THEN 
     
    9761014      CALL iom_rstput( 0, 0, inum, 'sometauy', vtau              )    ! j-wind stress 
    9771015      IF(  .NOT.ln_linssh  ) THEN              
    978          CALL iom_rstput( 0, 0, inum, 'vovvldep', gdept(:,:,:,Kmm)        )    !  T-cell depth  
    979          CALL iom_rstput( 0, 0, inum, 'vovvle3t', e3t(:,:,:,Kmm)          )    !  T-cell thickness   
     1016         CALL iom_rstput( 0, 0, inum, 'vovvldep', zgdept        )    !  T-cell depth  
     1017         CALL iom_rstput( 0, 0, inum, 'vovvle3t', ze3t          )    !  T-cell thickness   
    9801018      END IF 
    9811019      IF( ln_wave .AND. ln_sdw ) THEN 
Note: See TracChangeset for help on using the changeset viewer.