Groundbased

This module contains the OrbitalRadar class that runs the simulator for suborbital radar data. It is a subclass of the Simulator class.

Difference between ground-based and airborne radar geometry: - along-track coordinate from mean wind for ground based and mean flight vel. from airborne radar - no ground echo added to airborne radar - range grid of airborne radar assumed to be height above mean sea level and height above ground for groundbased - no attenuation correction for airborne radar - lat/lon coordinates included as input to airborne radar

class orbital_radar.suborbital.Suborbital(geometry, name, config_file, suborbital_radar, input_radar_format)[source]

Run the simulator for suborbital radar data.

Attributes:
summary

Prints short summary of simulator settings.

Methods

add_airborne_variables()

Add variables specific to airborne simulator to the dataset, i.e., the mean flight velocity.

add_attenuation(ds, da_gas_atten)

Add attenuation to dataset.

add_attributes()

Adds attributes to the variables of the dataset

add_ground_echo(ds)

Calculates artificial ground echo inside ground-based radar range grid.

add_groundbased_variables()

Add variables specific to groundbased simulator to the dataset, i.e., the mean horizontal wind.

add_noise(x, x_std, noise)

Equation to calculate the noise from values without noise, the uncertainty of the values, and random noise.

add_vmze_attrs(ds)

Adds attributes to Doppler velocity and radar reflectivity variables.

apply_detection_limit(var_ze, var_other)

Applies the detection limit of the spaceborne radar to the along-height convoluted data.

attenuation_correction(ds, ds_cloudnet)

Gas attenuation correction based on Cloudnet.

calculate_along_track_sat_bin_edges()

Calculate the bin edges of the along-track satellite grid.

calculate_height_sat_bin_edges()

Calculate the bin edges of the height satellite grid.

calculate_ms_flag()

Calculates the multiple scattering flag.

calculate_nubf()

Calculates the non-uniform beam filling from the standard deviation of Ze within the radar volume.

calculate_nubf_flag([threshold])

Calculate non-uniform beam filling flag.

calculate_signal_fraction()

Calculates the fraction of bins that contain a ze signal above the detection limit of the spaceborne radar.

calculate_vm_bias()

Calculate the satellite Doppler velocity bias between the estimate with and without satellite motion error.

calculate_vm_bias_flag([threshold])

Calculate the satellite Doppler velocity bias flag.

calculate_vm_noise()

Adds noise to satellite Doppler velocity based on the pre-defined lookup table with noise values for different radar reflectivity bins.

calculate_vm_std_nubf()

Calculate outstanding error in correcting Mean Doppler Velocity biases caused by non-uniform beam filling

calculate_ze_noise()

Adds noise to satellite radar reflectivity based on the pre-defined lookup table with noise values for different radar reflectivity bins.

check_input_dataset()

Check user input for consistency.

check_is_sea_level(ds)

Check if input radar range/height grid is defined with respect to ground level or sea level.

convert_frequency(ds)

Convert frequency from 35 to 94 GHz.

convolve_along_track()

Calculates the along-track convolution from the input suborbital data using the along-track weighting function of the spaceborne radar.

convolve_height()

Convolution of the along-track integrated data with the range weighting function of the spaceborne radar.

correct_dielectric_constant(ds)

Apply correction for dielectric constant assumed in Ze calculation of suborbital radar to match the dielectric constant of the spaceborne radar.

create_along_track(ds)

Creates along-track coordinates from time coordinates.

create_regular_along_track(ds)

Creates regular along-track coordinate for suborbital radar.

create_regular_height()

Creates regular height coordinate for suborbital radar.

fold_vm()

Doppler velocity folding correction.

integrate_along_track()

Integrates the along-track convoluted data to profiles, which represent the satellite's footprint.

interpolate_to_regular_grid(ds)

Interpolates radar data to regular grid in along-track and height.

plot(**kwds)

Along-track plot of the simulated radar reflectivity and Doppler velocity.

plot_histogram(**kwds)

Histogram plot of the simulated radar reflectivity and Doppler velocity.

plot_scatter(**kwds)

Scatter plot between satellite data and suborbital data.

prepare_dates(start_date, end_date)

Creates a date array from start to end date.

prepare_input_dataset()

Prepares input dataset for computations.

range_to_height(ds)

Convert range coordinate to height coordinate by adding the station height above mean sea level to the range coordinate.

run(start_date, end_date[, write_output])

Runs simulation for all days in the time frame.

run_date(date[, write_output])

Runs simulation for a single day.

to_netcdf(date)

Writes dataset to netcdf file.

transform(ds)

Runs the entire simulator.

vm_uncertainty_equation(vm_std_broad, ...)

Calculates the total Doppler velocity uncertainty based on the broadening Doppler velocity uncertainty and the non-uniform beam filling Doppler velocity uncertainty.

add_airborne_variables()[source]

Add variables specific to airborne simulator to the dataset, i.e., the mean flight velocity.

add_attenuation(ds, da_gas_atten)[source]

Add attenuation to dataset.

Parameters:
dsxarray.Dataset

Data with “ze” variable. Unit: mm6 m-3

da_gas_attenxarray.DataArray

Interpolated gas attenuation data on the same grid as ds. Unit: dBZ.

Returns:
dsxarray.Dataset

Data with added attenuation.

add_ground_echo(ds)[source]

Calculates artificial ground echo inside ground-based radar range grid. The values are chosen such that the final ground echo after the along-range convolution is equal to the ground echo of the satellite. The pulse length used here is not the same as the pulse length of the satellite.

Parameters:
dsxarray.Dataset

Data with “ze” and “vm” variables.

Returns:
dsxarray.Dataset

Data with added ground echo.

add_groundbased_variables()[source]

Add variables specific to groundbased simulator to the dataset, i.e., the mean horizontal wind.

add_vmze_attrs(ds)[source]

Adds attributes to Doppler velocity and radar reflectivity variables.

Parameters:
dsxarray.Dataset

Data with “ze” and “vm” variables.

Returns:
dsxarray.Dataset

Data with added attributes.

attenuation_correction(ds, ds_cloudnet)[source]

Gas attenuation correction based on Cloudnet.

Cloudnet contains gas attenuation (gas_atten) as a function of 137 ERA5 levels. The height of each level varies with time. Therefore, the height is first interpolated onto the height grid of the radar data. Then, the gas attenuation is interpolated onto the time and height grid of the radar data. Finally, the gas attenuation is added to the radar reflectivity.

There are major differences between the cloudnet_ecmwf and cloudnet_categorize attenuation products: - ecmwf height is time-dependent, categorize height is not - ecmwf height is wrt ground, categorize height is wrt mean sea level - ecmwf variable is named “radar_gas_atten”, categorize variable is named “gas_atten” - ecmwf is calculated for both frequencies, categorize only for 94 GHz

Parameters:
dsxarray.Dataset

Data with “ze” variable. Unit: mm6 m-3

ds_cloudnetxarray.Dataset

Cloudnet data with “gas_atten” variable. Unit: dBZ

Returns:
dsxarray.Dataset

Data with added attenuation.

check_is_sea_level(ds)[source]

Check if input radar range/height grid is defined with respect to ground level or sea level. The input to the simulator should be wrt. sea level.

Parameters:
dsxr.Dataset

Input radar data

convert_frequency(ds)[source]

Convert frequency from 35 to 94 GHz.

The conversion is based on Kollias et al. (2019) (doi: https://doi.org/10.5194/amt-12-4949-2019)

Parameters:
dsxarray.Dataset

Data with “ze” variable in mm6/mm3. Ze was measured at 35 GHz.

Returns:
dsxarray.Dataset

Data with converted “ze” variable. Ze is now transformed to 94 GHz.

correct_dielectric_constant(ds)[source]

Apply correction for dielectric constant assumed in Ze calculation of suborbital radar to match the dielectric constant of the spaceborne radar.

Correction equation with \(K_g\) and \(K_s\) as dielectric constants of the suborbital and spaceborne radar, respectively:

\[Z_e = 10 \log_{10} \left( \frac{K_s}{K_g} \right) + Z_e\]
Parameters:
dsxarray.Dataset

Data with “ze” variable.

create_along_track(ds)[source]

Creates along-track coordinates from time coordinates.

Parameters:
dsxarray.Dataset

Data with “time” and “height” coordinates.

Returns:
dsxarray.Dataset

Data with “along_track” coordinate.

create_regular_along_track(ds)[source]

Creates regular along-track coordinate for suborbital radar.

Parameters:
dsxarray.Dataset

Data with “along_track” coordinate.

Returns:
xarray.DataArray

Regular along-track coordinate for suborbital radar.

create_regular_height()[source]

Creates regular height coordinate for suborbital radar.

Returns:
xarray.DataArray

Regular height coordinate for suborbital radar.

interpolate_to_regular_grid(ds)[source]

Interpolates radar data to regular grid in along-track and height.

Parameters:
dsxarray.Dataset

Data with “time” and “height” coordinates.

Returns:
dsxarray.Dataset

Data with interpolated “along_track” and “height” coordinates.

static prepare_dates(start_date, end_date)[source]

Creates a date array from start to end date.

Parameters:
start_datenp.datetime64

Start date.

end_datenp.datetime64

End date.

Returns:
dates: pd.DatetimeIndex

Date range from start to end date.

range_to_height(ds)[source]

Convert range coordinate to height coordinate by adding the station height above mean sea level to the range coordinate.

The altitude is pre-defined for each station in the configuration file.

Parameters:
dsxarray.Dataset

Data with “range” coordinate.

run(start_date, end_date, write_output=True)[source]

Runs simulation for all days in the time frame.

Parameters:
start_datenp.datetime64

Start date.

end_datenp.datetime64

End date (inclusive).

run_date(date, write_output=True)[source]

Runs simulation for a single day.

Parameters:
datenp.datetime64

Date to simulate.

write_outputbool

If True, write output to netcdf file.

property summary

Prints short summary of simulator settings.

to_netcdf(date)[source]

Writes dataset to netcdf file. Note that not all variables are stored.

Parameters:
datenp.datetime64

Date of the simulation. Used to create the filename.