Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Update U10 to be resolved wind; add variable for U10+gusts #440

Merged
merged 3 commits into from
Mar 26, 2024
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
16 changes: 13 additions & 3 deletions cesm/flux_atmocn/shr_flux_mod.F90
Original file line number Diff line number Diff line change
Expand Up @@ -143,7 +143,8 @@ SUBROUTINE flux_atmOcn(logunit, nMax ,zbot ,ubot ,vbot ,thbot , &
& ocn_surface_flux_scheme, &
& add_gusts, &
& duu10n, &
& ugust_out, &
& ugust_out, &
& u10res, &
& ustar_sv ,re_sv ,ssq_sv, &
& missval)

Expand Down Expand Up @@ -194,6 +195,7 @@ SUBROUTINE flux_atmOcn(logunit, nMax ,zbot ,ubot ,vbot ,thbot , &
real(R8),intent(out) :: qref (nMax) ! diag: 2m ref humidity (kg/kg)
real(R8),intent(out) :: duu10n(nMax) ! diag: 10m wind speed squared (m/s)^2
real(R8),intent(out) :: ugust_out(nMax) ! diag: gustiness addition to U10 (m/s)
real(R8),intent(out) :: u10res(nMax) ! diag: gustiness addition to U10 (m/s)

real(R8),intent(out),optional :: ustar_sv(nMax) ! diag: ustar
real(R8),intent(out),optional :: re_sv (nMax) ! diag: sqrt of exchange coefficient (water)
Expand Down Expand Up @@ -243,6 +245,7 @@ SUBROUTINE flux_atmOcn(logunit, nMax ,zbot ,ubot ,vbot ,thbot , &
real(R8) :: cp ! specific heat of moist air
real(R8) :: fac ! vertical interpolation factor
real(R8) :: spval ! local missing value
real(R8) :: wind0 ! resolved large-scale 10m wind (no gust added)
!!++ COARE only
real(R8) :: zo,zot,zoq ! roughness lengths
real(R8) :: hsb,hlb ! sens & lat heat flxs at zbot
Expand Down Expand Up @@ -343,12 +346,15 @@ SUBROUTINE flux_atmOcn(logunit, nMax ,zbot ,ubot ,vbot ,thbot , &

!--- compute some needed quantities ---
if (add_gusts) then
vmag = max(seq_flux_atmocn_minwind, sqrt( (ubot(n)-us(n))**2 + (vbot(n)-vs(n))**2) + ugust(min(rainc(n),6.94444e-4_r8)) )
vmag = max(seq_flux_atmocn_minwind, sqrt( (ubot(n)-us(n))**2 + (vbot(n)-vs(n))**2 + (1.0_R8*ugust(min(rainc(n),6.94444e-4_r8))**2)) )

ugust_out(n) = ugust(min(rainc(n),6.94444e-4_r8))
else
vmag = max(seq_flux_atmocn_minwind, sqrt( (ubot(n)-us(n))**2 + (vbot(n)-vs(n))**2) )
ugust_out(n) = 0.0_r8
end if
wind0 = max(seq_flux_atmocn_minwind, sqrt( (ubot(n)-us(n))**2 + (vbot(n)-vs(n))**2) )


if (use_coldair_outbreak_mod) then
! Cold Air Outbreak Modification:
Expand Down Expand Up @@ -461,6 +467,8 @@ SUBROUTINE flux_atmOcn(logunit, nMax ,zbot ,ubot ,vbot ,thbot , &

duu10n(n) = u10n*u10n ! 10m wind speed squared

u10res(n) = u10n * (wind0/vmag) ! resolved 10m wind

!------------------------------------------------------------
! optional diagnostics, needed for water tracer fluxes (dcn)
!------------------------------------------------------------
Expand All @@ -472,6 +480,7 @@ SUBROUTINE flux_atmOcn(logunit, nMax ,zbot ,ubot ,vbot ,thbot , &
!------------------------------------------------------------
! no valid data here -- out of domain
!------------------------------------------------------------

sen (n) = spval ! sensible heat flux (W/m^2)
lat (n) = spval ! latent heat flux (W/m^2)
lwup (n) = spval ! long-wave upward heat flux (W/m^2)
Expand All @@ -484,7 +493,8 @@ SUBROUTINE flux_atmOcn(logunit, nMax ,zbot ,ubot ,vbot ,thbot , &
tref (n) = spval ! 2m reference height temperature (K)
qref (n) = spval ! 2m reference height humidity (kg/kg)
duu10n(n) = spval ! 10m wind speed squared (m/s)^2
ugust_out(n) = spval ! gustiness addition (m/s)
ugust_out(n) = spval ! gustiness addition (m/s)
u10res(n) = spval ! 10m resolved wind (no gusts) (m/s)

if (present(ustar_sv)) ustar_sv(n) = spval
if (present(re_sv )) re_sv (n) = spval
Expand Down
34 changes: 34 additions & 0 deletions mediator/esmFldsExchange_cesm_mod.F90
Original file line number Diff line number Diff line change
Expand Up @@ -1385,6 +1385,40 @@ subroutine esmFldsExchange_cesm(gcomp, phase, rc)
end if
end if

! ---------------------------------------------------------------------
! to atm: 10 m winds including/excluding gust component
! ---------------------------------------------------------------------
if (phase == 'advertise') then
call addfld_aoflux('So_u10withGust')
call addfld_to(compatm, 'So_u10withGust')
else
if ( fldchk(is_local%wrap%FBexp(compatm), 'So_u10withGust', rc=rc)) then
if (fldchk(is_local%wrap%FBMed_aoflux_o, 'So_u10withGust', rc=rc)) then
if (trim(is_local%wrap%aoflux_grid) == 'ogrid') then
call addmap_aoflux('So_u10withGust', compatm, mapconsf, 'ofrac', ocn2atm_map)
end if
call addmrg_to(compatm , 'So_u10withGust', &
mrg_from=compmed, mrg_fld='So_u10withGust', mrg_type='merge', mrg_fracname='ofrac')
end if
end if
end if

if (phase == 'advertise') then
call addfld_aoflux('So_u10res')
call addfld_to(compatm, 'So_u10res')
else
if ( fldchk(is_local%wrap%FBexp(compatm), 'So_u10res', rc=rc)) then
if (fldchk(is_local%wrap%FBMed_aoflux_o, 'So_u10res', rc=rc)) then
if (trim(is_local%wrap%aoflux_grid) == 'ogrid') then
call addmap_aoflux('So_u10res', compatm, mapconsf, 'ofrac', ocn2atm_map)
end if
call addmrg_to(compatm , 'So_u10res', &
mrg_from=compmed, mrg_fld='So_u10res', mrg_type='merge', mrg_fracname='ofrac')
end if
end if
end if


! ---------------------------------------------------------------------
! to atm: surface snow depth from ice (needed for cam)
! to atm: mean ice volume per unit area from ice
Expand Down
8 changes: 8 additions & 0 deletions mediator/fd_cesm.yaml
Original file line number Diff line number Diff line change
Expand Up @@ -491,6 +491,14 @@
canonical_units: m/s
description: atmosphere import
#
- standard_name: So_u10withGust
canonical_units: m/s
description: atmosphere import
#
- standard_name: So_u10res
canonical_units: m/s
description: atmosphere import
#
#-----------------------------------
# section: land-ice export
# Note that the fields sent from glc->med do NOT have elevation classes,
Expand Down
15 changes: 13 additions & 2 deletions mediator/med_phases_aofluxes_mod.F90
Original file line number Diff line number Diff line change
Expand Up @@ -149,6 +149,8 @@ module med_phases_aofluxes_mod
real(R8) , pointer :: u10 (:) => null() ! diagnostic: 10m wind speed
real(R8) , pointer :: duu10n (:) => null() ! diagnostic: 10m wind speed squared
real(R8) , pointer :: ugust_out (:) => null() ! diagnostic: gust wind added
real(R8) , pointer :: u10_withGust(:) => null() ! diagnostic: gust wind added
real(R8) , pointer :: u10res (:) => null() ! diagnostic: no gust wind added
real(R8) , pointer :: ustar (:) => null() ! saved ustar
real(R8) , pointer :: re (:) => null() ! saved re
real(R8) , pointer :: ssq (:) => null() ! saved sq
Expand Down Expand Up @@ -1075,6 +1077,7 @@ subroutine med_aofluxes_update(gcomp, aoflux_in, aoflux_out, rc)
add_gusts=add_gusts, &
duu10n=aoflux_out%duu10n, &
ugust_out = aoflux_out%ugust_out, &
u10res = aoflux_out%u10res, &
ustar_sv=aoflux_out%ustar, re_sv=aoflux_out%re, ssq_sv=aoflux_out%ssq, &
missval=0.0_r8)

Expand Down Expand Up @@ -1108,8 +1111,9 @@ subroutine med_aofluxes_update(gcomp, aoflux_in, aoflux_out, rc)
#endif

do n = 1,aoflux_in%lsize
if (aoflux_in%mask(n) /= 0) then
aoflux_out%u10(n) = sqrt(aoflux_out%duu10n(n))
if (aoflux_in%mask(n) /= 0) then
aoflux_out%u10(n) = aoflux_out%u10res(n)
aoflux_out%u10_withGust(n) = sqrt(aoflux_out%duu10n(n))
end if
enddo

Expand Down Expand Up @@ -1712,6 +1716,13 @@ subroutine set_aoflux_out_pointers(fldbun, lsize, aoflux_out, xgrid, rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
call fldbun_getfldptr(fldbun, 'So_duu10n', aoflux_out%duu10n, xgrid=xgrid, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return

call fldbun_getfldptr(fldbun, 'So_ugustOut', aoflux_out%ugust_out, xgrid=xgrid, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
call fldbun_getfldptr(fldbun, 'So_u10withGust', aoflux_out%u10_withGust, xgrid=xgrid, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
call fldbun_getfldptr(fldbun, 'So_u10res', aoflux_out%u10res, xgrid=xgrid, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
call fldbun_getfldptr(fldbun, 'Faox_taux', aoflux_out%taux, xgrid=xgrid, rc=rc)
if (chkerr(rc,__LINE__,u_FILE_u)) return
call fldbun_getfldptr(fldbun, 'Faox_tauy', aoflux_out%tauy, xgrid=xgrid, rc=rc)
Expand Down
Loading