MODULE sbcwave !!====================================================================== !! *** MODULE sbcwave *** !! Wave module !!====================================================================== !! History : 3.3 ! 2011-09 (Adani M) Original code: Drag Coefficient !! : 3.4 ! 2012-10 (Adani M) Stokes Drift !! 3.6 ! 2014-09 (Clementi E, Oddo P)New Stokes Drift Computation !! - ! 2016-12 (G. Madec, E. Clementi) update Stoke drift computation !! + add sbc_wave_ini routine !!---------------------------------------------------------------------- !!---------------------------------------------------------------------- !! sbc_stokes : calculate 3D Stokes-drift velocities !! sbc_wave : wave data from wave model in netcdf files !! sbc_wave_init : initialisation fo surface waves !!---------------------------------------------------------------------- USE oce ! ocean variables USE sbc_oce ! Surface boundary condition: ocean fields USE bdy_oce ! open boundary condition variables USE domvvl ! domain: variable volume layers ! USE iom ! I/O manager library USE in_out_manager ! I/O manager USE lib_mpp ! distribued memory computing library USE fldread ! read input fields USE wrk_nemo ! USE phycst ! physical constants IMPLICIT NONE PRIVATE PUBLIC sbc_stokes ! routine called in sbccpl PUBLIC sbc_stress ! routine called in sbcmod PUBLIC sbc_wave ! routine called in sbcmod PUBLIC sbc_wave_init ! routine called in sbcmod ! Variables checking if the wave parameters are coupled (if not, they are read from file) LOGICAL, PUBLIC :: cpl_hsig = .FALSE. LOGICAL, PUBLIC :: cpl_phioc = .FALSE. LOGICAL, PUBLIC :: cpl_sdrft = .FALSE. LOGICAL, PUBLIC :: cpl_wper = .FALSE. LOGICAL, PUBLIC :: cpl_wfreq = .FALSE. LOGICAL, PUBLIC :: cpl_wnum = .FALSE. LOGICAL, PUBLIC :: cpl_tauoc = .FALSE. LOGICAL, PUBLIC :: cpl_tauw = .FALSE. LOGICAL, PUBLIC :: cpl_wdrag = .FALSE. INTEGER :: nn_sdrift ! type of parameterization to calculate vertical Stokes drift INTEGER, PARAMETER :: jp_breivik = 0 ! Breivik 2015: v_z=v_0*[exp(2*k*z)/(1-8*k*z)] INTEGER, PARAMETER :: jp_phillips = 1 ! Phillips: v_z=v_o*[exp(2*k*z)-beta*sqrt(-2*k*pi*z)*erfc(sqrt(-2*k*z))] INTEGER, PARAMETER :: jp_peakph = 2 ! Phillips using the peak wave number read from wave model instead of the inverse depth scale INTEGER :: jpfld ! number of files to read for stokes drift INTEGER :: jp_usd ! index of stokes drift (i-component) (m/s) at T-point INTEGER :: jp_vsd ! index of stokes drift (j-component) (m/s) at T-point INTEGER :: jp_hsw ! index of significant wave hight (m) at T-point INTEGER :: jp_wmp ! index of mean wave period (s) at T-point INTEGER :: jp_wfr ! index of wave peak frequency (s^-1) at T-point TYPE(FLD), ALLOCATABLE, DIMENSION(:) :: sf_cd ! structure of input fields (file informations, fields read) Drag Coefficient TYPE(FLD), ALLOCATABLE, DIMENSION(:) :: sf_sd ! structure of input fields (file informations, fields read) Stokes Drift TYPE(FLD), ALLOCATABLE, DIMENSION(:) :: sf_wn ! structure of input fields (file informations, fields read) wave number for Qiao TYPE(FLD), ALLOCATABLE, DIMENSION(:) :: sf_tauoc ! structure of input fields (file informations, fields read) normalized wave stress into the ocean TYPE(FLD), ALLOCATABLE, DIMENSION(:) :: sf_tauw ! structure of input fields (file informations, fields read) ocean stress components from wave model TYPE(FLD), ALLOCATABLE, DIMENSION(:) :: sf_phioc ! structure of input fields (file informations, fields read) wave to ocean energy REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) :: cdn_wave !: REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) :: hsw, wmp, wnum !: REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) :: wfreq !: REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) :: rn_crban !: Craig and Banner constant for surface breaking waves mixing REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) :: tauoc_wave !: REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) :: tauw_x !: REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) :: tauw_y !: REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) :: tsd2d !: REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) :: div_sd !: barotropic stokes drift divergence REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:) :: ut0sd, vt0sd !: surface Stokes drift velocities at t-point REAL(wp), PUBLIC, ALLOCATABLE, DIMENSION(:,:,:) :: usd , vsd , wsd !: Stokes drift velocities at u-, v- & w-points, resp. # include "vectopt_loop_substitute.h90" !!---------------------------------------------------------------------- !! NEMO/OPA 4.0 , NEMO Consortium (2011) !! $Id$ !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) !!---------------------------------------------------------------------- CONTAINS SUBROUTINE sbc_stokes( ) !!--------------------------------------------------------------------- !! *** ROUTINE sbc_stokes *** !! !! ** Purpose : compute the 3d Stokes Drift according to Breivik et al., !! 2014 (DOI: 10.1175/JPO-D-14-0020.1) !! !! ** Method : - Calculate Stokes transport speed !! - Calculate horizontal divergence !! - Integrate the horizontal divergenze from the bottom !! ** action !!--------------------------------------------------------------------- INTEGER :: jj, ji, jk ! dummy loop argument INTEGER :: ik ! local integer REAL(wp) :: ztransp, zfac, ztemp REAL(wp) :: zdep_u, zdep_v, zkh_u, zkh_v, zda_u, zda_v REAL(wp), DIMENSION(:,:) , POINTER :: zk_t, zk_u, zk_v, zu0_sd, zv0_sd ! 2D workspace REAL(wp), DIMENSION(:,:,:), POINTER :: ze3divh ! 3D workspace !!--------------------------------------------------------------------- ! CALL wrk_alloc( jpi,jpj,jpk, ze3divh ) CALL wrk_alloc( jpi,jpj, zk_t, zk_u, zk_v, zu0_sd, zv0_sd ) ! ! select parameterization for the calculation of vertical Stokes drift ! exp. wave number at t-point IF( nn_sdrift==jp_breivik .OR. nn_sdrift==jp_phillips ) THEN ! (Eq. (19) in Breivick et al. (2014) ) zfac = 2.0_wp * rpi / 16.0_wp DO jj = 1, jpj DO ji = 1, jpi ! Stokes drift velocity estimated from Hs and Tmean ztransp = zfac * hsw(ji,jj)*hsw(ji,jj) / MAX( wmp(ji,jj), 0.0000001_wp ) ! Stokes surface speed tsd2d(ji,jj) = SQRT( ut0sd(ji,jj)*ut0sd(ji,jj) + vt0sd(ji,jj)*vt0sd(ji,jj)) ! Wavenumber scale zk_t(ji,jj) = ABS( tsd2d(ji,jj) ) / MAX( ABS( 5.97_wp*ztransp ), 0.0000001_wp ) END DO END DO DO jj = 1, jpjm1 ! exp. wave number & Stokes drift velocity at u- & v-points DO ji = 1, jpim1 zk_u(ji,jj) = 0.5_wp * ( zk_t(ji,jj) + zk_t(ji+1,jj) ) zk_v(ji,jj) = 0.5_wp * ( zk_t(ji,jj) + zk_t(ji,jj+1) ) ! zu0_sd(ji,jj) = 0.5_wp * ( ut0sd(ji,jj) + ut0sd(ji+1,jj) ) zv0_sd(ji,jj) = 0.5_wp * ( vt0sd(ji,jj) + vt0sd(ji,jj+1) ) END DO END DO ELSE IF( nn_sdrift==jp_peakph ) THEN ! peak wave number calculated from the peak frequency received by the wave model DO jj = 1, jpjm1 DO ji = 1, jpim1 zk_u(ji,jj) = 0.5_wp * ( wfreq(ji,jj)*wfreq(ji,jj) + wfreq(ji+1,jj)*wfreq(ji+1,jj) ) / grav zk_v(ji,jj) = 0.5_wp * ( wfreq(ji,jj)*wfreq(ji,jj) + wfreq(ji,jj+1)*wfreq(ji,jj+1) ) / grav ! zu0_sd(ji,jj) = 0.5_wp * ( ut0sd(ji,jj) + ut0sd(ji+1,jj) ) zv0_sd(ji,jj) = 0.5_wp * ( vt0sd(ji,jj) + vt0sd(ji,jj+1) ) END DO END DO ENDIF ! ! !== horizontal Stokes Drift 3D velocity ==! IF( nn_sdrift==jp_breivik ) THEN DO jk = 1, jpkm1 DO jj = 2, jpjm1 DO ji = 2, jpim1 zdep_u = 0.5_wp * ( gdept_n(ji,jj,jk) + gdept_n(ji+1,jj,jk) ) zdep_v = 0.5_wp * ( gdept_n(ji,jj,jk) + gdept_n(ji,jj+1,jk) ) ! zkh_u = zk_u(ji,jj) * zdep_u ! k * depth zkh_v = zk_v(ji,jj) * zdep_v ! ! Depth attenuation zda_u = EXP( -2.0_wp*zkh_u ) / ( 1.0_wp + 8.0_wp*zkh_u ) zda_v = EXP( -2.0_wp*zkh_v ) / ( 1.0_wp + 8.0_wp*zkh_v ) ! usd(ji,jj,jk) = zda_u * zu0_sd(ji,jj) * umask(ji,jj,jk) vsd(ji,jj,jk) = zda_v * zv0_sd(ji,jj) * vmask(ji,jj,jk) END DO END DO END DO ELSE IF( nn_sdrift==jp_phillips .OR. nn_sdrift==jp_peakph ) THEN DO jk = 1, jpkm1 DO jj = 2, jpjm1 DO ji = 2, jpim1 zdep_u = 0.5_wp * ( gdept_n(ji,jj,jk) + gdept_n(ji+1,jj,jk) ) zdep_v = 0.5_wp * ( gdept_n(ji,jj,jk) + gdept_n(ji,jj+1,jk) ) ! zkh_u = zk_u(ji,jj) * zdep_u ! k * depth zkh_v = zk_v(ji,jj) * zdep_v ! ! Depth attenuation zda_u = EXP( -2.0_wp*zkh_u ) - SQRT(2.0_wp*rpi*zkh_u) * ERFC(SQRT(2.0_wp*zkh_u)) zda_v = EXP( -2.0_wp*zkh_v ) - SQRT(2.0_wp*rpi*zkh_v) * ERFC(SQRT(2.0_wp*zkh_v)) ! usd(ji,jj,jk) = zda_u * zu0_sd(ji,jj) * umask(ji,jj,jk) vsd(ji,jj,jk) = zda_v * zv0_sd(ji,jj) * vmask(ji,jj,jk) END DO END DO END DO ENDIF CALL lbc_lnk( usd(:,:,:), 'U', vsd(:,:,:), 'V', -1. ) ! ! !== vertical Stokes Drift 3D velocity ==! ! DO jk = 1, jpkm1 ! Horizontal e3*divergence DO jj = 2, jpj DO ji = fs_2, jpi ze3divh(ji,jj,jk) = ( e2u(ji ,jj) * e3u_n(ji ,jj,jk) * usd(ji, jj,jk) & & - e2u(ji-1,jj) * e3u_n(ji-1,jj,jk) * usd(ji-1,jj,jk) & & + e1v(ji,jj ) * e3v_n(ji,jj ,jk) * vsd(ji,jj ,jk) & & - e1v(ji,jj-1) * e3v_n(ji,jj-1,jk) * vsd(ji,jj-1,jk) ) * r1_e12t(ji,jj) END DO END DO END DO ! IF( .NOT. AGRIF_Root() ) THEN IF( nbondi == 1 .OR. nbondi == 2 ) ze3divh(nlci-1, : ,:) = 0._wp ! east IF( nbondi == -1 .OR. nbondi == 2 ) ze3divh( 2 , : ,:) = 0._wp ! west IF( nbondj == 1 .OR. nbondj == 2 ) ze3divh( : ,nlcj-1,:) = 0._wp ! north IF( nbondj == -1 .OR. nbondj == 2 ) ze3divh( : , 2 ,:) = 0._wp ! south ENDIF ! CALL lbc_lnk( ze3divh, 'T', 1. ) ! IF( .NOT. lk_vvl ) THEN ; ik = 1 ! none zero velocity through the sea surface ELSE ; ik = 2 ! w=0 at the surface (set one for all in sbc_wave_init) ENDIF DO jk = jpkm1, ik, -1 ! integrate from the bottom the hor. divergence (NB: at k=jpk w is always zero) wsd(:,:,jk) = wsd(:,:,jk+1) - ze3divh(:,:,jk) END DO #if defined key_bdy IF( lk_bdy ) THEN DO jk = 1, jpkm1 wsd(:,:,jk) = wsd(:,:,jk) * bdytmask(:,:) END DO ENDIF #endif ! !== Horizontal divergence of barotropic Stokes transport ==! div_sd(:,:) = 0._wp DO jk = 1, jpkm1 ! div_sd(:,:) = div_sd(:,:) + ze3divh(:,:,jk) END DO ! CALL iom_put( "ustokes", usd ) CALL iom_put( "vstokes", vsd ) CALL iom_put( "wstokes", wsd ) ! CALL wrk_dealloc( jpi,jpj,jpk, ze3divh ) CALL wrk_dealloc( jpi,jpj, zk_t, zk_u, zk_v, zu0_sd, zv0_sd ) ! END SUBROUTINE sbc_stokes SUBROUTINE sbc_stress( ) !!--------------------------------------------------------------------- !! *** ROUTINE sbc_stress *** !! !! ** Purpose : Updates the ocean momentum modified by waves !! !! ** Method : - Calculate u,v components of stress depending on stress !! model !! - Calculate the stress module !! - The wind module is not modified by waves !! ** action !!--------------------------------------------------------------------- INTEGER :: jj, ji ! dummy loop argument ! IF( ln_tauoc ) THEN utau(:,:) = utau(:,:)*tauoc_wave(:,:) vtau(:,:) = vtau(:,:)*tauoc_wave(:,:) taum(:,:) = taum(:,:)*tauoc_wave(:,:) ENDIF ! IF( ln_tauw ) THEN DO jj = 1, jpjm1 DO ji = 1, jpim1 ! Stress components at u- & v-points utau(ji,jj) = 0.5_wp * ( tauw_x(ji,jj) + tauw_x(ji+1,jj) ) vtau(ji,jj) = 0.5_wp * ( tauw_y(ji,jj) + tauw_y(ji,jj+1) ) ! ! Stress module at t points taum(ji,jj) = SQRT( tauw_x(ji,jj)*tauw_x(ji,jj) + tauw_y(ji,jj)*tauw_y(ji,jj) ) END DO END DO CALL lbc_lnk_multi( utau(:,:), 'U', -1. , vtau(:,:), 'V', -1. , taum(:,: ), 'T', -1. ) ENDIF ! END SUBROUTINE sbc_stress SUBROUTINE sbc_wave( kt ) !!--------------------------------------------------------------------- !! *** ROUTINE sbc_wave *** !! !! ** Purpose : read wave parameters from wave model in netcdf files. !! !! ** Method : - Read namelist namsbc_wave !! - Read Cd_n10 fields in netcdf files !! - Read stokes drift 2d in netcdf files !! - Read wave number in netcdf files !! - Compute 3d stokes drift using Breivik et al.,2014 !! formulation !! ** action !!--------------------------------------------------------------------- INTEGER, INTENT(in ) :: kt ! ocean time step !!--------------------------------------------------------------------- ! IF( ln_cdgw .AND. .NOT. cpl_wdrag ) THEN !== Neutral drag coefficient ==! CALL fld_read( kt, nn_fsbc, sf_cd ) ! read from external forcing cdn_wave(:,:) = sf_cd(1)%fnow(:,:,1) ! check that the drag coefficient contains proper information even if ! the masks do not match - the momentum stress is not masked! WHERE( cdn_wave < 0.0 ) cdn_wave = 1.5e-3 WHERE( cdn_wave > 1.0 ) cdn_wave = 1.5e-3 ENDIF IF( ln_tauoc .AND. .NOT. cpl_tauoc ) THEN !== Wave induced stress ==! CALL fld_read( kt, nn_fsbc, sf_tauoc ) ! read wave norm stress from external forcing tauoc_wave(:,:) = sf_tauoc(1)%fnow(:,:,1) WHERE( tauoc_wave < 0.0 ) tauoc_wave = 1.0 WHERE( tauoc_wave > 100.0 ) tauoc_wave = 1.0 ENDIF IF( ln_tauw .AND. .NOT. cpl_tauw ) THEN !== Wave induced stress ==! CALL fld_read( kt, nn_fsbc, sf_tauw ) ! read ocean stress components from external forcing (T grid) tauw_x(:,:) = sf_tauw(1)%fnow(:,:,1) WHERE( tauw_x < -100.0 ) tauw_x = 0.0 WHERE( tauw_x > 100.0 ) tauw_x = 0.0 tauw_y(:,:) = sf_tauw(2)%fnow(:,:,1) WHERE( tauw_y < -100.0 ) tauw_y = 0.0 WHERE( tauw_y > 100.0 ) tauw_y = 0.0 ENDIF IF( ln_phioc .AND. .NOT. cpl_phioc ) THEN !== Wave to ocean energy ==! CALL fld_read( kt, nn_fsbc, sf_phioc ) ! read wave to ocean energy from external forcing rn_crban(:,:) = 29.0 * sf_phioc(1)%fnow(:,:,1) ! ! Alfa is phioc*sqrt(rau0/zrhoa) : rau0=water density, zhroa= air density WHERE( rn_crban > 1.e8 ) rn_crban = 0.0 !remove first mask mistmatch points, then cap values in case of low friction velocity WHERE( rn_crban < 0.0 ) rn_crban = 0.0 WHERE( rn_crban > 1000.0 ) rn_crban = 1000.0 ENDIF IF( ln_sdw .OR. ln_rough ) THEN !== Computation of the 3d Stokes Drift ==! ! IF( jpfld > 0 ) THEN ! Read from file only if the field is not coupled CALL fld_read( kt, nn_fsbc, sf_sd ) ! read wave parameters from external forcing IF( jp_hsw > 0 ) THEN hsw (:,:) = sf_sd(jp_hsw)%fnow(:,:,1) ! significant wave height WHERE( hsw > 100.0 ) hsw = 0.0 WHERE( hsw < 0.0 ) hsw = 0.0 ENDIF IF( jp_wmp > 0 ) THEN wmp (:,:) = sf_sd(jp_wmp)%fnow(:,:,1) ! wave mean period WHERE( wmp > 100.0 ) wmp = 0.0 WHERE( wmp < 0.0 ) wmp = 0.0 ENDIF IF( jp_wfr > 0 ) THEN wfreq(:,:) = sf_sd(jp_wfr)%fnow(:,:,1) ! Peak wave frequency WHERE( wfreq < 0.0 ) wfreq = 0.0 WHERE( wfreq > 100.0 ) wfreq = 0.0 ENDIF IF( jp_usd > 0 ) THEN ut0sd(:,:) = sf_sd(jp_usd)%fnow(:,:,1) ! 2D zonal Stokes Drift at T point WHERE( ut0sd < -100.0 ) ut0sd = 1.0 WHERE( ut0sd > 100.0 ) ut0sd = 1.0 ENDIF IF( jp_vsd > 0 ) THEN vt0sd(:,:) = sf_sd(jp_vsd)%fnow(:,:,1) ! 2D meridional Stokes Drift at T point WHERE( vt0sd < -100.0 ) vt0sd = 1.0 WHERE( vt0sd > 100.0 ) vt0sd = 1.0 ENDIF ENDIF ENDIF ! IF( ln_sdw ) THEN ! Read also wave number if needed, so that it is available in coupling routines IF( ln_zdfqiao .AND. .NOT.cpl_wnum ) THEN CALL fld_read( kt, nn_fsbc, sf_wn ) ! read wave parameters from external forcing wnum(:,:) = sf_wn(1)%fnow(:,:,1) ENDIF ! !== Computation of the 3d Stokes Drift ==! ! IF( ((nn_sdrift==jp_breivik .OR. nn_sdrift==jp_phillips) .AND. & jp_hsw>0 .AND. jp_wmp>0 .AND. jp_usd>0 .AND. jp_vsd>0) .OR. & (nn_sdrift==jp_peakph .AND. jp_wfr>0 .AND. jp_usd>0 .AND. jp_vsd>0) ) & CALL sbc_stokes() ! Calculate only if required fields are read ! ! In coupled wave model-NEMO case the call is done after coupling ! ENDIF ! END SUBROUTINE sbc_wave SUBROUTINE sbc_wave_init !!--------------------------------------------------------------------- !! *** ROUTINE sbc_wave_init *** !! !! ** Purpose : read wave parameters from wave model in netcdf files. !! !! ** Method : - Read namelist namsbc_wave !! - Read Cd_n10 fields in netcdf files !! - Read stokes drift 2d in netcdf files !! - Read wave number in netcdf files !! - Compute 3d stokes drift using Breivik et al.,2014 !! formulation !! ** action !!--------------------------------------------------------------------- INTEGER :: ierror, ios ! local integer INTEGER :: ifpr !! CHARACTER(len=100) :: cn_dir ! Root directory for location of drag coefficient files TYPE(FLD_N), ALLOCATABLE, DIMENSION(:) :: slf_i, slf_j ! array of namelist informations on the fields to read TYPE(FLD_N) :: sn_cdg, sn_usd, sn_vsd, sn_phioc, & & sn_hsw, sn_wmp, sn_wfr, sn_wnum , & & sn_tauoc, sn_tauwx, sn_tauwy ! information about the fields to be read ! NAMELIST/namsbc_wave/ sn_cdg, cn_dir, sn_usd, sn_vsd, sn_hsw, sn_wmp, sn_wfr, sn_wnum, sn_phioc, & sn_tauoc, sn_tauwx, sn_tauwy, & ln_cdgw, ln_sdw, ln_stcor, ln_phioc, ln_tauoc, ln_tauw, ln_zdfqiao, ln_rough, & nn_sdrift, nn_wmix !!--------------------------------------------------------------------- ! REWIND( numnam_ref ) ! Namelist namsbc_wave in reference namelist : File for drag coeff. from wave model READ ( numnam_ref, namsbc_wave, IOSTAT = ios, ERR = 901) 901 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namsbc_wave in reference namelist', lwp ) REWIND( numnam_cfg ) ! Namelist namsbc_wave in configuration namelist : File for drag coeff. from wave model READ ( numnam_cfg, namsbc_wave, IOSTAT = ios, ERR = 902 ) 902 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namsbc_wave in configuration namelist', lwp ) IF(lwm) WRITE ( numond, namsbc_wave ) ! IF(lwp) THEN !* Parameter print WRITE(numout,*) WRITE(numout,*) 'sbc_wave_init: wave physics' WRITE(numout,*) '~~~~~~~~' WRITE(numout,*) ' Namelist namsbc_wave : set wave physics parameters' WRITE(numout,*) ' Stokes drift corr. to vert. velocity ln_sdw = ', ln_sdw WRITE(numout,*) ' vertical parametrization nn_sdrift = ', nn_sdrift WRITE(numout,*) ' Stokes coriolis term ln_stcor = ', ln_stcor WRITE(numout,*) ' wave modified ocean stress ln_tauoc = ', ln_tauoc WRITE(numout,*) ' wave modified ocean stress components ln_tauw = ', ln_tauw WRITE(numout,*) ' wave to ocean energy ln_phioc = ', ln_phioc WRITE(numout,*) ' vertical mixing parametrization nn_wmix = ', nn_wmix WRITE(numout,*) ' neutral drag coefficient ln_cdgw = ', ln_cdgw WRITE(numout,*) ' wave roughness length modification ln_rough = ', ln_rough WRITE(numout,*) ' Qiao vertical mixing formulation ln_zdfqiao = ', ln_zdfqiao ENDIF IF ( ln_wave ) THEN ! Activated wave physics but no wave physics components activated IF ( .NOT.(ln_cdgw .OR. ln_sdw .OR. ln_tauoc .OR. ln_tauw .OR. ln_stcor .OR. ln_phioc & .OR. ln_rough .OR. ln_zdfqiao) ) THEN CALL ctl_warn( 'Ask for wave coupling but ln_cdgw=F, ln_sdw=F, ln_tauoc=F, ln_tauw=F, ln_stcor=F ', & 'ln_phioc=F, ln_rough=F, ln_zdfqiao=F' ) ELSE IF (ln_stcor .AND. .NOT. ln_sdw) & CALL ctl_stop( 'Stokes-Coriolis term calculated only if activated Stokes Drift ln_sdw=T') IF ( ln_cdgw .AND. .NOT.(nn_drag==jp_ukmo .OR. nn_drag==jp_std .OR. nn_drag==jp_const .OR. nn_drag==jp_mcore) ) & CALL ctl_stop( 'The chosen nn_drag for momentum calculation must be 0, 1, 2, or 3') IF ( ln_cdgw .AND. ln_blk_core .AND. nn_drag==0 ) & CALL ctl_stop( 'The chosen nn_drag for momentum calculation in core forcing must be 1, 2, or 3') IF ( ln_cdgw .AND. ln_flx .AND. nn_drag==3 ) & CALL ctl_stop( 'The chosen nn_drag for momentum calculation in direct forcing must be 0, 1, or 2') IF( ln_phioc .AND. .NOT.(nn_wmix==jp_craigbanner .OR. nn_wmix==jp_janssen) ) & CALL ctl_stop( 'The chosen nn_wmix for wave vertical mixing must be 0, or 1' ) IF( ln_sdw .AND. .NOT.(nn_sdrift==jp_breivik .OR. nn_sdrift==jp_phillips .OR. nn_sdrift==jp_peakph) ) & CALL ctl_stop( 'The chosen nn_sdrift for Stokes drift vertical velocity must be 0, 1, or 2' ) IF( ln_zdfqiao .AND. .NOT.ln_sdw ) & CALL ctl_stop( 'Qiao vertical mixing can not be used without Stokes drift (ln_sdw)' ) IF( ln_tauoc .AND. ln_tauw ) & CALL ctl_stop( 'More than one method for modifying the ocean stress has been selected ', & '(ln_tauoc=.true. and ln_tauw=.true.)' ) IF( ln_tauoc ) & CALL ctl_warn( 'You are subtracting the wave stress to the ocean (ln_tauoc=.true.)' ) IF( ln_tauw ) & CALL ctl_warn( 'The wave modified ocean stress components are used (ln_tauw=.true.) ', & 'This will override any other specification of the ocean stress' ) ENDIF ELSE IF ( ln_cdgw .OR. ln_sdw .OR. ln_tauoc .OR. ln_tauw .OR. ln_stcor .OR. ln_phioc .OR. ln_rough .OR. ln_zdfqiao ) & & CALL ctl_stop( 'Not Activated Wave Module (ln_wave=F) but asked coupling ', & & 'with drag coefficient (ln_cdgw =T) ' , & & 'or Stokes Drift (ln_sdw=T) ' , & & 'or Stokes-Coriolis term (ln_stcor=T)', & & 'or ocean stress modification due to waves (ln_tauoc=T) ', & & 'or ocean stress components from waves (ln_tauw=T) ', & & 'or wave to ocean energy modification (ln_phioc=T) ', & & 'or wave surface roughness (ln_rough=T) ', & & 'or Qiao vertical mixing formulation (ln_zdfqiao=T) ' ) ENDIF ! IF( ln_cdgw ) THEN IF( .NOT. cpl_wdrag ) THEN ALLOCATE( sf_cd(1), STAT=ierror ) !* allocate and fill sf_wave with sn_cdg IF( ierror > 0 ) CALL ctl_stop( 'STOP', 'sbc_wave_init: unable to allocate sf_cd structure' ) ! ALLOCATE( sf_cd(1)%fnow(jpi,jpj,1) ) IF( sn_cdg%ln_tint ) ALLOCATE( sf_cd(1)%fdta(jpi,jpj,1,2) ) CALL fld_fill( sf_cd, (/ sn_cdg /), cn_dir, 'sbc_wave_init', 'read wave input', 'namsbc_wave' ) ENDIF ALLOCATE( cdn_wave(jpi,jpj) ) ENDIF IF( ln_tauoc ) THEN IF( .NOT. cpl_tauoc ) THEN ALLOCATE( sf_tauoc(1), STAT=ierror ) !* allocate and fill sf_wave with sn_tauoc IF( ierror > 0 ) CALL ctl_stop( 'STOP', 'sbc_wave_init: unable to allocate sf_tauoc structure' ) ! ALLOCATE( sf_tauoc(1)%fnow(jpi,jpj,1) ) IF( sn_tauoc%ln_tint ) ALLOCATE( sf_tauoc(1)%fdta(jpi,jpj,1,2) ) CALL fld_fill( sf_tauoc, (/ sn_tauoc /), cn_dir, 'sbc_wave_init', 'read wave input', 'namsbc_wave' ) ENDIF ALLOCATE( tauoc_wave(jpi,jpj) ) ENDIF IF( ln_tauw ) THEN IF( .NOT. cpl_tauw ) THEN ALLOCATE( sf_tauw(2), STAT=ierror ) !* allocate and fill sf_wave with sn_tauwx/y IF( ierror > 0 ) CALL ctl_stop( 'STOP', 'sbc_wave_init: unable to allocate sf_tauw structure' ) ! ALLOCATE( slf_j(2) ) slf_j(1) = sn_tauwx slf_j(2) = sn_tauwy ALLOCATE( sf_tauw(1)%fnow(jpi,jpj,1) ) ALLOCATE( sf_tauw(2)%fnow(jpi,jpj,1) ) IF( slf_j(1)%ln_tint ) ALLOCATE( sf_tauw(1)%fdta(jpi,jpj,1,2) ) IF( slf_j(2)%ln_tint ) ALLOCATE( sf_tauw(2)%fdta(jpi,jpj,1,2) ) CALL fld_fill( sf_tauw, (/ slf_j /), cn_dir, 'sbc_wave_init', 'read wave input', 'namsbc_wave' ) ENDIF ALLOCATE( tauw_x(jpi,jpj) ) ALLOCATE( tauw_y(jpi,jpj) ) ENDIF IF( ln_phioc ) THEN IF( .NOT. cpl_phioc ) THEN ALLOCATE( sf_phioc(1), STAT=ierror ) !* allocate and fill sf_wave with sn_phioc IF( ierror > 0 ) CALL ctl_stop( 'STOP', 'sbc_wave_init: unable to allocate sf_phioc structure' ) ! ALLOCATE( sf_phioc(1)%fnow(jpi,jpj,1) ) IF( sn_phioc%ln_tint ) ALLOCATE( sf_phioc(1)%fdta(jpi,jpj,1,2) ) CALL fld_fill( sf_phioc, (/ sn_phioc /), cn_dir, 'sbc_wave_init', 'read wave input', 'namsbc_wave' ) ENDIF ALLOCATE( rn_crban(jpi,jpj) ) ENDIF ! Find out how many fields have to be read from file if not coupled jpfld=0 jp_usd=0 ; jp_vsd=0 ; jp_hsw=0 ; jp_wmp=0 ; jp_wfr=0 IF( ln_sdw ) THEN IF( .NOT. cpl_sdrft ) THEN jpfld = jpfld + 1 jp_usd = jpfld jpfld = jpfld + 1 jp_vsd = jpfld ENDIF IF( .NOT. cpl_hsig .AND. (nn_sdrift==jp_breivik .OR. nn_sdrift==jp_phillips) ) THEN jpfld = jpfld + 1 jp_hsw = jpfld ENDIF IF( .NOT. cpl_wper .AND. (nn_sdrift==jp_breivik .OR. nn_sdrift==jp_phillips) ) THEN jpfld = jpfld + 1 jp_wmp = jpfld ENDIF IF( .NOT. cpl_wfreq .AND. nn_sdrift==jp_peakph ) THEN jpfld = jpfld + 1 jp_wfr = jpfld ENDIF ENDIF IF( ln_rough .AND. .NOT. cpl_hsig .AND. jp_hsw==0 ) THEN jpfld = jpfld + 1 jp_hsw = jpfld ENDIF ! Read from file only the non-coupled fields IF( jpfld > 0 ) THEN ALLOCATE( slf_i(jpfld) ) IF( jp_usd > 0 ) slf_i(jp_usd) = sn_usd IF( jp_vsd > 0 ) slf_i(jp_vsd) = sn_vsd IF( jp_hsw > 0 ) slf_i(jp_hsw) = sn_hsw IF( jp_wmp > 0 ) slf_i(jp_wmp) = sn_wmp IF( jp_wfr > 0 ) slf_i(jp_wfr) = sn_wfr ALLOCATE( sf_sd(jpfld), STAT=ierror ) !* allocate and fill sf_sd with stokes drift IF( ierror > 0 ) CALL ctl_stop( 'STOP', 'sbc_wave_init: unable to allocate sf_sd structure' ) ! DO ifpr= 1, jpfld ALLOCATE( sf_sd(ifpr)%fnow(jpi,jpj,1) ) IF( slf_i(ifpr)%ln_tint ) ALLOCATE( sf_sd(ifpr)%fdta(jpi,jpj,1,2) ) END DO ! CALL fld_fill( sf_sd, slf_i, cn_dir, 'sbc_wave_init', 'read wave input', 'namsbc_wave' ) ENDIF IF( ln_sdw ) THEN ALLOCATE( usd (jpi,jpj,jpk), vsd (jpi,jpj,jpk), wsd(jpi,jpj,jpk) ) ALLOCATE( wmp (jpi,jpj) ) ALLOCATE( wfreq (jpi,jpj) ) ALLOCATE( ut0sd(jpi,jpj) , vt0sd(jpi,jpj) ) ALLOCATE( div_sd(jpi,jpj) ) ALLOCATE( tsd2d (jpi,jpj) ) usd(:,:,:) = 0._wp vsd(:,:,:) = 0._wp wsd(:,:,:) = 0._wp ! Wave number needed only if ln_zdfqiao=T IF( ln_zdfqiao .AND. .NOT.cpl_wnum ) THEN ALLOCATE( sf_wn(1), STAT=ierror ) !* allocate and fill sf_wave with sn_wnum IF( ierror > 0 ) CALL ctl_stop( 'STOP', 'sbc_wave_init: unable toallocate sf_wn structure' ) ALLOCATE( sf_wn(1)%fnow(jpi,jpj,1) ) IF( sn_wnum%ln_tint ) ALLOCATE( sf_wn(1)%fdta(jpi,jpj,1,2) ) CALL fld_fill( sf_wn, (/ sn_wnum /), cn_dir, 'sbc_wave', 'read wave input', 'namsbc_wave' ) ENDIF ALLOCATE( wnum(jpi,jpj) ) ENDIF IF( ln_sdw .OR. ln_rough ) THEN ALLOCATE( hsw (jpi,jpj) ) ENDIF ! END SUBROUTINE sbc_wave_init !!====================================================================== END MODULE sbcwave