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

Ignore:
Timestamp:
2020-07-03T11:12:53+02:00 (4 years ago)
Author:
smasson
Message:

trunk: Mid-year merge, merge back KERNEL-06_techene_e3

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

Legend:

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

    r13226 r13237  
    3939   !! * Substitutions 
    4040#  include "do_loop_substitute.h90" 
     41#  include "domzgr_substitute.h90" 
    4142   !!---------------------------------------------------------------------- 
    4243   !! NEMO/OCE 4.0 , NEMO Consortium (2018) 
     
    7677      ! 
    7778      REAL(wp), ALLOCATABLE, DIMENSION(:,:)     :: zarea_ssh , zbotpres       ! 2D workspace  
    78       REAL(wp), ALLOCATABLE, DIMENSION(:,:)     :: zpe, z2d                   ! 2D workspace  
    79       REAL(wp), ALLOCATABLE, DIMENSION(:,:,:)   :: zrhd , ztpot               ! 3D workspace 
     79      REAL(wp), ALLOCATABLE, DIMENSION(:,:)     :: z2d, zpe                   ! 2D workspace  
     80      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:)   :: z3d, zrhd, ztpot, zgdept   ! 3D workspace (zgdept: needed to use the substitute) 
    8081      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:) :: ztsn                       ! 4D workspace 
    8182 
     
    101102            zrhd(:,:,jk) = e1e2t(:,:) * e3t(:,:,jk,Kmm) * tmask(:,:,jk) 
    102103         END DO 
     104         DO jk = 1, jpk 
     105            z3d(:,:,jk) =  rho0 * e3t(:,:,jk,Kmm) * tmask(:,:,jk) 
     106         END DO  
    103107         CALL iom_put( 'volcello'  , zrhd(:,:,:)  )  ! WARNING not consistent with CMIP DR where volcello is at ca. 2000 
    104          CALL iom_put( 'masscello' , rho0 * e3t(:,:,:,Kmm) * tmask(:,:,:) )  ! ocean mass 
     108         CALL iom_put( 'masscello' , z3d (:,:,:) )   ! ocean mass 
    105109      ENDIF  
    106110      ! 
     
    128132         ztsn(:,:,:,jp_tem) = ts(:,:,:,jp_tem,Kmm)                    ! thermosteric ssh 
    129133         ztsn(:,:,:,jp_sal) = sn0(:,:,:) 
    130          CALL eos( ztsn, zrhd, gdept(:,:,:,Kmm) )                       ! now in situ density using initial salinity 
     134         ALLOCATE( zgdept(jpi,jpj,jpk) ) 
     135         DO jk = 1, jpk 
     136            zgdept(:,:,jk) = gdept(:,:,jk,Kmm) 
     137         END DO 
     138         CALL eos( ztsn, zrhd, zgdept)                       ! now in situ density using initial salinity 
    131139         ! 
    132140         zbotpres(:,:) = 0._wp                        ! no atmospheric surface pressure, levitating sea-ice 
     
    180188         CALL iom_put( 'botpres', zbotpres ) 
    181189         ! 
     190         DEALLOCATE( zgdept ) 
     191         ! 
    182192      ENDIF 
    183193 
  • NEMO/trunk/src/OCE/DIA/diacfl.F90

    r12489 r13237  
    3434   !! * Substitutions 
    3535#  include "do_loop_substitute.h90" 
     36#  include "domzgr_substitute.h90" 
    3637   !!---------------------------------------------------------------------- 
    3738   !! NEMO/OCE 4.0 , NEMO Consortium (2018) 
  • NEMO/trunk/src/OCE/DIA/diadct.F90

    r12489 r13237  
    1111   !!            3.4  ! 09/2011 (C Bricaud) 
    1212   !!---------------------------------------------------------------------- 
    13    !! does not work with agrif 
    1413#if ! defined key_agrif 
     14   !!                        ==>>  CAUTION: does not work with agrif 
    1515   !!---------------------------------------------------------------------- 
    1616   !!   dia_dct      :  Compute the transport through a sec. 
     
    6666   TYPE SECTION 
    6767      CHARACTER(len=60)                            :: name              ! name of the sec 
    68       LOGICAL                                      :: llstrpond         ! true if you want the computation of salt and 
    69                                                                        ! heat transports 
     68      LOGICAL                                      :: llstrpond         ! true if you want the computation of salt and heat transports 
    7069      LOGICAL                                      :: ll_ice_section    ! ice surface and ice volume computation 
    7170      LOGICAL                                      :: ll_date_line      ! = T if the section crosses the date-line 
     
    7473      INTEGER, DIMENSION(nb_point_max)             :: direction         ! vector direction of the point in the section 
    7574      CHARACTER(len=40),DIMENSION(nb_class_max)    :: classname         ! characteristics of the class 
    76       REAL(wp), DIMENSION(nb_class_max)            :: zsigi           ,&! in-situ   density classes    (99 if you don't want) 
    77                                                       zsigp           ,&! potential density classes    (99 if you don't want) 
    78                                                       zsal            ,&! salinity classes   (99 if you don't want) 
    79                                                       ztem            ,&! temperature classes(99 if you don't want) 
    80                                                       zlay              ! level classes      (99 if you don't want) 
     75      REAL(wp), DIMENSION(nb_class_max)            :: zsigi             ! in-situ   density classes    (99 if you don't want) 
     76      REAL(wp), DIMENSION(nb_class_max)            :: zsigp             ! potential density classes    (99 if you don't want) 
     77      REAL(wp), DIMENSION(nb_class_max)            :: zsal              ! salinity classes   (99 if you don't want) 
     78      REAL(wp), DIMENSION(nb_class_max)            :: ztem              ! temperature classes(99 if you don't want) 
     79      REAL(wp), DIMENSION(nb_class_max)            :: zlay              ! level classes      (99 if you don't want) 
    8180      REAL(wp), DIMENSION(nb_type_class,nb_class_max)  :: transport     ! transport output 
    8281      REAL(wp)                                         :: slopeSection  ! slope of the section 
     
    9089   REAL(wp), ALLOCATABLE, DIMENSION(:,:,:)   ::  transports_2d   
    9190 
     91 
     92   !! * Substitutions 
     93#  include "domzgr_substitute.h90" 
    9294   !!---------------------------------------------------------------------- 
    9395   !! NEMO/OCE 4.0 , NEMO Consortium (2018) 
     
    9597   !! Software governed by the CeCILL license (see ./LICENSE) 
    9698   !!---------------------------------------------------------------------- 
     99 
    97100CONTAINS 
    98101  
     
    11191122  !!    |               |                  |       interpolation between ptab(I,J,K) and ptab(I,J,K+1) 
    11201123  !!    |               |                  |       zbis =  
    1121   !!    |               |                  |      [ e3w(I+1,J,K)*ptab(I,J,K) + ( e3w(I,J,K) - e3w(I+1,J,K) ) * ptab(I,J,K-1) ] 
    1122   !!    |               |                  |      /[ e3w(I+1,J,K) + e3w(I,J,K) - e3w(I+1,J,K) ]  
     1124  !!    |               |                  |      [ e3w_n(I+1,J,K,NOW)*ptab(I,J,K) + ( e3w_n(I,J,K,NOW) - e3w_n(I+1,J,K,NOW) ) * ptab(I,J,K-1) ] 
     1125  !!    |               |                  |     /[ e3w_n(I+1,J,K,NOW)             +   e3w_n(I,J,K,NOW) - e3w_n(I+1,J,K,NOW) ]  
    11231126  !!    |               |                  |  
    11241127  !!    |               |                  |    2. Horizontal interpolation: compute value at U/V point 
     
    12131216 
    12141217     ze3t  = e3t(ii2,ij2,kk,Kmm) - e3t(ii1,ij1,kk,Kmm)  
    1215      zwgt1 = ( e3w(ii2,ij2,kk,Kmm) - e3w(ii1,ij1,kk,Kmm) ) / e3w(ii2,ij2,kk,Kmm) 
    1216      zwgt2 = ( e3w(ii1,ij1,kk,Kmm) - e3w(ii2,ij2,kk,Kmm) ) / e3w(ii1,ij1,kk,Kmm) 
     1218     zwgt1 = ( e3w(ii2,ij2,kk,Kmm) - e3w(ii1,ij1,kk,Kmm) )   & 
     1219        &    / e3w(ii2,ij2,kk,Kmm) 
     1220     zwgt2 = ( e3w(ii1,ij1,kk,Kmm) - e3w(ii2,ij2,kk,Kmm) )   & 
     1221        &    / e3w(ii1,ij1,kk,Kmm) 
    12171222 
    12181223     IF(kk .NE. 1)THEN 
  • NEMO/trunk/src/OCE/DIA/diahsb.F90

    r12489 r13237  
    5050   REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::   tmask_ini 
    5151 
     52   !! * Substitutions 
     53#  include "domzgr_substitute.h90" 
    5254   !!---------------------------------------------------------------------- 
    5355   !! NEMO/OCE 4.0 , NEMO Consortium (2018) 
     
    156158      ! 
    157159      DO jk = 1, jpkm1           ! volume variation (calculated with scale factors) 
    158          zwrk(:,:,jk) = surf(:,:)*e3t(:,:,jk,Kmm)*tmask(:,:,jk) - surf_ini(:,:)*e3t_ini(:,:,jk)*tmask_ini(:,:,jk) 
     160         zwrk(:,:,jk) =   surf    (:,:) * e3t    (:,:,jk,Kmm)*tmask    (:,:,jk)   & 
     161            &           - surf_ini(:,:) * e3t_ini(:,:,jk    )*tmask_ini(:,:,jk) 
    159162      END DO 
    160163      zdiff_v2 = glob_sum_full( 'diahsb', zwrk(:,:,:) )     ! glob_sum_full needed as tmask and tmask_ini could be different 
    161164      DO jk = 1, jpkm1           ! heat content variation 
    162          zwrk(:,:,jk) = ( surf(:,:)*e3t(:,:,jk,Kmm)*ts(:,:,jk,jp_tem,Kmm) - surf_ini(:,:)*hc_loc_ini(:,:,jk) ) 
     165         zwrk(:,:,jk) = ( surf    (:,:) * e3t(:,:,jk,Kmm)*ts(:,:,jk,jp_tem,Kmm)   & 
     166            &           - surf_ini(:,:) *         hc_loc_ini(:,:,jk) ) 
    163167      END DO 
    164168      zdiff_hc = glob_sum_full( 'diahsb', zwrk(:,:,:) ) 
    165169      DO jk = 1, jpkm1           ! salt content variation 
    166          zwrk(:,:,jk) = ( surf(:,:)*e3t(:,:,jk,Kmm)*ts(:,:,jk,jp_sal,Kmm) - surf_ini(:,:)*sc_loc_ini(:,:,jk) ) 
     170         zwrk(:,:,jk) = ( surf    (:,:) * e3t(:,:,jk,Kmm)*ts(:,:,jk,jp_sal,Kmm)   & 
     171            &           - surf_ini(:,:) *         sc_loc_ini(:,:,jk) ) 
    167172      END DO 
    168173      zdiff_sc = glob_sum_full( 'diahsb', zwrk(:,:,:) ) 
  • NEMO/trunk/src/OCE/DIA/diahth.F90

    r12489 r13237  
    4242   !! * Substitutions 
    4343#  include "do_loop_substitute.h90" 
     44#  include "domzgr_substitute.h90" 
    4445   !!---------------------------------------------------------------------- 
    4546   !! NEMO/OCE 4.0 , NEMO Consortium (2018) 
     
    361362         ik = ilevel(ji,jj) 
    362363         zthick(ji,jj) = pdep - zthick(ji,jj)   !   remaining thickness to reach depht pdep 
    363          phtc(ji,jj)   = phtc(ji,jj) + pt(ji,jj,ik+1) * MIN( e3t(ji,jj,ik+1,Kmm), zthick(ji,jj) ) & 
     364         phtc(ji,jj)   = phtc(ji,jj)    & 
     365            &           + pt (ji,jj,ik+1) * MIN( e3t(ji,jj,ik+1,Kmm), zthick(ji,jj) ) & 
    364366                                                       * tmask(ji,jj,ik+1) 
    365367      END_2D 
  • NEMO/trunk/src/OCE/DIA/diamlr.F90

    r12642 r13237  
    44   !! Management of the IOM context for multiple-linear-regression analysis 
    55   !!====================================================================== 
    6    !! History :       !  2019  (S. Mueller) 
     6   !! History :  4.0  !  2019  (S. Mueller)   Original code 
    77   !!---------------------------------------------------------------------- 
    88 
    99   USE par_oce        , ONLY :   wp, jpi, jpj 
    1010   USE phycst         , ONLY :   rpi 
     11   USE dom_oce        , ONLY :   adatrj 
     12   USE tide_mod 
     13   ! 
    1114   USE in_out_manager , ONLY :   lwp, numout, ln_timing 
    1215   USE iom            , ONLY :   iom_put, iom_use, iom_update_file_name 
    13    USE dom_oce        , ONLY :   adatrj 
    1416   USE timing         , ONLY :   timing_start, timing_stop 
    1517#if defined key_iomput 
    1618   USE xios 
    1719#endif 
    18    USE tide_mod 
    1920 
    2021   IMPLICIT NONE 
    2122   PRIVATE 
    2223 
    23    LOGICAL, PUBLIC ::   lk_diamlr = .FALSE. 
     24   LOGICAL, PUBLIC ::   lk_diamlr = .FALSE.   !:         ===>>>   NOT a DOCTOR norm name :  use l_diamlr 
     25   !                                                              lk_  is used only for logical controlled by a CPP key 
    2426 
    2527   PUBLIC ::   dia_mlr_init, dia_mlr_iom_init, dia_mlr 
     
    4244      !! 
    4345      !!---------------------------------------------------------------------- 
    44  
     46      ! 
    4547      lk_diamlr = .TRUE. 
    46  
     48      ! 
    4749      IF(lwp) THEN 
    4850         WRITE(numout, *) 
     
    5052         WRITE(numout, *) '~~~~~~~~~~~~   multiple-linear-regression analysis' 
    5153      END IF 
    52  
     54      ! 
    5355   END SUBROUTINE dia_mlr_init 
     56 
    5457 
    5558   SUBROUTINE dia_mlr_iom_init 
     
    397400   END SUBROUTINE dia_mlr_iom_init 
    398401 
     402 
    399403   SUBROUTINE dia_mlr 
    400404      !!---------------------------------------------------------------------- 
     
    404408      !! 
    405409      !!---------------------------------------------------------------------- 
    406  
    407410      REAL(wp), DIMENSION(jpi,jpj) ::   zadatrj2d 
     411      !!---------------------------------------------------------------------- 
    408412 
    409413      IF( ln_timing )   CALL timing_start('dia_mlr') 
     
    412416      ! (value of adatrj converted to time in units of seconds) 
    413417      ! 
    414       ! A 2-dimensional field of constant value is sent, and subsequently used 
    415       ! directly or transformed to a scalar or a constant 3-dimensional field as 
    416       ! required. 
     418      ! A 2-dimensional field of constant value is sent, and subsequently used directly  
     419      ! or transformed to a scalar or a constant 3-dimensional field as required. 
    417420      zadatrj2d(:,:) = adatrj*86400.0_wp 
    418421      IF ( iom_use('diamlr_time') ) CALL iom_put('diamlr_time', zadatrj2d) 
    419        
     422      ! 
    420423      IF( ln_timing )   CALL timing_stop('dia_mlr') 
    421  
     424      ! 
    422425   END SUBROUTINE dia_mlr 
    423426 
     427   !!====================================================================== 
    424428END MODULE diamlr 
  • NEMO/trunk/src/OCE/DIA/diaptr.F90

    r13226 r13237  
    6060 
    6161   LOGICAL ::   ll_init = .TRUE.        !: tracers  trend flag (set from namelist in trdini) 
     62    
    6263   !! * Substitutions 
    6364#  include "do_loop_substitute.h90" 
     65#  include "domzgr_substitute.h90" 
    6466   !!---------------------------------------------------------------------- 
    6567   !! NEMO/OCE 4.0 , NEMO Consortium (2018) 
  • NEMO/trunk/src/OCE/DIA/diawri.F90

    r13226 r13237  
    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 
    459487 
    460488 
     
    571599         DEALLOCATE(zw3d_abl) 
    572600         ENDIF 
     601         ! 
    573602 
    574603         ! Declare all the output fields as NETCDF variables 
     
    580609            &          jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout ) 
    581610         IF(  .NOT.ln_linssh  ) THEN 
    582             CALL histdef( nid_T, "vovvle3t", "Level thickness"                    , "m"      ,&  ! e3t(:,:,:,Kmm) 
     611            CALL histdef( nid_T, "vovvle3t", "Level thickness"                    , "m"      ,&  ! e3t n 
    583612            &             jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout ) 
    584             CALL histdef( nid_T, "vovvldep", "T point depth"                      , "m"      ,&  ! e3t(:,:,:,Kmm) 
     613            CALL histdef( nid_T, "vovvldep", "T point depth"                      , "m"      ,&  ! e3t n 
    585614            &             jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout ) 
    586             CALL histdef( nid_T, "vovvldef", "Squared level deformation"          , "%^2"    ,&  ! e3t(:,:,:,Kmm) 
     615            CALL histdef( nid_T, "vovvldef", "Squared level deformation"          , "%^2"    ,&  ! e3t n 
    587616            &             jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout ) 
    588617         ENDIF 
     
    768797 
    769798      IF( .NOT.ln_linssh ) THEN 
    770          CALL histwrite( nid_T, "votemper", it, ts(:,:,:,jp_tem,Kmm) * e3t(:,:,:,Kmm) , ndim_T , ndex_T  )   ! heat content 
    771          CALL histwrite( nid_T, "vosaline", it, ts(:,:,:,jp_sal,Kmm) * e3t(:,:,:,Kmm) , ndim_T , ndex_T  )   ! salt content 
    772          CALL histwrite( nid_T, "sosstsst", it, ts(:,:,1,jp_tem,Kmm) * e3t(:,:,1,Kmm) , ndim_hT, ndex_hT )   ! sea surface heat content 
    773          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 
    774803      ELSE 
    775804         CALL histwrite( nid_T, "votemper", it, ts(:,:,:,jp_tem,Kmm) , ndim_T , ndex_T  )   ! temperature 
     
    779808      ENDIF 
    780809      IF( .NOT.ln_linssh ) THEN 
    781          zw3d(:,:,:) = ( ( e3t(:,:,:,Kmm) - e3t_0(:,:,:) ) / e3t_0(:,:,:) * 100 * tmask(:,:,:) ) ** 2 
    782          CALL histwrite( nid_T, "vovvle3t", it, e3t (:,:,:,Kmm) , ndim_T , ndex_T  )   ! level thickness 
    783          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  
    784813         CALL histwrite( nid_T, "vovvldef", it, zw3d             , ndim_T , ndex_T  )   ! level thickness deformation 
    785814      ENDIF 
     
    920949      !! 
    921950      INTEGER :: inum, jk 
     951      REAL(wp), DIMENSION(jpi,jpj,jpk) :: ze3t, zgdept      ! 3D workspace !!st patch to use substitution 
    922952      !!---------------------------------------------------------------------- 
    923953      !  
    924       IF(lwp) WRITE(numout,*) 
    925       IF(lwp) WRITE(numout,*) 'dia_wri_state : single instantaneous ocean state' 
    926       IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~   and forcing fields file created ' 
    927       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 
    928965      ! 
    929966      CALL iom_open( TRIM(cdfile_name), inum, ldwrt = .TRUE. ) 
     
    940977      ENDIF 
    941978      CALL iom_rstput( 0, 0, inum, 'risfdep', risfdep            )    ! now k-velocity 
    942       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 
    943980      ! 
    944981      IF ( ln_isf ) THEN 
     
    9771014      CALL iom_rstput( 0, 0, inum, 'sometauy', vtau              )    ! j-wind stress 
    9781015      IF(  .NOT.ln_linssh  ) THEN              
    979          CALL iom_rstput( 0, 0, inum, 'vovvldep', gdept(:,:,:,Kmm)        )    !  T-cell depth  
    980          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   
    9811018      END IF 
    9821019      IF( ln_wave .AND. ln_sdw ) THEN 
Note: See TracChangeset for help on using the changeset viewer.