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 3625 for branches/2012/dev_NOC_2012_rev3555/NEMOGCM/NEMO/OPA_SRC/DIA – NEMO

Ignore:
Timestamp:
2012-11-21T14:19:18+01:00 (12 years ago)
Author:
acc
Message:

Branch dev_NOC_2012_r3555. #1006. Step 7. Check in code now merged with dev_r3385_NOCS04_HAMF

Location:
branches/2012/dev_NOC_2012_rev3555/NEMOGCM/NEMO/OPA_SRC/DIA
Files:
3 edited

Legend:

Unmodified
Added
Removed
  • branches/2012/dev_NOC_2012_rev3555/NEMOGCM/NEMO/OPA_SRC/DIA/diahsb.F90

    r3294 r3625  
    8383      z_frc_trd_s =           SUM( sbc_tsc(:,:,jp_sal) * surf(:,:) )     ! salt fluxes 
    8484      ! Add penetrative solar radiation 
    85       IF( ln_traqsr )   z_frc_trd_t = z_frc_trd_t + ro0cpr * SUM( qsr     (:,:) * surf(:,:) ) 
     85      IF( ln_traqsr )   z_frc_trd_t = z_frc_trd_t + r1_rau0_rcp * SUM( qsr     (:,:) * surf(:,:) ) 
    8686      ! Add geothermal heat flux 
    87       IF( ln_trabbc )   z_frc_trd_t = z_frc_trd_t + ro0cpr * SUM( qgh_trd0(:,:) * surf(:,:) ) 
     87      IF( ln_trabbc )   z_frc_trd_t = z_frc_trd_t + r1_rau0_rcp * SUM( qgh_trd0(:,:) * surf(:,:) ) 
    8888      IF( lk_mpp ) THEN 
    8989         CALL mpp_sum( z_frc_trd_v ) 
  • branches/2012/dev_NOC_2012_rev3555/NEMOGCM/NEMO/OPA_SRC/DIA/diawri.F90

    r3609 r3625  
    400400         CALL histdef( nid_T, "sossheig", "Sea Surface Height"                 , "m"      ,   &  ! ssh 
    401401            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
    402 !!$#if defined key_lim3 || defined key_lim2  
    403 !!$         ! sowaflup = sowaflep + sorunoff + sowafldp + a term associated to 
    404 !!$         !    internal damping to Levitus that can be diagnosed from others 
    405 !!$         ! sowaflcd = sowaflep + sorunoff + sowafldp + iowaflup 
    406 !!$         CALL histdef( nid_T, "iowaflup", "Ice=>ocean net freshwater"          , "kg/m2/s",   &  ! fsalt 
    407 !!$            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
    408 !!$         CALL histdef( nid_T, "sowaflep", "atmos=>ocean net freshwater"        , "kg/m2/s",   &  ! fmass 
    409 !!$            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
    410 !!$#endif 
    411402         CALL histdef( nid_T, "sowaflup", "Net Upward Water Flux"              , "Kg/m2/s",   &  ! (emp-rnf) 
    412403            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
    413 !!$         CALL histdef( nid_T, "sorunoff", "Runoffs"                            , "Kg/m2/s",   &  ! runoffs 
    414 !!$            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
    415          CALL histdef( nid_T, "sowaflcd", "concentration/dilution water flux"  , "kg/m2/s",   &  ! (emps-rnf) 
    416             &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
    417          CALL histdef( nid_T, "sosalflx", "Surface Salt Flux"                  , "Kg/m2/s",   &  ! (emps-rnf) * sn 
    418             &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
     404         CALL histdef( nid_T, "sosfldow", "downward salt flux"                 , "PSU/m2/s",  &  ! sfx 
     405            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
     406#if ! defined key_vvl 
     407         CALL histdef( nid_T, "sosst_cd", "Concentration/Dilution term on temperature"        &  ! emp * tsn(:,:,1,jp_tem) 
     408            &                                                                  , "KgC/m2/s",  &  ! sosst_cd 
     409            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
     410         CALL histdef( nid_T, "sosss_cd", "Concentration/Dilution term on salinity"           &  ! emp * tsn(:,:,1,jp_sal) 
     411            &                                                                  , "KgPSU/m2/s",&  ! sosss_cd 
     412            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
     413#endif 
    419414         CALL histdef( nid_T, "sohefldo", "Net Downward Heat Flux"             , "W/m2"   ,   &  ! qns + qsr 
    420415            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout ) 
     
    602597      CALL histwrite( nid_T, "sosaline", it, tsn(:,:,1,jp_sal), ndim_hT, ndex_hT )   ! sea surface salinity 
    603598      CALL histwrite( nid_T, "sossheig", it, sshn          , ndim_hT, ndex_hT )   ! sea surface height 
    604 !!$#if  defined key_lim3 || defined key_lim2  
    605 !!$      CALL histwrite( nid_T, "iowaflup", it, fsalt(:,:)    , ndim_hT, ndex_hT )   ! ice=>ocean water flux 
    606 !!$      CALL histwrite( nid_T, "sowaflep", it, fmass(:,:)    , ndim_hT, ndex_hT )   ! atmos=>ocean water flux 
    607 !!$#endif 
    608599      CALL histwrite( nid_T, "sowaflup", it, ( emp-rnf )   , ndim_hT, ndex_hT )   ! upward water flux 
    609 !!$      CALL histwrite( nid_T, "sorunoff", it, runoff        , ndim_hT, ndex_hT )   ! runoff 
    610       CALL histwrite( nid_T, "sowaflcd", it, ( emps-rnf )  , ndim_hT, ndex_hT )   ! c/d water flux 
    611       zw2d(:,:) = ( emps(:,:) - rnf(:,:) ) * tsn(:,:,1,jp_sal) * tmask(:,:,1) 
    612       CALL histwrite( nid_T, "sosalflx", it, zw2d          , ndim_hT, ndex_hT )   ! c/d salt flux 
     600      CALL histwrite( nid_T, "sosfldow", it, sfx           , ndim_hT, ndex_hT )   ! downward salt flux  
     601                                                                                  ! (includes virtual salt flux beneath ice  
     602                                                                                  ! in linear free surface case) 
     603#if ! defined key_vvl 
     604      zw2d(:,:) = emp (:,:) * tsn(:,:,1,jp_tem) 
     605      CALL histwrite( nid_T, "sosst_cd", it, zw2d, ndim_hT, ndex_hT )             ! c/d term on sst 
     606      zw2d(:,:) = emp (:,:) * tsn(:,:,1,jp_sal) 
     607      CALL histwrite( nid_T, "sosss_cd", it, zw2d, ndim_hT, ndex_hT )             ! c/d term on sss 
     608#endif 
    613609      CALL histwrite( nid_T, "sohefldo", it, qns + qsr     , ndim_hT, ndex_hT )   ! total heat flux 
    614610      CALL histwrite( nid_T, "soshfldo", it, qsr           , ndim_hT, ndex_hT )   ! solar heat flux 
     
    782778      !!---------------------------------------------------------------------- 
    783779      !  
    784       IF( nn_timing == 1 )   CALL timing_start('dia_wri_state') 
     780!     IF( nn_timing == 1 )   CALL timing_start('dia_wri_state') ! not sure this works for routines not called in first timestep 
    785781 
    786782      ! 0. Initialisation 
     
    879875#endif 
    880876        
    881       IF( nn_timing == 1 )   CALL timing_stop('dia_wri_state') 
     877!     IF( nn_timing == 1 )   CALL timing_stop('dia_wri_state') ! not sure this works for routines not called in first timestep 
    882878      !  
    883879 
  • branches/2012/dev_NOC_2012_rev3555/NEMOGCM/NEMO/OPA_SRC/DIA/diawri_dimg.h90

    r3294 r3625  
    5454    !!  level 14: qct(:,:)                 equivalent flux due to treshold SST 
    5555    !!  level 15: fbt(:,:)                 feedback term . 
    56     !!  level 16: ( emps(:,:) - rnf(:,:) ) concentration/dilution water flux 
     56    !!  level 16: ( emp * sss )            concentration/dilution term on salinity 
     57    !!  level 17: ( emp * sst )            concentration/dilution term on temperature 
    5758    !!  level 17: fsalt(:,:)               Ice=>ocean net freshwater 
    5859    !!  level 18: gps(:,:)                 the surface pressure (m). 
     
    107108 
    108109 
    109     inbsel = 17 
     110    inbsel = 18 
    110111 
    111112    IF( inbsel >  jpk ) THEN 
     
    172173       !        fsel(:,:,14) = fsel(:,:,14) + qct(:,:) 
    173174       !        fsel(:,:,15) = fsel(:,:,15) + fbt(:,:) 
    174        fsel(:,:,16) = fsel(:,:,16) + ( emps(:,:)-rnf(:,:) )  
     175       fsel(:,:,16) = fsel(:,:,16) + ( emp(:,:)*tsn(:,:,1,jp_sal) )  
     176       fsel(:,:,17) = fsel(:,:,17) + ( emp(:,:)*tsn(:,:,1,jp_tem) )  
    175177       ! 
    176178       ! Output of dynamics and tracer fields and selected fields 
     
    240242          !         fsel(:,:,14) =  qct(:,:) 
    241243          !         fsel(:,:,15) =  fbt(:,:) 
    242           fsel(:,:,16) = ( emps(:,:)-rnf(:,:) ) * tmask(:,:,1)  
     244          fsel(:,:,16) = ( emp(:,:)-tsn(:,:,1,jp_sal) ) * tmask(:,:,1)  
     245          fsel(:,:,17) = ( emp(:,:)-tsn(:,:,1,jp_tem) ) * tmask(:,:,1)  
    243246          ! 
    244247          !         qct(:,:) = 0._wp 
Note: See TracChangeset for help on using the changeset viewer.