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 7099 – NEMO

Changeset 7099


Ignore:
Timestamp:
2016-10-25T18:15:30+02:00 (8 years ago)
Author:
jcastill
Message:

Changes as in r7078 of the original ING branch: branches/2015/dev_r5936_INGV1_WAVE

Location:
branches/UKMO/r5518_INGV1_WAVE-coupling/NEMOGCM
Files:
2 added
11 edited

Legend:

Unmodified
Added
Removed
  • branches/UKMO/r5518_INGV1_WAVE-coupling/NEMOGCM/CONFIG/SHARED/namelist_ref

    r5501 r7099  
    267267   ln_cdgw = .false.       !  Neutral drag coefficient read from wave model (T => fill namsbc_wave) 
    268268   ln_sdw  = .false.       !  Computation of 3D stokes drift                (T => fill namsbc_wave) 
     269   ln_tauoc= .false.       !  Activate ocean stress modified by external wave induced stress (T => fill namsbc_wave) 
     270   ln_stcor= .false.       !  Activate Stokes Coriolis term (T => fill namsbc_wave) 
    269271   nn_lsm  = 0             !  =0 land/sea mask for input fields is not applied (keep empty land/sea mask filename field) , 
    270272                           !  =1:n number of iterations of land/sea mask application for input fields (fill land/sea mask filename field) 
     
    921923   ln_zdfexp   = .false.   !  time-stepping: split-explicit (T) or implicit (F) time stepping 
    922924   nn_zdfexp   =    3            !  number of sub-timestep for ln_zdfexp=T 
     925   ln_zdfqiao  = .false.   !  Enhanced wave vertical mixing Qiao (2010) 
    923926/ 
    924927!----------------------------------------------------------------------- 
     
    12771280!              !  file name  ! frequency (hours) ! variable     ! time interp. !  clim   ! 'yearly'/ ! weights  ! rotation ! land/sea mask ! 
    12781281!              !             !  (if <0  months)  !   name       !   (logical)  !  (T/F)  ! 'monthly' ! filename ! pairing  ! filename      ! 
    1279    sn_cdg      =  'cdg_wave' ,        1          , 'drag_coeff' ,     .true.   , .false. , 'daily'   ,  ''      , ''       , '' 
     1282   sn_cdg      =  'sdw_wave' ,        1          , 'drag_coeff' ,     .true.   , .false. , 'daily'   ,  ''      , ''       , '' 
    12801283   sn_usd      =  'sdw_wave' ,        1          , 'u_sd2d'     ,     .true.   , .false. , 'daily'   ,  ''      , ''       , '' 
    12811284   sn_vsd      =  'sdw_wave' ,        1          , 'v_sd2d'     ,     .true.   , .false. , 'daily'   ,  ''      , ''       , '' 
    1282    sn_wn       =  'sdw_wave' ,        1          , 'wave_num'   ,     .true.   , .false. , 'daily'   ,  ''      , ''       , '' 
    1283 ! 
    1284    cn_dir_cdg  = './'  !  root directory for the location of drag coefficient files 
     1285   sn_swh      =  'sdw_wave' ,        1          , 'hs'         ,     .true.   , .false. , 'daily'   ,  ''      , ''       , ''  
     1286   sn_wmp      =  'sdw_wave' ,        1          , 'wmp'        ,     .true.   , .false. , 'daily'   ,  ''      , ''       , ''  
     1287   sn_wnum     =  'sdw_wave' ,        1          , 'wave_num'   ,     .true.   , .false. , 'daily'   ,  ''      , ''       , ''  
     1288   sn_tauoc    -  'sdw_wave' ,        1          , 'wave_stress',     .true.   , .false. , 'daily'   ,  ''      , ''       , ''  
     1289! 
     1290   cn_dir  = './'  !  root directory for the location of drag coefficient files 
    12851291/ 
    12861292!----------------------------------------------------------------------- 
  • branches/UKMO/r5518_INGV1_WAVE-coupling/NEMOGCM/NEMO/OPA_SRC/DYN/dynspg_ts.F90

    r5436 r7099  
    4242   USE sbcapr          ! surface boundary condition: atmospheric pressure 
    4343   USE dynadv, ONLY: ln_dynadv_vec 
     44   USE sbcwave, ONLY: usd2d, vsd2d 
    4445#if defined key_agrif 
    4546   USE agrif_opa_interp ! agrif 
  • branches/UKMO/r5518_INGV1_WAVE-coupling/NEMOGCM/NEMO/OPA_SRC/DYN/sshwzv.F90

    r5467 r7099  
    99   !!             -   !  2010-09  (D.Storkey and E.O'Dea) bug fixes for BDY module 
    1010   !!            3.3  !  2011-10  (M. Leclair) split former ssh_wzv routine and remove all vvl related work 
     11   !!            3.6  !  2014-10  (E. Clementi, P. Oddo) add wave contribution to surface vertical velocity 
    1112   !!---------------------------------------------------------------------- 
    1213 
     
    3940   USE wrk_nemo        ! Memory Allocation 
    4041   USE timing          ! Timing 
     42   USE sbcwave,  ONLY:  zusd2dt, zvsd2dt,wsd3d 
    4143 
    4244   IMPLICIT NONE 
     
    164166      REAL(wp), POINTER, DIMENSION(:,:  ) ::  z2d 
    165167      REAL(wp), POINTER, DIMENSION(:,:,:) ::  z3d, zhdiv 
     168      REAL(wp), ALLOCATABLE, DIMENSION(:,:  ) ::  sshnu, sshnv, dsshnu, dsshnv 
    166169      ! 
    167170      INTEGER             ::   ji, jj, jk   ! dummy loop indices 
     
    216219      ENDIF 
    217220 
     221!  
     222!     In case ln_wave and ln_sdw the surface vertical velocity is modified  
     223!     accounting for Sokes Drift velocity  
     224      IF ( ln_wave .AND. ln_sdw )  THEN  
     225        ALLOCATE(sshnu(jpi,jpj),sshnv(jpi,jpj),dsshnu(jpi,jpj),dsshnv(jpi,jpj))  
     226        sshnu (:,:) = 0._wp  
     227        sshnv (:,:) = 0._wp  
     228        dsshnu(:,:) = 0._wp  
     229        dsshnv(:,:) = 0._wp  
     230        ! sshn interpolated on U-V grid  
     231        !---------------------------------  
     232        DO jj = 1 , jpjm1  
     233          DO ji = 1 , jpim1  
     234            sshnu(ji,jj) =  0.5 * ( 2. - umask(ji,jj,1) ) * &  
     235                         &        ( sshn(ji  ,jj) * tmask(ji  ,jj,1)      &  
     236                         &        + sshn(ji+1,jj) * tmask(ji+1,jj,1) )  
     237            sshnv(ji,jj) =  0.5 * ( 2. - vmask(ji,jj,1) ) * &  
     238                         &        ( sshn(ji,jj  ) * tmask(ji,jji ,1)      &  
     239                         &        + sshn(ji,jj+1) * tmask(ji,jj+1,1) )  
     240          ENDDO  
     241        ENDDO  
     242        ! Compute d(ssh)/dx  and d(ssh)/dy    
     243        !---------------------------------  
     244        DO jj = 1 , jpjm1  
     245          DO ji = 1 , jpim1  
     246            dsshnu(ji,jj) = ( sshnu(ji+1,jj) - sshnu(ji,jj) ) / e1u(ji,jj)  
     247            dsshnv(ji,jj) = ( sshnv(ji,jj+1) - sshnv(ji,jj) ) / e2v(ji,jj)  
     248          ENDDO  
     249        ENDDO  
     250        ! Compute the surface vertical velocity accounting for the Stokes Drift   
     251        !---------------------------------------------------------------------  
     252        wn(:,:,1) = wn(:,:,1) + zusd2dt(:,:) * dsshnu(:,:)     &  
     253                  &           + zvsd2dt(:,:) * dsshnv(:,:)      &  
     254                  &           - ( wsd3d (:,:,1) ) * tmask(:,:,1)  
     255      ENDIF  
     256      !  
    218257#if defined key_bdy 
    219258      IF (lk_bdy) THEN 
  • branches/UKMO/r5518_INGV1_WAVE-coupling/NEMOGCM/NEMO/OPA_SRC/SBC/sbc_oce.F90

    r5407 r7099  
    6565   LOGICAL , PUBLIC ::   ln_cdgw        !: true if neutral drag coefficient from wave model 
    6666   LOGICAL , PUBLIC ::   ln_sdw         !: true if 3d stokes drift from wave model 
     67   LOGICAL , PUBLIC ::   ln_tauoc       !: true if normalized stress from wave is used  
     68   LOGICAL , PUBLIC ::   ln_stcor       !: true if Stokes-Coriolis term is used  
    6769   ! 
    6870   LOGICAL , PUBLIC ::   ln_icebergs    !: Icebergs 
  • branches/UKMO/r5518_INGV1_WAVE-coupling/NEMOGCM/NEMO/OPA_SRC/SBC/sbcmod.F90

    r5501 r7099  
    8989         &             ln_blk_mfs, ln_apr_dyn, nn_ice, nn_ice_embd, ln_dm2dc   , ln_rnf   ,   & 
    9090         &             ln_ssr    , nn_isf    , nn_fwb, ln_cdgw    , ln_wave    , ln_sdw   ,   & 
    91          &             nn_lsm    , nn_limflx , nn_components, ln_cpl 
     91         &             ln_tauoc  , ln_stcor  , nn_lsm, nn_limflx  , nn_components, ln_cpl  
    9292      INTEGER  ::   ios 
    9393      INTEGER  ::   ierr, ierr0, ierr1, ierr2, ierr3, jpm 
     
    381381      END SELECT 
    382382 
     383      IF (ln_wave .AND. ln_tauoc) THEN  
     384            utau(:,:) = utau(:,:)*tauoc_wave(:,:)  
     385            vtau(:,:) = vtau(:,:)*tauoc_wave(:,:)  
     386            taum(:,:) = taum(:,:)*tauoc_wave(:,:)  
     387      !  
     388            SELECT CASE( nsbc )  
     389            CASE(  0,1,2,3,5,-1 )  ;  
     390                IF(lwp) WRITE(numout,*) 'WARNING: You are subtracting the wave stress to the ocean. &  
     391                        & If not requested select ln_tauoc=.false'  
     392            END SELECT  
     393      !  
     394      END IF  
    383395      IF( ln_mixcpl )        CALL sbc_cpl_rcv ( kt, nn_fsbc, nn_ice )   ! forced-coupled mixed formulation after forcing 
    384396 
  • branches/UKMO/r5518_INGV1_WAVE-coupling/NEMOGCM/NEMO/OPA_SRC/SBC/sbcwave.F90

    r5215 r7099  
    44   !! Wave module  
    55   !!====================================================================== 
    6    !! History :  3.3.1  !   2011-09  (Adani M)  Original code: Drag Coefficient  
    7    !!         :  3.4    !   2012-10  (Adani M)                 Stokes Drift  
    8    !!---------------------------------------------------------------------- 
    9    USE iom             ! I/O manager library 
    10    USE in_out_manager  ! I/O manager 
    11    USE lib_mpp         ! distribued memory computing library 
    12    USE fldread        ! read input fields 
    13    USE oce 
    14    USE sbc_oce        ! Surface boundary condition: ocean fields 
    15    USE domvvl 
    16  
    17     
    18    !!---------------------------------------------------------------------- 
    19    !!   sbc_wave       : read drag coefficient from wave model in netcdf files  
    20    !!---------------------------------------------------------------------- 
     6   !! History :  3.3  !   2011-09  (Adani M)  Original code: Drag Coefficient  
     7   !!         :  3.4  !   2012-10  (Adani M)                 Stokes Drift  
     8   !!            3.6  !   2014-09  (Clementi E, Oddo P)New Stokes Drift Computation 
     9   !!---------------------------------------------------------------------- 
     10 
     11   !!---------------------------------------------------------------------- 
     12   !!   sbc_wave      : wave data from wave model in netcdf files  
     13   !!---------------------------------------------------------------------- 
     14   USE oce            !  
     15   USE sbc_oce       ! Surface boundary condition: ocean fields 
     16   USE bdy_oce        ! 
     17   USE domvvl         ! 
     18   ! 
     19   USE iom            ! I/O manager library 
     20   USE in_out_manager ! I/O manager 
     21   USE lib_mpp        ! distribued memory computing library 
     22   USE fldread       ! read input fields 
     23   USE wrk_nemo       ! 
     24   USE phycst         ! physical constants  
    2125 
    2226   IMPLICIT NONE 
    2327   PRIVATE 
    2428 
    25    PUBLIC   sbc_wave    ! routine called in sbc_blk_core or sbc_blk_mfs 
     29   PUBLIC   sbc_wave    ! routine called in sbcmod 
    2630    
    27    INTEGER , PARAMETER ::   jpfld  = 3           ! maximum number of files to read for srokes drift 
     31   INTEGER , PARAMETER ::   jpfld  = 4           ! number of files to read for stokes drift 
    2832   INTEGER , PARAMETER ::   jp_usd = 1           ! index of stokes drift  (i-component) (m/s)    at T-point 
    2933   INTEGER , PARAMETER ::   jp_vsd = 2           ! index of stokes drift  (j-component) (m/s)    at T-point 
    30    INTEGER , PARAMETER ::   jp_wn  = 3           ! index of wave number                 (1/m)    at T-point 
    31    TYPE(FLD), ALLOCATABLE, DIMENSION(:)  :: sf_cd    ! structure of input fields (file informations, fields read) Drag Coefficient 
    32    TYPE(FLD), ALLOCATABLE, DIMENSION(:)  :: sf_sd    ! structure of input fields (file informations, fields read) Stokes Drift 
    33    REAL(wp),PUBLIC,ALLOCATABLE,DIMENSION (:,:)       :: cdn_wave  
    34    REAL(wp),ALLOCATABLE,DIMENSION (:,:)              :: usd2d,vsd2d,uwavenum,vwavenum  
    35    REAL(wp),PUBLIC,ALLOCATABLE,DIMENSION (:,:,:)     :: usd3d,vsd3d,wsd3d  
     34   INTEGER , PARAMETER ::   jp_swh = 3           ! index of significant wave hight      (m)      at T-point 
     35   INTEGER , PARAMETER ::   jp_wmp = 4           ! index of mean wave period            (s)      at T-point 
     36 
     37   TYPE(FLD), ALLOCATABLE, DIMENSION(:)  :: sf_cd    ! structure of input fields (file informations, fields read) Drag Coefficient 
     38   TYPE(FLD), ALLOCATABLE, DIMENSION(:)  :: sf_sd    ! structure of input fields (file informations, fields read) Stokes Drift 
     39   TYPE(FLD), ALLOCATABLE, DIMENSION(:)  :: sf_wn    ! structure of input fields (file informations, fields read) wave number for Qiao 
     40   TYPE(FLD), ALLOCATABLE, DIMENSION(:)  :: sf_tauoc ! structure of input fields (file informations, fields read) normalized wave stress into the ocean 
     41   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:)       :: cdn_wave  
     42   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:)       :: swh,wmp, wnum 
     43   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:)       :: tauoc_wave 
     44   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:)       :: tsd2d 
     45   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:)       :: usd2d, vsd2d 
     46   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:)     :: usd3d, vsd3d, wsd3d  
     47   REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:)       :: zusd2dt, zvsd2dt 
    3648 
    3749   !! * Substitutions 
    3850#  include "domzgr_substitute.h90" 
    39    !!---------------------------------------------------------------------- 
    40    !! NEMO/OPA 4.0 , NEMO Consortium (2011)  
     51#  include "vectopt_loop_substitute.h90" 
     52   !!---------------------------------------------------------------------- 
     53   !! NEMO/OPA 3.7 , NEMO Consortium (2014)  
    4154   !! $Id$ 
    4255   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt) 
     
    4659   SUBROUTINE sbc_wave( kt ) 
    4760      !!--------------------------------------------------------------------- 
    48       !!                     ***  ROUTINE sbc_apr  *** 
     61      !!                     ***  ROUTINE sbc_wave  *** 
    4962      !! 
    50       !! ** Purpose :   read drag coefficient from wave model  in netcdf files. 
     63      !! ** Purpose :   read wave parameters from wave model  in netcdf files. 
    5164      !! 
    5265      !! ** Method  : - Read namelist namsbc_wave 
    5366      !!              - Read Cd_n10 fields in netcdf files  
    5467      !!              - Read stokes drift 2d in netcdf files  
    55       !!              - Read wave number      in netcdf files  
    56       !!              - Compute 3d stokes drift using monochromatic 
    57       !! ** action  :    
    58       !!                
     68      !!              - Read wave number in netcdf files  
     69      !!              - Compute 3d stokes drift using Breivik et al.,2014 
     70      !!                formulation 
     71      !! ** action   
    5972      !!--------------------------------------------------------------------- 
    60       USE oce,  ONLY : un,vn,hdivn,rotn 
    61       USE divcur 
    62       USE wrk_nemo 
    63 #if defined key_bdy 
    64       USE bdy_oce, ONLY : bdytmask 
    65 #endif 
    66       INTEGER, INTENT( in  ) ::  kt       ! ocean time step 
    67       INTEGER                ::  ierror   ! return error code 
    68       INTEGER                ::  ifpr, jj,ji,jk  
    69       INTEGER                ::   ios     ! Local integer output status for namelist read 
    70       REAL(wp),DIMENSION(:,:,:),POINTER             ::  udummy,vdummy,hdivdummy,rotdummy 
    71       REAL                                          ::  z2dt,z1_2dt 
     73      USE zdf_oce,  ONLY : ln_zdfqiao 
     74 
     75      INTEGER, INTENT( in  ) ::   kt       ! ocean time step 
     76      ! 
     77      INTEGER                ::   ierror   ! return error code 
     78      INTEGER                ::   ifpr, jj,ji,jk  
     79      INTEGER                ::   ios      ! Local integer output status for namelist read 
     80      ! 
     81      CHARACTER(len=100)     ::  cn_dir                          ! Root directory for location of drag coefficient files 
     82      REAL(wp)                       ::  ztransp, zsp0, zk, zus, zvs 
     83      REAL(wp), DIMENSION(jpi,jpj)   ::  zfac  
     84      REAL(wp), DIMENSION(:,:,:), POINTER :: ze3hdiv   ! 3D workspace 
    7285      TYPE(FLD_N), DIMENSION(jpfld) ::   slf_i     ! array of namelist informations on the fields to read 
    73       CHARACTER(len=100)     ::  cn_dir                          ! Root directory for location of drag coefficient files 
    74       TYPE(FLD_N)            ::  sn_cdg, sn_usd, sn_vsd, sn_wn   ! informations about the fields to be read 
     86      TYPE(FLD_N)            ::  sn_cdg, sn_usd, sn_vsd,  & 
     87                             &   sn_swh, sn_wmp, sn_wnum, sn_tauoc      ! informations about the fields to be read 
     88      !! 
     89      NAMELIST/namsbc_wave/  sn_cdg, cn_dir, sn_usd, sn_vsd, sn_swh, sn_wmp, sn_wnum, sn_tauoc 
    7590      !!--------------------------------------------------------------------- 
    76       NAMELIST/namsbc_wave/  sn_cdg, cn_dir, sn_usd, sn_vsd, sn_wn 
    77       !!--------------------------------------------------------------------- 
    78  
    79       !!---------------------------------------------------------------------- 
    80       ! 
    8191      ! 
    8292      !                                         ! -------------------- ! 
     
    92102         IF(lwm) WRITE ( numond, namsbc_wave ) 
    93103         ! 
    94  
    95104         IF ( ln_cdgw ) THEN 
    96105            ALLOCATE( sf_cd(1), STAT=ierror )           !* allocate and fill sf_wave with sn_cdg 
     
    102111            ALLOCATE( cdn_wave(jpi,jpj) ) 
    103112            cdn_wave(:,:) = 0.0 
     113         ENDIF 
     114 
     115         IF ( ln_tauoc ) THEN 
     116            ALLOCATE( sf_tauoc(1), STAT=ierror )           !* allocate and fill sf_wave with sn_tauoc 
     117            IF( ierror > 0 )   CALL ctl_stop( 'STOP', 'sbc_wave: unable to allocate sf_wave structure' ) 
     118            ! 
     119                                   ALLOCATE( sf_tauoc(1)%fnow(jpi,jpj,1)   ) 
     120            IF( sn_tauoc%ln_tint )   ALLOCATE( sf_tauoc(1)%fdta(jpi,jpj,1,2) ) 
     121            CALL fld_fill( sf_tauoc, (/ sn_tauoc /), cn_dir, 'sbc_wave', 'Wave module', 'namsbc_wave' ) 
     122            ALLOCATE( tauoc_wave(jpi,jpj) ) 
     123            tauoc_wave(:,:) = 0.0 
    104124        ENDIF 
     125 
    105126         IF ( ln_sdw ) THEN 
    106             slf_i(jp_usd) = sn_usd ; slf_i(jp_vsd) = sn_vsd; slf_i(jp_wn) = sn_wn 
    107             ALLOCATE( sf_sd(3), STAT=ierror )           !* allocate and fill sf_wave with sn_cdg 
     127            slf_i(jp_usd) = sn_usd ; slf_i(jp_vsd) = sn_vsd; 
     128            slf_i(jp_swh) = sn_swh ; slf_i(jp_wmp) = sn_wmp; 
     129            ALLOCATE( sf_sd(jpfld), STAT=ierror )           !* allocate and fill sf_sd with stokes drift 
    108130            IF( ierror > 0 )   CALL ctl_stop( 'STOP', 'sbc_wave: unable to allocate sf_wave structure' ) 
    109131            ! 
     
    112134               IF( slf_i(ifpr)%ln_tint )   ALLOCATE( sf_sd(ifpr)%fdta(jpi,jpj,1,2) ) 
    113135            END DO 
     136 
    114137            CALL fld_fill( sf_sd, slf_i, cn_dir, 'sbc_wave', 'Wave module ', 'namsbc_wave' ) 
    115             ALLOCATE( usd2d(jpi,jpj),vsd2d(jpi,jpj),uwavenum(jpi,jpj),vwavenum(jpi,jpj) ) 
     138            ALLOCATE( usd2d(jpi,jpj),vsd2d(jpi,jpj) ) 
    116139            ALLOCATE( usd3d(jpi,jpj,jpk),vsd3d(jpi,jpj,jpk),wsd3d(jpi,jpj,jpk) ) 
    117             usd2d(:,:) = 0.0 ;  vsd2d(:,:) = 0.0 ; uwavenum(:,:) = 0.0 ; vwavenum(:,:) = 0.0 
    118             usd3d(:,:,:) = 0.0 ;vsd3d(:,:,:) = 0.0 ; wsd3d(:,:,:) = 0.0 
    119          ENDIF 
    120       ENDIF 
    121          ! 
    122          ! 
    123       IF ( ln_cdgw ) THEN 
    124          CALL fld_read( kt, nn_fsbc, sf_cd )      !* read drag coefficient from external forcing 
     140            ALLOCATE( swh(jpi,jpj), wmp(jpi,jpj) ) 
     141            ALLOCATE( zusd2dt(jpi,jpj), zvsd2dt(jpi,jpj) ) 
     142            usd3d(:,:,:) = 0._wp   ;   usd2d(:,:) = 0._wp   ;    
     143            vsd3d(:,:,:) = 0._wp   ;   vsd2d(:,:) = 0._wp   ;  
     144            wsd3d(:,:,:) = 0._wp   ; 
     145            swh  (:,:)   = 0._wp   ;    wmp (:,:) = 0._wp ; 
     146            IF ( ln_zdfqiao ) THEN     !==  Vertical mixing enhancement using Qiao,2010  ==! 
     147               ALLOCATE( sf_wn(1), STAT=ierror )           !* allocate and fill sf_wave with sn_wnum 
     148               IF( ierror > 0 )   CALL ctl_stop( 'STOP', 'sbc_wave: unable toallocate sf_wave structure' ) 
     149                                      ALLOCATE( sf_wn(1)%fnow(jpi,jpj,1)   ) 
     150               IF( sn_wnum%ln_tint )  ALLOCATE( sf_wn(1)%fdta(jpi,jpj,1,2) ) 
     151               CALL fld_fill( sf_wn, (/ sn_wnum /), cn_dir, 'sbc_wave', 'Wave module', 'namsbc_wave' ) 
     152               ALLOCATE( wnum(jpi,jpj),tsd2d(jpi,jpj) ) 
     153               wnum(:,:) = 0._wp ; tsd2d(:,:) = 0._wp 
     154            ENDIF 
     155         ENDIF 
     156      ENDIF 
     157      ! 
     158      IF ( ln_cdgw ) THEN              !==  Neutral drag coefficient  ==! 
     159         CALL fld_read( kt, nn_fsbc, sf_cd )      ! read from external forcing 
    125160         cdn_wave(:,:) = sf_cd(1)%fnow(:,:,1) 
    126161      ENDIF 
    127       IF ( ln_sdw )  THEN 
    128           CALL fld_read( kt, nn_fsbc, sf_sd )      !* read drag coefficient from external forcing 
    129  
    130          ! Interpolate wavenumber, stokes drift into the grid_V and grid_V 
    131          !------------------------------------------------- 
    132  
    133          DO jj = 1, jpjm1 
    134             DO ji = 1, jpim1 
    135                uwavenum(ji,jj)=0.5 * ( 2. - umask(ji,jj,1) ) * ( sf_sd(3)%fnow(ji,jj,1) * tmask(ji,jj,1) & 
    136                &                                + sf_sd(3)%fnow(ji+1,jj,1) * tmask(ji+1,jj,1) ) 
    137  
    138                vwavenum(ji,jj)=0.5 * ( 2. - vmask(ji,jj,1) ) * ( sf_sd(3)%fnow(ji,jj,1) * tmask(ji,jj,1) & 
    139                &                                + sf_sd(3)%fnow(ji,jj+1,1) * tmask(ji,jj+1,1) ) 
    140  
    141                usd2d(ji,jj) = 0.5 * ( 2. - umask(ji,jj,1) ) * ( sf_sd(1)%fnow(ji,jj,1) * tmask(ji,jj,1) & 
    142                &                                + sf_sd(1)%fnow(ji+1,jj,1) * tmask(ji+1,jj,1) ) 
    143  
    144                vsd2d(ji,jj) = 0.5 * ( 2. - vmask(ji,jj,1) ) * ( sf_sd(2)%fnow(ji,jj,1) * tmask(ji,jj,1) & 
    145                &                                + sf_sd(2)%fnow(ji,jj+1,1) * tmask(ji,jj+1,1) ) 
     162 
     163      IF ( ln_tauoc ) THEN             !==  Wave induced stress  ==! 
     164         CALL fld_read( kt, nn_fsbc, sf_tauoc )      !* read wave norm stress from external forcing 
     165         tauoc_wave(:,:) = sf_tauoc(1)%fnow(:,:,1) 
     166      ENDIF 
     167 
     168      IF ( ln_sdw )  THEN                         !==  Computation of the 3d Stokes Drift  ==!  
     169         ! 
     170         CALL fld_read( kt, nn_fsbc, sf_sd )      !* read wave parameters from external forcing 
     171         swh(:,:)     = sf_sd(jp_swh)%fnow(:,:,1)   ! significant wave height 
     172         wmp(:,:)     = sf_sd(jp_wmp)%fnow(:,:,1)   ! wave mean period 
     173         zusd2dt(:,:) = sf_sd(jp_usd)%fnow(:,:,1)   ! 2D zonal Stokes Drift at T point 
     174         zvsd2dt(:,:) = sf_sd(jp_vsd)%fnow(:,:,1)   ! 2D meridional Stokes Drift at T point 
     175         ! 
     176         !==  Computation of the 3d Stokes Drift according to Breivik et al.,2014 
     177         !(DOI: 10.1175/JPO-D-14-0020.1)==!  
     178         ! 
     179         CALL wrk_alloc( jpi,jpj,jpk, ze3hdiv ) 
     180         DO jk = 1, jpk 
     181            DO jj = 1, jpj 
     182               DO ji = 1, jpi 
     183               ! On T grid 
     184               ! Stokes transport speed estimated from Hs and Tmean 
     185               ztransp = 2.0_wp*rpi*swh(ji,jj)**2.0_wp/(16.0_wp*MAX(wmp(ji,jj),0.0000001_wp)) 
     186               ! Stokes surface speed 
     187               zsp0 = SQRT( sf_sd(jp_usd)%fnow(ji,jj,1)**2 +  sf_sd(jp_vsd)%fnow(ji,jj,1)**2) 
     188               ! Wavenumber scale 
     189               zk = ABS(zsp0)/MAX(ABS(5.97_wp*ztransp),0.0000001_wp) 
     190               ! Depth attenuation 
     191               zfac(ji,jj) = EXP(-2.0_wp*zk*fsdept(ji,jj,jk))/(1.0_wp+8.0_wp*zk*fsdept(ji,jj,jk)) 
     192               END DO 
     193            END DO 
     194 
     195            DO jj = 1, jpjm1 
     196               DO ji = 1, jpim1 
     197                  ! Into the U and V Grid  
     198                  zus = 0.5 * ( 2. - umask(ji,jj,1) ) * ( zfac(ji,jj) * tmask(ji,jj,1) & 
     199                 &                                + zfac(ji+1,jj) * tmask(ji+1,jj,1) ) 
     200 
     201                  zvs = 0.5 * ( 2. - vmask(ji,jj,1) ) * ( zfac(ji,jj) * tmask(ji,jj,1) & 
     202                 &                                + zfac(ji,jj+1) * tmask(ji,jj+1,1) ) 
     203 
     204                  usd2d(ji,jj) = 0.5 * ( 2. - umask(ji,jj,1) ) * ( zusd2dt(ji,jj) * tmask(ji,jj,1) & 
     205                 &                                             +  zusd2dt(ji+1,jj) * tmask(ji+1,jj,1) ) 
     206 
     207                  vsd2d(ji,jj) = 0.5 * ( 2. - vmask(ji,jj,1) ) * ( zvsd2dt(ji,jj) * tmask(ji,jj,1) & 
     208                 &                                              + zvsd2dt(ji,jj+1) * tmask(ji,jj+1,1) ) 
     209 
     210                  usd3d(ji,jj,jk) = usd2d(ji,jj) * zus 
     211                  vsd3d(ji,jj,jk) = vsd2d(ji,jj) * zvs 
     212               END DO 
    146213            END DO 
    147214         END DO 
    148  
    149           !Computation of the 3d Stokes Drift 
    150           DO jk = 1, jpk 
    151              DO jj = 1, jpj-1 
    152                 DO ji = 1, jpi-1 
    153                    usd3d(ji,jj,jk) = usd2d(ji,jj)*exp(2.0*uwavenum(ji,jj)*(-MIN( gdept_0(ji,jj,jk) , gdept_0(ji+1,jj  ,jk)))) 
    154                    vsd3d(ji,jj,jk) = vsd2d(ji,jj)*exp(2.0*vwavenum(ji,jj)*(-MIN( gdept_0(ji,jj,jk) , gdept_0(ji  ,jj+1,jk)))) 
    155                 END DO 
    156              END DO 
    157              usd3d(jpi,:,jk) = usd2d(jpi,:)*exp( 2.0*uwavenum(jpi,:)*(-gdept_0(jpi,:,jk)) ) 
    158              vsd3d(:,jpj,jk) = vsd2d(:,jpj)*exp( 2.0*vwavenum(:,jpj)*(-gdept_0(:,jpj,jk)) ) 
    159           END DO 
    160  
    161           CALL wrk_alloc( jpi,jpj,jpk,udummy,vdummy,hdivdummy,rotdummy) 
    162            
    163           udummy(:,:,:)=un(:,:,:) 
    164           vdummy(:,:,:)=vn(:,:,:) 
    165           hdivdummy(:,:,:)=hdivn(:,:,:) 
    166           rotdummy(:,:,:)=rotn(:,:,:) 
    167           un(:,:,:)=usd3d(:,:,:) 
    168           vn(:,:,:)=vsd3d(:,:,:) 
    169           CALL div_cur(kt) 
    170       !                                           !------------------------------! 
    171       !                                           !     Now Vertical Velocity    ! 
    172       !                                           !------------------------------! 
    173           z2dt = 2._wp * rdt                              ! set time step size (Euler/Leapfrog) 
    174  
    175           z1_2dt = 1.e0 / z2dt 
    176           DO jk = jpkm1, 1, -1                             ! integrate from the bottom the hor. divergence 
    177              ! - ML - need 3 lines here because replacement of fse3t by its expression yields too long lines otherwise 
    178              wsd3d(:,:,jk) = wsd3d(:,:,jk+1) -   fse3t_n(:,:,jk) * hdivn(:,:,jk)        & 
    179                 &                      - ( fse3t_a(:,:,jk) - fse3t_b(:,:,jk) )    & 
    180                 &                         * tmask(:,:,jk) * z1_2dt 
     215         ! 
     216         CALL lbc_lnk( usd3d(:,:,:), 'U', -1. ) 
     217         CALL lbc_lnk( vsd3d(:,:,:), 'V', -1. ) 
     218         ! 
     219         DO jk = 1, jpkm1                       ! e3t * Horizontal divergence 
     220            DO jj = 2, jpjm1 
     221               DO ji = fs_2, fs_jpim1   ! vector opt. 
     222                  ze3hdiv(ji,jj,jk) = (  e2u(ji  ,jj) * fse3u_n(ji  ,jj,jk) * usd3d(ji  ,jj,jk)     & 
     223                     &                 - e2u(ji-1,jj) * fse3u_n(ji-1,jj,jk) * usd3d(ji-1,jj,jk)     & 
     224                     &                 + e1v(ji,jj  ) * fse3v_n(ji,jj  ,jk) * vsd3d(ji,jj  ,jk)     & 
     225                     &                 - e1v(ji,jj-1) * fse3v_n(ji,jj-1,jk) * vsd3d(ji,jj-1,jk)   ) * r1_e1e2t(ji,jj) 
     226               END DO   
     227            END DO   
     228            IF( .NOT. AGRIF_Root() ) THEN 
     229               IF( nbondi ==  1 .OR. nbondi == 2 )   ze3hdiv(nlci-1,   :  ,jk) = 0._wp      ! east 
     230               IF( nbondi == -1 .OR. nbondi == 2 )   ze3hdiv(  2   ,   :  ,jk) = 0._wp      ! west 
     231               IF( nbondj ==  1 .OR. nbondj == 2 )   ze3hdiv(  :   ,nlcj-1,jk) = 0._wp      ! north 
     232               IF( nbondj == -1 .OR. nbondj == 2 )   ze3hdiv(  :   ,  2   ,jk) = 0._wp      ! south 
     233            ENDIF 
     234         END DO 
     235         CALL lbc_lnk( ze3hdiv, 'T', 1. )  
     236         ! 
     237         DO jk = jpkm1, 1, -1                   !* integrate from the bottom the e3t * hor. divergence 
     238            wsd3d(:,:,jk) = wsd3d(:,:,jk+1) - ze3hdiv(:,:,jk) 
     239         END DO 
    181240#if defined key_bdy 
    182              wsd3d(:,:,jk) = wsd3d(:,:,jk) * bdytmask(:,:) 
     241         IF( lk_bdy ) THEN 
     242            DO jk = 1, jpkm1 
     243               wsd3d(:,:,jk) = wsd3d(:,:,jk) * bdytmask(:,:) 
     244            END DO 
     245         ENDIF 
    183246#endif 
    184           END DO 
    185           hdivn(:,:,:)=hdivdummy(:,:,:) 
    186           rotn(:,:,:)=rotdummy(:,:,:) 
    187           vn(:,:,:)=vdummy(:,:,:) 
    188           un(:,:,:)=udummy(:,:,:) 
    189           CALL wrk_dealloc( jpi,jpj,jpk,udummy,vdummy,hdivdummy,rotdummy) 
    190       ENDIF 
     247         CALL wrk_dealloc( jpi,jpj,jpk, ze3hdiv ) 
     248 
     249         IF ( ln_zdfqiao ) THEN 
     250            wnum(:,:) = sf_wn(1)%fnow(:,:,1) 
     251            
     252            ! Calculate the module of the stokes drift on T grid 
     253            !------------------------------------------------- 
     254            DO jj = 1, jpj 
     255               DO ji = 1, jpi 
     256                   tsd2d(ji,jj) = ((sf_sd(jp_usd)%fnow(ji,jj,1) * tmask(ji,jj,1))**2.0  +     & 
     257                &               (sf_sd(jp_vsd)%fnow(ji,jj,1) * tmask(ji,jj,1))**2.0)**0.5 
     258               END DO 
     259            END DO 
     260         ENDIF 
     261         ! 
     262      ENDIF 
     263      ! 
    191264   END SUBROUTINE sbc_wave 
    192265       
  • branches/UKMO/r5518_INGV1_WAVE-coupling/NEMOGCM/NEMO/OPA_SRC/TRA/traadv.F90

    r5147 r7099  
    77   !!            3.3  !  2010-09  (C. Ethe, G. Madec)  merge TRC-TRA + switch from velocity to transport 
    88   !!            4.0  !  2011-06  (G. Madec)  Addition of Mixed Layer Eddy parameterisation 
     9   !!            3.6  !  2015-06  (E. Clementi) Addition of Stokes drift in case of wave coupling  
    910   !!---------------------------------------------------------------------- 
    1011 
     
    3435   USE timing          ! Timing 
    3536   USE sbc_oce 
    36    USE diaptr          ! Poleward heat transport  
    37  
     37   USE sbcwave        ! wave module  
     38   USE sbc_oce        ! surface boundary condition: ocean  
     39    
     40   USE diaptr         ! Poleward heat transport   
     41   USE sbcwave        ! wave module  
     42   USE sbc_oce        ! surface boundary condition: ocean  
    3843 
    3944   IMPLICIT NONE 
     
    8590      CALL wrk_alloc( jpi, jpj, jpk, zun, zvn, zwn ) 
    8691      !                                          ! set time step 
     92      zun(:,:,:) = 0.0  
     93      zvn(:,:,:) = 0.0  
     94      zwn(:,:,:) = 0.0  
     95      !      
    8796      IF( neuler == 0 .AND. kt == nit000 ) THEN     ! at nit000 
    8897         r2dtra(:) =  rdttra(:)                          ! = rdtra (restarting with Euler time stepping) 
     
    94103      ! 
    95104      !                                               !==  effective transport  ==! 
     105      IF (ln_wave .AND. ln_sdw) THEN  
     106         DO jk = 1, jpkm1  
     107            zun(:,:,jk) = e2u(:,:) * fse3u(:,:,jk) *      &  
     108                        &  ( un(:,:,jk) + usd3d(:,:,jk) )                    !eulerian transport + Stokes Drift  
     109            zvn(:,:,jk) = e1v(:,:) * fse3v(:,:,jk) *      &  
     110                   &  ( vn(:,:,jk) + vsd3d(:,:,jk) )  
     111            zwn(:,:,jk) = e1e2t(:,:) *                    &  
     112                   &  ( wn(:,:,jk) + wsd3d(:,:,jk) )  
     113         END DO  
     114      ELSE  
    96115      DO jk = 1, jpkm1 
    97116         zun(:,:,jk) = e2u(:,:) * fse3u(:,:,jk) * un(:,:,jk)                  ! eulerian transport only 
     
    99118         zwn(:,:,jk) = e1t(:,:) * e2t(:,:)      * wn(:,:,jk) 
    100119      END DO 
     120      ENDIF 
    101121      ! 
    102122      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN 
  • branches/UKMO/r5518_INGV1_WAVE-coupling/NEMOGCM/NEMO/OPA_SRC/ZDF/zdf_oce.F90

    r7096 r7099  
    3535   INTEGER , PUBLIC ::   nn_npc      !: non penetrative convective scheme call  frequency 
    3636   INTEGER , PUBLIC ::   nn_npcp     !: non penetrative convective scheme print frequency 
     37   LOGICAL , PUBLIC ::   ln_zdfqiao  !: Enhanced wave vertical mixing Qiao(2010) formulation flag  
    3738 
    3839 
  • branches/UKMO/r5518_INGV1_WAVE-coupling/NEMOGCM/NEMO/OPA_SRC/ZDF/zdfini.F90

    r5386 r7099  
    5353      INTEGER ::   ios 
    5454      !! 
    55       NAMELIST/namzdf/ rn_avm0, rn_avt0, nn_avb, nn_havtb, ln_zdfexp, nn_zdfexp,   & 
    56          &              ln_zdfevd, nn_evdm, rn_avevd, ln_zdfnpc, nn_npc, nn_npcp 
     55      NAMELIST/namzdf/ rn_avm0, rn_avt0, nn_avb, nn_havtb, ln_zdfexp, nn_zdfexp,  &  
     56         &        ln_zdfevd, nn_evdm, rn_avevd, ln_zdfnpc, nn_npc, nn_npcp,       &  
     57         &        ln_zdfqiao  
    5758      !!---------------------------------------------------------------------- 
    5859 
     
    8384         WRITE(numout,*) '      npc call  frequency                 nn_npc    = ', nn_npc 
    8485         WRITE(numout,*) '      npc print frequency                 nn_npcp   = ', nn_npcp 
     86         WRITE(numout,*) '      Qiao formulation flag ln_zdfqiao=', ln_zdfqiao  
    8587      ENDIF 
    8688 
  • branches/UKMO/r5518_INGV1_WAVE-coupling/NEMOGCM/NEMO/OPA_SRC/step.F90

    r5510 r7099  
    2626   !!                 !  2012-07  (J. Simeon, G. Madec, C. Ethe) Online coarsening of outputs 
    2727   !!            3.7  !  2014-04  (F. Roquet, G. Madec) New equations of state 
     28   !!            3.6  !  2014-10  (E. Clementi, P. Oddo) Add Qiao vertical mixing in case of waves  
    2829   !!---------------------------------------------------------------------- 
    2930 
     
    7273      !!              -8- Outputs and diagnostics 
    7374      !!---------------------------------------------------------------------- 
    74       INTEGER ::   jk      ! dummy loop indice 
     75      INTEGER ::   ji,jj,jk ! dummy loop indice 
    7576      INTEGER ::   indic    ! error indicator if < 0 
    7677      INTEGER ::   kcall    ! optional integer argument (dom_vvl_sf_nxt) 
     
    128129      IF( lk_zdfgls  )   CALL zdf_gls( kstp )            ! GLS closure scheme for Kz 
    129130      IF( lk_zdfkpp  )   CALL zdf_kpp( kstp )            ! KPP closure scheme for Kz 
     131      IF( ln_zdfqiao )   THEN  
     132        CALL zdf_qiao(kstp )                             ! Qiao vertical mixing   
     133        DO jk = 1, jpkm1  
     134          DO jj = 1, jpj  
     135            DO ji = 1, jpi  
     136              avmu(ji,jj,jk) = (avmu(ji,jj,jk) + QBvu(ji,jj,jk)) * umask(ji,jj,jk)  
     137              avmv(ji,jj,jk) = (avmv(ji,jj,jk) + QBvv(ji,jj,jk)) * vmask(ji,jj,jk)  
     138              avt( ji,jj,jk) = (avt( ji,jj,jk) + QBv(ji,jj,jk))  * tmask(ji,jj,jk)  
     139            END DO  
     140          END DO  
     141        END DO  
     142      ENDIF  
     143      !  
    130144      IF( lk_zdfcst  ) THEN                              ! Constant Kz (reset avt, avm[uv] to the background value) 
    131145         avt (:,:,:) = rn_avt0 * wmask (:,:,:) 
     
    207221                                  CALL dyn_vor      ( kstp )   ! vorticity term including Coriolis 
    208222                                  CALL dyn_ldf      ( kstp )   ! lateral mixing 
     223          IF( ln_stcor    )       CALL dyn_stcor    ( kstp )   ! Stokes-Coriolis forcing 
    209224          IF( ln_neptsimp )       CALL dyn_nept_cor ( kstp )   ! add Neptune velocities (simplified) 
    210225#if defined key_agrif 
  • branches/UKMO/r5518_INGV1_WAVE-coupling/NEMOGCM/NEMO/OPA_SRC/step_oce.F90

    r5501 r7099  
    2828   USE sbctide          ! Tide initialisation 
    2929   USE sbcapr           ! surface boundary condition: ssh_ib required by bdydta  
     30   USE sbcwave          ! Wave intialisation 
    3031 
    3132   USE traqsr           ! solar radiation penetration      (tra_qsr routine) 
     
    5253   USE dynspg           ! surface pressure gradient        (dyn_spg routine) 
    5354   USE dynnept          ! simp. form of Neptune effect(dyn_nept_cor routine) 
     55   USE dynstcor         ! simp. form of Stokes-Coriolis 
    5456 
    5557   USE dynnxt           ! time-stepping                    (dyn_nxt routine) 
     
    8486   USE zdfric           ! Richardson vertical mixing       (zdf_ric routine) 
    8587   USE zdfmxl           ! Mixed-layer depth                (zdf_mxl routine) 
     88   USE zdfqiao          !Qiao module wave induced mixing   (zdf_qiao routine) 
    8689 
    8790   USE zpshde           ! partial step: hor. derivative     (zps_hde routine) 
Note: See TracChangeset for help on using the changeset viewer.