climkern

climkern.calc_RH_feedback(ctrl_q, ctrl_ta, ctrl_ps, pert_q, pert_ta, pert_ps, pert_trop=None, kern='GFDL', sky='all-sky', method=1)

Calculate the TOA radiative perturbations from changes in relative humidity following Held & Shell (2012). Horizontal resolution is kept at input data’s resolution.

Parameters:
  • ctrl_q (xarray DataArray) – DataArray containing specific humidity on pressure levels in the control simulation. Must be 4D with coords of time, lat, lon, and plev with units of either “1”, “kg/kg”, or “g/kg”.

  • ctrl_ta (xarray DataArray) – DataArray containing air temperature on pressure levels in the control simulation. Must be 4D with coords of time, lat, lon, and plev with units “K”.

  • ctrl_ps (xarray DataArray) – DataArray containing surface pressure in the control simulation. 3D with coordinates of time, latitude, and longitude and units of Pa or hPa.

  • pert_q (xarray DataArray) – DataArray containing specific humidity on pressure levels in the perturbed simulation. 4D with coordinates of time, latitude, longitude, and pressure level and units of “1”, “kg/kg”, or “g/kg”.

  • pert_ps (xarray DataArray) – DataArray containing surface pressure in the perturbed simulation. 3D with coordinates of time, latitude, and longitude and units of Pa or hPa.

  • pert_trop (xarray DataArray, optional) – DataArray containing tropopause pressure in the perturbed simulation. 3D with coordinates of time, latitude, and longitude and units of Pa or hPa. If not provided, ClimKern will assume a tropopause height of 100 hPa at the equator, linearly descending with the cosine of latitude to 300 hPa at the poles.

  • kern (string, optional) – String specifying the name of the desired kernel. Defaults to “GFDL”.

  • sky (string, optional) – String, either “all-sky” or “clear-sky”, specifying whether to calculate the all-sky or clear-sky feedbacks. Defaults to “all-sky”.

  • method (int, optional) – Specifies the method to use to calculate the specific humidity feedback. Options 1, 2, and 3 use the change in the natural logarithm of specific humidity, while 4 uses the linear change. The options are: 1 – Uses the actual logarithm for both the specific humidity response and the normalization factor. 2 – Uses the fractional change approximation of logarithms only in the normalization factor, with the actual logarithm used in the specific humidity response. 3 – Uses the fractional change approximation of logarithms in the specific humidity response & normalization factor. 4 – Uses the linear change in specific humidity for both.

Returns:

RH_feedback – 3D DataArray containing the vertically integrated radiative perturbations from changes in relative humidity (LW+SW). Has coordinates of time, lat, and lon.

Return type:

xarray DataArray

climkern.calc_T_feedbacks(ctrl_ta, ctrl_ts, ctrl_ps, pert_ta, pert_ts, pert_ps, pert_trop=None, kern='GFDL', sky='all-sky', fixRH=False)

Calculate the raditive pertubations (W/m^2) at the TOA from changes in surface skin and air temperature using user-specified kernel. Horizontal resolution is kept at input data’s resolution.

Parameters:
  • ctrl_ta (xarray DataArray) – DataArray containing air temperature on pressure levels in the control simulation. 4D with coordinates of time, latitude, longitude, and pressure level and units of K.

  • ctrl_ts (xarray DataArray) – DataArray containing surface skin temperature in the control simulation. 3D with coordinates of time, latitude, and longitude and units of K.

  • ctrl_ps (xarray DataArray) – DataArray containing surface pressure in the control simulation. 3D with coordinates of time, latitude, and longitude and units of Pa or hPa.

  • pert_ta (xarray DataArray) – DataArray containing air temperature on pressure levels in the perturbed simulation. 4D with coordinates of time, latitude, longitude, and pressure level and units of K.

  • pert_ts (xarray DataArray) – DataArray containing surface skin temperature in the perturbed simulation. 3D with coordinates of time, latitude, and longitude and units of K.

  • pert_ps (xarray DataArray) – DataArray containing surface pressure in the perturbed simulation. 3D with coordinates of time, latitude, and longitude and units of Pa or hPa.

  • pert_trop (xarray DataArray, optional) – DataArray containing tropopause pressure in the perturbed simulation. 3D with coordinates of time, latitude, and longitude and units of Pa or hPa. If not provided, ClimKern will assume a tropopause height of 100 hPa at the equator, linearly descending with the cosine of latitude to 300 hPa at the poles.

  • kern (string, optional) – String specifying the name of the desired kernel. Defaults to “GFDL”.

  • sky (string, optional) – String, either “all-sky” or “clear-sky”, specifying whether to calculate the all-sky or clear-sky feedbacks. Defaults to “all-sky”.

  • fixRH (boolean, optional) – Specifies whether to calculate alternative Planck and lapse rate feedbacks using relative humidity as a state variable, as outlined in Held & Shell (2012).

Returns:

  • lr_feedback (xarray DataArray) – 3D DataArray containing radiative perturbations at TOA caused by changes in tropospheric lapse rate with coordinates of time, latitude, and longitude.

  • planck_feedback (xarray DataArray) – 3D DataArray containing radiative perturbations at TOA caused by vertically-uniform temperature change with coordinates of time, latitude, and longitude.

climkern.calc_alb_feedback(ctrl_rsus: xarray.DataArray, ctrl_rsds: xarray.DataArray, pert_rsus: xarray.DataArray, pert_rsds: xarray.DataArray, kern: str = 'GFDL', sky: str = 'all-sky') xarray.DataArray

Calculate the radiative perturbation (W/m^2) from changes in surface albedo using user-specified radiative kernel. Horizontal resolution is kept at input data’s resolution.

Parameters:
  • ctrl_rsus (xarray DataArray) – DataArray containing the upwelling SW radiation at the surface in the control simulation. Must be 3D with coordinates of time, latitude, and longitude with units of W/m^2.

  • ctrl_rsds (xarray DataArray) – DataArray containing the downwelling SW radiation at the surface in the control simulation. Must be 3D with coordinates of time, latitude, and longitude with units of W/m^2.

  • pert_rsus (xarray DataArray) – DataArray containing the upwelling SW radiation at the surface in the perturbed simulation. Must be 3D with coordinates of time, latitude, and longitude with units of W/m^2.

  • pert_rsds (xarray DataArray) – DataArray containing the downwelling SW radiation at the surface in the perturbed simulation. Must be 3D with coordinates of time, latitude, and longitude with units of W/m^2.

  • kern (string, optional) – String specifying the name of the desired kernel. Defaults to “GFDL”.

  • sky (string, optional) – String, either “all-sky” or “clear-sky”, specifying whether to calculate the all-sky or clear-sky feedbacks. Defaults to “all-sky”.

Returns:

alb_feedback – 3D DataArray containing radiative perturbations at TOA caused by changes in surface albedo with coordinates of time, latitude, and longitude.

Return type:

xarray DataArray

climkern.calc_cloud_LW(t_as, t_cs, q_lwas, q_lwcs, dCRE_lw, rf_lwas, rf_lwcs)

Calculate the radiative perturbation from the longwave cloud feedback using the adjustment method outlined in Soden et al. (2008).

Parameters:
  • t_as (xarray DataArray) – DataArray containing the vertically integrated all-sky radiative perturbation at the TOA from the total temperature feedback. The total temperature feedback is the sum of the Planck and lapse rate feedbacks. Should have dims of lat, lon, and time.

  • t_cs (xarray DataArray) – DataArray containing the vertically integrated clear-sky radiative perturbation at the TOA from the total temperature feedback. The total temperature feedback is the sum of the Planck and lapse rate feedbacks. Should have dims of lat, lon, and time.

  • q_lwas (xarray DataArray) – DataArray containing the vertically integrated LW all-sky radiative perturbation at the TOA from the water vapor feedback. Should have coords of lat, lon, and time.

  • q_lwcs (xarray DataArray) – DataArray containing the vertically integrated LW clear-sky radiative perturbation at the TOA from the water vapor feedback. Should have coords of lat, lon, and time.

  • dCRE_lw (xarray DataArray) – DataArray containing the change in LW cloud radiative effect at the TOA with coords of time, lat, and lon and units of Wm^-2. positive = downwards.

  • rf_lwas (xarray DataArray) – DataArray containing the LW all-sky radiative forcing in units of Wm^-2 with coords of lat, lon, and time.

  • rf_lwcs (xarray DataArray) – DataArray containing the LW clear-sky radiative forcing in units of Wm^-2 with coords of lat, lon, and time.

Returns:

lw_cld_feedback – Three-dimensional DataArray containing the TOA radiative perturbation from the longwave cloud feedback.

Return type:

xarray DataArray

climkern.calc_cloud_LW_res(ctrl_FLNT, pert_FLNT, RF_lw, t_lw, q_lw)

Calculate the radiative perturbation from the shortwave cloud feedback using the residual method outlined in Soden & Held (2006).

Parameters:
  • ctrl_FLNT (xarray DataArray) – Three-dimensional DataArray containing the all-sky net longwave flux at the top-of-atmosphere in the control simulation with coords of time, lat, and lon and units of Wm^-2. It should be oriented such that positive = downwards.

  • pert_FLNT (xarray DataArray) – Three-dimensional DataArray containing the all-sky net longwave flux at the top-of-atmosphere in the perturbed simulation with coords of time, lat, and lon and units of Wm^-2. It should be oriented such that positive = downwards.

  • RF_lw (xarray DataArray) – The longwave all-sky radiative forcing in units of Wm^-2 with coords of lat, lon, and time.

  • t_lw (xarray DataArray) – DataArray containing the vertically integrated all-sky radiative perturbation at the TOA from the total temperature feedback. The total temperature feedback is the sum of the Planck and lapse rate feedbacks. Should have dims of lat, lon, and time.

  • q_lw (xarray DataArray) – DataArray containing the vertically integrated LW all-sky radiative perturbation at the TOA from the water vapor feedback. Should have coords of lat, lon, and time.

Returns:

lw_cld_feedback – Three-dimensional DataArray containing the TOA radiative perturbation from the longwave cloud feedback.

Return type:

xarray DataArray

climkern.calc_cloud_SW(alb_as, alb_cs, q_swas, q_swcs, dCRE_sw, RF_swas, RF_swcs)

Calculate the radiative perturbation from the shortwave cloud feedback using the adjustment method outlined in Soden et al. (2008).

Parameters:
  • alb_as (xarray DataArray) – DataArray containing the all-sky radiative perturbation at the TOA from the surface albedo feedback. Should have coords of lat, lon, and time.

  • alb_cs (xarray DataArray) – DataArray containing the clear-sky radiative perturbation at the TOA from the surface albedo feedback. Should have coords of lat, lon, and time.

  • q_swas (xarray DataArray) – DataArray containing the vertically integrated all-sky LW radiative perturbation at the TOA from the shortwave water vapor feedback. Should have coords of lat, lon, and time.

  • q_swcs (xarray DataArray) – DataArray containing the vertically integrated clear-sky LW radiative perturbation at the TOA from the shortwave water vapor feedback. Should have coords of lat, lon, and time.

  • dCRE_sw (xarray DataArray) – Three-dimensional DataArray containing the change in shortwave cloud radiative effect at the top-of-atmosphere with coords of time, lat, and lon and units of Wm^-2. positive = downwards.

  • RF_swas (xarray DataArray) – The shortwave all-sky radiative forcing in units of Wm^-2 with coords of lat, lon, and time.

  • RF_swcs (xarray DataArray) – The shortwave clear-sky radiative forcing in units of Wm^-2 with coords of lat, lon, and time.

Returns:

sw_cld_feedback – Three-dimensional DataArray containing the TOA radiative perturbation from the shortwave cloud feedback.

Return type:

xarray DataArray

climkern.calc_cloud_SW_res(ctrl_FSNT, pert_FSNT, RF_sw, q_sw, alb_sw)

Calculate the radiative perturbation from the shortwave cloud feedback using the residual method outlined in Soden & Held (2006).

Parameters:
  • ctrl_FSNT (xarray DataArray) – Three-dimensional DataArray containing the all-sky net shortwave flux at the top-of-atmosphere in the control simulation with coords of time, lat, and lon and units of Wm^-2. It should be oriented such that positive = downwards/incoming.

  • pert_FSNT (xarray DataArray) – Three-dimensional DataArray containing the all-sky net shortwave flux at the top-of-atmosphere in the perturbed simulation with coords of time, lat, and lon and units of Wm^-2. It should be oriented such that positive = downwards/incoming.

  • RF_sw (xarray DataArray) – The shortwave all-sky radiative forcing in units of Wm^-2 with coords of lat, lon, and time. This is usually the stratosphere- adjusted RF.

  • q_sw (xarray DataArray) – DataArray containing the vertically integrated all-sky LW radiative perturbation at the TOA from the shortwave water vapor feedback. Should have coords of lat, lon, and time.

  • alb_sw (xarray DataArray) – DataArray containing the all-sky radiative perturbation at the TOA from the surface albedo feedback. Should have coords of lat, lon, and time.

Returns:

sw_cld_feedback – Three-dimensional DataArray containing the TOA radiative perturbation from the shortwave cloud feedback.

Return type:

xarray DataArray

climkern.calc_dCRE_LW(ctrl_FLNT, pert_FLNT, ctrl_FLNTC, pert_FLNTC)

Calculate the change in the LW cloud radiative effect at the TOA.

Parameters:
  • ctrl_FLNT (xarray DataArray) – Three-dimensional DataArray containing the all-sky net longwave flux at the top-of-atmosphere in the control simulation with coords of time, lat, and lon and units of Wm^-2. It should be oriented such that positive = downwards.

  • pert_FLNT (xarray DataArray) – Three-dimensional DataArray containing the all-sky net longwave flux at the top-of-atmosphere in the perturbed simulation with coords of time, lat, and lon and units of Wm^-2. It should be oriented such that positive = downwards.

  • ctrl_FLNTC (xarray DataArray) – Three-dimensional DataArray containing the clear-sky net longwave flux at the top-of-atmosphere in the control simulation with coords of time, lat, and lon and units of Wm^-2. It should be oriented such that positive = downwards.

  • pert_FLNTC (xarray DataArray) – Three-dimensional DataArray containing the clear-sky net longwave flux at the top-of-atmosphere in the perturbed simulation with coords of time, lat, and lon and units of Wm^-2. It should be oriented such that positive = downwards.

Returns:

dCRE_LW – Three-dimensional DataArray containing the change in longwave cloud radiative effect at the top-of-atmosphere with coords of time, lat, and lon and units of Wm^-2. positive = downwards.

Return type:

xarray DataArray

climkern.calc_dCRE_SW(ctrl_FSNT, pert_FSNT, ctrl_FSNTC, pert_FSNTC)

Calculate the change in the SW cloud radiative effect at the TOA.

Parameters:
  • ctrl_FSNT (xarray DataArray) – Three-dimensional DataArray containing the all-sky net shortwave flux at the top-of-atmosphere in the control simulation with coords of time, lat, and lon and units of Wm^-2. It should be oriented such that positive = downwards/incoming.

  • pert_FSNT (xarray DataArray) – Three-dimensional DataArray containing the all-sky net shortwave flux at the top-of-atmosphere in the perturbed simulation with coords of time, lat, and lon and units of Wm^-2. It should be oriented such that positive = downwards/incoming.

  • ctrl_FSNTC (xarray DataArray) – Three-dimensional DataArray containing the clear-sky net shortwave flux at the top-of-atmosphere in the control simulation with coords of time, lat, and lon and units of Wm^-2. It should be oriented such that positive = downwards/incoming.

  • pert_FSNTC (xarray DataArray) – Three-dimensional DataArray containing the clear-sky net shortwave flux at the top-of-atmosphere in the perturbed simulation with coords of time, lat, and lon and units of Wm^-2. It should be oriented such that positive = downwards/incoming.

Returns:

dCRE_SW – Three-dimensional DataArray containing the change in shortwave cloud radiative effect at the top-of-atmosphere with coords of time, lat, and lon and units of Wm^-2. positive = downwards.

Return type:

xarray DataArray

climkern.calc_q_feedbacks(ctrl_q, ctrl_ta, ctrl_ps, pert_q, pert_ps, pert_trop=None, kern='GFDL', sky='all-sky', method=1)

Calculate the raditive pertubations (W/m^2), LW & SW, at the TOA from changes in specific humidity using user-specified kernel. Horizontal resolution is kept at input data’s resolution.

Parameters:
  • ctrl_q (xarray DataArray) – DataArray containing specific humidity on pressure levels in the control simulation. 4D with coordinates of time, latitude, longitude, and pressure level and units of “1”, “kg/kg”, or “g/kg”.

  • ctrl_ta (xarray DataArray) – DataArray containing air temperature on pressure levels in the control simulation. 4D with coordinates of time, latitude, longitude, and pressure level and units of K.

  • ctrl_ps (xarray DataArray) – DataArray containing surface pressure in the control simulation. 3D with coordinates of time, latitude, and longitude and units of Pa or hPa.

  • pert_q (xarray DataArray) – DataArray containing specific humidity on pressure levels in the perturbed simulation. 4D with coordinates of time, latitude, longitude, and pressure level and units of “1”, “kg/kg”, or “g/kg”.

  • pert_ps (xarray DataArray) – DataArray containing surface pressure in the perturbed simulation. 3D with coordinates of time, latitude, and longitude and units of Pa or hPa.

  • pert_trop (xarray DataArray, optional) – DataArray containing tropopause pressure in the perturbed simulation. 3D with coordinates of time, latitude, and longitude and units of Pa or hPa. If not provided, ClimKern will assume a tropopause height of 100 hPa at the equator, linearly descending with the cosine of latitude to 300 hPa at the poles.

  • kern (string, optional) – String specifying the name of the desired kernel. Defaults to “GFDL”.

  • sky (string, optional) – String, either “all-sky” or “clear-sky”, specifying whether to calculate the all-sky or clear-sky feedbacks. Defaults to “all-sky”.

  • method (int, optional) – Specifies the method to use to calculate the specific humidity feedback. Options 1, 2, and 3 use the change in the natural logarithm of specific humidity, while 4 uses the linear change. The options are: 1 – Uses the actual logarithm for both the specific humidity response and the normalization factor. 2 – Uses the fractional change approximation of logarithms only in the normalization factor, with the actual logarithm used in the specific humidity response. 3 – Uses the fractional change approximation of logarithms in the specific humidity response & normalization factor. 4 – Uses the linear change in specific humidity for both.

Returns:

  • lw_q_feedback (xarray DataArray) – 3D DataArray containing LW radiative perturbations at TOA caused by changes in specific humidity with coordinates of time, latitude, and longitude.

  • sw_q_feedback (xarray DataArray) – 3D DataArray containing SW radiative perturbations at TOA caused by changes in specific humidity with coordinates of time, latitude, and longitude.

climkern.calc_strato_T(ctrl_ta, pert_ta, pert_ps, pert_trop=None, kern='GFDL', sky='all-sky')

Calculate the raditive pertubations (W/m^2) at the TOA from changes in statosphere air temperature using user-specified kernel. Horizontal resolution is kept at input data’s resolution.

Parameters:
  • ctrl_ta (xarray DataArray) – DataArray containing air temperature on pressure levels in the control simulation. 4D with coordinates of time, latitude, longitude, and pressure level and units of K.

  • pert_ta (xarray DataArray) – DataArray containing air temperature on pressure levels in the perturbed simulation. 4D with coordinates of time, latitude, longitude, and pressure level and units of K.

  • pert_ps (xarray DataArray) – DataArray containing surface pressure in the perturbed simulation. 3D with coordinates of time, latitude, and longitude and units of Pa or hPa.

  • pert_trop (xarray DataArray, optional) – DataArray containing tropopause pressure in the perturbed simulation. 3D with coordinates of time, latitude, and longitude and units of Pa or hPa. If not provided, ClimKern will assume a tropopause height of 100 hPa at the equator, linearly descending with the cosine of latitude to 300 hPa at the poles.

  • kern (string, optional) – String specifying the name of the desired kernel. Defaults to “GFDL”.

  • sky (string, optional) – String, either “all-sky” or “clear-sky”, specifying whether to calculate the all-sky or clear-sky feedbacks. Defaults to “all-sky”.

Returns:

T_feedback – 3D DataArray containing the vertically integrated radiative perturbations caused by changes in stratospheric temperature. Has coordinates of time, lat, and lon.

Return type:

xarray DataArray

climkern.calc_strato_q(ctrl_q, ctrl_ta, pert_q, pert_ps, pert_trop=None, kern='GFDL', sky='all-sky', method=1)

Calculate the raditive pertubations (W/m^2), LW & SW, at the TOA from changes in stratospheric specific humidity using user-specified kernel. Horizontal resolution is kept at input data’s resolution.

Parameters:
  • ctrl_q (xarray DataArray) – DataArray containing specific humidity on pressure levels in the control simulation. 4D with coordinates of time, latitude, longitude, and pressure level and units of “1”, “kg/kg”, or “g/kg”.

  • ctrl_ta (xarray DataArray) – DataArray containing air temperature on pressure levels in the control simulation. 4D with coordinates of time, latitude, longitude, and pressure level and units of K.

  • pert_q (xarray DataArray) – DataArray containing specific humidity on pressure levels in the perturbed simulation. 4D with coordinates of time, latitude, longitude, and pressure level and units of “1”, “kg/kg”, or “g/kg”.

  • pert_ps (xarray DataArray) – DataArray containing surface pressure in the perturbed simulation. 3D with coordinates of time, latitude, and longitude and units of Pa or hPa.

  • pert_trop (xarray DataArray, optional) – DataArray containing tropopause pressure in the perturbed simulation. 3D with coordinates of time, latitude, and longitude and units of Pa or hPa. If not provided, ClimKern will assume a tropopause height of 100 hPa at the equator, linearly descending with the cosine of latitude to 300 hPa at the poles.

  • kern (string, optional) – String specifying the name of the desired kernel. Defaults to “GFDL”.

  • sky (string, optional) – String, either “all-sky” or “clear-sky”, specifying whether to calculate the all-sky or clear-sky feedbacks. Defaults to “all-sky”.

  • method (int, optional) – Specifies the method to use to calculate the specific humidity feedback. Options 1, 2, and 3 use the change in the natural logarithm of specific humidity, while 4 uses the linear change. The options are: 1 – Uses the actual logarithm for both the specific humidity response and the normalization factor. 2 – Uses the fractional change approximation of logarithms only in the normalization factor, with the actual logarithm used in the specific humidity response. 3 – Uses the fractional change approximation of logarithms in the specific humidity response & normalization factor. 4 – Uses the linear change in specific humidity for both.

Returns:

  • lw_q_feedback (xarray DataArray) – 3D DataArray containing the vertically integrated radiative perturbations from changes in specific humidity in the stratosphere (longwave). Has coordinates of time, lat, and lon.

  • sw_q_feedback (xarray DataArray) – 3D DataArray containing the vertically integrated radiative perturbations from changes in specific humidity in the stratosphere (shortwave). Has coordinates of time, lat, and lon.

climkern.spat_avg(data, lat_bound_s=-90, lat_bound_n=90)

Compute the spatial average while weighting for cos(latitude), optionally specifying latitudinal boundaries.

Parameters:
  • data (Xarray DataArray) – 3D input data over which to compute the spatial average.

  • lat_bound_s (float, optional) – Latitude of the southern boundary of the area over which to average. Defaults to -90 to compute a global average.

  • lat_bound_n (float, optional) – Latitude of the northern boundary of the area over which to average. Defaults to 90 to compute a global average.

Returns:

avg – New spatially averaged DataArray with latitude and longitude coordinates now removed.

Return type:

Xarray DataArray

climkern.tutorial_data(label)

Retrieve tutorial data which should be located in the package’s “data” folder.

Parameters:

label (string) –

Specifies which data to access. Choices are “ctrl”, “pert”, “IRF”, “adjRF”,

or “ERF” for the 1xCO2, 2xCO2, instantaneous radiative forcing,

stratosphere-adjusted radiative forcing, and effective radiative forcing, respectively.

Returns:

data – Requested tutorial dataset.

Return type:

xarray Dataset