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/LIM_SRC_3/limwri.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/LIM_SRC_3/limwri.F90

    r4765 r5581  
    2424   USE lib_mpp         ! MPP library 
    2525   USE wrk_nemo        ! work arrays 
    26    USE par_ice 
    2726   USE iom 
    2827   USE timing          ! Timing 
     
    3534   PUBLIC lim_wri_state  ! called by dia_wri_state  
    3635 
    37    REAL(wp)  ::   epsi06 = 1.e-6_wp 
    3836   !!---------------------------------------------------------------------- 
    3937   !! NEMO/LIM3 4.0 , UCL - NEMO Consortium (2011) 
     
    5957      INTEGER, INTENT(in) ::   kindic   ! if kindic < 0 there has been an error somewhere 
    6058      ! 
    61       INTEGER ::  ji, jj, jk, jl  ! dummy loop indices 
    62       REAL(wp) ::  zinda, zindb, z1_365 
    63       REAL(wp), POINTER, DIMENSION(:,:,:) ::   zoi, zei 
    64       REAL(wp), POINTER, DIMENSION(:,:)   :: z2d, z2da, z2db, zind    ! 2D workspace 
     59      INTEGER  ::  ji, jj, jk, jl  ! dummy loop indices 
     60      REAL(wp) ::  z1_365 
     61      REAL(wp) ::  ztmp 
     62      REAL(wp), POINTER, DIMENSION(:,:,:) ::  zoi, zei, zt_i, zt_s 
     63      REAL(wp), POINTER, DIMENSION(:,:)   ::  z2d, z2da, z2db, zswi    ! 2D workspace 
    6564      !!------------------------------------------------------------------- 
    6665 
    6766      IF( nn_timing == 1 )  CALL timing_start('limwri') 
    6867 
    69       CALL wrk_alloc( jpi, jpj, jpl, zoi, zei ) 
    70       CALL wrk_alloc( jpi, jpj     , z2d, z2da, z2db, zind ) 
     68      CALL wrk_alloc( jpi, jpj, jpl, zoi, zei, zt_i, zt_s ) 
     69      CALL wrk_alloc( jpi, jpj     , z2d, z2da, z2db, zswi ) 
    7170 
    7271      !----------------------------- 
    7372      ! Mean category values 
    7473      !----------------------------- 
     74      z1_365 = 1._wp / 365._wp 
    7575 
    7676      CALL lim_var_icetm      ! mean sea ice temperature 
     
    8080      DO jj = 1, jpj          ! presence indicator of ice 
    8181         DO ji = 1, jpi 
    82             zind(ji,jj)  = MAX( 0._wp , SIGN( 1._wp , at_i(ji,jj) - epsi06 ) ) 
     82            zswi(ji,jj)  = MAX( 0._wp , SIGN( 1._wp , at_i(ji,jj) - epsi06 ) ) 
    8383         END DO 
    8484      END DO 
     
    8989         DO jj = 1, jpj  
    9090            DO ji = 1, jpi 
    91                z2d(ji,jj)  = vt_i(ji,jj) / MAX( at_i(ji,jj), epsi06 ) * zind(ji,jj) 
     91               z2d(ji,jj)  = vt_i(ji,jj) / MAX( at_i(ji,jj), epsi06 ) * zswi(ji,jj) 
    9292            END DO 
    9393         END DO 
     
    9898         DO jj = 1, jpj                                             
    9999            DO ji = 1, jpi 
    100                z2d(ji,jj)  = vt_s(ji,jj) / MAX( at_i(ji,jj), epsi06 ) * zind(ji,jj) 
     100               z2d(ji,jj)  = vt_s(ji,jj) / MAX( at_i(ji,jj), epsi06 ) * zswi(ji,jj) 
    101101            END DO 
    102102         END DO 
     
    107107         DO jj = 2 , jpjm1 
    108108            DO ji = 2 , jpim1 
    109                z2da(ji,jj)  = (  u_ice(ji,jj) * tmu(ji,jj) + u_ice(ji-1,jj) * tmu(ji-1,jj) ) * 0.5_wp 
    110                z2db(ji,jj)  = (  v_ice(ji,jj) * tmv(ji,jj) + v_ice(ji,jj-1) * tmv(ji,jj-1) ) * 0.5_wp 
     109               z2da(ji,jj)  = (  u_ice(ji,jj) * umask(ji,jj,1) + u_ice(ji-1,jj) * umask(ji-1,jj,1) ) * 0.5_wp 
     110               z2db(ji,jj)  = (  v_ice(ji,jj) * vmask(ji,jj,1) + v_ice(ji,jj-1) * vmask(ji,jj-1,1) ) * 0.5_wp 
    111111           END DO 
    112112         END DO 
    113113         CALL lbc_lnk( z2da, 'T', -1. ) 
    114114         CALL lbc_lnk( z2db, 'T', -1. ) 
    115          CALL iom_put( "uice_ipa"     , z2da                )       ! ice velocity u component 
    116          CALL iom_put( "vice_ipa"     , z2db                )       ! ice velocity v component 
     115         CALL iom_put( "uice_ipa"     , z2da             )       ! ice velocity u component 
     116         CALL iom_put( "vice_ipa"     , z2db             )       ! ice velocity v component 
    117117         DO jj = 1, jpj                                  
    118118            DO ji = 1, jpi 
     
    120120            END DO 
    121121         END DO 
    122          CALL iom_put( "icevel"       , z2d                 )       ! ice velocity module 
     122         CALL iom_put( "icevel"       , z2d              )       ! ice velocity module 
    123123      ENDIF 
    124124      ! 
     
    128128            DO jj = 1, jpj 
    129129               DO ji = 1, jpi 
    130                   z2d(ji,jj) = z2d(ji,jj) + zind(ji,jj) * oa_i(ji,jj,jl) 
     130                  rswitch    = MAX( 0._wp , SIGN( 1._wp , at_i(ji,jj) - 0.1 ) ) 
     131                  z2d(ji,jj) = z2d(ji,jj) + rswitch * oa_i(ji,jj,jl) / MAX( at_i(ji,jj), 0.1 ) 
    131132               END DO 
    132133            END DO 
    133134         END DO 
    134          z1_365 = 1._wp / 365._wp 
    135          CALL iom_put( "miceage"     , z2d * z1_365         )        ! mean ice age 
     135         CALL iom_put( "miceage"     , z2d * z1_365      )        ! mean ice age 
    136136      ENDIF 
    137137 
     
    139139         DO jj = 1, jpj 
    140140            DO ji = 1, jpi 
    141                z2d(ji,jj) = ( tm_i(ji,jj) - rtt ) * zind(ji,jj) 
    142             END DO 
    143          END DO 
    144          CALL iom_put( "micet"       , z2d                  )        ! mean ice temperature 
     141               z2d(ji,jj) = ( tm_i(ji,jj) - rt0 ) * zswi(ji,jj) 
     142            END DO 
     143         END DO 
     144         CALL iom_put( "micet"       , z2d               )        ! mean ice temperature 
    145145      ENDIF 
    146146      ! 
     
    150150            DO jj = 1, jpj 
    151151               DO ji = 1, jpi 
    152                   z2d(ji,jj) = z2d(ji,jj) + zind(ji,jj) * ( t_su(ji,jj,jl) - rtt ) * a_i(ji,jj,jl) / MAX( at_i(ji,jj) , epsi06 ) 
     152                  z2d(ji,jj) = z2d(ji,jj) + zswi(ji,jj) * ( t_su(ji,jj,jl) - rt0 ) * a_i(ji,jj,jl) / MAX( at_i(ji,jj) , epsi06 ) 
    153153               END DO 
    154154            END DO 
    155155         END DO 
    156          CALL iom_put( "icest"       , z2d                 )        ! ice surface temperature 
     156         CALL iom_put( "icest"       , z2d              )        ! ice surface temperature 
    157157      ENDIF 
    158158 
     
    160160         DO jj = 1, jpj 
    161161            DO ji = 1, jpi 
    162                zindb  = MAX( 0._wp , SIGN( 1._wp , at_i(ji,jj) ) ) 
    163                z2d(ji,jj) = hicol(ji,jj) * zindb 
    164             END DO 
    165          END DO 
    166          CALL iom_put( "icecolf"     , z2d                 )        ! frazil ice collection thickness 
     162               rswitch  = MAX( 0._wp , SIGN( 1._wp , at_i(ji,jj) ) ) 
     163               z2d(ji,jj) = hicol(ji,jj) * rswitch 
     164            END DO 
     165         END DO 
     166         CALL iom_put( "icecolf"     , z2d              )        ! frazil ice collection thickness 
    167167      ENDIF 
    168168 
     
    176176      CALL iom_put( "utau_ice"    , utau_ice            )        ! wind stress over ice along i-axis at I-point 
    177177      CALL iom_put( "vtau_ice"    , vtau_ice            )        ! wind stress over ice along j-axis at I-point 
    178       CALL iom_put( "snowpre"     , sprecip             )        ! snow precipitation  
     178      CALL iom_put( "snowpre"     , sprecip * 86400.    )        ! snow precipitation  
    179179      CALL iom_put( "micesalt"    , smt_i               )        ! mean ice salinity 
    180180 
     
    186186      CALL iom_put( "icetrp"      , diag_trp_vi * rday  )        ! ice volume transport 
    187187      CALL iom_put( "snwtrp"      , diag_trp_vs * rday  )        ! snw volume transport 
     188      CALL iom_put( "saltrp"      , diag_trp_smv * rday * rhoic ) ! salt content transport 
    188189      CALL iom_put( "deitrp"      , diag_trp_ei         )        ! advected ice enthalpy (W/m2) 
    189190      CALL iom_put( "destrp"      , diag_trp_es         )        ! advected snw enthalpy (W/m2) 
     
    199200      CALL iom_put( "sfx"         , sfx     * rday      )        ! total salt flux 
    200201 
    201       CALL iom_put( "vfxres"     , wfx_res * rday / rhoic  )        ! daily prod./melting due to limupdate  
    202       CALL iom_put( "vfxopw"     , wfx_opw * rday / rhoic  )        ! daily lateral thermodynamic ice production 
    203       CALL iom_put( "vfxsni"     , wfx_sni * rday / rhoic  )        ! daily snowice ice production 
    204       CALL iom_put( "vfxbog"     , wfx_bog * rday / rhoic  )       ! daily bottom thermodynamic ice production 
    205       CALL iom_put( "vfxdyn"     , wfx_dyn * rday / rhoic  )       ! daily dynamic ice production (rid/raft) 
    206       CALL iom_put( "vfxsum"     , wfx_sum * rday / rhoic  )        ! surface melt  
    207       CALL iom_put( "vfxbom"     , wfx_bom * rday / rhoic  )        ! bottom melt  
    208       CALL iom_put( "vfxice"     , wfx_ice * rday / rhoic  )        ! total ice growth/melt  
    209       CALL iom_put( "vfxsnw"     , wfx_snw * rday / rhoic  )        ! total snw growth/melt  
    210       CALL iom_put( "vfxsub"     , wfx_sub * rday / rhoic  )        ! sublimation (snow)  
    211       CALL iom_put( "vfxspr"     , wfx_spr * rday / rhoic  )        ! precip (snow)  
    212  
    213       CALL iom_put ('hfxthd', hfx_thd(:,:) )   !   
    214       CALL iom_put ('hfxdyn', hfx_dyn(:,:) )   !   
    215       CALL iom_put ('hfxres', hfx_res(:,:) )   !   
    216       CALL iom_put ('hfxout', hfx_out(:,:) )   !   
    217       CALL iom_put ('hfxin' , hfx_in(:,:) )   !   
    218       CALL iom_put ('hfxsnw', hfx_snw(:,:) )   !   
    219       CALL iom_put ('hfxsub', hfx_sub(:,:) )   !   
    220       CALL iom_put ('hfxerr', hfx_err(:,:) )   !   
    221       CALL iom_put ('hfxerr_rem', hfx_err_rem(:,:) )   !   
    222        
    223       CALL iom_put ('hfxsum', hfx_sum(:,:) )   !   
    224       CALL iom_put ('hfxbom', hfx_bom(:,:) )   !   
    225       CALL iom_put ('hfxbog', hfx_bog(:,:) )   !   
    226       CALL iom_put ('hfxdif', hfx_dif(:,:) )   !   
    227       CALL iom_put ('hfxopw', hfx_opw(:,:) )   !   
    228       CALL iom_put ('hfxtur', fhtur(:,:) * at_i(:,:) )   ! turbulent heat flux at ice base  
    229       CALL iom_put ('hfxdhc', diag_heat_dhc(:,:) )          ! Heat content variation in snow and ice  
    230       CALL iom_put ('hfxspr', hfx_spr(:,:) )          ! Heat content of snow precip  
     202      ztmp = rday / rhoic 
     203      CALL iom_put( "vfxres"     , wfx_res * ztmp       )        ! daily prod./melting due to limupdate  
     204      CALL iom_put( "vfxopw"     , wfx_opw * ztmp       )        ! daily lateral thermodynamic ice production 
     205      CALL iom_put( "vfxsni"     , wfx_sni * ztmp       )        ! daily snowice ice production 
     206      CALL iom_put( "vfxbog"     , wfx_bog * ztmp       )        ! daily bottom thermodynamic ice production 
     207      CALL iom_put( "vfxdyn"     , wfx_dyn * ztmp       )        ! daily dynamic ice production (rid/raft) 
     208      CALL iom_put( "vfxsum"     , wfx_sum * ztmp       )        ! surface melt  
     209      CALL iom_put( "vfxbom"     , wfx_bom * ztmp       )        ! bottom melt  
     210      CALL iom_put( "vfxice"     , wfx_ice * ztmp       )        ! total ice growth/melt  
     211      CALL iom_put( "vfxsnw"     , wfx_snw * ztmp       )        ! total snw growth/melt  
     212      CALL iom_put( "vfxsub"     , wfx_sub * ztmp       )        ! sublimation (snow)  
     213      CALL iom_put( "vfxspr"     , wfx_spr * ztmp       )        ! precip (snow) 
     214       
     215      CALL iom_put( "afxtot"     , afx_tot * rday       )        ! concentration tendency (total) 
     216      CALL iom_put( "afxdyn"     , afx_dyn * rday       )        ! concentration tendency (dynamics) 
     217      CALL iom_put( "afxthd"     , afx_thd * rday       )        ! concentration tendency (thermo) 
     218 
     219      CALL iom_put ('hfxthd'     , hfx_thd(:,:)         )   !   
     220      CALL iom_put ('hfxdyn'     , hfx_dyn(:,:)         )   !   
     221      CALL iom_put ('hfxres'     , hfx_res(:,:)         )   !   
     222      CALL iom_put ('hfxout'     , hfx_out(:,:)         )   !   
     223      CALL iom_put ('hfxin'      , hfx_in(:,:)          )   !   
     224      CALL iom_put ('hfxsnw'     , hfx_snw(:,:)         )   !   
     225      CALL iom_put ('hfxsub'     , hfx_sub(:,:)         )   !   
     226      CALL iom_put ('hfxerr'     , hfx_err(:,:)         )   !   
     227      CALL iom_put ('hfxerr_rem' , hfx_err_rem(:,:)     )   !   
     228       
     229      CALL iom_put ('hfxsum'     , hfx_sum(:,:)         )   !   
     230      CALL iom_put ('hfxbom'     , hfx_bom(:,:)         )   !   
     231      CALL iom_put ('hfxbog'     , hfx_bog(:,:)         )   !   
     232      CALL iom_put ('hfxdif'     , hfx_dif(:,:)         )   !   
     233      CALL iom_put ('hfxopw'     , hfx_opw(:,:)         )   !   
     234      CALL iom_put ('hfxtur'     , fhtur(:,:) * SUM(a_i_b(:,:,:), dim=3) ) ! turbulent heat flux at ice base  
     235      CALL iom_put ('hfxdhc'     , diag_heat(:,:)       )   ! Heat content variation in snow and ice  
     236      CALL iom_put ('hfxspr'     , hfx_spr(:,:)         )   ! Heat content of snow precip  
    231237       
    232238      !-------------------------------- 
     
    238244      CALL iom_put( "salinity_cat"     , sm_i        )        ! salinity for categories 
    239245 
     246      ! ice temperature 
     247      IF ( iom_use( "icetemp_cat" ) ) THEN  
     248         zt_i(:,:,:) = SUM( t_i(:,:,:,:), dim=3 ) * r1_nlay_i 
     249         CALL iom_put( "icetemp_cat"   , zt_i - rt0  ) 
     250      ENDIF 
     251       
     252      ! snow temperature 
     253      IF ( iom_use( "snwtemp_cat" ) ) THEN  
     254         zt_s(:,:,:) = SUM( t_s(:,:,:,:), dim=3 ) * r1_nlay_s 
     255         CALL iom_put( "snwtemp_cat"   , zt_s - rt0  ) 
     256      ENDIF 
     257 
    240258      ! Compute ice age 
    241259      IF ( iom_use( "iceage_cat" ) ) THEN  
     
    243261            DO jj = 1, jpj 
    244262               DO ji = 1, jpi 
    245                   zinda = MAX( 0._wp , SIGN( 1._wp , a_i(ji,jj,jl) - epsi06 ) ) 
    246                   zoi(ji,jj,jl) = oa_i(ji,jj,jl)  / MAX( a_i(ji,jj,jl) , epsi06 ) * zinda 
     263                  rswitch = MAX( 0._wp , SIGN( 1._wp , at_i(ji,jj) - 0.1 ) ) 
     264                  rswitch = rswitch * MAX( 0._wp , SIGN( 1._wp , a_i(ji,jj,jl) - 0.1 ) ) 
     265                  zoi(ji,jj,jl) = oa_i(ji,jj,jl)  / MAX( a_i(ji,jj,jl) , 0.1 ) * rswitch 
    247266               END DO 
    248267            END DO 
    249268         END DO 
    250          CALL iom_put( "iceage_cat"     , zoi        )        ! ice age for categories 
     269         CALL iom_put( "iceage_cat"   , zoi * z1_365 )        ! ice age for categories 
    251270      ENDIF 
    252271 
     
    258277               DO jj = 1, jpj 
    259278                  DO ji = 1, jpi 
    260                      zinda = MAX( 0._wp , SIGN( 1._wp , a_i(ji,jj,jl) - epsi06 ) ) 
    261                      zei(ji,jj,jl) = zei(ji,jj,jl) + 100.0* & 
    262                         ( - tmut * s_i(ji,jj,jk,jl) / MIN( ( t_i(ji,jj,jk,jl) - rtt ), - epsi06 ) ) * & 
    263                         zinda / nlay_i 
     279                     rswitch = MAX( 0._wp , SIGN( 1._wp , a_i(ji,jj,jl) - epsi06 ) ) 
     280                     zei(ji,jj,jl) = zei(ji,jj,jl) + 100.0 * & 
     281                        ( - tmut * s_i(ji,jj,jk,jl) / MIN( ( t_i(ji,jj,jk,jl) - rt0 ), - epsi06 ) ) * & 
     282                        rswitch * r1_nlay_i 
    264283                  END DO 
    265284               END DO 
    266285            END DO 
    267286         END DO 
    268          CALL iom_put( "brinevol_cat"     , zei         )        ! brine volume for categories 
     287         CALL iom_put( "brinevol_cat"     , zei      )        ! brine volume for categories 
    269288      ENDIF 
    270289 
     
    273292      !     not yet implemented 
    274293       
    275       CALL wrk_dealloc( jpi, jpj, jpl, zoi, zei ) 
    276       CALL wrk_dealloc( jpi, jpj     , z2d, zind, z2da, z2db ) 
     294      CALL wrk_dealloc( jpi, jpj, jpl, zoi, zei, zt_i, zt_s ) 
     295      CALL wrk_dealloc( jpi, jpj     , z2d, zswi, z2da, z2db ) 
    277296 
    278297      IF( nn_timing == 1 )  CALL timing_stop('limwri') 
     
    347366      CALL histwrite( kid, "iicethic", kt, icethi        , jpi*jpj, (/1/) )     
    348367      CALL histwrite( kid, "iiceconc", kt, at_i          , jpi*jpj, (/1/) ) 
    349       CALL histwrite( kid, "iicetemp", kt, tm_i - rtt    , jpi*jpj, (/1/) ) 
     368      CALL histwrite( kid, "iicetemp", kt, tm_i - rt0    , jpi*jpj, (/1/) ) 
    350369      CALL histwrite( kid, "iicevelu", kt, u_ice          , jpi*jpj, (/1/) ) 
    351370      CALL histwrite( kid, "iicevelv", kt, v_ice          , jpi*jpj, (/1/) ) 
Note: See TracChangeset for help on using the changeset viewer.