Skip to content

Commit

Permalink
Merge branch 'jackreeveseyre/cime/alternate_surface_flux' (PR #2972)
Browse files Browse the repository at this point in the history
Add UA and COARE ocean surface flux schemes

This branch adds two algorithms to calculate ocean surface fluxes
(wind stress, latent and sensible heat fluxes) in CIME. The new
algorithms are selected by a new CIME namelist variable:

ocn_surface_flux_scheme = 0 [default] -- original scheme as used in E3SMv1;
= 1 -- COARE v3.0 scheme (from Fairall et al. 2003);
= 2 -- University of Arizona (UA) scheme (from Zeng et al., 1998).

The desired option should be added to user_nl_cpl.

[BFB]
[FCC]
[NML]

EW-108, EW-127

* jackreeveseyre/cime/alternate_surface_flux:
  Reverts change to parameter specification in Large&Pond scheme.
  Adds COAREv3.0 ocean surface flux scheme to CIME.
  Fixes bug in calling routine for default ocean surface flux scheme.
  Adds UA ocean surface flux scheme to CIME.
  • Loading branch information
singhbalwinder committed Jun 20, 2019
2 parents fcea3e2 + 5b0241c commit 1aa18d3
Show file tree
Hide file tree
Showing 3 changed files with 77 additions and 5 deletions.
17 changes: 16 additions & 1 deletion cime_config/namelist_definition_drv.xml
Original file line number Diff line number Diff line change
Expand Up @@ -696,7 +696,22 @@
<values>
<value>.false.</value>
</values>
</entry>
</entry>

<entry id="ocn_surface_flux_scheme">
<type>integer</type>
<category>control</category>
<group>seq_infodata_inparm</group>
<desc>
0: standard surface flux calculation as in E3SMv1
1: COAREv3.0 flux computation (Fairall et al., 2003)
2: University of Arizona algorithm (Zeng et al., 1998)
default: 0
</desc>
<values>
<value>0</value>
</values>
</entry>

<entry id="flux_convergence">
<type>real</type>
Expand Down
48 changes: 45 additions & 3 deletions main/seq_flux_mct.F90
Original file line number Diff line number Diff line change
Expand Up @@ -2,7 +2,7 @@ module seq_flux_mct

use shr_kind_mod, only: r8 => shr_kind_r8, in=>shr_kind_in
use shr_sys_mod, only: shr_sys_abort
use shr_flux_mod, only: shr_flux_atmocn, shr_flux_atmocn_diurnal, shr_flux_adjust_constants
use shr_flux_mod, only: shr_flux_atmocn, shr_flux_atmocn_ua, shr_flux_atmocn_diurnal, shr_flux_adjust_constants
use shr_orb_mod, only: shr_orb_params, shr_orb_cosz, shr_orb_decl
use shr_mct_mod, only: shr_mct_queryConfigFile, shr_mct_sMatReaddnc

Expand Down Expand Up @@ -54,6 +54,7 @@ module seq_flux_mct
real(r8), allocatable :: roce_18O (:) ! ocn H218O ratio
real(r8), allocatable :: dens (:) ! atm density
real(r8), allocatable :: tbot (:) ! atm bottom surface T
real(r8), allocatable :: pslv (:) ! sea level pressure (Pa)
real(r8), allocatable :: sen (:) ! heat flux: sensible
real(r8), allocatable :: lat (:) ! heat flux: latent
real(r8), allocatable :: lwup (:) ! lwup over ocean
Expand Down Expand Up @@ -118,6 +119,7 @@ module seq_flux_mct
integer :: index_a2x_Sa_shum_HDO
integer :: index_a2x_Sa_shum_18O
integer :: index_a2x_Sa_dens
integer :: index_a2x_Sa_pslv
integer :: index_a2x_Faxa_swndr
integer :: index_a2x_Faxa_swndf
integer :: index_a2x_Faxa_swvdr
Expand Down Expand Up @@ -248,6 +250,9 @@ subroutine seq_flux_init_mct(comp, fractions)
allocate(tbot(nloc),stat=ier)
if(ier/=0) call mct_die(subName,'allocate tbot',ier)
tbot = 0.0_r8
allocate(pslv(nloc),stat=ier)
if(ier/=0) call mct_die(subName,'allocate pslv',ier)
pslv = 0.0_r8
allocate(ustar(nloc),stat=ier)
if(ier/=0) call mct_die(subName,'allocate ustar',ier)
ustar = 0.0_r8
Expand Down Expand Up @@ -605,6 +610,8 @@ subroutine seq_flux_initexch_mct(atm, ocn, mpicom_cplid, cplid)
if(ier/=0) call mct_die(subName,'allocate dens',ier)
allocate(tbot(nloc_a2o),stat=ier)
if(ier/=0) call mct_die(subName,'allocate tbot',ier)
allocate(pslv(nloc_a2o),stat=ier)
if(ier/=0) call mct_die(subName,'allocate pslv',ier)
allocate(ustar(nloc_a2o),stat=ier)
if(ier/=0) call mct_die(subName,'allocate ustar',ier)
allocate(re(nloc_a2o), stat=ier)
Expand Down Expand Up @@ -906,6 +913,7 @@ subroutine seq_flux_atmocnexch_mct( infodata, atm, ocn, fractions_a, fractions_o
logical :: read_restart ! .true. => model starting from restart
logical :: ocn_prognostic ! .true. => ocn is prognostic
logical :: flux_diurnal ! .true. => turn on diurnal cycle in atm/ocn fluxes
integer :: ocn_surface_flux_scheme ! 0: E3SMv1 1: COARE 2: UA
logical :: cold_start ! .true. to initialize internal fields in shr_flux diurnal
character(len=256) :: fldlist ! subset of xao fields
!
Expand All @@ -931,7 +939,8 @@ subroutine seq_flux_atmocnexch_mct( infodata, atm, ocn, fractions_a, fractions_o
atm_nx=atm_nx, atm_ny=atm_ny, &
ocn_nx=ocn_nx, ocn_ny=ocn_ny, &
ocn_prognostic=ocn_prognostic, &
flux_diurnal=flux_diurnal)
flux_diurnal=flux_diurnal, &
ocn_surface_flux_scheme=ocn_surface_flux_scheme)

cold_start = .false. ! use restart data or data from last timestep

Expand All @@ -958,6 +967,7 @@ subroutine seq_flux_atmocnexch_mct( infodata, atm, ocn, fractions_a, fractions_o
roce_18O(n) = 1.0_r8 ! H218O ratio ~ mol/mol
dens(n) = 1.0_r8 ! atm density ~ kg/m^3
tbot(n) = 300.0_r8 ! atm temperature ~ Kelvin
pslv(n) = 101300.0_r8 ! sea level pressure ~ Pa
enddo
else

Expand Down Expand Up @@ -990,6 +1000,7 @@ subroutine seq_flux_atmocnexch_mct( infodata, atm, ocn, fractions_a, fractions_o
shum_18O(n) = a2x_e%rAttr(index_a2x_Sa_shum_18O,ia)
dens(n) = a2x_e%rAttr(index_a2x_Sa_dens,ia)
tbot(n) = a2x_e%rAttr(index_a2x_Sa_tbot,ia)
pslv(n) = a2x_e%rAttr(index_a2x_Sa_pslv,ia)
tocn(n) = o2x_e%rAttr(index_o2x_So_t ,io)
uocn(n) = o2x_e%rAttr(index_o2x_So_u ,io)
vocn(n) = o2x_e%rAttr(index_o2x_So_v ,io)
Expand All @@ -1002,26 +1013,39 @@ subroutine seq_flux_atmocnexch_mct( infodata, atm, ocn, fractions_a, fractions_o
end if

if (flux_diurnal) then
if (ocn_surface_flux_scheme.eq.2) then
call shr_sys_abort(trim(subname)//' ERROR cannot use flux_diurnal with UA flux scheme')
endif
call shr_flux_atmocn_diurnal (nloc_a2o , zbot , ubot, vbot, thbot, &
shum , shum_16O , shum_HDO, shum_18O, dens , tbot, uocn, vocn , &
tocn , emask, sen , lat , lwup , &
roce_16O, roce_HDO, roce_18O, &
evap , evap_16O, evap_HDO, evap_18O, taux , tauy, tref, qref , &
uGust, lwdn , swdn , swup, prec, &
fswpen, ocnsal, ocn_prognostic, flux_diurnal, &
ocn_surface_flux_scheme, &
lats , lons , warm , salt , speed, regime, &
warmMax, windMax, qSolAvg, windAvg, &
warmMaxInc, windMaxInc, qSolInc, windInc, nInc, &
tbulk, tskin, tskin_day, tskin_night, &
cskin, cskin_night, tod, dt, &
duu10n,ustar, re , ssq , missval = 0.0_r8, &
cold_start=cold_start)
else if (ocn_surface_flux_scheme.eq.2) then
call shr_flux_atmOcn_UA(nloc_a2o , zbot , ubot, vbot, thbot, &
shum , shum_16O , shum_HDO, shum_18O, dens , tbot, pslv, &
uocn, vocn , tocn , emask, sen , lat , lwup , &
roce_16O, roce_HDO, roce_18O, &
evap , evap_16O, evap_HDO, evap_18O, taux, tauy, tref, qref , &
duu10n,ustar, re , ssq , missval = 0.0_r8 )
else

call shr_flux_atmocn (nloc_a2o , zbot , ubot, vbot, thbot, &
shum , shum_16O , shum_HDO, shum_18O, dens , tbot, uocn, vocn , &
tocn , emask, sen , lat , lwup , &
roce_16O, roce_HDO, roce_18O, &
evap , evap_16O, evap_HDO, evap_18O, taux, tauy, tref, qref , &
ocn_surface_flux_scheme, &
duu10n,ustar, re , ssq , missval = 0.0_r8 )
endif

Expand Down Expand Up @@ -1187,6 +1211,7 @@ subroutine seq_flux_atmocn_mct(infodata, tod, dt, a2x, o2x, xao)
logical :: read_restart ! .true. => continue run
logical :: ocn_prognostic ! .true. => ocn is prognostic
logical :: flux_diurnal ! .true. => turn on diurnal cycle in atm/ocn fluxes
integer :: ocn_surface_flux_scheme ! 0: E3SMv1 1: COARE 2: UA
real(r8) :: flux_convergence ! convergence criteria for imlicit flux computation
integer(in) :: flux_max_iteration ! maximum number of iterations for convergence
logical :: coldair_outbreak_mod ! cold air outbreak adjustment (Mahrt & Sun 1995,MWR)
Expand All @@ -1202,7 +1227,8 @@ subroutine seq_flux_atmocn_mct(infodata, tod, dt, a2x, o2x, xao)
flux_albav=flux_albav, &
dead_comps=dead_comps, &
ocn_prognostic=ocn_prognostic, &
flux_diurnal=flux_diurnal)
flux_diurnal=flux_diurnal, &
ocn_surface_flux_scheme=ocn_surface_flux_scheme)

cold_start = .false. ! use restart data or data from last timestep

Expand Down Expand Up @@ -1256,6 +1282,7 @@ subroutine seq_flux_atmocn_mct(infodata, tod, dt, a2x, o2x, xao)
index_a2x_Sa_u = mct_aVect_indexRA(a2x,'Sa_u')
index_a2x_Sa_v = mct_aVect_indexRA(a2x,'Sa_v')
index_a2x_Sa_tbot = mct_aVect_indexRA(a2x,'Sa_tbot')
index_a2x_Sa_pslv = mct_aVect_indexRA(a2x,'Sa_pslv')
index_a2x_Sa_ptem = mct_aVect_indexRA(a2x,'Sa_ptem')
index_a2x_Sa_shum = mct_aVect_indexRA(a2x,'Sa_shum')
index_a2x_Sa_shum_16O = mct_aVect_indexRA(a2x,'Sa_shum_16O', perrWith='quiet')
Expand Down Expand Up @@ -1318,6 +1345,7 @@ subroutine seq_flux_atmocn_mct(infodata, tod, dt, a2x, o2x, xao)
roce_18O(n) = 1.0_r8 ! H218O surface ratio ~ mol/mol
dens(n) = 1.0_r8 ! atm density ~ kg/m^3
tbot(n) = 300.0_r8 ! atm temperature ~ Kelvin
pslv(n) = 101300.0_r8 ! sea level pressure ~ Pa
uGust(n)= 0.0_r8
lwdn(n) = 0.0_r8
prec(n) = 0.0_r8
Expand Down Expand Up @@ -1360,6 +1388,7 @@ subroutine seq_flux_atmocn_mct(infodata, tod, dt, a2x, o2x, xao)
if ( index_a2x_Sa_shum_18O /= 0 ) shum_18O(n) = a2x%rAttr(index_a2x_Sa_shum_18O,n)
dens(n) = a2x%rAttr(index_a2x_Sa_dens,n)
tbot(n) = a2x%rAttr(index_a2x_Sa_tbot,n)
pslv(n) = a2x%rAttr(index_a2x_Sa_pslv,n)
tocn(n) = o2x%rAttr(index_o2x_So_t ,n)
uocn(n) = o2x%rAttr(index_o2x_So_u ,n)
vocn(n) = o2x%rAttr(index_o2x_So_v ,n)
Expand Down Expand Up @@ -1408,13 +1437,18 @@ subroutine seq_flux_atmocn_mct(infodata, tod, dt, a2x, o2x, xao)
end if

if (flux_diurnal) then
if (ocn_surface_flux_scheme.eq.2) then
call shr_sys_abort(trim(subname)//' ERROR cannot use flux_diurnal with UA flux scheme')
endif

call shr_flux_atmocn_diurnal (nloc , zbot , ubot, vbot, thbot, &
shum , shum_16O , shum_HDO, shum_18O, dens , tbot, uocn, vocn , &
tocn , emask, sen , lat , lwup , &
roce_16O, roce_HDO, roce_18O, &
evap , evap_16O, evap_HDO, evap_18O, taux , tauy, tref, qref , &
uGust, lwdn , swdn , swup, prec, &
fswpen, ocnsal, ocn_prognostic, flux_diurnal, &
ocn_surface_flux_scheme, &
lats, lons , warm , salt , speed, regime, &
warmMax, windMax, qSolAvg, windAvg, &
warmMaxInc, windMaxInc, qSolInc, windInc, nInc, &
Expand All @@ -1425,12 +1459,20 @@ subroutine seq_flux_atmocn_mct(infodata, tod, dt, a2x, o2x, xao)
!consistent with mrgx2a fraction
!duu10n,ustar, re , ssq, missval = 0.0_r8 )
cold_start=cold_start)
else if (ocn_surface_flux_scheme.eq.2) then
call shr_flux_atmOcn_UA(nloc , zbot , ubot, vbot, thbot, &
shum , shum_16O , shum_HDO, shum_18O, dens , tbot, pslv, &
uocn, vocn , tocn , emask, sen , lat , lwup , &
roce_16O, roce_HDO, roce_18O, &
evap , evap_16O, evap_HDO, evap_18O, taux , tauy, tref, qref , &
duu10n,ustar, re , ssq)
else
call shr_flux_atmocn (nloc , zbot , ubot, vbot, thbot, &
shum , shum_16O , shum_HDO, shum_18O, dens , tbot, uocn, vocn , &
tocn , emask, sen , lat , lwup , &
roce_16O, roce_HDO, roce_18O, &
evap , evap_16O, evap_HDO, evap_18O, taux , tauy, tref, qref , &
ocn_surface_flux_scheme, &
duu10n,ustar, re , ssq)
!missval should not be needed if flux calc
!consistent with mrgx2a fraction
Expand Down
17 changes: 16 additions & 1 deletion shr/seq_infodata_mod.F90
Original file line number Diff line number Diff line change
Expand Up @@ -105,6 +105,7 @@ MODULE seq_infodata_mod
character(SHR_KIND_CL) :: flux_epbal ! selects E,P,R adjustment technique
logical :: flux_albav ! T => no diurnal cycle in ocn albedos
logical :: flux_diurnal ! T => diurnal cycle in atm/ocn fluxes
integer :: ocn_surface_flux_scheme ! 0: E3SMv1 1: COARE 2: UA
logical :: coldair_outbreak_mod ! (Mahrt & Sun 1995,MWR)
real(SHR_KIND_R8) :: flux_convergence ! atmocn flux calc convergence value
integer :: flux_max_iteration ! max number of iterations of atmocn flux loop
Expand Down Expand Up @@ -351,6 +352,7 @@ SUBROUTINE seq_infodata_Init( infodata, nmlfile, ID, pioid, cpl_tag)
character(SHR_KIND_CL) :: flux_epbal ! selects E,P,R adjustment technique
logical :: flux_albav ! T => no diurnal cycle in ocn albedos
logical :: flux_diurnal ! T => diurnal cycle in atm/ocn fluxes
integer :: ocn_surface_flux_scheme ! 0: E3SMv1 1: COARE 2: UA
logical :: coldair_outbreak_mod ! (Mahrt & Sun 1995,MWR)
real(SHR_KIND_R8) :: flux_convergence ! atmocn flux calc convergence value
integer :: flux_max_iteration ! max number of iterations of atmocn flux loop
Expand Down Expand Up @@ -424,6 +426,7 @@ SUBROUTINE seq_infodata_Init( infodata, nmlfile, ID, pioid, cpl_tag)
restart_pfile, restart_file, run_barriers, &
single_column, scmlat, force_stop_at, &
scmlon, logFilePostFix, outPathRoot, flux_diurnal,&
ocn_surface_flux_scheme, &
coldair_outbreak_mod, &
flux_convergence, flux_max_iteration, &
perpetual, perpetual_ymd, flux_epbal, flux_albav, &
Expand Down Expand Up @@ -504,6 +507,7 @@ SUBROUTINE seq_infodata_Init( infodata, nmlfile, ID, pioid, cpl_tag)
flux_epbal = 'off'
flux_albav = .false.
flux_diurnal = .false.
ocn_surface_flux_scheme = 0
coldair_outbreak_mod = .false.
flux_convergence = 0.0_SHR_KIND_R8
flux_max_iteration = 2
Expand Down Expand Up @@ -631,6 +635,7 @@ SUBROUTINE seq_infodata_Init( infodata, nmlfile, ID, pioid, cpl_tag)
infodata%flux_epbal = flux_epbal
infodata%flux_albav = flux_albav
infodata%flux_diurnal = flux_diurnal
infodata%ocn_surface_flux_scheme = ocn_surface_flux_scheme
infodata%flux_convergence = flux_convergence
infodata%coldair_outbreak_mod = coldair_outbreak_mod
infodata%flux_max_iteration = flux_max_iteration
Expand Down Expand Up @@ -966,6 +971,7 @@ SUBROUTINE seq_infodata_GetData_explicit( infodata, cime_model, case_name, case_
nextsw_cday, precip_fact, flux_epbal, flux_albav, &
glc_g2lupdate, atm_aero, run_barriers, esmf_map_flag, &
do_budgets, do_histinit, drv_threading, flux_diurnal, &
ocn_surface_flux_scheme, &
coldair_outbreak_mod, &
flux_convergence, flux_max_iteration, &
budget_inst, budget_daily, budget_month, wall_time_limit, &
Expand Down Expand Up @@ -1036,6 +1042,7 @@ SUBROUTINE seq_infodata_GetData_explicit( infodata, cime_model, case_name, case_
character(len=*), optional, intent(OUT) :: flux_epbal ! selects E,P,R adjustment technique
logical, optional, intent(OUT) :: flux_albav ! T => no diurnal cycle in ocn albedos
logical, optional, intent(OUT) :: flux_diurnal ! T => diurnal cycle in atm/ocn flux
integer, optional, intent(OUT) :: ocn_surface_flux_scheme ! 0: E3SMv1 1: COARE 2: UA
real(SHR_KIND_R8), optional, intent(out) :: flux_convergence ! atmocn flux calc convergence value
logical, optional, intent(out) :: coldair_outbreak_mod ! (Mahrt & Sun 1995, MWR)
integer, optional, intent(OUT) :: flux_max_iteration ! max number of iterations of atmocn flux loop
Expand Down Expand Up @@ -1213,6 +1220,8 @@ SUBROUTINE seq_infodata_GetData_explicit( infodata, cime_model, case_name, case_
if ( present(flux_epbal) ) flux_epbal = infodata%flux_epbal
if ( present(flux_albav) ) flux_albav = infodata%flux_albav
if ( present(flux_diurnal) ) flux_diurnal = infodata%flux_diurnal
if ( present(ocn_surface_flux_scheme) ) ocn_surface_flux_scheme = &
infodata%ocn_surface_flux_scheme
if ( present(coldair_outbreak_mod)) coldair_outbreak_mod = infodata%coldair_outbreak_mod
if ( present(flux_convergence)) flux_convergence = infodata%flux_convergence
if ( present(flux_max_iteration)) flux_max_iteration = infodata%flux_max_iteration
Expand Down Expand Up @@ -1497,6 +1506,7 @@ SUBROUTINE seq_infodata_PutData_explicit( infodata, cime_model, case_name, case_
nextsw_cday, precip_fact, flux_epbal, flux_albav, &
glc_g2lupdate, atm_aero, esmf_map_flag, wall_time_limit, &
do_budgets, do_histinit, drv_threading, flux_diurnal, &
ocn_surface_flux_scheme, &
coldair_outbreak_mod, &
flux_convergence, flux_max_iteration, &
budget_inst, budget_daily, budget_month, force_stop_at, &
Expand Down Expand Up @@ -1566,6 +1576,7 @@ SUBROUTINE seq_infodata_PutData_explicit( infodata, cime_model, case_name, case_
character(len=*), optional, intent(IN) :: flux_epbal ! selects E,P,R adjustment technique
logical, optional, intent(IN) :: flux_albav ! T => no diurnal cycle in ocn albedos
logical, optional, intent(IN) :: flux_diurnal ! T => diurnal cycle in atm/ocn flux
integer, optional, intent(IN) :: ocn_surface_flux_scheme ! 0: E3SMv1 1: COARE 2:UA
logical, optional, intent(in) :: coldair_outbreak_mod
real(SHR_KIND_R8), optional, intent(IN) :: flux_convergence ! atmocn flux calc convergence value
integer, optional, intent(IN) :: flux_max_iteration ! max number of iterations of atmocn flux loop
Expand Down Expand Up @@ -1740,6 +1751,8 @@ SUBROUTINE seq_infodata_PutData_explicit( infodata, cime_model, case_name, case_
if ( present(flux_epbal) ) infodata%flux_epbal = flux_epbal
if ( present(flux_albav) ) infodata%flux_albav = flux_albav
if ( present(flux_diurnal) ) infodata%flux_diurnal = flux_diurnal
if ( present(ocn_surface_flux_scheme) ) infodata%ocn_surface_flux_scheme = &
ocn_surface_flux_scheme
if ( present(coldair_outbreak_mod) ) infodata%coldair_outbreak_mod = coldair_outbreak_mod
if ( present(flux_convergence)) infodata%flux_convergence = flux_convergence
if ( present(flux_max_iteration)) infodata%flux_max_iteration = flux_max_iteration
Expand Down Expand Up @@ -2039,7 +2052,8 @@ subroutine seq_infodata_bcast(infodata,mpicom)
call shr_mpi_bcast(infodata%flux_epbal, mpicom)
call shr_mpi_bcast(infodata%flux_albav, mpicom)
call shr_mpi_bcast(infodata%flux_diurnal, mpicom)
call shr_mpi_bcast(infodata%coldair_outbreak_mod, mpicom)
call shr_mpi_bcast(infodata%ocn_surface_flux_scheme, mpicom)
call shr_mpi_bcast(infodata%coldair_outbreak_mod, mpicom)
call shr_mpi_bcast(infodata%flux_convergence, mpicom)
call shr_mpi_bcast(infodata%flux_max_iteration, mpicom)
call shr_mpi_bcast(infodata%glc_renormalize_smb, mpicom)
Expand Down Expand Up @@ -2735,6 +2749,7 @@ SUBROUTINE seq_infodata_print( infodata )
write(logunit,F0A) subname,'flux_epbal = ', trim(infodata%flux_epbal)
write(logunit,F0L) subname,'flux_albav = ', infodata%flux_albav
write(logunit,F0L) subname,'flux_diurnal = ', infodata%flux_diurnal
write(logunit,F0L) subname,'ocn_surface_flux_scheme = ', infodata%ocn_surface_flux_scheme
write(logunit,F0L) subname,'coldair_outbreak_mod = ', infodata%coldair_outbreak_mod
write(logunit,F0R) subname,'flux_convergence = ', infodata%flux_convergence
write(logunit,F0I) subname,'flux_max_iteration = ', infodata%flux_max_iteration
Expand Down

0 comments on commit 1aa18d3

Please sign in to comment.