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 9019 for branches/2017/dev_merge_2017/NEMOGCM/NEMO/OPA_SRC/DIA/diawri.F90 – NEMO

Ignore:
Timestamp:
2017-12-13T15:58:53+01:00 (6 years ago)
Author:
timgraham
Message:

Merge of dev_CNRS_2017 into branch

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/2017/dev_merge_2017/NEMOGCM/NEMO/OPA_SRC/DIA/diawri.F90

    r8465 r9019  
    2525   !!   dia_wri_state : create an output NetCDF file for a single instantaeous ocean state and forcing fields 
    2626   !!---------------------------------------------------------------------- 
    27    USE oce             ! ocean dynamics and tracers  
    28    USE dom_oce         ! ocean space and time domain 
    29    USE dynadv, ONLY: ln_dynadv_vec 
    30    USE zdf_oce         ! ocean vertical physics 
    31    USE ldftra          ! lateral physics: eddy diffusivity coef. 
    32    USE ldfdyn          ! lateral physics: eddy viscosity   coef. 
    33    USE sbc_oce         ! Surface boundary condition: ocean fields 
    34    USE sbc_ice         ! Surface boundary condition: ice fields 
    35    USE icb_oce         ! Icebergs 
    36    USE icbdia          ! Iceberg budgets 
    37    USE sbcssr          ! restoring term toward SST/SSS climatology 
    38    USE phycst          ! physical constants 
    39    USE zdfmxl          ! mixed layer 
    40    USE dianam          ! build name of file (routine) 
    41    USE zdfddm          ! vertical  physics: double diffusion 
    42    USE diahth          ! thermocline diagnostics 
    43    USE wet_dry         ! wetting and drying 
    44    USE sbcwave         ! wave parameters 
     27   USE oce            ! ocean dynamics and tracers  
     28   USE dom_oce        ! ocean space and time domain 
     29   USE phycst         ! physical constants 
     30   USE dianam         ! build name of file (routine) 
     31   USE diahth         ! thermocline diagnostics 
     32   USE dynadv   , ONLY: ln_dynadv_vec 
     33   USE icb_oce        ! Icebergs 
     34   USE icbdia         ! Iceberg budgets 
     35   USE ldftra         ! lateral physics: eddy diffusivity coef. 
     36   USE ldfdyn         ! lateral physics: eddy viscosity   coef. 
     37   USE sbc_oce        ! Surface boundary condition: ocean fields 
     38   USE sbc_ice        ! Surface boundary condition: ice fields 
     39   USE sbcssr         ! restoring term toward SST/SSS climatology 
     40   USE sbcwave        ! wave parameters 
     41   USE wet_dry        ! wetting and drying 
     42   USE zdf_oce        ! ocean vertical physics 
     43   USE zdfdrg         ! ocean vertical physics: top/bottom friction 
     44   USE zdfmxl         ! mixed layer 
    4545   ! 
    46    USE lbclnk          ! ocean lateral boundary conditions (or mpp link) 
    47    USE in_out_manager  ! I/O manager 
    48    USE diatmb          ! Top,middle,bottom output 
    49    USE dia25h          ! 25h Mean output 
    50    USE iom 
    51    USE ioipsl 
    52  
    53 #if defined key_lim2 
    54    USE limwri_2  
    55 #elif defined key_lim3 
    56    USE limwri  
     46   USE lbclnk         ! ocean lateral boundary conditions (or mpp link) 
     47   USE in_out_manager ! I/O manager 
     48   USE diatmb         ! Top,middle,bottom output 
     49   USE dia25h         ! 25h Mean output 
     50   USE iom            !  
     51   USE ioipsl         !  
     52 
     53#if defined key_lim3 
     54   USE icewri  
    5755#endif 
    5856   USE lib_mpp         ! MPP library 
     
    6058   USE diurnal_bulk    ! diurnal warm layer 
    6159   USE cool_skin       ! Cool skin 
    62    USE wrk_nemo        ! working array 
    6360 
    6461   IMPLICIT NONE 
     
    8077 
    8178   !! * Substitutions 
    82 #  include "zdfddm_substitute.h90" 
    8379#  include "vectopt_loop_substitute.h90" 
    8480   !!---------------------------------------------------------------------- 
     
    120116      !! ** Method  :  use iom_put 
    121117      !!---------------------------------------------------------------------- 
    122       !! 
    123118      INTEGER, INTENT( in ) ::   kt      ! ocean time-step index 
    124119      !! 
    125       INTEGER                      ::   ji, jj, jk              ! dummy loop indices 
    126       INTEGER                      ::   jkbot                   ! 
    127       REAL(wp)                     ::   zztmp, zztmpx, zztmpy   !  
    128       !! 
    129       REAL(wp), POINTER, DIMENSION(:,:)   :: z2d      ! 2D workspace 
    130       REAL(wp), POINTER, DIMENSION(:,:,:) :: z3d      ! 3D workspace 
     120      INTEGER ::   ji, jj, jk       ! dummy loop indices 
     121      INTEGER ::   ikbot            ! local integer 
     122      REAL(wp)::   zztmp , zztmpx   ! local scalar 
     123      REAL(wp)::   zztmp2, zztmpy   !   -      - 
     124      REAL(wp), DIMENSION(jpi,jpj)     ::   z2d   ! 2D workspace 
     125      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   z3d   ! 3D workspace 
    131126      !!---------------------------------------------------------------------- 
    132127      !  
    133128      IF( nn_timing == 1 )   CALL timing_start('dia_wri') 
    134129      !  
    135       CALL wrk_alloc( jpi , jpj      , z2d ) 
    136       CALL wrk_alloc( jpi , jpj, jpk , z3d ) 
    137       ! 
    138130      ! Output the initial state and forcings 
    139131      IF( ninist == 1 ) THEN                        
     
    163155         DO jj = 1, jpj 
    164156            DO ji = 1, jpi 
    165                jkbot = mbkt(ji,jj) 
    166                z2d(ji,jj) = tsn(ji,jj,jkbot,jp_tem) 
     157               ikbot = mbkt(ji,jj) 
     158               z2d(ji,jj) = tsn(ji,jj,ikbot,jp_tem) 
    167159            END DO 
    168160         END DO 
     
    175167         DO jj = 1, jpj 
    176168            DO ji = 1, jpi 
    177                jkbot = mbkt(ji,jj) 
    178                z2d(ji,jj) = tsn(ji,jj,jkbot,jp_sal) 
     169               ikbot = mbkt(ji,jj) 
     170               z2d(ji,jj) = tsn(ji,jj,ikbot,jp_sal) 
    179171            END DO 
    180172         END DO 
     
    183175 
    184176      IF ( iom_use("taubot") ) THEN                ! bottom stress 
     177         zztmp = rau0 * 0.25 
    185178         z2d(:,:) = 0._wp 
    186179         DO jj = 2, jpjm1 
    187180            DO ji = fs_2, fs_jpim1   ! vector opt. 
    188                zztmpx = (  bfrua(ji  ,jj) * un(ji  ,jj,mbku(ji  ,jj))  & 
    189                       &  + bfrua(ji-1,jj) * un(ji-1,jj,mbku(ji-1,jj))  )       
    190                zztmpy = (  bfrva(ji,  jj) * vn(ji,jj  ,mbkv(ji,jj  ))  & 
    191                       &  + bfrva(ji,jj-1) * vn(ji,jj-1,mbkv(ji,jj-1))  )  
    192                z2d(ji,jj) = rau0 * SQRT( zztmpx * zztmpx + zztmpy * zztmpy ) * tmask(ji,jj,1)  
     181               zztmp2 = (  ( rCdU_bot(ji+1,jj)+rCdU_bot(ji  ,jj) ) * un(ji  ,jj,mbku(ji  ,jj))  )**2   & 
     182                  &   + (  ( rCdU_bot(ji  ,jj)+rCdU_bot(ji-1,jj) ) * un(ji-1,jj,mbku(ji-1,jj))  )**2   & 
     183                  &   + (  ( rCdU_bot(ji,jj+1)+rCdU_bot(ji,jj  ) ) * vn(ji,jj  ,mbkv(ji,jj  ))  )**2   & 
     184                  &   + (  ( rCdU_bot(ji,jj  )+rCdU_bot(ji,jj-1) ) * vn(ji,jj-1,mbkv(ji,jj-1))  )**2 
     185               z2d(ji,jj) = zztmp * SQRT( zztmp2 ) * tmask(ji,jj,1)  
    193186               ! 
    194             ENDDO 
    195          ENDDO 
     187            END DO 
     188         END DO 
    196189         CALL lbc_lnk( z2d, 'T', 1. ) 
    197190         CALL iom_put( "taubot", z2d )            
    198191      ENDIF 
    199192          
    200       CALL iom_put( "uoce", un(:,:,:)         )    ! 3D i-current 
    201       CALL iom_put(  "ssu", un(:,:,1)         )    ! surface i-current 
     193      CALL iom_put( "uoce", un(:,:,:) )            ! 3D i-current 
     194      CALL iom_put(  "ssu", un(:,:,1) )            ! surface i-current 
    202195      IF ( iom_use("sbu") ) THEN 
    203196         DO jj = 1, jpj 
    204197            DO ji = 1, jpi 
    205                jkbot = mbku(ji,jj) 
    206                z2d(ji,jj) = un(ji,jj,jkbot) 
     198               ikbot = mbku(ji,jj) 
     199               z2d(ji,jj) = un(ji,jj,ikbot) 
    207200            END DO 
    208201         END DO 
     
    210203      ENDIF 
    211204       
    212       CALL iom_put( "voce", vn(:,:,:)         )    ! 3D j-current 
    213       CALL iom_put(  "ssv", vn(:,:,1)         )    ! surface j-current 
     205      CALL iom_put( "voce", vn(:,:,:) )            ! 3D j-current 
     206      CALL iom_put(  "ssv", vn(:,:,1) )            ! surface j-current 
    214207      IF ( iom_use("sbv") ) THEN 
    215208         DO jj = 1, jpj 
    216209            DO ji = 1, jpi 
    217                jkbot = mbkv(ji,jj) 
    218                z2d(ji,jj) = vn(ji,jj,jkbot) 
     210               ikbot = mbkv(ji,jj) 
     211               z2d(ji,jj) = vn(ji,jj,ikbot) 
    219212            END DO 
    220213         END DO 
     
    233226      ENDIF 
    234227 
    235       CALL iom_put( "avt" , avt                        )    ! T vert. eddy diff. coef. 
    236       CALL iom_put( "avm" , avmu                       )    ! T vert. eddy visc. coef. 
    237       CALL iom_put( "avs" , fsavs(:,:,:)               )    ! S vert. eddy diff. coef. (useful only with key_zdfddm) 
    238  
    239       IF( iom_use('logavt') )   CALL iom_put( "logavt", LOG( MAX( 1.e-20_wp, avt  (:,:,:) ) ) ) 
    240       IF( iom_use('logavs') )   CALL iom_put( "logavs", LOG( MAX( 1.e-20_wp, fsavs(:,:,:) ) ) ) 
     228      CALL iom_put( "avt" , avt )                  ! T vert. eddy diff. coef. 
     229      CALL iom_put( "avs" , avs )                  ! S vert. eddy diff. coef. 
     230      CALL iom_put( "avm" , avm )                  ! T vert. eddy visc. coef. 
     231 
     232      IF( iom_use('logavt') )   CALL iom_put( "logavt", LOG( MAX( 1.e-20_wp, avt(:,:,:) ) ) ) 
     233      IF( iom_use('logavs') )   CALL iom_put( "logavs", LOG( MAX( 1.e-20_wp, avs(:,:,:) ) ) ) 
    241234 
    242235      IF ( iom_use("sstgrad") .OR. iom_use("sstgrad2") ) THEN 
     
    251244         END DO 
    252245         CALL lbc_lnk( z2d, 'T', 1. ) 
    253          CALL iom_put( "sstgrad2",  z2d               )    ! square of module of sst gradient 
     246         CALL iom_put( "sstgrad2",  z2d )          ! square of module of sst gradient 
    254247         z2d(:,:) = SQRT( z2d(:,:) ) 
    255          CALL iom_put( "sstgrad" ,  z2d               )    ! module of sst gradient 
     248         CALL iom_put( "sstgrad" ,  z2d )          ! module of sst gradient 
    256249      ENDIF 
    257250          
    258       ! clem: heat and salt content 
     251      ! heat and salt contents 
    259252      IF( iom_use("heatc") ) THEN 
    260253         z2d(:,:)  = 0._wp  
     
    266259            END DO 
    267260         END DO 
    268          CALL iom_put( "heatc", (rau0 * rcp) * z2d )    ! vertically integrated heat content (J/m2) 
     261         CALL iom_put( "heatc", rau0_rcp * z2d )   ! vertically integrated heat content (J/m2) 
    269262      ENDIF 
    270263 
     
    278271            END DO 
    279272         END DO 
    280          CALL iom_put( "saltc", rau0 * z2d )   ! vertically integrated salt content (PSU*kg/m2) 
     273         CALL iom_put( "saltc", rau0 * z2d )          ! vertically integrated salt content (PSU*kg/m2) 
    281274      ENDIF 
    282275      ! 
    283276      IF ( iom_use("eken") ) THEN 
    284          rke(:,:,jpk) = 0._wp                               !      kinetic energy  
     277         z3d(:,:,jk) = 0._wp  
    285278         DO jk = 1, jpkm1 
    286279            DO jj = 2, jpjm1 
    287280               DO ji = fs_2, fs_jpim1   ! vector opt. 
    288                   zztmp   = r1_e1e2t(ji,jj) / e3t_n(ji,jj,jk) 
    289                   zztmpx  = 0.5 * (  un(ji-1,jj,jk) * un(ji-1,jj,jk) * e1e2u(ji-1,jj) * e3u_n(ji-1,jj,jk)    & 
    290                      &             + un(ji  ,jj,jk) * un(ji  ,jj,jk) * e1e2u(ji  ,jj) * e3u_n(ji  ,jj,jk) )  & 
    291                      &          *  zztmp  
    292                   ! 
    293                   zztmpy  = 0.5 * (  vn(ji,jj-1,jk) * vn(ji,jj-1,jk) * e1e2v(ji,jj-1) * e3v_n(ji,jj-1,jk)    & 
    294                      &             + vn(ji,jj  ,jk) * vn(ji,jj  ,jk) * e1e2v(ji,jj  ) * e3v_n(ji,jj  ,jk) )  & 
    295                      &          *  zztmp  
    296                   ! 
    297                   rke(ji,jj,jk) = 0.5_wp * ( zztmpx + zztmpy ) 
    298                   ! 
    299                ENDDO 
    300             ENDDO 
    301          ENDDO 
    302          CALL lbc_lnk( rke, 'T', 1. ) 
    303          CALL iom_put( "eken", rke )            
     281                  zztmp  = 0.25_wp * r1_e1e2t(ji,jj) / e3t_n(ji,jj,jk) 
     282                  z3d(ji,jj,jk) = zztmp * (  un(ji-1,jj,jk)**2 * e2u(ji-1,jj) * e3u_n(ji-1,jj,jk)   & 
     283                     &                     + un(ji  ,jj,jk)**2 * e2u(ji  ,jj) * e3u_n(ji  ,jj,jk)   & 
     284                     &                     + vn(ji,jj-1,jk)**2 * e1v(ji,jj-1) * e3v_n(ji,jj-1,jk)   & 
     285                     &                     + vn(ji,jj  ,jk)**2 * e1v(ji,jj  ) * e3v_n(ji,jj  ,jk)   ) 
     286               END DO 
     287            END DO 
     288         END DO 
     289         CALL lbc_lnk( z3d, 'T', 1. ) 
     290         CALL iom_put( "eken", z3d )                 ! kinetic energy 
    304291      ENDIF 
    305292      ! 
     
    313300            z2d(:,:) = z2d(:,:) + z3d(:,:,jk) 
    314301         END DO 
    315          CALL iom_put( "u_masstr", z3d )                  ! mass transport in i-direction 
    316          CALL iom_put( "u_masstr_vint", z2d )             ! mass transport in i-direction vertical sum 
     302         CALL iom_put( "u_masstr"     , z3d )         ! mass transport in i-direction 
     303         CALL iom_put( "u_masstr_vint", z2d )         ! mass transport in i-direction vertical sum 
    317304      ENDIF 
    318305       
    319306      IF( iom_use("u_heattr") ) THEN 
    320          z2d(:,:) = 0.e0  
     307         z2d(:,:) = 0._wp  
    321308         DO jk = 1, jpkm1 
    322309            DO jj = 2, jpjm1 
     
    327314         END DO 
    328315         CALL lbc_lnk( z2d, 'U', -1. ) 
    329          CALL iom_put( "u_heattr", (0.5 * rcp) * z2d )    ! heat transport in i-direction 
     316         CALL iom_put( "u_heattr", 0.5*rcp * z2d )    ! heat transport in i-direction 
    330317      ENDIF 
    331318 
     
    340327         END DO 
    341328         CALL lbc_lnk( z2d, 'U', -1. ) 
    342          CALL iom_put( "u_salttr", 0.5 * z2d )            ! heat transport in i-direction 
     329         CALL iom_put( "u_salttr", 0.5 * z2d )        ! heat transport in i-direction 
    343330      ENDIF 
    344331 
     
    349336            z3d(:,:,jk) = rau0 * vn(:,:,jk) * e1v(:,:) * e3v_n(:,:,jk) * vmask(:,:,jk) 
    350337         END DO 
    351          CALL iom_put( "v_masstr", z3d )                  ! mass transport in j-direction 
     338         CALL iom_put( "v_masstr", z3d )              ! mass transport in j-direction 
    352339      ENDIF 
    353340       
     
    362349         END DO 
    363350         CALL lbc_lnk( z2d, 'V', -1. ) 
    364          CALL iom_put( "v_heattr", (0.5 * rcp) * z2d )    !  heat transport in j-direction 
     351         CALL iom_put( "v_heattr", 0.5*rcp * z2d )    !  heat transport in j-direction 
    365352      ENDIF 
    366353 
    367354      IF( iom_use("v_salttr") ) THEN 
    368          z2d(:,:) = 0.e0  
     355         z2d(:,:) = 0._wp  
    369356         DO jk = 1, jpkm1 
    370357            DO jj = 2, jpjm1 
     
    375362         END DO 
    376363         CALL lbc_lnk( z2d, 'V', -1. ) 
    377          CALL iom_put( "v_salttr", 0.5 * z2d )            !  heat transport in j-direction 
    378       ENDIF 
    379  
    380       ! Vertical integral of temperature 
     364         CALL iom_put( "v_salttr", 0.5 * z2d )        !  heat transport in j-direction 
     365      ENDIF 
     366 
    381367      IF( iom_use("tosmint") ) THEN 
    382          z2d(:,:)=0._wp 
     368         z2d(:,:) = 0._wp 
    383369         DO jk = 1, jpkm1 
    384370            DO jj = 2, jpjm1 
    385371               DO ji = fs_2, fs_jpim1   ! vector opt. 
    386                   z2d(ji,jj) = z2d(ji,jj) + rau0 * e3t_n(ji,jj,jk) *  tsn(ji,jj,jk,jp_tem) 
     372                  z2d(ji,jj) = z2d(ji,jj) + e3t_n(ji,jj,jk) *  tsn(ji,jj,jk,jp_tem) 
    387373               END DO 
    388374            END DO 
    389375         END DO 
    390376         CALL lbc_lnk( z2d, 'T', -1. ) 
    391          CALL iom_put( "tosmint", z2d )  
    392       ENDIF 
    393  
    394       ! Vertical integral of salinity 
     377         CALL iom_put( "tosmint", rau0 * z2d )        ! Vertical integral of temperature 
     378      ENDIF 
    395379      IF( iom_use("somint") ) THEN 
    396380         z2d(:,:)=0._wp 
     
    398382            DO jj = 2, jpjm1 
    399383               DO ji = fs_2, fs_jpim1   ! vector opt. 
    400                   z2d(ji,jj) = z2d(ji,jj) + rau0 * e3t_n(ji,jj,jk) * tsn(ji,jj,jk,jp_sal) 
     384                  z2d(ji,jj) = z2d(ji,jj) + e3t_n(ji,jj,jk) * tsn(ji,jj,jk,jp_sal) 
    401385               END DO 
    402386            END DO 
    403387         END DO 
    404388         CALL lbc_lnk( z2d, 'T', -1. ) 
    405          CALL iom_put( "somint", z2d )  
    406       ENDIF 
    407  
    408       CALL iom_put( "bn2", rn2 )  !Brunt-Vaisala buoyancy frequency (N^2) 
    409       ! 
    410       CALL wrk_dealloc( jpi , jpj      , z2d ) 
    411       CALL wrk_dealloc( jpi , jpj, jpk , z3d ) 
    412       ! 
    413       ! If we want tmb values  
    414  
    415       IF (ln_diatmb) THEN 
    416          CALL dia_tmb  
    417       ENDIF  
    418       IF (ln_dia25h) THEN 
    419          CALL dia_25h( kt ) 
    420       ENDIF  
     389         CALL iom_put( "somint", rau0 * z2d )         ! Vertical integral of salinity 
     390      ENDIF 
     391 
     392      CALL iom_put( "bn2", rn2 )                      ! Brunt-Vaisala buoyancy frequency (N^2) 
     393      ! 
     394 
     395      IF (ln_diatmb)   CALL dia_tmb                   ! tmb values  
     396           
     397      IF (ln_dia25h)   CALL dia_25h( kt )             ! 25h averaging 
    421398 
    422399      IF( nn_timing == 1 )   CALL timing_stop('dia_wri') 
     
    452429      REAL(wp) ::   zsto, zout, zmax, zjulian                ! local scalars 
    453430      ! 
    454       REAL(wp), POINTER, DIMENSION(:,:)   :: zw2d       ! 2D workspace 
    455       REAL(wp), POINTER, DIMENSION(:,:,:) :: zw3d       ! 3D workspace 
     431      REAL(wp), DIMENSION(jpi,jpj)   :: zw2d       ! 2D workspace 
     432      REAL(wp), DIMENSION(jpi,jpj,jpk) :: zw3d       ! 3D workspace 
    456433      !!---------------------------------------------------------------------- 
    457434      !  
    458435      IF( nn_timing == 1 )   CALL timing_start('dia_wri') 
    459436      ! 
    460                              CALL wrk_alloc( jpi,jpj      , zw2d ) 
    461       IF( .NOT.ln_linssh )   CALL wrk_alloc( jpi,jpj,jpk  , zw3d ) 
    462       ! 
    463       ! Output the initial state and forcings 
    464       IF( ninist == 1 ) THEN                        
     437      IF( ninist == 1 ) THEN     !==  Output the initial state and forcings  ==! 
    465438         CALL dia_wri_state( 'output.init', kt ) 
    466439         ninist = 0 
     
    470443      ! ----------------- 
    471444 
    472       ! local variable for debugging 
    473       ll_print = .FALSE. 
     445      ll_print = .FALSE.                  ! local variable for debugging 
    474446      ll_print = ll_print .AND. lwp 
    475447 
     
    707679#endif 
    708680 
    709          IF( ln_cpl .AND. nn_ice == 2 ) THEN 
    710             CALL histdef( nid_T,"soicetem" , "Ice Surface Temperature"            , "K"      ,   &  ! tn_ice 
    711                &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
    712             CALL histdef( nid_T,"soicealb" , "Ice Albedo"                         , "[0,1]"  ,   &  ! alb_ice 
    713                &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
    714          ENDIF 
    715  
    716681         CALL histend( nid_T, snc4chunks=snc4set ) 
    717682 
     
    747712         CALL histdef( nid_W, "votkeavt", "Vertical Eddy Diffusivity"          , "m2/s"   ,   &  ! avt 
    748713            &          jpi, jpj, nh_W, ipk, 1, ipk, nz_W, 32, clop, zsto, zout ) 
    749          CALL histdef( nid_W, "votkeavm", "Vertical Eddy Viscosity"             , "m2/s"  ,   &  ! avmu 
     714         CALL histdef( nid_W, "votkeavm", "Vertical Eddy Viscosity"             , "m2/s"  ,   &  ! avm 
    750715            &          jpi, jpj, nh_W, ipk, 1, ipk, nz_W, 32, clop, zsto, zout ) 
    751716 
    752          IF( lk_zdfddm ) THEN 
     717         IF( ln_zdfddm ) THEN 
    753718            CALL histdef( nid_W,"voddmavs","Salt Vertical Eddy Diffusivity"    , "m2/s"   ,   &  ! avs 
    754719               &          jpi, jpj, nh_W, ipk, 1, ipk, nz_W, 32, clop, zsto, zout ) 
     
    861826#endif 
    862827 
    863       IF( ln_cpl .AND. nn_ice == 2 ) THEN 
    864          CALL histwrite( nid_T, "soicetem", it, tn_ice(:,:,1) , ndim_hT, ndex_hT )   ! surf. ice temperature 
    865          CALL histwrite( nid_T, "soicealb", it, alb_ice(:,:,1), ndim_hT, ndex_hT )   ! ice albedo 
    866       ENDIF 
    867  
    868828      CALL histwrite( nid_U, "vozocrtx", it, un            , ndim_U , ndex_U )    ! i-current 
    869829      CALL histwrite( nid_U, "sozotaux", it, utau          , ndim_hU, ndex_hU )   ! i-wind stress 
     
    874834      CALL histwrite( nid_W, "vovecrtz", it, wn             , ndim_T, ndex_T )    ! vert. current 
    875835      CALL histwrite( nid_W, "votkeavt", it, avt            , ndim_T, ndex_T )    ! T vert. eddy diff. coef. 
    876       CALL histwrite( nid_W, "votkeavm", it, avmu           , ndim_T, ndex_T )    ! T vert. eddy visc. coef. 
    877       IF( lk_zdfddm ) THEN 
    878          CALL histwrite( nid_W, "voddmavs", it, fsavs(:,:,:), ndim_T, ndex_T )    ! S vert. eddy diff. coef. 
     836      CALL histwrite( nid_W, "votkeavm", it, avm            , ndim_T, ndex_T )    ! T vert. eddy visc. coef. 
     837      IF( ln_zdfddm ) THEN 
     838         CALL histwrite( nid_W, "voddmavs", it, avs         , ndim_T, ndex_T )    ! S vert. eddy diff. coef. 
    879839      ENDIF 
    880840 
    881841      IF( ln_wave .AND. ln_sdw ) THEN 
    882          CALL histwrite( nid_U, "sdzocrtx", it, usd           , ndim_U , ndex_U )    ! i-StokesDrift-current 
    883          CALL histwrite( nid_V, "sdmecrty", it, vsd           , ndim_V , ndex_V )    ! j-StokesDrift-current 
    884          CALL histwrite( nid_W, "sdvecrtz", it, wsd           , ndim_T , ndex_T )    ! StokesDrift vert. current 
     842         CALL histwrite( nid_U, "sdzocrtx", it, usd         , ndim_U , ndex_U )    ! i-StokesDrift-current 
     843         CALL histwrite( nid_V, "sdmecrty", it, vsd         , ndim_V , ndex_V )    ! j-StokesDrift-current 
     844         CALL histwrite( nid_W, "sdvecrtz", it, wsd         , ndim_T , ndex_T )    ! StokesDrift vert. current 
    885845      ENDIF 
    886846 
     
    893853         CALL histclo( nid_W ) 
    894854      ENDIF 
    895       ! 
    896                              CALL wrk_dealloc( jpi , jpj        , zw2d ) 
    897       IF( .NOT.ln_linssh )   CALL wrk_dealloc( jpi , jpj , jpk  , zw3d ) 
    898855      ! 
    899856      IF( nn_timing == 1 )   CALL timing_stop('dia_wri') 
     
    1009966      ENDIF 
    1010967 
    1011 #if defined key_lim2 
    1012       CALL lim_wri_state_2( kt, id_i, nh_i ) 
    1013 #elif defined key_lim3 
    1014       CALL lim_wri_state( kt, id_i, nh_i ) 
     968#if defined key_lim3 
     969      IF( nn_ice == 2 ) THEN   ! condition needed in case agrif + ice-model but no-ice in child grid 
     970         CALL ice_wri_state( kt, id_i, nh_i ) 
     971      ENDIF 
    1015972#else 
    1016973      CALL histend( id_i, snc4chunks=snc4set ) 
Note: See TracChangeset for help on using the changeset viewer.