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

Ignore:
Timestamp:
2015-07-10T13:28:53+02:00 (9 years ago)
Author:
timgraham
Message:

Merged head of trunk into branch

File:
1 edited

Legend:

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

    r4761 r5581  
    4444   USE in_out_manager  ! I/O manager 
    4545   USE diadimg         ! dimg direct access file format output 
    46    USE diaar5, ONLY :   lk_diaar5 
    4746   USE iom 
    4847   USE ioipsl 
     48   USE dynspg_oce, ONLY: un_adv, vn_adv ! barotropic velocities      
     49 
    4950#if defined key_lim2 
    5051   USE limwri_2  
     
    7980   !!---------------------------------------------------------------------- 
    8081   !! NEMO/OPA 3.3 , NEMO Consortium (2010) 
    81    !! $Id $ 
     82   !! $Id$ 
    8283   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt) 
    8384   !!---------------------------------------------------------------------- 
     
    8889      INTEGER, DIMENSION(2) :: ierr 
    8990      !!---------------------------------------------------------------------- 
    90       ! 
    9191      ierr = 0 
    92       ! 
    9392      ALLOCATE( ndex_hT(jpi*jpj) , ndex_T(jpi*jpj*jpk) ,     & 
    9493         &      ndex_hU(jpi*jpj) , ndex_U(jpi*jpj*jpk) ,     & 
     
    128127      !! 
    129128      INTEGER                      ::   ji, jj, jk              ! dummy loop indices 
     129      INTEGER                      ::   jkbot                   ! 
    130130      REAL(wp)                     ::   zztmp, zztmpx, zztmpy   !  
    131131      !! 
    132132      REAL(wp), POINTER, DIMENSION(:,:)   :: z2d      ! 2D workspace 
    133       REAL(wp), POINTER, DIMENSION(:,:)   :: z2ds     ! 2D workspace 
    134133      REAL(wp), POINTER, DIMENSION(:,:,:) :: z3d      ! 3D workspace 
    135134      !!---------------------------------------------------------------------- 
     
    137136      IF( nn_timing == 1 )   CALL timing_start('dia_wri') 
    138137      !  
    139       CALL wrk_alloc( jpi , jpj      , z2d , z2ds ) 
     138      CALL wrk_alloc( jpi , jpj      , z2d ) 
    140139      CALL wrk_alloc( jpi , jpj, jpk , z3d ) 
    141140      ! 
     
    146145      ENDIF 
    147146 
    148       IF( lk_vvl ) THEN 
    149          z3d(:,:,:) = tsn(:,:,:,jp_tem) * fse3t_n(:,:,:) 
    150          CALL iom_put( "toce" , z3d                        )   ! heat content 
    151          CALL iom_put( "sst"  , z3d(:,:,1)                 )   ! sea surface heat content 
    152          z3d(:,:,1) = tsn(:,:,1,jp_tem) * z3d(:,:,1) 
    153          CALL iom_put( "sst2" , z3d(:,:,1)                 )   ! sea surface content of squared temperature 
    154          z3d(:,:,:) = tsn(:,:,:,jp_sal) * fse3t_n(:,:,:)             
    155          CALL iom_put( "soce" , z3d                        )   ! salinity content 
    156          CALL iom_put( "sss"  , z3d(:,:,1)                 )   ! sea surface salinity content 
    157          z3d(:,:,1) = tsn(:,:,1,jp_sal) * z3d(:,:,1) 
    158          CALL iom_put( "sss2" , z3d(:,:,1)                 )   ! sea surface content of squared salinity 
    159       ELSE 
    160          CALL iom_put( "toce" , tsn(:,:,:,jp_tem)          )   ! temperature 
    161          CALL iom_put( "sst"  , tsn(:,:,1,jp_tem)          )   ! sea surface temperature 
    162          CALL iom_put( "sst2" , tsn(:,:,1,jp_tem) * tsn(:,:,1,jp_tem) ) ! square of sea surface temperature 
    163          CALL iom_put( "soce" , tsn(:,:,:,jp_sal)          )   ! salinity 
    164          CALL iom_put( "sss"  , tsn(:,:,1,jp_sal)          )   ! sea surface salinity 
    165          CALL iom_put( "sss2" , tsn(:,:,1,jp_sal) * tsn(:,:,1,jp_sal) ) ! square of sea surface salinity 
    166       END IF 
    167       IF( lk_vvl .AND. (.NOT. ln_dynadv_vec) ) THEN 
    168          CALL iom_put( "uoce" , un(:,:,:) * fse3u_n(:,:,:) )    ! i-transport 
    169          CALL iom_put( "voce" , vn(:,:,:) * fse3v_n(:,:,:) )    ! j-transport 
    170       ELSE 
    171          CALL iom_put( "uoce" , un                         )    ! i-current 
    172          CALL iom_put( "voce" , vn                         )    ! j-current 
    173       END IF 
    174       CALL iom_put(    "avt"  , avt                        )    ! T vert. eddy diff. coef. 
    175       CALL iom_put(    "avm"  , avmu                       )    ! T vert. eddy visc. coef. 
    176       IF( lk_zdfddm ) THEN 
    177          CALL iom_put( "avs" , fsavs(:,:,:)                          )    ! S vert. eddy diff. coef. 
    178       ENDIF 
    179  
    180       DO jj = 2, jpjm1                                    ! sst gradient 
    181          DO ji = fs_2, fs_jpim1   ! vector opt. 
    182             zztmp      = tsn(ji,jj,1,jp_tem) 
    183             zztmpx     = ( tsn(ji+1,jj  ,1,jp_tem) - zztmp ) / e1u(ji,jj) + ( zztmp - tsn(ji-1,jj  ,1,jp_tem) ) / e1u(ji-1,jj  ) 
    184             zztmpy     = ( tsn(ji  ,jj+1,1,jp_tem) - zztmp ) / e2v(ji,jj) + ( zztmp - tsn(ji  ,jj-1,1,jp_tem) ) / e2v(ji  ,jj-1) 
    185             z2d(ji,jj) = 0.25 * ( zztmpx * zztmpx + zztmpy * zztmpy )   & 
    186                &              * umask(ji,jj,1) * umask(ji-1,jj,1) * vmask(ji,jj,1) * umask(ji,jj-1,1) 
    187          END DO 
    188       END DO 
    189       CALL lbc_lnk( z2d, 'T', 1. ) 
    190       CALL iom_put( "sstgrad2",  z2d               )    ! square of module of sst gradient 
    191 !CDIR NOVERRCHK 
    192       z2d(:,:) = SQRT( z2d(:,:) ) 
    193       CALL iom_put( "sstgrad" ,  z2d               )    ! module of sst gradient 
    194  
    195       ! clem: heat and salt content 
    196       z2d(:,:)  = 0._wp  
    197       z2ds(:,:) = 0._wp  
    198       DO jk = 1, jpkm1 
     147      IF( .NOT.lk_vvl ) THEN 
     148         CALL iom_put( "e3t" , fse3t_n(:,:,:) ) 
     149         CALL iom_put( "e3u" , fse3u_n(:,:,:) ) 
     150         CALL iom_put( "e3v" , fse3v_n(:,:,:) ) 
     151         CALL iom_put( "e3w" , fse3w_n(:,:,:) ) 
     152      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 
     156       
     157      CALL iom_put( "toce", tsn(:,:,:,jp_tem) )    ! 3D temperature 
     158      CALL iom_put(  "sst", tsn(:,:,1,jp_tem) )    ! surface temperature 
     159      IF ( iom_use("sbt") ) THEN 
     160         DO jj = 1, jpj 
     161            DO ji = 1, jpi 
     162               jkbot = mbkt(ji,jj) 
     163               z2d(ji,jj) = tsn(ji,jj,jkbot,jp_tem) 
     164            END DO 
     165         END DO 
     166         CALL iom_put( "sbt", z2d )                ! bottom temperature 
     167      ENDIF 
     168       
     169      CALL iom_put( "soce", tsn(:,:,:,jp_sal) )    ! 3D salinity 
     170      CALL iom_put(  "sss", tsn(:,:,1,jp_sal) )    ! surface salinity 
     171      IF ( iom_use("sbs") ) THEN 
     172         DO jj = 1, jpj 
     173            DO ji = 1, jpi 
     174               jkbot = mbkt(ji,jj) 
     175               z2d(ji,jj) = tsn(ji,jj,jkbot,jp_sal) 
     176            END DO 
     177         END DO 
     178         CALL iom_put( "sbs", z2d )                ! bottom salinity 
     179      ENDIF 
     180 
     181      IF ( iom_use("taubot") ) THEN                ! bottom stress 
     182         z2d(:,:) = 0._wp 
    199183         DO jj = 2, jpjm1 
    200184            DO ji = fs_2, fs_jpim1   ! vector opt. 
    201                z2d(ji,jj) = z2d(ji,jj) + rau0 * rcp * fse3t(ji,jj,jk) * tsn(ji,jj,jk,jp_tem) * tmask(ji,jj,jk) 
    202                z2ds(ji,jj) = z2ds(ji,jj) + rau0 * fse3t(ji,jj,jk) * tsn(ji,jj,jk,jp_sal) * tmask(ji,jj,jk) 
    203             END DO 
    204          END DO 
    205       END DO 
    206       CALL lbc_lnk( z2d, 'T', 1. ) 
    207       CALL lbc_lnk( z2ds, 'T', 1. ) 
    208       CALL iom_put( "heatc", z2d )    ! vertically integrated heat content (J/m2) 
    209       CALL iom_put( "saltc", z2ds )   ! vertically integrated salt content (PSU*kg/m2) 
     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 )            
     195      ENDIF 
     196          
     197      CALL iom_put( "uoce", un(:,:,:)         )    ! 3D i-current 
     198      CALL iom_put(  "ssu", un(:,:,1)         )    ! surface i-current 
     199      IF ( iom_use("sbu") ) THEN 
     200         DO jj = 1, jpj 
     201            DO ji = 1, jpi 
     202               jkbot = mbku(ji,jj) 
     203               z2d(ji,jj) = un(ji,jj,jkbot) 
     204            END DO 
     205         END DO 
     206         CALL iom_put( "sbu", z2d )                ! bottom i-current 
     207      ENDIF 
     208#if defined key_dynspg_ts 
     209      CALL iom_put(  "ubar", un_adv(:,:)      )    ! barotropic i-current 
     210#else 
     211      CALL iom_put(  "ubar", un_b(:,:)        )    ! barotropic i-current 
     212#endif 
    210213       
    211  
    212       IF( lk_diaar5 ) THEN 
     214      CALL iom_put( "voce", vn(:,:,:)         )    ! 3D j-current 
     215      CALL iom_put(  "ssv", vn(:,:,1)         )    ! surface j-current 
     216      IF ( iom_use("sbv") ) THEN 
     217         DO jj = 1, jpj 
     218            DO ji = 1, jpi 
     219               jkbot = mbkv(ji,jj) 
     220               z2d(ji,jj) = vn(ji,jj,jkbot) 
     221            END DO 
     222         END DO 
     223         CALL iom_put( "sbv", z2d )                ! bottom j-current 
     224      ENDIF 
     225#if defined key_dynspg_ts 
     226      CALL iom_put(  "vbar", vn_adv(:,:)      )    ! barotropic j-current 
     227#else 
     228      CALL iom_put(  "vbar", vn_b(:,:)        )    ! barotropic j-current 
     229#endif 
     230 
     231      CALL iom_put( "woce", wn )                   ! vertical velocity 
     232      IF( iom_use('w_masstr') .OR. iom_use('w_masstr2') ) THEN   ! vertical mass transport & its square value 
     233         ! Caution: in the VVL case, it only correponds to the baroclinic mass transport. 
     234         z2d(:,:) = rau0 * e12t(:,:) 
     235         DO jk = 1, jpk 
     236            z3d(:,:,jk) = wn(:,:,jk) * z2d(:,:) 
     237         END DO 
     238         CALL iom_put( "w_masstr" , z3d )   
     239         IF( iom_use('w_masstr2') )   CALL iom_put( "w_masstr2", z3d(:,:,:) * z3d(:,:,:) ) 
     240      ENDIF 
     241 
     242      CALL iom_put( "avt" , avt                        )    ! T vert. eddy diff. coef. 
     243      CALL iom_put( "avm" , avmu                       )    ! T vert. eddy visc. coef. 
     244      CALL iom_put( "avs" , fsavs(:,:,:)               )    ! S vert. eddy diff. coef. (useful only with key_zdfddm) 
     245 
     246      IF ( iom_use("sstgrad") .OR. iom_use("sstgrad2") ) THEN 
     247         DO jj = 2, jpjm1                                    ! sst gradient 
     248            DO ji = fs_2, fs_jpim1   ! vector opt. 
     249               zztmp      = tsn(ji,jj,1,jp_tem) 
     250               zztmpx     = ( tsn(ji+1,jj  ,1,jp_tem) - zztmp ) / e1u(ji,jj) + ( zztmp - tsn(ji-1,jj  ,1,jp_tem) ) / e1u(ji-1,jj  ) 
     251               zztmpy     = ( tsn(ji  ,jj+1,1,jp_tem) - zztmp ) / e2v(ji,jj) + ( zztmp - tsn(ji  ,jj-1,1,jp_tem) ) / e2v(ji  ,jj-1) 
     252               z2d(ji,jj) = 0.25 * ( zztmpx * zztmpx + zztmpy * zztmpy )   & 
     253                  &              * umask(ji,jj,1) * umask(ji-1,jj,1) * vmask(ji,jj,1) * umask(ji,jj-1,1) 
     254            END DO 
     255         END DO 
     256         CALL lbc_lnk( z2d, 'T', 1. ) 
     257         CALL iom_put( "sstgrad2",  z2d               )    ! square of module of sst gradient 
     258         z2d(:,:) = SQRT( z2d(:,:) ) 
     259         CALL iom_put( "sstgrad" ,  z2d               )    ! module of sst gradient 
     260      ENDIF 
     261          
     262      ! clem: heat and salt content 
     263      IF( iom_use("heatc") ) THEN 
     264         z2d(:,:)  = 0._wp  
     265         DO jk = 1, jpkm1 
     266            DO jj = 1, jpj 
     267               DO ji = 1, jpi 
     268                  z2d(ji,jj) = z2d(ji,jj) + fse3t(ji,jj,jk) * tsn(ji,jj,jk,jp_tem) * tmask(ji,jj,jk) 
     269               END DO 
     270            END DO 
     271         END DO 
     272         CALL iom_put( "heatc", (rau0 * rcp) * z2d )    ! vertically integrated heat content (J/m2) 
     273      ENDIF 
     274 
     275      IF( iom_use("saltc") ) THEN 
     276         z2d(:,:)  = 0._wp  
     277         DO jk = 1, jpkm1 
     278            DO jj = 1, jpj 
     279               DO ji = 1, jpi 
     280                  z2d(ji,jj) = z2d(ji,jj) + fse3t(ji,jj,jk) * tsn(ji,jj,jk,jp_sal) * tmask(ji,jj,jk) 
     281               END DO 
     282            END DO 
     283         END DO 
     284         CALL iom_put( "saltc", rau0 * z2d )   ! vertically integrated salt content (PSU*kg/m2) 
     285      ENDIF 
     286      ! 
     287      IF ( iom_use("eken") ) THEN 
     288         rke(:,:,jk) = 0._wp                               !      kinetic energy  
     289         DO jk = 1, jpkm1 
     290            DO jj = 2, jpjm1 
     291               DO ji = fs_2, fs_jpim1   ! vector opt. 
     292                  zztmp   = 1._wp / ( e1e2t(ji,jj) * fse3t(ji,jj,jk) ) 
     293                  zztmpx  = 0.5 * (  un(ji-1,jj,jk) * un(ji-1,jj,jk) * e2u(ji-1,jj) * fse3u(ji-1,jj,jk)    & 
     294                     &             + un(ji  ,jj,jk) * un(ji  ,jj,jk) * e2u(ji  ,jj) * fse3u(ji  ,jj,jk) )  & 
     295                     &          *  zztmp  
     296                  ! 
     297                  zztmpy  = 0.5 * (  vn(ji,jj-1,jk) * vn(ji,jj-1,jk) * e1v(ji,jj-1) * fse3v(ji,jj-1,jk)    & 
     298                     &             + vn(ji,jj  ,jk) * vn(ji,jj  ,jk) * e1v(ji,jj  ) * fse3v(ji,jj  ,jk) )  & 
     299                     &          *  zztmp  
     300                  ! 
     301                  rke(ji,jj,jk) = 0.5_wp * ( zztmpx + zztmpy ) 
     302                  ! 
     303               ENDDO 
     304            ENDDO 
     305         ENDDO 
     306         CALL lbc_lnk( rke, 'T', 1. ) 
     307         CALL iom_put( "eken", rke )            
     308      ENDIF 
     309          
     310      IF( iom_use("u_masstr") .OR. iom_use("u_heattr") .OR. iom_use("u_salttr") ) THEN 
    213311         z3d(:,:,jpk) = 0.e0 
    214312         DO jk = 1, jpkm1 
     
    216314         END DO 
    217315         CALL iom_put( "u_masstr", z3d )                  ! mass transport in i-direction 
    218  
    219          zztmp = 0.5 * rcp 
     316      ENDIF 
     317       
     318      IF( iom_use("u_heattr") ) THEN 
    220319         z2d(:,:) = 0.e0  
    221          z2ds(:,:) = 0.e0  
    222320         DO jk = 1, jpkm1 
    223321            DO jj = 2, jpjm1 
    224322               DO ji = fs_2, fs_jpim1   ! vector opt. 
    225                   z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * zztmp * ( tsn(ji,jj,jk,jp_tem) + tsn(ji+1,jj,jk,jp_tem) ) 
    226                   z2ds(ji,jj) = z2ds(ji,jj) + z3d(ji,jj,jk) * 0.5_wp * ( tsn(ji,jj,jk,jp_sal) + tsn(ji+1,jj,jk,jp_sal) ) 
     323                  z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( tsn(ji,jj,jk,jp_tem) + tsn(ji+1,jj,jk,jp_tem) ) 
    227324               END DO 
    228325            END DO 
    229326         END DO 
    230327         CALL lbc_lnk( z2d, 'U', -1. ) 
    231          CALL lbc_lnk( z2ds, 'U', -1. ) 
    232          CALL iom_put( "u_heattr", z2d )                  ! heat transport in i-direction 
    233          CALL iom_put( "u_salttr", z2ds )                 ! salt transport in i-direction 
    234  
     328         CALL iom_put( "u_heattr", (0.5 * rcp) * z2d )    ! heat transport in i-direction 
     329      ENDIF 
     330 
     331      IF( iom_use("u_salttr") ) THEN 
     332         z2d(:,:) = 0.e0  
     333         DO jk = 1, jpkm1 
     334            DO jj = 2, jpjm1 
     335               DO ji = fs_2, fs_jpim1   ! vector opt. 
     336                  z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( tsn(ji,jj,jk,jp_sal) + tsn(ji+1,jj,jk,jp_sal) ) 
     337               END DO 
     338            END DO 
     339         END DO 
     340         CALL lbc_lnk( z2d, 'U', -1. ) 
     341         CALL iom_put( "u_salttr", 0.5 * z2d )            ! heat transport in i-direction 
     342      ENDIF 
     343 
     344       
     345      IF( iom_use("v_masstr") .OR. iom_use("v_heattr") .OR. iom_use("v_salttr") ) THEN 
    235346         z3d(:,:,jpk) = 0.e0 
    236347         DO jk = 1, jpkm1 
     
    238349         END DO 
    239350         CALL iom_put( "v_masstr", z3d )                  ! mass transport in j-direction 
    240  
     351      ENDIF 
     352       
     353      IF( iom_use("v_heattr") ) THEN 
    241354         z2d(:,:) = 0.e0  
    242          z2ds(:,:) = 0.e0  
    243355         DO jk = 1, jpkm1 
    244356            DO jj = 2, jpjm1 
    245357               DO ji = fs_2, fs_jpim1   ! vector opt. 
    246                   z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * zztmp * ( tsn(ji,jj,jk,jp_tem) + tsn(ji,jj+1,jk,jp_tem) ) 
    247                   z2ds(ji,jj) = z2ds(ji,jj) + z3d(ji,jj,jk) * 0.5_wp * ( tsn(ji,jj,jk,jp_sal) + tsn(ji,jj+1,jk,jp_sal) ) 
     358                  z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( tsn(ji,jj,jk,jp_tem) + tsn(ji,jj+1,jk,jp_tem) ) 
    248359               END DO 
    249360            END DO 
    250361         END DO 
    251362         CALL lbc_lnk( z2d, 'V', -1. ) 
    252          CALL lbc_lnk( z2ds, 'V', -1. ) 
    253          CALL iom_put( "v_heattr", z2d )                  !  heat transport in j-direction 
    254          CALL iom_put( "v_salttr", z2ds )                 !  salt transport in j-direction 
    255       ENDIF 
    256       ! 
    257       CALL wrk_dealloc( jpi , jpj      , z2d , z2ds ) 
     363         CALL iom_put( "v_heattr", (0.5 * rcp) * z2d )    !  heat transport in j-direction 
     364      ENDIF 
     365 
     366      IF( iom_use("v_salttr") ) THEN 
     367         z2d(:,:) = 0.e0  
     368         DO jk = 1, jpkm1 
     369            DO jj = 2, jpjm1 
     370               DO ji = fs_2, fs_jpim1   ! vector opt. 
     371                  z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( tsn(ji,jj,jk,jp_sal) + tsn(ji,jj+1,jk,jp_sal) ) 
     372               END DO 
     373            END DO 
     374         END DO 
     375         CALL lbc_lnk( z2d, 'V', -1. ) 
     376         CALL iom_put( "v_salttr", 0.5 * z2d )            !  heat transport in j-direction 
     377      ENDIF 
     378      ! 
     379      CALL wrk_dealloc( jpi , jpj      , z2d ) 
    258380      CALL wrk_dealloc( jpi , jpj, jpk , z3d ) 
    259381      ! 
     
    316438      zdt = rdt 
    317439      IF( nacc == 1 ) zdt = rdtmin 
    318       IF( ln_mskland )   THEN   ;   clop = "only(x)"   ! put 1.e+20 on land (very expensive!!) 
    319       ELSE                      ;   clop = "x"         ! no use of the mask value (require less cpu time) 
    320       ENDIF 
     440      clop = "x"         ! no use of the mask value (require less cpu time and otherwise the model crashes) 
    321441#if defined key_diainstant 
    322442      zsto = nwrite * zdt 
     
    518638         ENDIF 
    519639 
    520 #if ! defined key_coupled  
    521          CALL histdef( nid_T, "sohefldp", "Surface Heat Flux: Damping"         , "W/m2"   ,   &  ! qrp 
    522             &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
    523          CALL histdef( nid_T, "sowafldp", "Surface Water Flux: Damping"        , "Kg/m2/s",   &  ! erp 
    524             &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
    525          CALL histdef( nid_T, "sosafldp", "Surface salt flux: damping"         , "Kg/m2/s",   &  ! erp * sn 
    526             &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
    527 #endif 
    528  
    529  
    530  
    531 #if ( defined key_coupled && ! defined key_lim3 && ! defined key_lim2 )  
    532          CALL histdef( nid_T, "sohefldp", "Surface Heat Flux: Damping"         , "W/m2"   ,   &  ! qrp 
    533             &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
    534          CALL histdef( nid_T, "sowafldp", "Surface Water Flux: Damping"        , "Kg/m2/s",   &  ! erp 
    535             &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
    536          CALL histdef( nid_T, "sosafldp", "Surface salt flux: Damping"         , "Kg/m2/s",   &  ! erp * sn 
    537             &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
    538 #endif 
     640         IF( .NOT. ln_cpl ) THEN 
     641            CALL histdef( nid_T, "sohefldp", "Surface Heat Flux: Damping"         , "W/m2"   ,   &  ! qrp 
     642               &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
     643            CALL histdef( nid_T, "sowafldp", "Surface Water Flux: Damping"        , "Kg/m2/s",   &  ! erp 
     644               &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
     645            CALL histdef( nid_T, "sosafldp", "Surface salt flux: damping"         , "Kg/m2/s",   &  ! erp * sn 
     646               &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
     647         ENDIF 
     648 
     649         IF( ln_cpl .AND. nn_ice <= 1 ) THEN 
     650            CALL histdef( nid_T, "sohefldp", "Surface Heat Flux: Damping"         , "W/m2"   ,   &  ! qrp 
     651               &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
     652            CALL histdef( nid_T, "sowafldp", "Surface Water Flux: Damping"        , "Kg/m2/s",   &  ! erp 
     653               &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
     654            CALL histdef( nid_T, "sosafldp", "Surface salt flux: Damping"         , "Kg/m2/s",   &  ! erp * sn 
     655               &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
     656         ENDIF 
     657          
    539658         clmx ="l_max(only(x))"    ! max index on a period 
    540659         CALL histdef( nid_T, "sobowlin", "Bowl Index"                         , "W-point",   &  ! bowl INDEX  
     
    551670#endif 
    552671 
    553 #if defined key_coupled  
    554 # if defined key_lim3 
    555          Must be adapted to LIM3 
    556 # endif  
    557 # if defined key_lim2 
    558          CALL histdef( nid_T,"soicetem" , "Ice Surface Temperature"            , "K"      ,   &  ! tn_ice 
    559             &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
    560          CALL histdef( nid_T,"soicealb" , "Ice Albedo"                         , "[0,1]"  ,   &  ! alb_ice 
    561             &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
    562 # endif  
    563 #endif  
     672         IF( ln_cpl .AND. nn_ice == 2 ) THEN 
     673            CALL histdef( nid_T,"soicetem" , "Ice Surface Temperature"            , "K"      ,   &  ! tn_ice 
     674               &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
     675            CALL histdef( nid_T,"soicealb" , "Ice Albedo"                         , "[0,1]"  ,   &  ! alb_ice 
     676               &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
     677         ENDIF 
    564678 
    565679         CALL histend( nid_T, snc4chunks=snc4set ) 
     
    652766      ENDIF 
    653767 
    654       ! Write fields on T grid 
    655768      IF( lk_vvl ) THEN 
    656769         CALL histwrite( nid_T, "votemper", it, tsn(:,:,:,jp_tem) * fse3t_n(:,:,:) , ndim_T , ndex_T  )   ! heat content 
     
    663776         CALL histwrite( nid_T, "sosstsst", it, tsn(:,:,1,jp_tem) , ndim_hT, ndex_hT )   ! sea surface temperature 
    664777         CALL histwrite( nid_T, "sosaline", it, tsn(:,:,1,jp_sal) , ndim_hT, ndex_hT )   ! sea surface salinity 
    665  
    666778      ENDIF 
    667779      IF( lk_vvl ) THEN 
     
    713825      ENDIF 
    714826 
    715 #if ! defined key_coupled 
    716       CALL histwrite( nid_T, "sohefldp", it, qrp           , ndim_hT, ndex_hT )   ! heat flux damping 
    717       CALL histwrite( nid_T, "sowafldp", it, erp           , ndim_hT, ndex_hT )   ! freshwater flux damping 
    718       IF( ln_ssr ) zw2d(:,:) = erp(:,:) * tsn(:,:,1,jp_sal) * tmask(:,:,1) 
    719       CALL histwrite( nid_T, "sosafldp", it, zw2d          , ndim_hT, ndex_hT )   ! salt flux damping 
    720 #endif 
    721 #if ( defined key_coupled && ! defined key_lim3 && ! defined key_lim2 )  
    722       CALL histwrite( nid_T, "sohefldp", it, qrp           , ndim_hT, ndex_hT )   ! heat flux damping 
    723       CALL histwrite( nid_T, "sowafldp", it, erp           , ndim_hT, ndex_hT )   ! freshwater flux damping 
     827      IF( .NOT. ln_cpl ) THEN 
     828         CALL histwrite( nid_T, "sohefldp", it, qrp           , ndim_hT, ndex_hT )   ! heat flux damping 
     829         CALL histwrite( nid_T, "sowafldp", it, erp           , ndim_hT, ndex_hT )   ! freshwater flux damping 
    724830         IF( ln_ssr ) zw2d(:,:) = erp(:,:) * tsn(:,:,1,jp_sal) * tmask(:,:,1) 
    725       CALL histwrite( nid_T, "sosafldp", it, zw2d          , ndim_hT, ndex_hT )   ! salt flux damping 
    726 #endif 
    727       zw2d(:,:) = FLOAT( nmln(:,:) ) * tmask(:,:,1) 
    728       CALL histwrite( nid_T, "sobowlin", it, zw2d          , ndim_hT, ndex_hT )   ! ??? 
     831         CALL histwrite( nid_T, "sosafldp", it, zw2d          , ndim_hT, ndex_hT )   ! salt flux damping 
     832      ENDIF 
     833      IF( ln_cpl .AND. nn_ice <= 1 ) THEN 
     834         CALL histwrite( nid_T, "sohefldp", it, qrp           , ndim_hT, ndex_hT )   ! heat flux damping 
     835         CALL histwrite( nid_T, "sowafldp", it, erp           , ndim_hT, ndex_hT )   ! freshwater flux damping 
     836         IF( ln_ssr ) zw2d(:,:) = erp(:,:) * tsn(:,:,1,jp_sal) * tmask(:,:,1) 
     837         CALL histwrite( nid_T, "sosafldp", it, zw2d          , ndim_hT, ndex_hT )   ! salt flux damping 
     838      ENDIF 
     839!      zw2d(:,:) = FLOAT( nmln(:,:) ) * tmask(:,:,1) 
     840!      CALL histwrite( nid_T, "sobowlin", it, zw2d          , ndim_hT, ndex_hT )   ! ??? 
    729841 
    730842#if defined key_diahth 
     
    735847#endif 
    736848 
    737 #if defined key_coupled  
    738 # if defined key_lim3 
    739       Must be adapted for LIM3 
    740       CALL histwrite( nid_T, "soicetem", it, tn_ice        , ndim_hT, ndex_hT )   ! surf. ice temperature 
    741       CALL histwrite( nid_T, "soicealb", it, alb_ice       , ndim_hT, ndex_hT )   ! ice albedo 
    742 # endif 
    743 # if defined key_lim2 
    744       CALL histwrite( nid_T, "soicetem", it, tn_ice(:,:,1) , ndim_hT, ndex_hT )   ! surf. ice temperature 
    745       CALL histwrite( nid_T, "soicealb", it, alb_ice(:,:,1), ndim_hT, ndex_hT )   ! ice albedo 
    746 # endif 
    747 #endif 
    748          ! Write fields on U grid 
     849      IF( ln_cpl .AND. nn_ice == 2 ) THEN 
     850         CALL histwrite( nid_T, "soicetem", it, tn_ice(:,:,1) , ndim_hT, ndex_hT )   ! surf. ice temperature 
     851         CALL histwrite( nid_T, "soicealb", it, alb_ice(:,:,1), ndim_hT, ndex_hT )   ! ice albedo 
     852      ENDIF 
     853 
    749854      CALL histwrite( nid_U, "vozocrtx", it, un            , ndim_U , ndex_U )    ! i-current 
    750855      IF( ln_traldf_gdia ) THEN 
     
    768873      CALL histwrite( nid_U, "sozotaux", it, utau          , ndim_hU, ndex_hU )   ! i-wind stress 
    769874 
    770          ! Write fields on V grid 
    771875      CALL histwrite( nid_V, "vomecrty", it, vn            , ndim_V , ndex_V  )   ! j-current 
    772876      IF( ln_traldf_gdia ) THEN 
     
    783887      CALL histwrite( nid_V, "sometauy", it, vtau          , ndim_hV, ndex_hV )   ! j-wind stress 
    784888 
    785          ! Write fields on W grid 
    786889      CALL histwrite( nid_W, "vovecrtz", it, wn             , ndim_T, ndex_T )    ! vert. current 
    787890      IF( ln_traldf_gdia ) THEN 
Note: See TracChangeset for help on using the changeset viewer.