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 1705 for trunk/NEMO/OPA_SRC/SBC – NEMO

Ignore:
Timestamp:
2009-11-03T17:37:00+01:00 (15 years ago)
Author:
smasson
Message:

impact of HF winds in TKE, see ticket:585

Location:
trunk/NEMO/OPA_SRC/SBC
Files:
4 edited

Legend:

Unmodified
Added
Removed
  • trunk/NEMO/OPA_SRC/SBC/sbc_oce.F90

    r1695 r1705  
    2929   !                                             !:  = 1 global mean of e-p-r set to zero at each nn_fsbc time step 
    3030   !                                             !:  = 2 annual global mean of e-p-r set to zero 
    31   INTEGER , PUBLIC ::   nn_ico_cpl  = 0          !: ice-ocean coupling indicator 
     31   INTEGER , PUBLIC ::   nn_ico_cpl  = 0          !: ice-ocean coupling indicator 
    3232   !                                             !:  = 0   LIM-3 old case 
    3333   !                                             !:  = 1   stresses computed using now ocean velocity 
     
    3737   !!              Ocean Surface Boundary Condition fields 
    3838   !!---------------------------------------------------------------------- 
     39   LOGICAL , PUBLIC ::   lhftau = .FALSE.              !: HF tau contribution: mean of stress module - module of the mean stress 
    3940   REAL(wp), PUBLIC, DIMENSION(jpi,jpj) ::   utau      !: sea surface i-stress (ocean referential)     [N/m2] 
    4041   REAL(wp), PUBLIC, DIMENSION(jpi,jpj) ::   vtau      !: sea surface j-stress (ocean referential)     [N/m2] 
    4142   REAL(wp), PUBLIC, DIMENSION(jpi,jpj) ::   taum      !: module of sea surface stress (at T-point)    [N/m2]  
    42    REAL(wp), PUBLIC, DIMENSION(jpi,jpj) ::   wndm      !: wind speed module at T-point (=|U10m-Uoce|)  [m/s] Used only in PISCES 
     43   !! wndm is used only in PISCES to compute gases exchanges at the surface of the free ocean or in the leads in sea-ice parts 
     44   REAL(wp), PUBLIC, DIMENSION(jpi,jpj) ::   wndm      !: wind speed module at T-point (=|U10m-Uoce|)  [m/s]  
    4345   REAL(wp), PUBLIC, DIMENSION(jpi,jpj) ::   qsr       !: sea heat flux:     solar                     [W/m2] 
    4446   REAL(wp), PUBLIC, DIMENSION(jpi,jpj) ::   qns       !: sea heat flux: non solar                     [W/m2] 
  • trunk/NEMO/OPA_SRC/SBC/sbcblk_core.F90

    r1695 r1705  
    4343   PUBLIC   blk_ice_core       ! routine called in sbc_ice_lim module 
    4444       
    45    INTEGER , PARAMETER ::   jpfld   = 8           ! maximum number of files to read  
     45   INTEGER , PARAMETER ::   jpfld   = 9           ! maximum number of files to read  
    4646   INTEGER , PARAMETER ::   jp_wndi = 1           ! index of 10m wind velocity (i-component) (m/s)    at T-point 
    4747   INTEGER , PARAMETER ::   jp_wndj = 2           ! index of 10m wind velocity (j-component) (m/s)    at T-point 
     
    5252   INTEGER , PARAMETER ::   jp_prec = 7           ! index of total precipitation (rain+snow) (Kg/m2/s) 
    5353   INTEGER , PARAMETER ::   jp_snow = 8           ! index of snow (solid prcipitation)       (kg/m2/s) 
     54   INTEGER , PARAMETER ::   jp_tdif = 9           ! index of tau diff associated to HF tau   (N/m2)   at T-point 
    5455   TYPE(FLD), ALLOCATABLE, DIMENSION(:) ::   sf   ! structure of input fields (file informations, fields read) 
    5556          
     
    6364 
    6465   !                                !!* Namelist namsbc_core : CORE bulk parameters 
    65    LOGICAL  ::   ln_2m = .FALSE.     ! logical flag for height of air temp. and hum 
    66    REAL(wp) ::   rn_pfac = 1.        ! multiplication factor for precipitation 
     66   LOGICAL  ::   ln_2m     = .FALSE.     ! logical flag for height of air temp. and hum 
     67   LOGICAL  ::   ln_taudif = .FALSE.     ! logical flag to use the "mean of stress module - module of mean stress" data 
     68   REAL(wp) ::   rn_pfac   = 1.          ! multiplication factor for precipitation 
    6769 
    6870   !! * Substitutions 
     
    9395      !!      the total precipitation (rain+snow) (Kg/m2/s) 
    9496      !!      the snow (solid prcipitation)       (kg/m2/s) 
     97      !!   OPTIONAL parameter (see ln_taudif namelist flag): 
     98      !!      the tau diff associated to HF tau   (N/m2)   at T-point  
    9599      !!              (2) CALL blk_oce_core 
    96100      !! 
     
    110114      INTEGER  ::   ierror   ! return error code 
    111115      INTEGER  ::   ifpr     ! dummy loop indice 
     116      INTEGER  ::   jfld     ! dummy loop arguments 
    112117      !! 
    113118      CHARACTER(len=100) ::  cn_dir   !   Root directory for location of core files 
     
    115120      TYPE(FLD_N) ::   sn_wndi, sn_wndj, sn_humi, sn_qsr       ! informations about the fields to be read 
    116121      TYPE(FLD_N) ::   sn_qlw , sn_tair, sn_prec, sn_snow      !   "                                 " 
    117       NAMELIST/namsbc_core/ cn_dir, ln_2m, rn_pfac, sn_wndi, sn_wndj, sn_humi, sn_qsr,   & 
    118          &                                                sn_qlw , sn_tair, sn_prec, sn_snow 
     122      TYPE(FLD_N) ::   sn_tdif                                 !   "                                 " 
     123      NAMELIST/namsbc_core/ cn_dir , ln_2m  , ln_taudif, rn_pfac,           & 
     124         &                  sn_wndi, sn_wndj, sn_humi  , sn_qsr ,           & 
     125         &                  sn_qlw , sn_tair, sn_prec  , sn_snow, sn_tdif 
    119126      !!--------------------------------------------------------------------- 
    120127 
     
    136143         sn_prec = FLD_N( 'precip'  ,    -1.    ,  'precip'  ,  .true.    , .false. ,   'yearly'  , ''       , ''         ) 
    137144         sn_snow = FLD_N( 'snow'    ,    -1.    ,  'snow'    ,  .true.    , .false. ,   'yearly'  , ''       , ''         ) 
     145         sn_tdif = FLD_N( 'taudif'  ,    24.    ,  'taudif'  ,  .true.    , .false. ,   'yearly'  , ''       , ''         ) 
    138146         ! 
    139147         REWIND( numnam )                    ! ... read in namlist namsbc_core 
     
    145153         slf_i(jp_tair) = sn_tair   ;   slf_i(jp_humi) = sn_humi 
    146154         slf_i(jp_prec) = sn_prec   ;   slf_i(jp_snow) = sn_snow 
     155         slf_i(jp_tdif) = sn_tdif 
     156         ! 
     157         ! do we use HF tau information? 
     158         lhftau = ln_taudif 
     159         jfld = jpfld - COUNT( (/.NOT. lhftau/) ) 
    147160         ! 
    148161         ! set sf structure 
    149          ALLOCATE( sf(jpfld), STAT=ierror ) 
     162         ALLOCATE( sf(jfld), STAT=ierror ) 
    150163         IF( ierror > 0 ) THEN 
    151164            CALL ctl_stop( 'sbc_blk_core: unable to allocate sf structure' )   ;   RETURN 
    152165         ENDIF 
    153          DO ifpr= 1, jpfld 
     166         DO ifpr= 1, jfld 
    154167            ALLOCATE( sf(ifpr)%fnow(jpi,jpj) ) 
    155168            ALLOCATE( sf(ifpr)%fdta(jpi,jpj,2) ) 
     
    291304         END DO 
    292305      END DO 
     306 
     307      ! ... add the HF tau contribution to the wind stress module? 
     308      IF( lhftau ) THEN  
     309!CDIR COLLAPSE 
     310         taum(:,:) = taum(:,:) + sf(jp_tdif)%fnow(:,:) 
     311      ENDIF 
     312      CALL iom_put( "taum_oce", taum )   ! output wind stress module 
     313 
    293314      ! ... utau, vtau at U- and V_points, resp. 
    294315      !     Note the use of 0.5*(2-umask) in order to unmask the stress along coastlines 
  • trunk/NEMO/OPA_SRC/SBC/sbccpl.F90

    r1698 r1705  
    419419      !                                                      ! ------------------------- ! 
    420420      srcv(jpr_taum)%clname = 'O_TauMod'   ;   IF( TRIM(cn_rcv_taumod) == 'coupled' )   srcv(jpr_taum)%laction = .TRUE. 
     421      lhftau = srcv(jpr_taum)%laction 
    421422 
    422423#if defined key_cpl_carbon_cycle 
     
    693694         taum(:,:) = frcv(:,:,jpr_taum) 
    694695         wndm(:,:) = frcv(:,:,jpr_w10m) 
     696         CALL iom_put( "taum_oce", taum )   ! output wind stress module 
    695697         !   
    696698      ENDIF 
  • trunk/NEMO/OPA_SRC/SBC/sbcmod.F90

    r1649 r1705  
    253253      CALL iom_put( "utau", utau )   ! i-wind stress   (stress can be updated at  
    254254      CALL iom_put( "vtau", vtau )   ! j-wind stress    each time step in sea-ice) 
    255       CALL iom_put( "wspd", wndm )   ! wind speed module  
     255      CALL iom_put( "taum", taum )   ! wind stress module  
     256      CALL iom_put( "wspd", wndm )   ! wind speed  module  
    256257      ! 
    257258      IF(ln_ctl) THEN         ! print mean trends (used for debugging) 
Note: See TracChangeset for help on using the changeset viewer.