Sionna
  • Installation
    • Using pip
    • From source
    • Testing
    • Documentation
    • Developing
  • Ray Tracing (RT)
    • Tutorials
      • Introduction to Sionna RT
        • Imports
        • Loading and Visualizing Scenes
        • Inspecting SceneObjects and Editing of Scenes
        • Ray tracing of Propagation Paths
        • From Paths to Channel Impulse and Frequency Responses
        • Radio Maps
        • Summary
      • Tutorial on Mobility
        • Background Information
        • GPU Configuration and Imports
        • Controlling Position and Orientation of Scene Objects
        • Time Evolution of Channels Via Doppler Shift
          • Example: Delay-Doppler Spectrum
        • Comparison of Doppler- vs Position-based Time Evolution
        • Summary
      • Tutorial on Radio Maps
        • Imports
        • Understanding radio maps
          • Metrics
          • Multiple transmitters
          • User association
          • Sampling of random user positions
          • Directional antennas and precoding vectors
        • Radio map for a realistic scene
          • Channel impulse responses for random user locations
        • Radio maps on meshes
        • Conclusions
      • Tutorial on Scattering
        • Imports
        • Scattering Basics
        • Scattering Patterns
        • Validation Against the “Far”-Wall Approximation
        • Radio Maps With Scattering
        • Impact on Channel Impulse Response
        • Summary
        • References
      • Tutorial on Loading and Editing of Scenes
        • Imports
        • Loading Scenes and Merging Objects
        • Editing Scenes
        • Path Computation with the Edited Scene
        • Summary
    • Primer on Electromagnetics
      • Coordinate system, rotations, and vector fields
      • Planar Time-Harmonic Waves
      • Far Field of a Transmitting Antenna
      • Modelling of a Receiving Antenna
      • General Propagation Path
      • Frequency & Impulse Response
      • Reflection and Refraction
        • Single-layer Slab
      • Diffraction
      • Scattering
      • Reconfigurable Intelligent Surfaces (RIS)
    • API Documentation
      • Antenna Arrays
        • AntennaArray
          • AntennaArray.antenna_pattern
          • AntennaArray.array_size
          • AntennaArray.normalized_positions
          • AntennaArray.num_ant
          • AntennaArray.positions()
          • AntennaArray.rotate()
        • PlanarArray
          • PlanarArray.antenna_pattern
          • PlanarArray.array_size
          • PlanarArray.normalized_positions
          • PlanarArray.num_ant
          • PlanarArray.positions()
          • PlanarArray.rotate()
          • PlanarArray.show()
      • Antenna Patterns
        • AntennaPattern
          • AntennaPattern.compute_gain()
          • AntennaPattern.patterns
          • AntennaPattern.show()
        • PolarizedAntennaPattern
          • PolarizedAntennaPattern.compute_gain()
          • PolarizedAntennaPattern.patterns
          • PolarizedAntennaPattern.show()
        • Vertically Polarized Antenna Pattern Functions
          • v_iso_pattern()
          • v_dipole_pattern()
          • v_hw_dipole_pattern()
          • v_tr38901_pattern()
        • Polarization Models
          • polarization_model_tr38901_1()
          • polarization_model_tr38901_2()
        • Utility Functions
          • antenna_pattern_to_world_implicit()
          • complex2real_antenna_pattern()
          • register_antenna_pattern()
          • register_polarization()
          • register_polarization_model()
      • Cameras
        • Camera
          • Camera.look_at()
          • Camera.orientation
          • Camera.position
      • Paths
        • Paths
          • Paths.a
          • Paths.cfr()
          • Paths.cir()
          • Paths.doppler
          • Paths.interactions
          • Paths.num_rx
          • Paths.num_tx
          • Paths.objects
          • Paths.phi_r
          • Paths.phi_t
          • Paths.primitives
          • Paths.rx_array
          • Paths.sources
          • Paths.synthetic_array
          • Paths.taps()
          • Paths.targets
          • Paths.tau
          • Paths.theta_r
          • Paths.theta_t
          • Paths.tx_array
          • Paths.valid
          • Paths.vertices
        • Constants
          • InteractionType
          • INVALID_SHAPE
          • INVALID_PRIMITIVE
      • Path Solvers
        • PathSolver
          • PathSolver.__call__()
          • PathSolver.loop_mode
      • Radio Devices
        • RadioDevice
          • RadioDevice.color
          • RadioDevice.display_radius
          • RadioDevice.look_at()
          • RadioDevice.name
          • RadioDevice.orientation
          • RadioDevice.position
          • RadioDevice.velocity
        • Receiver
        • Transmitter
          • Transmitter.power
          • Transmitter.power_dbm
      • Radio Maps
        • PlanarRadioMap
          • PlanarRadioMap.cell_centers
          • PlanarRadioMap.cell_size
          • PlanarRadioMap.cells_count
          • PlanarRadioMap.cells_per_dim
          • PlanarRadioMap.center
          • PlanarRadioMap.measurement_surface
          • PlanarRadioMap.orientation
          • PlanarRadioMap.path_gain
          • PlanarRadioMap.rx_cell_indices
          • PlanarRadioMap.sample_positions()
          • PlanarRadioMap.show()
          • PlanarRadioMap.show_association()
          • PlanarRadioMap.size
          • PlanarRadioMap.tx_association()
          • PlanarRadioMap.tx_cell_indices
        • MeshRadioMap
          • MeshRadioMap.cell_centers
          • MeshRadioMap.cells_count
          • MeshRadioMap.measurement_surface
          • MeshRadioMap.path_gain
          • MeshRadioMap.sample_positions()
        • RadioMap
          • RadioMap.cdf()
          • RadioMap.cell_centers
          • RadioMap.cells_count
          • RadioMap.measurement_surface
          • RadioMap.num_rx
          • RadioMap.num_tx
          • RadioMap.path_gain
          • RadioMap.rss
          • RadioMap.sample_cells()
          • RadioMap.sinr
          • RadioMap.transmitter_radio_map()
          • RadioMap.tx_association()
      • Radio Map Solvers
        • RadioMapSolver
          • RadioMapSolver.__call__()
          • RadioMapSolver.loop_mode
      • Radio Materials
        • RadioMaterialBase
          • RadioMaterialBase.color
          • RadioMaterialBase.eval()
          • RadioMaterialBase.name
          • RadioMaterialBase.pdf()
          • RadioMaterialBase.sample()
          • RadioMaterialBase.scene
          • RadioMaterialBase.to_string()
          • RadioMaterialBase.traverse()
        • RadioMaterial
          • RadioMaterial.conductivity
          • RadioMaterial.frequency_update()
          • RadioMaterial.frequency_update_callback
          • RadioMaterial.relative_permittivity
          • RadioMaterial.scattering_coefficient
          • RadioMaterial.scattering_pattern
          • RadioMaterial.thickness
          • RadioMaterial.xpd_coefficient
        • ITURadioMaterial
          • ITURadioMaterial.itu_type
          • ITURadioMaterial.to_string()
        • HolderMaterial
          • HolderMaterial.radio_material
          • HolderMaterial.velocity
        • Scattering Patterns
          • ScatteringPattern
          • LambertianPattern
          • DirectivePattern
          • BackscatteringPattern
          • register_scattering_pattern()
      • Scenes
        • load_scene()
        • Scene
          • Scene.add()
          • Scene.all_set()
          • Scene.angular_frequency
          • Scene.bandwidth
          • Scene.edit()
          • Scene.frequency
          • Scene.get()
          • Scene.mi_scene
          • Scene.mi_scene_params
          • Scene.objects
          • Scene.paths_solver
          • Scene.preview()
          • Scene.radio_materials
          • Scene.receivers
          • Scene.remove()
          • Scene.render()
          • Scene.render_to_file()
          • Scene.rx_array
          • Scene.scene_geometry_updated()
          • Scene.sources()
          • Scene.targets()
          • Scene.temperature
          • Scene.thermal_noise_power
          • Scene.transmitters
          • Scene.tx_array
          • Scene.wavelength
        • Built-in scenes
          • box
          • box_one_screen
          • box_two_screens
          • double_reflector
          • etoile
          • floor_wall
          • florence
          • munich
          • san_francisco
          • simple_reflector
          • simple_street_canyon
          • simple_street_canyon_with_cars
          • simple_wedge
          • triple_reflector
        • Built-in meshes
          • low_poly_car
          • sphere
      • Scene Objects
        • SceneObject
          • SceneObject.clone()
          • SceneObject.look_at()
          • SceneObject.mi_mesh
          • SceneObject.name
          • SceneObject.object_id
          • SceneObject.orientation
          • SceneObject.position
          • SceneObject.radio_material
          • SceneObject.scaling
          • SceneObject.scene
          • SceneObject.velocity
      • Utility Functions
        • Complex-valued tensors
          • cpx_abs()
          • cpx_abs_square()
          • cpx_add()
          • cpx_convert()
          • cpx_div()
          • cpx_exp()
          • cpx_mul()
          • cpx_sqrt()
          • cpx_sub()
        • Electromagnetics
          • complex_relative_permittivity()
          • fresnel_reflection_coefficients_simplified()
          • itu_coefficients_single_layer_slab()
        • Geometry
          • phi_hat()
          • theta_hat()
          • r_hat()
          • theta_phi_from_unit_vec()
          • rotation_matrix()
        • Jones calculus
          • implicit_basis_vector()
          • jones_matrix_rotator()
          • jones_matrix_rotator_flip_forward()
          • to_world_jones_rotator()
          • jones_matrix_to_world_implicit()
          • jones_vec_dot()
        • Meshes
          • load_mesh()
          • transform_mesh()
        • Miscellaneous
          • complex_sqrt()
          • dbm_to_watt()
          • isclose()
          • log10()
          • sigmoid()
          • sinc()
          • subcarrier_frequencies()
          • watt_to_dbm()
        • Ray tracing
          • fibonacci_lattice()
          • spawn_ray_from_sources()
          • offset_p()
          • spawn_ray_towards()
          • spawn_ray_to()
    • Developer Guides
      • Compatibility with other Frameworks
        • Type conversions
        • Gradients
        • Training-Loop in PyTorch
      • Understanding the Paths Object
      • Custom Antenna Patterns
        • Gradient-based Optimization
      • Understanding Radio Materials
        • Modifying Parameters of Radio Materials
        • Calibrating Material Parameters Through Gradient Descent
        • Custom Radio Materials
          • Representation of Jones vector and Matrices
          • Implicit Basis
          • The Local Interaction Basis
          • Mandatory Subclass Methods
          • Implementation of a Simple Radio Material Model
          • A More Complex Material Model
      • Custom Scattering Patterns
        • Differentiable Parameters
    • Technical Report
  • Physical Layer (PHY)
    • Tutorials
      • Beginners
        • “Hello, world!”
        • Part 1: Getting Started with Sionna
          • Imports & Basics
          • A note on random number generation
          • Sionna Data-flow and Design Paradigms
          • Hello, Sionna!
          • Communication Systems as Sionna Blocks
          • Forward Error Correction (FEC)
          • Eager vs Graph Mode
          • Exercise
        • Part 2: Differentiable Communication Systems
          • Imports
          • Gradient Computation Through End-to-end Systems
          • Creating Custom Layers
          • Setting up Training Loops
        • Part 3: Advanced Link-level Simulations
          • Imports
          • OFDM Resource Grid and Stream Management
          • Antenna Arrays
          • Channel Model
          • Uplink Transmission in the Frequency Domain
        • Part 4: Toward Learned Receivers
          • Imports
          • Simulation Parameters
          • Implemention of an Advanced Neural Receiver
          • Training the Neural Receiver
          • Benchmarking the Neural Receiver
          • Conclusion
          • References
        • Basic MIMO Simulations
          • Table of Contents
          • Simple uncoded transmission
          • Extension to channel coding
        • Pulse-shaping Basics
          • Table of Contents
          • GPU Configuration and Imports
          • Pulse-shaping of a sequence of QAM symbols
          • Recovering the QAM symbols through matched filtering and downsampling
          • Investigating the ACLR
          • Windowing
        • Optical Channel with Lumped Amplification
          • Table of Contents
          • Setup
          • Impulse Generation
          • Attenuation
          • Amplified Spontaneous Emission Noise
          • Chromatic Dispersion
          • Kerr Nonlinearity
          • Split-Step Fourier Method
          • References
      • Experts
        • 5G Channel Coding and Rate-Matching: Polar vs. LDPC Codes
          • Table of Contents
          • GPU Configuration and Imports
          • BER Performance of 5G Coding Schemes
          • A Deeper Look into the Polar Code Module
          • Rate-Matching and Rate-Recovery
          • Throughput and Decoding Complexity
          • References
        • 5G NR PUSCH Tutorial
          • Table of Contents
          • GPU Configuration and Imports
          • A Hello World Example
          • Carrier Configuration
          • Understanding the DMRS Configuration
          • Transport Blocks and MCS
          • Looking into the PUSCHTransmitter
          • Components of the PUSCHReceiver
          • End-to-end PUSCH Simulations
        • Bit-Interleaved Coded Modulation (BICM)
          • Table of Contents
          • System Block Diagram
          • GPU Configuration and Imports
          • A Simple BICM System
          • All-zero Codeword Simulations
          • EXIT Charts
          • Mismatched Demapping and the Advantages of Min-sum Decoding
          • References
        • MIMO OFDM Transmissions over the CDL Channel Model
          • Table of Contents
          • System Setup
          • Simulations
        • Neural Receiver for OFDM SIMO Systems
          • GPU Configuration and Imports
          • Simulation Parameters
          • Neural Receiver
          • End-to-end System
          • End-to-end System as a Sionna Block
          • Evaluation of the Baselines
          • Training the Neural Receiver
          • Evaluation of the Neural Receiver
          • Pre-computed Results
          • References
        • Realistic Multiuser MIMO OFDM Simulations
          • Table of Contents
          • GPU Configuration and Imports
          • System Setup
          • Uplink Transmissions in the Frequency Domain
        • OFDM MIMO Channel Estimation and Detection
          • Table of Contents
          • GPU Configuration and Imports
          • Simulation parameters
          • Estimation of the channel time, frequency, and spatial covariance matrices
          • Loading the channel covariance matrices
          • Comparison of OFDM estimators
          • Comparison of MIMO detectors
        • Introduction to Iterative Detection and Decoding
          • Iterative Detection and Decoding
          • Table of contents
          • GPU Configuration and Imports
          • Simulation Parameters
          • Setting-up an End-to-end Block
          • Non-IDD versus IDD Benchmarks
          • Discussion-Optimizing IDD with Machine Learning
          • Comments
          • List of References
        • End-to-end Learning with Autoencoders
          • GPU Configuration and Imports
          • Simulation Parameters
          • Neural Demapper
          • Trainable End-to-end System: Conventional Training
          • Trainable End-to-end System: RL-based Training
          • Evaluation
          • Visualizing the Learned Constellations
          • References
        • Weighted Belief Propagation Decoding
          • Table of Contents
          • GPU Configuration and Imports
          • Weighted BP for BCH Codes
          • Further Experiments
          • References
        • Channel Models from Datasets
          • GPU Configuration and Imports
          • Simulation Parameters
          • Creating a Simple Dataset
          • Generators
          • Use the Channel Model for OFDM Transmissions
        • Using the DeepMIMO Dataset with Sionna
          • Table of Contents
          • GPU Configuration and Imports
          • Configuration of DeepMIMO
          • Using DeepMIMO with Sionna
          • Link-level Simulations using Sionna and DeepMIMO
          • DeepMIMO License and Citation
        • Link-level simulations with Sionna RT
          • Background Information
          • Imports
          • Setting up the Ray Tracer
          • Creating a CIR Dataset
          • PUSCH Link-Level Simulations
    • API Documentation
      • Configuration
        • Config
          • Config.np_cdtype
          • Config.np_rdtype
          • Config.np_rng
          • Config.precision
          • Config.py_rng
          • Config.seed
          • Config.tf_cdtype
          • Config.tf_rdtype
          • Config.tf_rng
      • Forward Error Correction (FEC)
        • Linear Codes
          • LinearEncoder
          • OSDecoder
        • Low-Density Parity-Check (LDPC)
          • LDPC5GEncoder
          • LDPCBPDecoder
          • LDPC5GDecoder
          • Node Update Functions
          • Decoder Callbacks
        • Polar Codes
          • Polar5GEncoder
          • PolarEncoder
          • Polar5GDecoder
          • PolarSCDecoder
          • PolarSCLDecoder
          • PolarBPDecoder
          • Utility Functions
        • Convolutional Codes
          • ConvEncoder
          • ViterbiDecoder
          • BCJRDecoder
          • Utility Functions
        • Turbo Codes
          • TurboEncoder
          • TurboDecoder
          • Utility Functions
        • Cyclic Redundancy Check (CRC)
          • CRCEncoder
          • CRCDecoder
        • Interleaving
          • RowColumnInterleaver
          • RandomInterleaver
          • Turbo3GPPInterleaver
          • Deinterleaver
        • Scrambling
          • Scrambler
          • TB5GScrambler
          • Descrambler
        • Utility Functions
          • (Binary) Linear Codes
          • EXIT Analysis
          • Miscellaneous
      • Mapping
        • Constellation
          • Constellation.center
          • Constellation.constellation_type
          • Constellation.normalize
          • Constellation.num_bits_per_symbol
          • Constellation.num_points
          • Constellation.points
          • Constellation.show()
        • Mapper
          • Mapper.constellation
        • Demapper
          • Demapper.constellation
        • SymbolDemapper
        • Utility Functions
          • BinarySource
          • LLRs2SymbolLogits
          • PAMSource
          • PAM2QAM
          • QAMSource
          • QAM2PAM
          • SymbolInds2Bits
          • SymbolLogits2LLRs
          • SymbolLogits2Moments
          • SymbolSource
      • Wireless Channel Models
        • AWGN
        • Flat-fading channel
          • FlatFadingChannel
          • GenerateFlatFadingChannel
          • ApplyFlatFadingChannel
          • SpatialCorrelation
          • KroneckerModel
          • PerColumnModel
        • Channel model interface
          • ChannelModel
        • Time domain channel
          • TimeChannel
          • GenerateTimeChannel
          • ApplyTimeChannel
          • cir_to_time_channel()
          • time_to_ofdm_channel()
        • Channel with OFDM waveform
          • OFDMChannel
          • GenerateOFDMChannel
          • ApplyOFDMChannel
          • cir_to_ofdm_channel()
        • Rayleigh block fading
          • RayleighBlockFading
        • 3GPP 38.901 channel models
          • PanelArray
          • Antenna
          • AntennaArray
          • TDL
          • CDL
          • UMi
          • UMa
          • RMa
        • External datasets
          • CIRDataset
        • Utility functions
          • subcarrier_frequencies()
          • time_lag_discrete_time_channel()
          • deg_2_rad()
          • rad_2_deg()
          • wrap_angle_0_360()
          • drop_uts_in_sector()
          • relocate_uts()
          • set_3gpp_scenario_parameters()
          • gen_single_sector_topology()
          • gen_single_sector_topology_interferers()
          • exp_corr_mat()
          • one_ring_corr_mat()
      • Optical Channel Models
        • SSFM
        • EDFA
        • Utility functions
          • time_frequency_vector()
      • Discrete Channel Models
        • BinaryErasureChannel
        • BinaryMemorylessChannel
          • BinaryMemorylessChannel.llr_max
          • BinaryMemorylessChannel.temperature
        • BinarySymmetricChannel
        • BinaryZChannel
      • Orthogonal Frequency-Division Multiplexing (OFDM)
        • Resource Grid
          • ResourceGrid
          • ResourceGridMapper
          • ResourceGridDemapper
          • RemoveNulledSubcarriers
        • Modulation & Demodulation
          • OFDMModulator
          • OFDMDemodulator
        • Pilot Pattern
          • PilotPattern
          • EmptyPilotPattern
          • KroneckerPilotPattern
        • Channel Estimation
          • BaseChannelEstimator
          • BaseChannelInterpolator
          • LSChannelEstimator
          • LinearInterpolator
          • LMMSEInterpolator
          • NearestNeighborInterpolator
          • tdl_time_cov_mat()
          • tdl_freq_cov_mat()
        • Precoding
          • RZFPrecoder
          • PrecodedChannel
          • CBFPrecodedChannel
          • EyePrecodedChannel
          • RZFPrecodedChannel
        • Equalization
          • OFDMEqualizer
          • LMMSEEqualizer
          • MFEqualizer
          • ZFEqualizer
          • PostEqualizationSINR
          • LMMSEPostEqualizationSINR
        • Detection
          • OFDMDetector
          • OFDMDetectorWithPrior
          • EPDetector
          • KBestDetector
          • LinearDetector
          • MaximumLikelihoodDetector
          • MaximumLikelihoodDetectorWithPrior
          • MMSEPICDetector
      • Multiple-Input Multiple-Output (MIMO)
        • Stream Management
          • StreamManagement
        • Precoding
          • cbf_precoding_matrix()
          • rzf_precoding_matrix()
          • rzf_precoder()
          • grid_of_beams_dft_ula()
          • grid_of_beams_dft()
          • flatten_precoding_mat()
          • normalize_precoding_power()
        • Equalization
          • lmmse_matrix()
          • lmmse_equalizer()
          • mf_equalizer()
          • zf_equalizer()
        • Detection
          • EPDetector
          • KBestDetector
          • LinearDetector
          • MaximumLikelihoodDetector
          • MMSEPICDetector
        • Utility Functions
          • List2LLR
          • List2LLRSimple
          • complex2real_vector()
          • real2complex_vector()
          • complex2real_matrix()
          • real2complex_matrix()
          • complex2real_covariance()
          • real2complex_covariance()
          • complex2real_channel()
          • real2complex_channel()
          • whiten_channel()
      • 5G NR
        • Carrier
          • CarrierConfig
        • Layer Mapping
          • LayerMapper
          • LayerDemapper
        • PUSCH
          • PUSCHConfig
          • PUSCHDMRSConfig
          • PUSCHLSChannelEstimator
          • PUSCHPilotPattern
          • PUSCHPrecoder
          • PUSCHReceiver
          • PUSCHTransmitter
        • Transport Block
          • TBConfig
          • TBEncoder
          • TBDecoder
        • Utils
          • calculate_tb_size()
          • generate_prng_seq()
          • decode_mcs_index()
          • calculate_num_coded_bits()
          • TransportBlockNR
          • CodedAWGNChannelNR
          • MCSDecoderNR
      • Signal
        • Filters
          • SincFilter
          • RaisedCosineFilter
          • RootRaisedCosineFilter
          • CustomFilter
          • Filter
        • Window functions
          • HannWindow
          • HammingWindow
          • BlackmanWindow
          • CustomWindow
          • Window
        • Utility Functions
          • convolve()
          • fft()
          • ifft()
          • Upsampling
          • Downsampling
          • empirical_psd()
          • empirical_aclr()
      • Utility Functions
        • Linear Algebra
          • inv_cholesky()
          • matrix_pinv()
        • Metrics
          • compute_ber()
          • compute_bler()
          • compute_ser()
          • count_block_errors()
          • count_errors()
        • Miscellaneous
          • dbm_to_watt()
          • db_to_lin
          • DeepUpdateDict
          • dict_keys_to_int()
          • ebnodb2no()
          • complex_normal()
          • hard_decisions()
          • Interpolate
          • lin_to_db
          • log2()
          • log10()
          • MCSDecoder
          • scalar_to_shaped_tensor()
          • sim_ber()
          • SingleLinkChannel
          • SplineGriddataInterpolation
          • to_list()
          • TransportBlock
          • watt_to_dbm()
        • Numerics
          • bisection_method()
        • Plotting
          • plot_ber()
          • PlotBER
        • Tensors
          • expand_to_rank()
          • flatten_dims()
          • flatten_last_dims()
          • insert_dims()
          • split_dim()
          • diag_part_axis()
          • flatten_multi_index()
          • gather_from_batched_indices()
          • tensor_values_are_in_set()
          • enumerate_indices()
      • For Developers
        • Object
          • Object.cdtype
          • Object.precision
          • Object.rdtype
        • Block
          • Block.build()
          • Block.built
          • Block.call()
    • Developer Guides
      • Matrix inversion
        • Solving linear systems
        • Correlated random vectors
      • Random number generation
      • Sionna Block and Object
        • Understanding Sionna Blocks
  • System Level (SYS)
    • Tutorials
      • Beginners
        • Physical Layer Abstraction
          • Imports
          • Instantiate a PHYAbstraction object
          • Retrieve BLER values from interpolated tables
          • Generate a new BLER table
          • Bypass physical layer computations
          • Effective SINR
          • Conclusions
          • References
        • Link Adaptation
          • Imports
          • Simulation parameters
          • The principles of link adaptation
          • Inner-Loop Link Adaptation (ILLA)
          • Outer-Loop Link Adaption (OLLA)
          • Conclusions
          • References
        • Proportional Fairness Scheduler
          • Imports
          • The main principle
          • Basic scenario
          • Schedule users
          • Evaluate performance
          • Conclusions
          • References
        • Hexagonal Grid Topology
          • Imports
          • Generate a multicell topology
          • Drop users
          • Set up a 3GPP multicell scenario
          • Per-stream SINR computation
          • Conclusions
        • Power Control
          • Imports
          • Multicell scenario
          • Uplink power control
          • Downlink power control
          • Conclusions
          • References
      • Experts
        • Sionna SYS meets Sionna RT
          • Imports
          • Simulation parameters
          • Scene creation
          • Channel generation via Sionna RT
          • System-level simulation
          • Conclusions
        • System-Level Simulations
          • Imports
          • Utils
          • Simulation
          • Performance metric analysis
          • Conclusions
    • API Documentation
      • PHY Abstraction
        • EffectiveSINR
          • EffectiveSINR.calibrate()
        • EESM
          • EESM.beta_table
          • EESM.beta_table_filenames
          • EESM.beta_tensor
          • EESM.validate_beta_table()
        • PHYAbstraction
          • PHYAbstraction.bler_interp_delta
          • PHYAbstraction.bler_table
          • PHYAbstraction.bler_table_filenames
          • PHYAbstraction.bler_table_interp
          • PHYAbstraction.cbs_interp_min_max_delta
          • PHYAbstraction.get_bler()
          • PHYAbstraction.get_idx_from_grid()
          • PHYAbstraction.load_table()
          • PHYAbstraction.new_bler_table()
          • PHYAbstraction.plot()
          • PHYAbstraction.snr_db_interp_min_max_delta
          • PHYAbstraction.snr_table_interp
          • PHYAbstraction.validate_bler_table()
      • Link Adaptation
        • InnerLoopLinkAdaptation
          • InnerLoopLinkAdaptation.bler_target
        • OuterLoopLinkAdaptation
          • OuterLoopLinkAdaptation.bler_target
          • OuterLoopLinkAdaptation.delta_down
          • OuterLoopLinkAdaptation.delta_up
          • OuterLoopLinkAdaptation.offset
          • OuterLoopLinkAdaptation.offset_max
          • OuterLoopLinkAdaptation.offset_min
          • OuterLoopLinkAdaptation.reset()
          • OuterLoopLinkAdaptation.sinr_eff_db_last
        • Note
      • Power Control
        • Uplink
          • open_loop_uplink_power_control()
        • Downlink
          • downlink_fair_power_control()
      • Scheduling
        • PFSchedulerSUMIMO
          • PFSchedulerSUMIMO.beta
          • PFSchedulerSUMIMO.pf_metric
          • PFSchedulerSUMIMO.rate_achieved_past
      • Multicell Topology
        • Hexagon
          • Hexagon.coord_axial
          • Hexagon.coord_dict()
          • Hexagon.coord_euclid
          • Hexagon.coord_offset
          • Hexagon.corners()
          • Hexagon.neighbor()
          • Hexagon.radius
        • HexGrid
          • HexGrid.cell_height
          • HexGrid.cell_loc
          • HexGrid.cell_radius
          • HexGrid.center_loc
          • HexGrid.grid
          • HexGrid.isd
          • HexGrid.mirror_cell_loc
          • HexGrid.num_cells
          • HexGrid.num_rings
          • HexGrid.show()
        • gen_hexgrid_topology()
        • get_num_hex_in_grid()
        • convert_hex_coord()
      • Utils
        • get_pathloss()
        • is_scheduled_in_slot()
        • spread_across_subcarriers()
  • Research Kit (RK)
    • Quickstart
      • Hardware Requirements
      • Step 1: Jetson Setup
      • Step 2: USRP Setup
      • Step 3: UE Setup
      • Step 4: Deploy 5G Stack
      • Have Your First Call
    • Setup
      • Platform Preparation
        • Bill of Materials
          • Computing Platform
          • RF Components
          • User Equipment
          • Optional Components
        • Jetson Setup
          • OS Installation
          • Post-Installation Setup
          • Version Information
        • Custom Linux Kernel
          • Prerequisites
          • Source Code
          • Kernel Configuration
          • Building the Kernel
          • Installing Kernel Image and Modules
          • Configure Boot Sequence
        • Performance Tweaks
          • Power Management
          • Default System Mode
          • Real-Time Scheduling
          • Verifying Settings
        • USRP Driver Installation
          • Key Resources
          • Prerequisites
          • Building UHD
          • Post-build Configuration
          • Testing Installation
        • SIM Card Programming
          • Prerequisites
          • UICC Software Setup
          • Programming SIM Card
          • Additional Resources
        • Quectel Modem Setup
          • Basic Configuration
          • AT Command Reference
          • Additional Resources
      • Software Configuration
        • OpenAirInterface Setup
          • Required Components
          • TLDR
          • Getting the Source Code
          • Patching
          • Create configuration files
          • Building OAI RAN Images
          • Building OAI 5G Core Images
          • Building Components Manually
        • 5G System Configuration
          • Environment Variables
          • Additional Resources
        • Using RF Simulator Mode
          • Basic Configuration
          • Dynamic Re-configuration
        • Installing Sionna
      • Your First Call
        • Connect & Test Performance
      • Documentation of Scripts
        • Setup Scripts
          • quickstart-oai.sh
          • quickstart-cn5g.sh
          • configure-system.sh
          • install-usrp.sh
        • Linux Kernel Customization
          • build-custom-kernel.sh
          • install-custom-kernel.sh
        • Configuration Files
          • start-system.sh
          • stop-system.sh
          • generate-configs.sh
        • Building Docker Images
          • build-cn5g-images.sh
          • build-oai-images.sh
        • Development / Patch tracking
          • get-config-changes.sh
          • get-oai-changed-files.sh
          • get-oai-cn5g-changed-files.sh
          • get-oai-commit-versions.sh
          • get-oai-cn5g-commit-versions.sh
    • Tutorials
      • Running the Tutorials
        • Command Cheat-sheet
        • GPU-Accelerated LDPC
        • Demapper Plugin
        • Data Capture Plugin
        • TensorRT Neural Demapper
      • GPU-Accelerated LDPC Decoding
        • Part 1: Background & Python Reference Implementation
          • Background: Channel Coding in 5G
          • Overview Decoder Implementation
          • Memory Layout
          • CN Update Function
          • VN Update Function
          • Hard-decision
          • Main Decoding Function
          • Run and Test the Decoder
        • Part 2: CUDA Implementation
          • Overview
          • CUDA Integration in OAI
          • Running the Decoder
          • Implementation Aspects
          • Outlook - Weighted Belief Propagation
        • References
      • Plugins & Data Acquisition
        • Part 1: Create a Plugin
          • Architecture Overview
          • Select Functions
          • Define Module Interface
          • Loading the Module
          • Use Module Functions
          • Module Implementation
          • Compiling
          • Incremental Builds
          • Container Changes
        • Part 2: Capture Data
          • Adding New Plugin Variant
          • Using New Plugin Variant
          • Access Captured Files
          • Capture Format
          • Add Timestamps
          • Account for Multi-Threading
          • Final Source Code
        • References
      • Integration of a Neural Demapper
        • Part 1: Neural Demapper Training and TensorRT Export
          • Python Imports
          • Background: APP-Demapping
          • Understanding the OAI Data Structure
          • Neural Demapper
          • Blind Demapping & Training with Synthetic Data
          • Improved Demapping & Training with Captured Data
          • Export TensorRT Engine
        • Part 2: GPU-Accelerated Inference
          • Demapper Implementation Overview
          • Setting up the TensorRT Inference Engine
          • Running Batched Inference
          • Converting Data Types between Host and Device
          • Demapper Integration in OAI
          • Running the Demapper
          • Implementation Aspects
          • Unit tests
          • Outlook
        • References
      • Software-defined End-to-End 5G Network
        • Run the gNB
        • Run the UE
        • Test performance
      • Debugging & Troubleshooting
        • Attaching a debugger (gdb and VS code)
        • Inspecting and debugging inside a container interactively
        • Running memcheck within an interactive Docker compose session
        • Profiling with NVIDIA Nsight Systems
        • Fixing missing linker error messages in docker build ran-build
    • Get the Code
  • “Made with Sionna”
  • Discussions
  • Report an Issue
  • Contribute
  • Citation

Older Versions

  • v0.19.2
Sionna
  • Physical Layer (PHY)
  • Tutorials
  • Link-level simulations with Sionna RT

Colab logo Run in Google Colab View on GitHub Download notebook

Link-level simulations with Sionna RT

In this notebook, you will use ray-traced channels for link-level simulations instead of stochastic channel models

Background Information

Ray tracing is a technique to simulate environment-specific and physically accurate channel realizations for a given scene and user position. Sionna RT is a ray tracing extension for radio propagation modeling which is built on top of Mitsuba 3. Like all of Sionna’s components, it is differentiable.

For an introduction about how to use Sionna RT, please see the corresponding tutorials. The EM Primer provides further details on the theoretical background of ray tracing of wireless channels.

In this notebook, we will use Sionna RT for site-specific link-level simulations. For this, we evaluate the BER performance for a MU-MIMO 5G NR system in the uplink direction based on ray traced CIRs for random user positions.

We use the 5G NR PUSCH transmitter and receiver from the 5G NR PUSCH Tutorial notebook. Note that also the systems from the MIMO OFDM Transmissions over the CDL Channel Model or the Neural Receiver for OFDM SIMO Systems tutorials could be used instead.

There are different ways to implement uplink scenarios in Sionna RT. In this example, we configure the basestation as transmitter and the user equipments (UEs) as receivers which simplifies the ray tracing. Due to channel reciprocity, one can reverse the direction of the ray traced channels afterwards. For the ray tracer itself, the direction (uplink/downlink) does not change the simulated paths.

Imports

[1]:
import os # Configure which GPU
if os.getenv("CUDA_VISIBLE_DEVICES") is None:
    gpu_num = 0 # Use "" to use the CPU
    os.environ["CUDA_VISIBLE_DEVICES"] = f"{gpu_num}"
os.environ['TF_CPP_MIN_LOG_LEVEL'] = '3'

# Import or install Sionna
try:
    import sionna.phy
    import sionna.rt
except ImportError as e:
    import sys
    if 'google.colab' in sys.modules:
       # Install Sionna in Google Colab
       print("Installing Sionna and restarting the runtime. Please run the cell again.")
       os.system("pip install sionna")
       os.kill(os.getpid(), 5)
    else:
       raise e

# Configure the notebook to use only a single GPU and allocate only as much memory as needed
# For more details, see https://www.tensorflow.org/guide/gpu
import tensorflow as tf
gpus = tf.config.list_physical_devices('GPU')
if gpus:
    try:
        tf.config.experimental.set_memory_growth(gpus[0], True)
    except RuntimeError as e:
        print(e)

# Avoid warnings from TensorFlow
tf.get_logger().setLevel('ERROR')

import numpy as np

# For link-level simulations
from sionna.phy.channel import OFDMChannel, CIRDataset
from sionna.phy.nr import PUSCHConfig, PUSCHTransmitter, PUSCHReceiver
from sionna.phy.utils import ebnodb2no, PlotBER
from sionna.phy.ofdm import KBestDetector, LinearDetector
from sionna.phy.mimo import StreamManagement

# Import Sionna RT components
from sionna.rt import load_scene, Camera, Transmitter, Receiver, PlanarArray,\
                      PathSolver, RadioMapSolver

no_preview = True # Toggle to False to use the preview widget
                  # instead of rendering for scene visualization

Setting up the Ray Tracer

Let’s start by defining some constants that control the system we want to simulate.

[2]:
# System parameters
subcarrier_spacing = 30e3 # Hz
num_time_steps = 14 # Total number of ofdm symbols per slot

num_tx = 4 # Number of users
num_rx = 1 # Only one receiver considered
num_tx_ant = 4 # Each user has 4 antennas
num_rx_ant = 16 # The receiver is equipped with 16 antennas

# batch_size for CIR generation
batch_size_cir = 1000

We then set up the radio propagation environment. We start by loading a scene and then add a transmitter that acts as a base station. We will later use channel reciprocity to simulate the uplink direction.

[3]:
# Load an integrated scene.
# You can try other scenes, such as `sionna.rt.scene.etoile`. Note that this would require
# updating the position of the transmitter (see below in this cell).
scene = load_scene(sionna.rt.scene.munich)

# Transmitter (=basestation) has an antenna pattern from 3GPP 38.901
scene.tx_array = PlanarArray(num_rows=1,
                             num_cols=num_rx_ant//2, # We want to transmitter to be equiped with 16 antennas
                             vertical_spacing=0.5,
                             horizontal_spacing=0.5,
                             pattern="tr38901",
                             polarization="cross")

# Create transmitter
tx = Transmitter(name="tx",
                 position=[8.5,21,27],
                 look_at=[45,90,1.5], # optional, defines view direction
                 display_radius=3.) # optinal, radius of the sphere for visualizing the device
scene.add(tx)

# Create new camera
bird_cam = Camera(position=[0,80,500], orientation=np.array([0,np.pi/2,-np.pi/2]))

We then compute a radio map for the instantiated transmitter.

[4]:
max_depth = 5

# Radio map solver
rm_solver = RadioMapSolver()

# Compute the radio map
rm = rm_solver(scene,
               max_depth=5,
               cell_size=(1., 1.),
               samples_per_tx=10**7)

Let’s visualize the computed radio map.

[5]:
if no_preview:
    # Render an image
    scene.render(camera=bird_cam,
                 radio_map=rm,
                 rm_vmin=-110,
                 clip_at=12.); # Clip the scene at rendering for visualizing the refracted field
else:
    # Show preview
    scene.preview(radio_map=rm,
                  rm_vmin=-110,
                  clip_at=12.); # Clip the scene at rendering for visualizing the refracted field
../../_images/phy_tutorials_Link_Level_Simulations_with_RT_14_0.png

The function RadioMap.sample_positions() allows sampling of random user positions from a radio map. It ensures that only positions that have a path gain of at least min_gain_dB dB and at most max_gain_dB dB are sampled, i.e., it ignores positions without a connection to the transmitter. Further, one can set the distances min_dist and max_dist to sample only points within a certain distance range from the transmitter.

[6]:
min_gain_db = -130 # in dB; ignore any position with less than -130 dB path gain
max_gain_db = 0 # in dB; ignore strong paths

# Sample points in a 5-400m range around the receiver
min_dist = 5 # in m
max_dist = 400 # in m

# Sample batch_size random user positions from the radio map
ue_pos, _ = rm.sample_positions(num_pos=batch_size_cir,
                                metric="path_gain",
                                min_val_db=min_gain_db,
                                max_val_db=max_gain_db,
                                min_dist=min_dist,
                                max_dist=max_dist)

We now add new receivers (=UEs) at the sampled positions.

Remark: This is an example for 5G NR PUSCH (uplink direction), we will reverse the direction of the channel for later BER simulations.

[7]:
# Configure antenna array for all receivers (=UEs)
scene.rx_array = PlanarArray(num_rows=1,
                             num_cols=num_tx_ant//2, # Each receiver is equipped with 4 antennas
                             vertical_spacing=0.5,
                             horizontal_spacing=0.5,
                             pattern="iso",
                             polarization="cross")

# Create batch_size receivers
for i in range(batch_size_cir):
    scene.remove(f"rx-{i}") # Remove old receiver if any
    rx = Receiver(name=f"rx-{i}",
                  position=ue_pos[0][i], # Position sampled from radio map
                  velocity=(3.,3.,0),
                  display_radius=1., # optional, radius of the sphere for visualizing the device
                  color=(1,0,0) # optional, color for visualizing the device
                  )
    scene.add(rx)

# And visualize the scene
if no_preview:
    # Render an image
    scene.render(camera=bird_cam,
                 radio_map=rm,
                 rm_vmin=-110,
                 clip_at=12.); # Clip the scene at rendering for visualizing the refracted field
else:
    # Show preview
    scene.preview(radio_map=rm,
                  rm_vmin=-110,
                  clip_at=12.); # Clip the scene at rendering for visualizing the refracted field
../../_images/phy_tutorials_Link_Level_Simulations_with_RT_18_0.png

Each dot represents a receiver position drawn from the random sampling function of the radio map. This allows to efficiently sample batches of random channel realizations even in complex scenarios.

Creating a CIR Dataset

We can now simulate the CIRs for many different positions which will be used later on as a dataset for link-level simulations.

Remark: Running the cells below can take some time depending on the requested number of CIRs.

[8]:
target_num_cirs = 5000 # Defines how many different CIRs are generated.
# Remark: some path are removed if no path was found for this position

max_depth = 5
min_gain_db = -130 # in dB / ignore any position with less than -130 dB path gain
max_gain_db = 0 # in dB / ignore any position with more than 0 dB path gain

# Sample points within a 10-400m radius around the transmitter
min_dist = 10 # in m
max_dist = 400 # in m

# List of channel impulse reponses
a_list = []
tau_list = []

# Maximum number of paths over all batches of CIRs.
# This is used later to concatenate all CIRs.
max_num_paths = 0

# Path solver
p_solver = PathSolver()

# Each simulation returns batch_size_cir results
num_runs = int(np.ceil(target_num_cirs/batch_size_cir))
for idx in range(num_runs):
    print(f"Progress: {idx+1}/{num_runs}", end="\r")

    # Sample random user positions
    ue_pos, _ = rm.sample_positions(
                        num_pos=batch_size_cir,
                        metric="path_gain",
                        min_val_db=min_gain_db,
                        max_val_db=max_gain_db,
                        min_dist=min_dist,
                        max_dist=max_dist,
                        seed=idx) # Change the seed from one run to the next to avoid sampling the same positions

    # Update all receiver positions
    for rx in range(batch_size_cir):
        scene.receivers[f"rx-{rx}"].position = ue_pos[0][rx]

    # Simulate CIR
    paths = p_solver(scene, max_depth=max_depth, max_num_paths_per_src=10**7)

    # Transform paths into channel impulse responses
    a, tau = paths.cir(sampling_frequency=subcarrier_spacing,
                         num_time_steps=14,
                         out_type='numpy')
    a_list.append(a)
    tau_list.append(tau)

    # Update maximum number of paths over all batches of CIRs
    num_paths = a.shape[-2]
    if num_paths > max_num_paths:
        max_num_paths = num_paths

# Concatenate all the CIRs into a single tensor along the num_rx dimension.
# First, we need to pad the CIRs to ensure they all have the same number of paths.
a = []
tau = []
for a_,tau_ in zip(a_list, tau_list):
    num_paths = a_.shape[-2]
    a.append(np.pad(a_, [[0,0],[0,0],[0,0],[0,0],[0,max_num_paths-num_paths],[0,0]], constant_values=0))
    tau.append(np.pad(tau_, [[0,0],[0,0],[0,max_num_paths-num_paths]], constant_values=0))
a = np.concatenate(a, axis=0) # Concatenate along the num_rx dimension
tau = np.concatenate(tau, axis=0)

# Let's now convert to uplink direction, by switing the receiver and transmitter
# dimensions
a = np.transpose(a, (2,3,0,1,4,5))
tau = np.transpose(tau, (1,0,2))

# Add a batch_size dimension
a = np.expand_dims(a, axis=0)
tau = np.expand_dims(tau, axis=0)

# Exchange the num_tx and batchsize dimensions
a = np.transpose(a, [3, 1, 2, 0, 4, 5, 6])
tau = np.transpose(tau, [2, 1, 0, 3])

# Remove CIRs that have no active link (i.e., a is all-zero)
p_link = np.sum(np.abs(a)**2, axis=(1,2,3,4,5,6))
a = a[p_link>0.,...]
tau = tau[p_link>0.,...]

print("Shape of a:", a.shape)
print("Shape of tau: ", tau.shape)
Shape of a: (4727, 1, 16, 1, 4, 39, 14)
Shape of tau:  (4727, 1, 1, 39)

Note that transmitters and receivers have been reversed, i.e., the transmitter now denotes the UE (with 4 antennas each) and the receiver is the base station (with 16 antennas).

Remark: We have removed all positions for which the resulting CIR had zero gain, i.e., there was no path between the transmitter and the receiver. This comes from the fact that the RadioMap.sample_positions() function samples from a radio map subdivided into cells and randomizes the position within the cells. Therefore, randomly sampled positions may have no paths connecting them to the transmitter.

Let us now define a data generator that samples random UEs from the dataset and yields the previously simulated CIRs.

[9]:
class CIRGenerator:
    """Creates a generator from a given dataset of channel impulse responses

    The generator samples ``num_tx`` different transmitters from the given path
    coefficients `a` and path delays `tau` and stacks the CIRs into a single tensor.

    Note that the generator internally samples ``num_tx`` random transmitters
    from the dataset. For this, the inputs ``a`` and ``tau`` must be given for
    a single transmitter (i.e., ``num_tx`` =1) which will then be stacked
    internally.

    Parameters
    ----------
    a : [batch size, num_rx, num_rx_ant, 1, num_tx_ant, num_paths, num_time_steps], complex
        Path coefficients per transmitter

    tau : [batch size, num_rx, 1, num_paths], float
        Path delays [s] per transmitter

    num_tx : int
        Number of transmitters

    Output
    -------
    a : [batch size, num_rx, num_rx_ant, num_tx, num_tx_ant, num_paths, num_time_steps], tf.complex
        Path coefficients

    tau : [batch size, num_rx, num_tx, num_paths], tf.float
        Path delays [s]
    """

    def __init__(self,
                 a,
                 tau,
                 num_tx):

        # Copy to tensorflow
        self._a = tf.constant(a, tf.complex64)
        self._tau = tf.constant(tau, tf.float32)
        self._dataset_size = self._a.shape[0]

        self._num_tx = num_tx

    def __call__(self):

        # Generator implements an infinite loop that yields new random samples
        while True:
            # Sample 4 random users and stack them together
            idx,_,_ = tf.random.uniform_candidate_sampler(
                            tf.expand_dims(tf.range(self._dataset_size, dtype=tf.int64), axis=0),
                            num_true=self._dataset_size,
                            num_sampled=self._num_tx,
                            unique=True,
                            range_max=self._dataset_size)

            a = tf.gather(self._a, idx)
            tau = tf.gather(self._tau, idx)

            # Transpose to remove batch dimension
            a = tf.transpose(a, (3,1,2,0,4,5,6))
            tau = tf.transpose(tau, (2,1,0,3))

            # And remove batch-dimension
            a = tf.squeeze(a, axis=0)
            tau = tf.squeeze(tau, axis=0)

            yield a, tau

We use Sionna’s built-in CIRDataset to initialize a channel model that can be directly used in Sionna’s OFDMChannel layer.

[10]:
batch_size = 20 # Must be the same for the BER simulations as CIRDataset returns fixed batch_size

# Init CIR generator
cir_generator = CIRGenerator(a,
                             tau,
                             num_tx)
# Initialises a channel model that can be directly used by OFDMChannel layer
channel_model = CIRDataset(cir_generator,
                           batch_size,
                           num_rx,
                           num_rx_ant,
                           num_tx,
                           num_tx_ant,
                           max_num_paths,
                           num_time_steps)

PUSCH Link-Level Simulations

Let’s now define an end-to-end model that simulates a PUSCH transmission.

[11]:
class Model:
    """Simulate PUSCH transmissions

    This model runs BER simulations for a multi-user MIMO uplink channel
    compliant with the 5G NR PUSCH specifications.
    You can pick different scenarios, i.e., channel models, perfect or
    estimated CSI, as well as different MIMO detectors (LMMSE or KBest).

    Parameters
    ----------
    channel_model : :class:`~sionna.channel.ChannelModel` object
        An instance of a :class:`~sionna.channel.ChannelModel` object, such as
        :class:`~sionna.channel.RayleighBlockFading` or
        :class:`~sionna.channel.tr38901.UMi` or
        :class:`~sionna.channel.CIRDataset`.

    perfect_csi : bool
        Determines if perfect CSI is assumed or if the CSI is estimated

    detector : str, one of ["lmmse", "kbest"]
        MIMO detector to be used. Note that each detector has additional
        parameters that can be configured in the source code of the _init_ call.

    Input
    -----
    batch_size : int
        Number of simultaneously simulated slots

    ebno_db : float
        Signal-to-noise-ratio

    Output
    ------
    b : [batch_size, num_tx, tb_size], tf.float
        Transmitted information bits

    b_hat : [batch_size, num_tx, tb_size], tf.float
        Decoded information bits
    """
    def __init__(self,
                 channel_model,
                 perfect_csi, # bool
                 detector,    # "lmmse", "kbest"
                ):
        super().__init__()

        self._channel_model = channel_model
        self._perfect_csi = perfect_csi

        # System configuration
        self._num_prb = 16
        self._mcs_index = 14
        self._num_layers = 1
        self._mcs_table = 1
        self._domain = "freq"

        # Below parameters must equal the Path2CIR parameters
        self._num_tx_ant = 4
        self._num_tx = 4
        self._subcarrier_spacing = 30e3 # must be the same as used for Path2CIR

        # PUSCHConfig for the first transmitter
        pusch_config = PUSCHConfig()
        pusch_config.carrier.subcarrier_spacing = self._subcarrier_spacing/1000
        pusch_config.carrier.n_size_grid = self._num_prb
        pusch_config.num_antenna_ports = self._num_tx_ant
        pusch_config.num_layers = self._num_layers
        pusch_config.precoding = "codebook"
        pusch_config.tpmi = 1
        pusch_config.dmrs.dmrs_port_set = list(range(self._num_layers))
        pusch_config.dmrs.config_type = 1
        pusch_config.dmrs.length = 1
        pusch_config.dmrs.additional_position = 1
        pusch_config.dmrs.num_cdm_groups_without_data = 2
        pusch_config.tb.mcs_index = self._mcs_index
        pusch_config.tb.mcs_table = self._mcs_table

        # Create PUSCHConfigs for the other transmitters by cloning of the first PUSCHConfig
        # and modifying the used DMRS ports.
        pusch_configs = [pusch_config]
        for i in range(1, self._num_tx):
            pc = pusch_config.clone()
            pc.dmrs.dmrs_port_set = list(range(i*self._num_layers, (i+1)*self._num_layers))
            pusch_configs.append(pc)

        # Create PUSCHTransmitter
        self._pusch_transmitter = PUSCHTransmitter(pusch_configs, output_domain=self._domain)

        # Create PUSCHReceiver
        rx_tx_association = np.ones([1, self._num_tx], bool)
        stream_management = StreamManagement(rx_tx_association,
                                             self._num_layers)

        assert detector in["lmmse", "kbest"], "Unsupported MIMO detector"
        if detector=="lmmse":
            detector = LinearDetector(equalizer="lmmse",
                                      output="bit",
                                      demapping_method="maxlog",
                                      resource_grid=self._pusch_transmitter.resource_grid,
                                      stream_management=stream_management,
                                      constellation_type="qam",
                                      num_bits_per_symbol=pusch_config.tb.num_bits_per_symbol)
        elif detector=="kbest":
            detector = KBestDetector(output="bit",
                                     num_streams=self._num_tx*self._num_layers,
                                     k=64,
                                     resource_grid=self._pusch_transmitter.resource_grid,
                                     stream_management=stream_management,
                                     constellation_type="qam",
                                     num_bits_per_symbol=pusch_config.tb.num_bits_per_symbol)

        if self._perfect_csi:
            self._pusch_receiver = PUSCHReceiver(self._pusch_transmitter,
                                                 mimo_detector=detector,
                                                 input_domain=self._domain,
                                                 channel_estimator="perfect")
        else:
            self._pusch_receiver = PUSCHReceiver(self._pusch_transmitter,
                                                 mimo_detector=detector,
                                                 input_domain=self._domain)


        # Configure the actual channel
        self._channel = OFDMChannel(
                            self._channel_model,
                            self._pusch_transmitter.resource_grid,
                            normalize_channel=True,
                            return_channel=True)

    # XLA currently not supported by the CIRDataset function
    @tf.function(jit_compile=False)
    def __call__(self, batch_size, ebno_db):

        x, b = self._pusch_transmitter(batch_size)
        no = ebnodb2no(ebno_db,
                       self._pusch_transmitter._num_bits_per_symbol,
                       self._pusch_transmitter._target_coderate,
                       self._pusch_transmitter.resource_grid)
        y, h = self._channel(x, no)
        if self._perfect_csi:
            b_hat = self._pusch_receiver(y, no, h)
        else:
            b_hat = self._pusch_receiver(y, no)
        return b, b_hat

We now initialize the end-to-end model that uses the CIRDataset.

[12]:
ebno_db = 10.
e2e_model = Model(channel_model,
                  perfect_csi=False, # bool
                  detector="lmmse")  # "lmmse", "kbest"

# We can draw samples from the end-2-end link-level simulations
b, b_hat = e2e_model(batch_size, ebno_db)

Now, let’s run the BER evaluation for different system configurations.

Remark: Running the cell below can take some time.

[13]:
ebno_db = np.arange(-3, 18, 2) # sim SNR range
ber_plot = PlotBER(f"Site-Specific MU-MIMO 5G NR PUSCH")

for detector in ["lmmse", "kbest"]:
    for perf_csi in [True, False]:
        e2e_model = Model(channel_model,
                          perfect_csi=perf_csi,
                          detector=detector)
        # define legend
        csi = "Perf. CSI" if perf_csi else "Imperf. CSI"
        det = "K-Best" if detector=="kbest" else "LMMSE"
        l = det + " " + csi
        ber_plot.simulate(
                    e2e_model,
                    ebno_dbs=ebno_db, # SNR to simulate
                    legend=l, # legend string for plotting
                    max_mc_iter=500,
                    num_target_block_errors=2000,
                    batch_size=batch_size, # batch-size per Monte Carlo run
                    soft_estimates=False, # the model returns hard-estimates
                    early_stop=True,
                    show_fig=False,
                    add_bler=True,
                    forward_keyboard_interrupt=True);
EbNo [dB] |        BER |       BLER |  bit errors |    num bits | block errors |  num blocks | runtime [s] |    status
---------------------------------------------------------------------------------------------------------------------------------------
     -3.0 | 1.3210e-01 | 9.5185e-01 |     1424408 |    10782720 |         2056 |        2160 |        23.8 |reached target block errors
     -1.0 | 6.5100e-02 | 5.8488e-01 |     1117934 |    17172480 |         2012 |        3440 |         5.5 |reached target block errors
      1.0 | 3.2190e-02 | 2.7486e-01 |     1169833 |    36341760 |         2001 |        7280 |        11.6 |reached target block errors
      3.0 | 1.8199e-02 | 1.5847e-01 |     1148305 |    63098880 |         2003 |       12640 |        20.1 |reached target block errors
      5.0 | 9.6161e-03 | 8.6979e-02 |     1106003 |   115015680 |         2004 |       23040 |        36.5 |reached target block errors
      7.0 | 4.7417e-03 | 4.3675e-02 |      946815 |   199680000 |         1747 |       40000 |        63.4 |reached max iterations
      9.0 | 2.7623e-03 | 2.4625e-02 |      551580 |   199680000 |          985 |       40000 |        63.3 |reached max iterations
     11.0 | 1.2288e-03 | 1.1350e-02 |      245357 |   199680000 |          454 |       40000 |        63.2 |reached max iterations
     13.0 | 7.0918e-04 | 6.4000e-03 |      141609 |   199680000 |          256 |       40000 |        63.3 |reached max iterations
     15.0 | 4.6758e-04 | 3.9750e-03 |       93367 |   199680000 |          159 |       40000 |        63.3 |reached max iterations
     17.0 | 2.0871e-04 | 1.6000e-03 |       41675 |   199680000 |           64 |       40000 |        63.3 |reached max iterations
EbNo [dB] |        BER |       BLER |  bit errors |    num bits | block errors |  num blocks | runtime [s] |    status
---------------------------------------------------------------------------------------------------------------------------------------
     -3.0 | 2.3236e-01 | 1.0000e+00 |     2319874 |     9984000 |         2000 |        2000 |        23.0 |reached target block errors
     -1.0 | 1.8175e-01 | 9.9952e-01 |     1887149 |    10383360 |         2079 |        2080 |         3.3 |reached target block errors
      1.0 | 1.2062e-01 | 9.2593e-01 |     1300640 |    10782720 |         2000 |        2160 |         3.4 |reached target block errors
      3.0 | 6.7142e-02 | 5.3245e-01 |     1260249 |    18769920 |         2002 |        3760 |         5.9 |reached target block errors
      5.0 | 3.7267e-02 | 2.8594e-01 |     1309694 |    35143680 |         2013 |        7040 |        11.1 |reached target block errors
      7.0 | 2.1437e-02 | 1.6896e-01 |     1275606 |    59504640 |         2014 |       11920 |        18.9 |reached target block errors
      9.0 | 1.0997e-02 | 8.9235e-02 |     1234137 |   112220160 |         2006 |       22480 |        35.6 |reached target block errors
     11.0 | 6.3717e-03 | 5.3613e-02 |     1188324 |   186501120 |         2003 |       37360 |        59.1 |reached target block errors
     13.0 | 3.1225e-03 | 2.6450e-02 |      623502 |   199680000 |         1058 |       40000 |        63.5 |reached max iterations
     15.0 | 1.7128e-03 | 1.3050e-02 |      342016 |   199680000 |          522 |       40000 |        63.5 |reached max iterations
     17.0 | 9.6857e-04 | 7.4750e-03 |      193404 |   199680000 |          299 |       40000 |        63.4 |reached max iterations
EbNo [dB] |        BER |       BLER |  bit errors |    num bits | block errors |  num blocks | runtime [s] |    status
---------------------------------------------------------------------------------------------------------------------------------------
     -3.0 | 1.4725e-01 | 9.9904e-01 |     1528976 |    10383360 |         2078 |        2080 |        52.4 |reached target block errors
     -1.0 | 7.6806e-02 | 8.5125e-01 |      920195 |    11980800 |         2043 |        2400 |        33.2 |reached target block errors
      1.0 | 3.3255e-02 | 3.3273e-01 |     1009341 |    30351360 |         2023 |        6080 |        83.7 |reached target block errors
      3.0 | 1.5603e-02 | 1.4101e-01 |     1109188 |    71086080 |         2008 |       14240 |       196.6 |reached target block errors
      5.0 | 6.1264e-03 | 5.6948e-02 |     1074081 |   175319040 |         2000 |       35120 |       484.8 |reached target block errors
      7.0 | 1.8097e-03 | 1.9000e-02 |      361368 |   199680000 |          760 |       40000 |       552.2 |reached max iterations
      9.0 | 4.2337e-04 | 4.6250e-03 |       84538 |   199680000 |          185 |       40000 |       552.6 |reached max iterations
     11.0 | 1.4656e-04 | 1.3750e-03 |       29265 |   199680000 |           55 |       40000 |       552.5 |reached max iterations
     13.0 | 1.0936e-04 | 6.7500e-04 |       21838 |   199680000 |           27 |       40000 |       553.1 |reached max iterations
     15.0 | 4.9920e-05 | 3.2500e-04 |        9968 |   199680000 |           13 |       40000 |       552.6 |reached max iterations
     17.0 | 5.0691e-05 | 2.0000e-04 |       10122 |   199680000 |            8 |       40000 |       552.9 |reached max iterations
EbNo [dB] |        BER |       BLER |  bit errors |    num bits | block errors |  num blocks | runtime [s] |    status
---------------------------------------------------------------------------------------------------------------------------------------
     -3.0 | 2.4172e-01 | 1.0000e+00 |     2413315 |     9984000 |         2000 |        2000 |        48.9 |reached target block errors
     -1.0 | 1.9193e-01 | 1.0000e+00 |     1916234 |     9984000 |         2000 |        2000 |        27.7 |reached target block errors
      1.0 | 1.4189e-01 | 9.8317e-01 |     1473333 |    10383360 |         2045 |        2080 |        28.8 |reached target block errors
      3.0 | 7.5092e-02 | 7.0521e-01 |     1079597 |    14376960 |         2031 |        2880 |        40.0 |reached target block errors
      5.0 | 3.4595e-02 | 2.9412e-01 |     1174339 |    33945600 |         2000 |        6800 |        94.0 |reached target block errors
      7.0 | 1.7845e-02 | 1.4357e-01 |     1247159 |    69888000 |         2010 |       14000 |       193.9 |reached target block errors
      9.0 | 7.9984e-03 | 6.4832e-02 |     1232980 |   154152960 |         2002 |       30880 |       427.4 |reached target block errors
     11.0 | 2.7478e-03 | 2.4400e-02 |      548684 |   199680000 |          976 |       40000 |       554.0 |reached max iterations
     13.0 | 8.9915e-04 | 7.7500e-03 |      179543 |   199680000 |          310 |       40000 |       554.0 |reached max iterations
     15.0 | 2.8446e-04 | 2.1000e-03 |       56801 |   199680000 |           84 |       40000 |       554.4 |reached max iterations
     17.0 | 1.0446e-04 | 7.2500e-04 |       20858 |   199680000 |           29 |       40000 |       554.0 |reached max iterations
[14]:
ber_plot(show_ber=False)
../../_images/phy_tutorials_Link_Level_Simulations_with_RT_31_0.png
Previous Next

© Copyright 2021-2025 NVIDIA CORPORATION.

Built with Sphinx using a theme provided by Read the Docs.