CCPP Scientific Documentation
v4.0
CIRES Unified Gravity Wave Physics Module

The physics of NGWs in the UGWP framework (Yudin et al. 2018 [179]) is represented by four GW-solvers, which have been introduced in Lindzen (1981) [108], Hines (1997) [81], Alexander and Dunkerton (1999) [2], and Scinocca (2003) [155]. These GW solvers have been modified by considring the background dissipation of temperature and winds. This feature, which is important in the mesosphere and thermosphere for WAM applications, considers appropriate scale-dependent dissipation of waves near the model top lid providing momentum and energy conservation in the vertical column physics (Shaw and Shepherd 2009 [156]). In the UGWP v0, a modification of Scinocca (2003) [155] scheme for NGWs with non-hydrostatic and rotational effects for GW propagations and background dissipation is represented by the subroutine fv3_ugwp_solv2_v0(). In the next release of UGWP, additional GW-solvers will be implemented along with physics-based triggering of waves and stochastic approaches for selection of GW modes characterized by horizontal phase velocities, azimuthal directions and magnitude of the vertical momentum flux (VMF). More...

Detailed Description

In UGWP v0, the specification for the VMF function is adopted from the Goddard Earth Observing System Version 5 (GEOS-5) global atmosphere model of the National Aeronautic and Space Administration (NASA) Goddard Space Flight Center (GSFC) Global Modeling and Assimilation Office (GMAO), as described in Molod et al. (2015) [126] and employed in the Modern-Era Retrospective analysis for Research and Applications (MERRA)-2 reanalysis (Gelaro et al., 2017 [62]). The Fortran subroutine slat_geos5_tamp() describes the latitudinal shape of VMF-function as displayed in Figure 3 of Molod et al. (2015) [126]. It shows that the enhanced values of VMF in the equatorial region gives opportunity to simulate the QBO-like oscillations in the equatorial zonal winds and lead to more realistic simulations of the equatorial dynamics in GEOS-5 operational and MERRA-2 reanalysis products.

The UGWP v0 code has been enhanced by scientists from NOAA's EMC to modulate the zonal mean NGW forcing by three-dimensional distributions of the total precipitation (as a proxy for the excitation of NGWs by convection) and vertically-integrated (surface-tropopause) turbulent kinetic energy. The vertically extented configuration of the UFS weather model is being tuned using reanalysis products from MERRA-2 and the European Centre for Medium-Range Weather Forecasts (ERA-5), along with temperature, ozone, and water vapor observations of middle atmosphere satellite missions. The verification scores with updated NGW forcing, as reported elsewhere by EMC researchers, display noticeable improvements into the mesosphere.

Argument Table

cires_ugwp_run argument table
local_namestandard_namelong_nameunitstypedimensionskindintentoptional
do_ugwpdo_ugwpflag to activate CIRES UGWPflaglogical()inFalse
mempi_rankMPI rank of current processindexinteger()inFalse
mastermpi_rootMPI rank of master processindexinteger()inFalse
imhorizontal_loop_extenthorizontalcountinteger()inFalse
levsvertical_dimensionnumber of vertical levelscountinteger()inFalse
ntracnumber_of_tracersnumber of tracerscountinteger()inFalse
dtptime_step_for_physicsphysics timestepsreal()kind_physinFalse
kdtindex_of_time_stepcurrent forecast iterationindexinteger()inFalse
lonrnumber_of_equatorial_longitude_pointsnumber of global points in x-dir (i) along the equatorcountinteger()inFalse
oroorographyorographymreal(horizontal_dimension)kind_physinFalse
oro_uforography_unfilteredunfiltered orographymreal(horizontal_dimension)kind_physinFalse
hprimestandard_deviation_of_subgrid_orographystandard deviation of subgrid orographymreal(horizontal_dimension)kind_physinFalse
nmtvrnumber_of_statistical_measures_of_subgrid_orographynumber of topographic variables in GWDcountinteger()inFalse
occonvexity_of_subgrid_orographyconvexity of subgrid orographynonereal(horizontal_dimension)kind_physinFalse
oa4asymmetry_of_subgrid_orographyasymmetry of subgrid orographynonereal(horizontal_dimension, 4)kind_physinFalse
clxfraction_of_grid_box_with_subgrid_orography_higher_than_critical_heighthorizontal fraction of grid box covered by subgrid orography higher than critical heightfracreal(horizontal_dimension, 4)kind_physinFalse
thetaangle_from_east_of_maximum_subgrid_orographic_variationsangle with_respect to east of maximum subgrid orographic variationsdegreesreal(horizontal_dimension)kind_physinFalse
sigmaslope_of_subgrid_orographyslope of subgrid orographynonereal(horizontal_dimension)kind_physinFalse
gammaanisotropy_of_subgrid_orographyanisotropy of subgrid orographynonereal(horizontal_dimension)kind_physinFalse
elvmaxmaximum_subgrid_orographymaximum of subgrid orographymreal(horizontal_dimension)kind_physinFalse
do_tofdturb_oro_form_drag_flagflag for turbulent orographic form dragflaglogical()inFalse
ldiag_ugwpdiag_ugwp_flagflag for CIRES UGWP Diagnosticsflaglogical()inFalse
cdmbgwdmultiplication_factors_for_mountain_blocking_and_orographic_gravity_wave_dragmultiplication factors for cdmb and gwdnonereal(4)kind_physinFalse
xlatlatitudegrid latitude in radiansradiansreal(horizontal_dimension)kind_physinFalse
xlat_dlatitude_degreelatitude in degreesdegreereal(horizontal_dimension)kind_physinFalse
sinlatsine_of_latitudesine of the grid latitudenonereal(horizontal_dimension)kind_physinFalse
coslatcosine_of_latitudecosine of the grid latitudenonereal(horizontal_dimension)kind_physinFalse
areacell_areaarea of the grid cellm2real(horizontal_dimension)kind_physinFalse
ugrsx_windzonal windm s-1real(horizontal_dimension, vertical_dimension)kind_physinFalse
vgrsy_windmeridional windm s-1real(horizontal_dimension, vertical_dimension)kind_physinFalse
tgrsair_temperaturemodel layer mean temperatureKreal(horizontal_dimension, vertical_dimension)kind_physinFalse
qgrstracer_concentrationmodel layer mean tracer concentrationkg kg-1real(horizontal_dimension, vertical_dimension, number_of_tracers)kind_physinFalse
prsiair_pressure_at_interfaceair pressure at model layer interfacesPareal(horizontal_dimension, vertical_dimension_plus_one)kind_physinFalse
prslair_pressuremean layer pressurePareal(horizontal_dimension, vertical_dimension)kind_physinFalse
prslkdimensionless_exner_function_at_model_layersdimensionless Exner function at model layer centersnonereal(horizontal_dimension, vertical_dimension)kind_physinFalse
phiigeopotential_at_interfacegeopotential at model layer interfacesm2 s-2real(horizontal_dimension, vertical_dimension_plus_one)kind_physinFalse
philgeopotentialgeopotential at model layer centersm2 s-2real(horizontal_dimension, vertical_dimension)kind_physinFalse
delair_pressure_difference_between_midlayersair pressure difference between midlayersPareal(horizontal_dimension, vertical_dimension)kind_physinFalse
kpblvertical_index_at_top_of_atmosphere_boundary_layervertical index at top atmospheric boundary layerindexinteger(horizontal_dimension)inFalse
dusfcginstantaneous_x_stress_due_to_gravity_wave_dragzonal surface stress due to orographic gravity wave dragPareal(horizontal_dimension)kind_physoutFalse
dvsfcginstantaneous_y_stress_due_to_gravity_wave_dragmeridional surface stress due to orographic gravity wave dragPareal(horizontal_dimension)kind_physoutFalse
gw_dudttendency_of_x_wind_due_to_ugwpzonal wind tendency due to UGWPm s-2real(horizontal_dimension, vertical_dimension)kind_physoutFalse
gw_dvdttendency_of_y_wind_due_to_ugwpmeridional wind tendency due to UGWPm s-2real(horizontal_dimension, vertical_dimension)kind_physoutFalse
gw_dtdttendency_of_air_temperature_due_to_ugwpair temperature tendency due to UGWPK s-1real(horizontal_dimension, vertical_dimension)kind_physoutFalse
gw_kdiseddy_mixing_due_to_ugwpeddy mixing due to UGWPm2 s-1real(horizontal_dimension, vertical_dimension)kind_physoutFalse
tau_tofdinstantaneous_momentum_flux_due_to_turbulent_orographic_form_dragmomentum flux or stress due to TOFDPareal(horizontal_dimension)kind_physoutFalse
tau_mtbinstantaneous_momentum_flux_due_to_mountain_blocking_dragmomentum flux or stress due to mountain blocking dragPareal(horizontal_dimension)kind_physoutFalse
tau_ogwinstantaneous_momentum_flux_due_to_orographic_gravity_wave_dragmomentum flux or stress due to orographic gravity wave dragPareal(horizontal_dimension)kind_physoutFalse
tau_ngwinstantaneous_momentum_flux_due_to_nonstationary_gravity_wavemomentum flux or stress due to nonstationary gravity wavesPareal(horizontal_dimension)kind_physoutFalse
zmtbheight_of_mountain_blockingheight of mountain blocking dragmreal(horizontal_dimension)kind_physoutFalse
zlwbheight_of_low_level_wave_breakingheight of low level wave breakingmreal(horizontal_dimension)kind_physoutFalse
zogwheight_of_launch_level_of_orographic_gravity_waveheight of launch level of orographic gravity wavemreal(horizontal_dimension)kind_physoutFalse
dudt_mtbinstantaneous_change_in_x_wind_due_to_mountain_blocking_draginstantaneous change in x wind due to mountain blocking dragm s-2real(horizontal_dimension, vertical_dimension)kind_physoutFalse
dudt_ogwinstantaneous_change_in_x_wind_due_to_orographic_gravity_wave_draginstantaneous change in x wind due to orographic gw dragm s-2real(horizontal_dimension, vertical_dimension)kind_physoutFalse
dudt_tmsinstantaneous_change_in_x_wind_due_to_turbulent_orographic_form_draginstantaneous change in x wind due to TOFDm s-2real(horizontal_dimension, vertical_dimension)kind_physoutFalse
du3dt_mtbtime_integral_of_change_in_x_wind_due_to_mountain_blocking_dragtime integral of change in x wind due to mountain blocking dragm s-2real(horizontal_dimension, vertical_dimension)kind_physinoutFalse
du3dt_ogwtime_integral_of_change_in_x_wind_due_to_orographic_gravity_wave_dragtime integral of change in x wind due to orographic gw dragm s-2real(horizontal_dimension, vertical_dimension)kind_physinoutFalse
du3dt_tmstime_integral_of_change_in_x_wind_due_to_turbulent_orographic_form_dragtime integral of change in x wind due to TOFDm s-2real(horizontal_dimension, vertical_dimension)kind_physinoutFalse
dudttendency_of_x_wind_due_to_model_physicszonal wind tendency due to model physicsm s-2real(horizontal_dimension, vertical_dimension)kind_physinoutFalse
dvdttendency_of_y_wind_due_to_model_physicsmeridional wind tendency due to model physicsm s-2real(horizontal_dimension, vertical_dimension)kind_physinoutFalse
dtdttendency_of_air_temperature_due_to_model_physicsair temperature tendency due to model physicsK s-1real(horizontal_dimension, vertical_dimension)kind_physinoutFalse
rdxzblevel_of_dividing_streamlinelevel of the dividing streamlinenonereal(horizontal_dimension)kind_physoutFalse
con_ggravitational_accelerationgravitational accelerationm s-2real()kind_physinFalse
con_pipiratio of a circle's circumference to its diameterradiansreal()kind_physinFalse
con_cpspecific_heat_of_dry_air_at_constant_pressurespecific heat !of dry air at constant pressureJ kg-1 K-1real()kind_physinFalse
con_rdgas_constant_dry_airideal gas constant for dry airJ kg-1 K-1real()kind_physinFalse
con_rvgas_constant_water_vaporideal gas constant for water vaporJ kg-1 K-1real()kind_physinFalse
con_fvirtratio_of_vapor_to_dry_air_gas_constants_minus_onerv/rd - 1 (rv = ideal gas constant for water vapor)nonereal()kind_physinFalse
rainlwe_thickness_of_precipitation_amount_on_dynamics_timesteptotal rain at this time stepmreal(horizontal_dimension)kind_physinFalse
ntkeindex_for_turbulent_kinetic_energytracer index for turbulent kinetic energyindexinteger()inFalse
q_tketurbulent_kinetic_energyturbulent kinetic energyJreal(horizontal_dimension, vertical_dimension)kind_physinFalse
dqdt_tketendency_of_turbulent_kinetic_energy_due_to_model_physicsturbulent kinetic energy tendency due to model physicsJ s-1real(horizontal_dimension, vertical_dimension)kind_physinFalse
lprntflag_printcontrol flag for diagnostic print outflaglogical()inFalse
iprhorizontal_index_of_printed_columnhorizontal index of printed columnindexinteger()inFalse
errmsgccpp_error_messageerror message for error handling in CCPPnonecharacter()len=*outFalse
errflgccpp_error_flagerror flag for error handling in CCPPflaginteger()outFalse

CIRES UGWP Scheme General Algorithm

Functions/Subroutines

subroutine, public cires_ugwp::cires_ugwp_init (me, master, nlunit, input_nml_file, logunit, fn_nml2, lonr, latr, levs, ak, bk, dtp, cdmbgwd, cgwf, pa_rf_in, tau_rf_in, con_p0, do_ugwp, errmsg, errflg)
 The subroutine initializes the CIRES UGWP. More...
 
subroutine, public cires_ugwp::cires_ugwp_run (do_ugwp, me, master, im, levs, ntrac, dtp, kdt, lonr, oro, oro_uf, hprime, nmtvr, oc, theta, sigma, gamma, elvmax, clx, oa4, do_tofd, ldiag_ugwp, cdmbgwd, xlat, xlat_d, sinlat, coslat, area, ugrs, vgrs, tgrs, qgrs, prsi, prsl, prslk, phii, phil, del, kpbl, dusfcg, dvsfcg, gw_dudt, gw_dvdt, gw_dtdt, gw_kdis, tau_tofd, tau_mtb, tau_ogw, tau_ngw, zmtb, zlwb, zogw, dudt_mtb, dudt_ogw, dudt_tms, du3dt_mtb, du3dt_ogw, du3dt_tms, dudt, dvdt, dtdt, rdxzb, con_g, con_pi, con_cp, con_rd, con_rv, con_fvirt, rain, ntke, q_tke, dqdt_tke, lprnt, ipr, errmsg, errflg)
 
subroutine gwdps_v0 (IM, km, imx, do_tofd, Pdvdt, Pdudt, Pdtdt, Pkdis, U1, V1, T1, Q1, KPBL, PRSI, DEL, PRSL, PRSLK, PHII, PHIL, DTP, KDT, sgh30, HPRIME, OC, OA4, CLX4, THETA, vSIGMA, vGAMMA, ELVMAXD, DUSFC, DVSFC, xlatd, sinlat, coslat, sparea, cdmbgwd, me, master, rdxzb, zmtb, zogw, tau_mtb, tau_ogw, tau_tofd, dudt_mtb, dudt_ogw, dudt_tms)
 Note for the sub-grid scale orography scheme in UGWP v0: Due to degraded forecast scores of simulations with revised schemes for subgrid-scale orography effects in FV3GFS, EMC reinstalled the original gwdps-code with updated efficiency factors for the mountain blocking and OGW drag. The GFS OGW is described in the separate section (GFS Orographic Gravity Wave Drag and Mountain Blocking Scheme Module) and its “call” moved into UGWP-driver subroutine. This combination of NGW and OGW schemes was tested in the FV3GFS-L127 medium-range forecasts (15-30 days) for C96, C192, C384 and C768 resolutions and work in progress to introduce the optimal choice for the scale-aware representations of the efficiency factors that will reflect the better simulations of GW activity by FV3 dynamical core at higher horizontal resolutions. With the MERRA-2 VMF function for NGWs (slat_geos5_tamp) and operational OGW drag scheme (gwdps), FV3GFS simulations can successfully forecast the recent major mid-winter sudden stratospheric warming (SSW) events of 2018-02-12 and 2018-12-31 (10-14 days before the SSW onset; Yudin et al. 2019 [180]). The first multi-year (2015-2018) FV3GFS simulations with UGWP v0 also produce the equatorial QBO-like oscillations in the zonal wind and temperature anomalies. More...
 
subroutine fv3_ugwp_solv2_v0 (klon, klev, dtime, tm1, um1, vm1, qm1, prsl, prsi, philg, xlatd, sinlat, coslat, pdudt, pdvdt, pdtdt, dked, tau_ngw, mpi_id, master, kdt)
 
subroutine edmix_ugwp_v0 (im, levs, dtp, t1, u1, v1, q1, del, prsl, prsi, phil, prslk, pdudt, pdvdt, pdTdt, pkdis, ed_dudt, ed_dvdt, ed_dTdt, me, master, kdt)
 Part-3 of UGWP v1 Dissipative (eddy) effects of UGWP it will be activated after tests of OGW (new revision) and NGW with MERRA-2 forcing. More...
 
subroutine diff_1d_wtend (levs, dt, F, F1, Km, rdp, rdpm, S, S1)
 explicit diffusion solver More...
 
subroutine diff_1d_ptend (levs, dt, F, Km, rdp, rdpm, S)
 explicit "eddy" smoother for tendencies. More...
 
subroutine ugwp_triggers
 
subroutine get_xy_pt (V, Vx, Vy, nx, ny, dlam1, dlat)
 
subroutine get_xyd_wind (V, Vx, Vy, Vyd, nx, ny, dlam1, dlat, divJp, divJm)
 
subroutine trig3d_fjets (nx, ny, nz, U, V, T, Q, P3D, PS, delp, delz, lon, lat, pmid, trig3d_fgf)
 
subroutine trig3d_okubo (nx, ny, nz, U, V, T, Q, P3d, PS, delp, delz, lon, lat, pmid, trig3d_okw)
 
subroutine trig3d_dconv (nx, ny, nz, U, V, T, Q, P3d, PS, delp, delz, lon, lat, pmid, trig3d_conv, dcheat3d, precip2d, cld_klevs2d, scheat3d)
 
subroutine cires_3d_triggers (nx, ny, nz, lon, lat, pmid, U, V, W, T, Q, delp, delz, p3d, PS, HS, Hyam, Hybm, Hyai, Hybi, trig3d_okw, trig3d_fgf, trig3d_conv, dcheat3d, precip2d, cld_klevs2d, scheat3d)
 
subroutine get_spectra_tau_convgw (nw, im, levs, dcheat, scheat, precip, icld, xlatd, sinlat, coslat, taub, klev, if_src, nf_src)
 
subroutine get_spectra_tau_nstgw (nw, im, levs, trig_fgf, xlatd, sinlat, coslat, taub, klev, if_src, nf_src)
 
subroutine get_spectra_tau_okw (nw, im, levs, trig_okw, xlatd, sinlat, coslat, taub, klev, if_src, nf_src)
 
subroutine slat_geos5_tamp (im, tau_amp, xlatdeg, tau_gw)
 This subroutine describes the latitudinal shape of the VMF-function as displayed in Figure 3 of Molod et al. (2015). The enhanced values of VMF in the equatorial region result in QBO-like oscillations in the equatorial zonal winds and more realistic simulations of the equatorial dynamics in the GEOS-5 operational and MERRA-2 reanalysis products. More...
 
subroutine slat_geos5 (im, xlatdeg, tau_gw)
 
subroutine init_nazdir (naz, xaz, yaz)
 
subroutine cires_ugwp_module::cires_ugwp_mod_init (me, master, nlunit, input_nml_file, logunit, fn_nml, lonr, latr, levs, ak, bk, pref, dtp, cdmvgwd, cgwf, pa_rf_in, tau_rf_in)
 
subroutine cires_ugwp_module::cires_ugwp_driver (im, levs, dtp, kdt, me, lprnt, lonr, pa_rf, tau_rf, cdmbgwd, xlat, xlatd, sinlat, coslat, ugrs, vgrs, tgrs, qgrs, prsi, prsl, prslk, phii, phil, delp, orostat, kpbl, dusfc, dvsfc, dudt, dvdt, dtdt, kdis, axtot, axo, axc, axf, aytot, ayo, ayc, ayf, eps_tot, ekdis, trig_okw, trig_fgf, dcheat, precip, cld_klevs, zmtb, scheat, dlength, cldf, taus_sso, taus_ogw, tauf_ogw, tauf_ngw, ugw_zmtb, ugw_zlwb, ugw_zogw, ugw_axmtb, ugw_axlwb, ugw_axtms)
 
subroutine cires_ugwp_module::cires_ugwp_advance
 
subroutine cires_ugwp_module::cires_ugwp_mod_finalize
 
subroutine um_flow (nz, klow, ktop, up, vp, tp, qp, dp, zpm, zpi, pmid, pint, bn2, uhm, vhm, bn2hm, rhohm)
 
subroutine mflow_tauz (levs, up, vp, tp, qp, dp, zpm, zpi, pmid, pint, rho, ui, vi, ti, bn2i, bvi, rhoi)
 
subroutine get_unit_vector (u, v, u_n, v_n, mag)
 
subroutine ugwp_wmsdis_naz (levs, ksrc, nw, naz, kxw, taub_lat, ch, xaz, yaz, fcor, c2f2, dp, zmid, zint, pmid, pint, rho, ui, vi, ti, kvg, ktg, krad, kion, bn2i, bvi, rhoi, ax, ay, eps, ked, tau1)
 
subroutine ugwp_wmsdis_az1 (levs, ksrc, nw, kxw, ch, dch, taub_sp, tau_bulk, spnorm, fcor, c2f2, zm, zi, rho, um, tm, bn2, bn, rhoi, bnrho, dzirho, dzpi, kvg, ktg, krad, kion, kmol, eps, ked, tau)
 
subroutine fvs93_ugwps (nw, ch, dch, taub_sp, spnorm, nslope, bn2, bn, bnrhos)
 
subroutine ugwp_drag_mtb (iemax, nz, elvpd, elvp, hprime, sigma, theta, oc, oa4, clx4, gam, zpbl, up, vp, tp, qp, dp, zpm, zpi, pmid, pint, idxzb, drmtb, taumtb)
 
subroutine ugwp_taub_oro (levs, izb, kxw, tau_izb, fcor, hprime, sigma, theta, oc, oa4, clx4, gamm, elvp, up, vp, tp, qp, dp, zpm, zpi, pmid, pint, xn, yn, umag, tautot, tauogw, taulee, drlee, tau_src, kxridge, kdswj, krefj, kotr)
 
subroutine ugwp_oro_lsatdis (krefj, levs, tauogw, tautot, tau_src, kxw, fcor, kxridge, up, vp, tp, qp, dp, zpm, zpi, pmid, pint, xn, yn, umag, drtau, kdis)
 
subroutine ugwp_tofd (im, levs, sigflt, elvmax, zpbl, u, v, zmid, utofd, vtofd, epstofd, krf_tofd)
 
subroutine ugwp_tofd1d (levs, sigflt, elvmax, zsurf, zpbl, u, v, zmid, utofd, vtofd, epstofd, krf_tofd)
 
subroutine ugwp_lsatdis_naz (levs, ksrc, nw, naz, kxw, taub_spect, ch, xaz, yaz, fcor, c2f2, dp, zmid, zint, pmid, pint, rho, ui, vi, ti, kvg, ktg, krad, kion, bn2i, bvi, rhoi, ax, ay, eps, ked, tau1)
 
subroutine ugwp_lsatdis_az1 (levs, ksrc, nw, kxw, ch, taub_sp, fcor, c2f2, zm, zi, rho, um, tm, bn2, bn, rhoi, dzirho, dzpi, kvg, ktg, krad, kion, kmol, eps, ked, tau)
 
subroutine ugwp_limit_1d (ax, ay, eps, ked, levs)