Quadriga-Lib
C++/MEX/Python Utility library for radio channel modelling and simulations
C++ API Documentation for Quadriga-Lib v0.11.5


Markdown version of this page:
General usage notes
Overview
Array antenna class
Array antenna functions
Channel class
Channel functions
Channel generation functions
Channel statistics
Math functions
Site-specific simulation tools




Array antenna class



arrayant - Class for storing and manipulating array antenna models
Attributes:

Attribute Size Description
arma::Cube<dtype> e_theta_re [n_elevation, n_azimuth, n_elements] E-theta (vertical) field, real part
arma::Cube<dtype> e_theta_im [n_elevation, n_azimuth, n_elements] E-theta (vertical) field, imaginary part
arma::Cube<dtype> e_phi_re [n_elevation, n_azimuth, n_elements] E-phi (horizontal) field, real part
arma::Cube<dtype> e_phi_im [n_elevation, n_azimuth, n_elements] E-phi (horizontal) field, imaginary part
arma::Col<dtype> azimuth_grid [n_azimuth] Azimuth angles in rad, in [-pi, pi], sorted
arma::Col<dtype> elevation_grid [n_elevation] Elevation angles in rad, in [-pi/2, pi/2], sorted
arma::Mat<dtype> element_pos [3, n_elements] or empty Element positions in local Cartesian coords
arma::Mat<dtype> coupling_re [n_elements, n_ports] Coupling matrix, real part
arma::Mat<dtype> coupling_im [n_elements, n_ports] Coupling matrix, imaginary part
dtype center_frequency scalar Center frequency

Simple member functions:

Function Description
.n_elevation() Number of elevation angles
.n_azimuth() Number of azimuth angles
.n_elements() Number of antenna elements
.n_ports() Number of ports (columns of coupling matrix)
.copy() Returns a deep copy of the arrayant object
.reset() Clears all data, resetting size to zero
.is_valid() Returns "" if valid, or an error message string

Complex member functions:

Function Description
.append Append elements of another arrayant
.calc_directivity_dBi Calculate per-element directivity in dBi
.combine_pattern Compute effective patterns from elements, positions, and coupling
.copy_element Copy a single element to one or more destination slots
.export_obj_file Export pattern geometry to Wavefront OBJ
.interpolate Interpolate field patterns at given azimuth/elevation angles
.qdant_write Write arrayant to QDANT file
.remove_zeros Remove zero-valued entries from pattern data
.rotate_pattern Rotate pattern and/or polarization via Euler angles
.set_size Resize the arrayant to new dimensions
.is_valid Validate arrayant integrity



.append - Append elements of another arrayant to the current one
Declaration:
quadriga_lib::arrayant<dtype> quadriga_lib::arrayant<dtype>::append(
    const arrayant<dtype> *new_arrayant) const;
Inputs:
Returns:
See also:


.calc_beamwidth_deg - Calculate the beam width of an antenna element in degree
Declaration:
void calc_beamwidth_deg(arma::uword i_element,
    dtype threshold_dB = 3.0,
    dtype *beamwidth_az = nullptr,
    dtype *beamwidth_el = nullptr,
    dtype *z_point_ang = nullptr,
    dtype *el_point_ang = nullptr) const;
Inputs:
Outputs:
See also:


.calc_directivity_dBi - Calculate the directivity in dBi of a single array element
Declaration:
dtype quadriga_lib::arrayant<dtype>::calc_directivity_dBi(arma::uword i_element) const;
Inputs:
Returns:
See also:


.combine_pattern - Combine element patterns, positions, and coupling weights into effective radiation patterns
Declaration:
quadriga_lib::arrayant<dtype> quadriga_lib::arrayant<dtype>::combine_pattern(
    const arma::Col<dtype> *azimuth_grid_new = nullptr,
    const arma::Col<dtype> *elevation_grid_new = nullptr) const;
Inputs:
Returns:
See also:


.copy_element - Copy a single antenna element to one or more destination slots
Declaration:
void quadriga_lib::arrayant<dtype>::copy_element(arma::uword source, arma::uword destination);
void quadriga_lib::arrayant<dtype>::copy_element(arma::uword source, arma::uvec destination);
Inputs:
# See also:


.export_obj_file - Export antenna pattern geometry to a Wavefront OBJ file for 3D visualization
Declaration:
void quadriga_lib::arrayant<dtype>::export_obj_file(
    std::string fn,
    dtype directivity_range = 30.0,
    std::string colormap = "jet",
    dtype object_radius = 1.0,
    arma::uword icosphere_n_div = 4,
    arma::uvec i_element = {}) const;
Inputs:
See also:


.interpolate - Interpolate polarimetric antenna field patterns for given azimuth/elevation angles
Declaration:
void quadriga_lib::arrayant<dtype>::interpolate(
    const arma::Mat<dtype> *azimuth,
    const arma::Mat<dtype> *elevation,
    arma::Mat<dtype> *V_re, arma::Mat<dtype> *V_im,
    arma::Mat<dtype> *H_re, arma::Mat<dtype> *H_im,
    arma::uvec i_element = {},
    const arma::Cube<dtype> *orientation = nullptr,
    const arma::Mat<dtype> *element_pos_i = nullptr,
    arma::Mat<dtype> *dist = nullptr,
    arma::Mat<dtype> *azimuth_loc = nullptr,
    arma::Mat<dtype> *elevation_loc = nullptr,
    arma::Mat<dtype> *gamma = nullptr) const;
Inputs:
Outputs:
Example:
auto ant = quadriga_lib::generate_arrayant_custom<double>(90.0, 90.0, 0.0);
arma::mat azimuth = {0.0, 0.5*pi, -0.5*pi, pi};
arma::mat elevation(1, azimuth.n_elem);  // zeros
arma::mat V_re, V_im, H_re, H_im;
ant.interpolate(&azimuth, &elevation, &V_re, &V_im, &H_re, &H_im);
See also:


.is_valid - Validate the integrity of an arrayant object

Declaration:
std::string quadriga_lib::arrayant<dtype>::is_valid(bool quick_check = true) const;
Inputs:
Returns:
See also:


.qdant_write - Write arrayant data to a QDANT (XML) file
Declaration:
unsigned quadriga_lib::arrayant<dtype>::qdant_write(
    std::string fn,
    unsigned id = 0,
    arma::u32_mat layout = {}) const;
Inputs:
Returns:
See also:


.remove_zeros - Remove zero-valued entries from antenna pattern data, reducing its size
Declaration:
void quadriga_lib::arrayant<dtype>::remove_zeros(arrayant<dtype> *output = nullptr);
Inputs:


.rotate_pattern - Rotate antenna radiation patterns around the principal axes using Euler rotations
Declaration:
void quadriga_lib::arrayant<dtype>::rotate_pattern(
    dtype x_deg = 0.0,
    dtype y_deg = 0.0,
    dtype z_deg = 0.0,
    unsigned usage = 0,
    unsigned element = -1,
    arrayant<dtype> *output = nullptr);
Inputs:
See also:


.set_size - Resize an arrayant object to new dimensions
Declaration:
void quadriga_lib::arrayant<dtype>::set_size(
    arma::uword n_elevation,
    arma::uword n_azimuth,
    arma::uword n_elements,
    arma::uword n_ports);
Inputs:




Array antenna functions



arrayant_concat_multi - Concatenate two multi-frequency arrayant vectors into a single multi-element model
Declaration:
std::vector<quadriga_lib::arrayant<dtype>> quadriga_lib::arrayant_concat_multi(
        const std::vector<quadriga_lib::arrayant<dtype>> &arrayant_vec1,
        const std::vector<quadriga_lib::arrayant<dtype>> &arrayant_vec2);
Inputs:
Returns:
See also:


arrayant_copy_element_multi - Copy an antenna element to one or more destinations across all entries in a multi-frequency arrayant vector
Declaration:
void quadriga_lib::arrayant_copy_element_multi(
        std::vector<quadriga_lib::arrayant<dtype>> &arrayant_vec,
        arma::uword source,
        arma::uvec destination);

void quadriga_lib::arrayant_copy_element_multi(
        std::vector<quadriga_lib::arrayant<dtype>> &arrayant_vec,
        arma::uword source,
        arma::uword destination);
Inputs:
Example:
arma::vec freqs = {500.0, 1000.0, 2000.0, 5000.0};
auto driver = quadriga_lib::generate_speaker<double>(
    "piston", 0.05, 80.0, 12000.0, 12.0, 12.0, 85.0, "hemisphere",
    0.0, 0.0, 0.0, 0.15, 0.25, freqs, 10.0);
quadriga_lib::arrayant_copy_element_multi(driver, 0, arma::uvec{1, 2, 3});
See also:


arrayant_interpolate_multi - Interpolate multi-frequency arrayant patterns at arbitrary angles and frequencies
Declaration:
void quadriga_lib::arrayant_interpolate_multi(
        const std::vector<quadriga_lib::arrayant<dtype>> &arrayant_vec,
        const arma::Mat<dtype> *azimuth,
        const arma::Mat<dtype> *elevation,
        const arma::Col<dtype> *frequency,
        arma::Cube<dtype> *V_re,
        arma::Cube<dtype> *V_im,
        arma::Cube<dtype> *H_re,
        arma::Cube<dtype> *H_im,
        arma::uvec i_element = {},
        const arma::Cube<dtype> *orientation = nullptr,
        const arma::Mat<dtype> *element_pos_i = nullptr,
        bool validate_input = true);
Inputs:
Outputs:
Example:
auto speaker = quadriga_lib::arrayant_concat_multi(woofer, tweeter);
arma::mat az = {0.0, 1.5708, -1.5708, 3.14159};
arma::mat el(1, 4, arma::fill::zeros);
arma::vec qf = {250.0, 1500.0, 8000.0};
arma::cube V_re, V_im, H_re, H_im;
quadriga_lib::arrayant_interpolate_multi(speaker, &az, &el, &qf, &V_re, &V_im, &H_re, &H_im);
See also:


arrayant_is_valid_multi - Validate a vector of arrayant objects for multi-frequency consistency
Declaration:
std::string quadriga_lib::arrayant_is_valid_multi(
        const std::vector<quadriga_lib::arrayant<dtype>> &arrayant_vec,
        bool quick_check = true);
Inputs:
Returns:
See also:


arrayant_rotate_pattern_multi - Apply Euler rotations to all entries in a multi-frequency arrayant vector
Declaration:
void quadriga_lib::arrayant_rotate_pattern_multi(
        std::vector<quadriga_lib::arrayant<dtype>> &arrayant_vec,
        dtype x_deg = 0.0,
        dtype y_deg = 0.0,
        dtype z_deg = 0.0,
        unsigned usage = 0,
        arma::uvec i_element = arma::uvec())
Inputs:
See also:


arrayant_set_element_pos_multi - Set element positions for all entries in a multi-frequency arrayant vector
Declaration:
void quadriga_lib::arrayant_set_element_pos_multi(
        std::vector<quadriga_lib::arrayant<dtype>> &arrayant_vec,
        const arma::Mat<dtype> &element_pos,
        arma::uvec i_element = arma::uvec());
Inputs:
See also:


generate_arrayant_3GPP - Generate a 3GPP-NR compliant antenna array model
Declaration:
quadriga_lib::arrayant<dtype> quadriga_lib::generate_arrayant_3GPP(
    arma::uword M = 1, 
    arma::uword N = 1, 
    dtype center_freq = 299792458.0,
    unsigned pol = 1, 
    dtype tilt = 0.0, 
    dtype spacing = 0.5, 
    arma::uword Mg = 1,
    arma::uword Ng = 1, 
    dtype dgv = 0.5, 
    dtype dgh = 0.5,
    const quadriga_lib::arrayant<dtype> *pattern = nullptr, 
    dtype res = 1.0);
Inputs:
Returns:


generate_arrayant_custom - Generate an antenna with custom 3dB beamwidth
Declaration:
quadriga_lib::arrayant<dtype> quadriga_lib::generate_arrayant_custom(
    dtype az_3dB = 90.0,
    dtype el_3dB = 90.0, 
    dtype rear_gain_lin = 0.0, 
    dtype res = 1.0);
Inputs:
Returns:


generate_arrayant_dipole - Generate a short dipole antenna with vertical polarization
Declaration:
quadriga_lib::arrayant<dtype> quadriga_lib::generate_arrayant_dipole(dtype res = 1.0);
Inputs:
Returns:


generate_arrayant_half_wave_dipole - Generate a half-wave dipole antenna with vertical polarization
Declaration:
quadriga_lib::arrayant<dtype> quadriga_lib::generate_arrayant_half_wave_dipole(dtype res = 1.0);
Inputs:
Returns:


generate_arrayant_multibeam - Generate a planar multi-element antenna array with multiple beam directions
Declaration:
quadriga_lib::arrayant<dtype> quadriga_lib::generate_arrayant_multibeam(
    arma::uword M = 1,
    arma::uword N = 1,
    arma::Col<dtype> az = {0.0},
    arma::Col<dtype> el = {0.0},
    arma::Col<dtype> weight = {1.0},
    dtype center_freq = 299792458.0,
    unsigned pol = 1,
    dtype spacing = 0.5,
    dtype az_3dB = 120.0,
    dtype el_3dB = 120.0,
    dtype rear_gain_lin = 0.0,
    dtype res = 1.0,
    bool separate_beams = false,
    bool apply_weights = false);
Inputs:
Returns:


generate_arrayant_omni - Generate an isotropic radiator with vertical polarization
Declaration:
quadriga_lib::arrayant<dtype> quadriga_lib::generate_arrayant_omni(dtype res = 1.0);
Inputs:
Returns:


generate_arrayant_ula - Generate a uniform linear array (ULA)
Declaration:
quadriga_lib::arrayant<dtype> quadriga_lib::generate_arrayant_ula(
    arma::uword N = 1, 
    dtype center_freq = 299792458.0, 
    dtype spacing = 0.5,
    const quadriga_lib::arrayant<dtype> *pattern = nullptr, 
    dtype res = 1.0);
Inputs:
Returns:


generate_arrayant_xpol - Generate a cross-polarized isotropic radiator
Declaration:
quadriga_lib::arrayant<dtype> quadriga_lib::generate_arrayant_xpol(dtype res = 1.0);
Inputs:
Returns:


generate_speaker - Generate a parametric frequency-dependent loudspeaker directivity model
Declaration:
std::vector<quadriga_lib::arrayant<dtype>> quadriga_lib::generate_speaker(
        std::string driver_type = "piston",
        dtype radius = 0.05,
        dtype lower_cutoff = 80.0,
        dtype upper_cutoff = 12000.0,
        dtype lower_rolloff_slope = 12.0,
        dtype upper_rolloff_slope = 12.0,
        dtype sensitivity = 85.0,
        std::string radiation_type = "hemisphere",
        dtype hor_coverage = 0.0,
        dtype ver_coverage = 0.0,
        dtype horn_control_freq = 0.0,
        dtype baffle_width = 0.15,
        dtype baffle_height = 0.25,
        arma::Col<dtype> frequencies = arma::Col<dtype>(),
        dtype angular_resolution = 5.0);
Inputs:
Returns:
Example:
arma::vec freqs = {100.0, 500.0, 1000.0, 5000.0, 10000.0};
auto spk = quadriga_lib::generate_speaker<double>("piston", 0.05, 80.0, 12000.0,
               12.0, 12.0, 85.0, "hemisphere", 0.0, 0.0, 0.0, 0.15, 0.25, freqs, 5.0);
auto horn = quadriga_lib::generate_speaker<double>("horn");
auto sub = quadriga_lib::generate_speaker<double>("omni", 0.13, 30.0, 200.0,
               12.0, 24.0, 92.0, "monopole", 0.0, 0.0, 0.0, 0.15, 0.25, {30.,50.,80.,120.,200.}, 10.0);


qdant_read - Read an arrayant object from a QDANT file
Declaration:
quadriga_lib::arrayant<dtype> quadriga_lib::qdant_read(
        std::string fn,
        unsigned id = 1,
        arma::u32_mat *layout = nullptr);
Inputs:
Returns:
See also:


qdant_read_multi - Read all arrayant objects from a QDANT file into a vector
Declaration:
std::vector<quadriga_lib::arrayant<dtype>> quadriga_lib::qdant_read_multi(
        const std::string &fn,
        arma::u32_mat *layout = nullptr);
Inputs:
Returns:
See also:


qdant_write_multi - Write a vector of arrayant objects to a single QDANT file
Declaration:
void quadriga_lib::qdant_write_multi(
        const std::string &fn,
        const std::vector<quadriga_lib::arrayant<dtype>> &arrayant_vec);
Inputs:
See also:




Channel class



channel - Class for storing and managing MIMO channel data and metadata across multiple snapshots
Attributes:

Attribute Size Description
std::string name Name of the channel object
arma::Col<dtype> center_frequency [1]; [n_snap], or [] Center frequency
arma::Mat<dtype> tx_pos [3, n_snap] or [3, 1] = static Transmitter positions
arma::Mat<dtype> rx_pos [3, n_snap] or [3, 1] = static Receiver positions
arma::Mat<dtype> tx_orientation [3, n_snap]; [3, 1] = static, or [] = no rotation Transmitter orientation (Euler angles)
arma::Mat<dtype> rx_orientation [3, n_snap]; [3, 1] = static, or [] = no rotation Receiver orientation (Euler angles)
std::vector<arma::Cube<dtype>> coeff_re per snap [n_rx, n_tx, n_path] Channel coefficients, real part
std::vector<arma::Cube<dtype>> coeff_im per snap [n_rx, n_tx, n_path] Channel coefficients, imaginary part
std::vector<arma::Cube<dtype>> delay per snap [n_rx, n_tx, n_path] or [1, 1, n_path] = broadcast Path delays in seconds
std::vector<arma::Col<dtype>> path_gain per snap [n_path] Path gains before antenna pattern
std::vector<arma::Col<dtype>> path_length per snap [n_path] Path lengths TX to RX
std::vector<arma::Mat<dtype>> path_polarization per snap [8, n_path] Interleaved polarization transfer matrices
std::vector<arma::Mat<dtype>> path_angles per snap [n_path, 4] Angles {AOD, EOD, AOA, EOA} in rad
std::vector<arma::Mat<dtype>> path_fbs_pos per snap [3, n_path] First-bounce scatterer positions
std::vector<arma::Mat<dtype>> path_lbs_pos per snap [3, n_path] Last-bounce scatterer positions
std::vector<arma::Col<unsigned>> no_interact per snap [n_path] Number of interactions per path
std::vector<arma::Mat<dtype>> interact_coord per snap [3, sum(no_interact)] Interaction point coordinates
std::vector<std::string> par_names Names of unstructured metadata fields
std::vector<std::any> par_data Unstructured metadata values (string, scalar, matrix, etc.)
int initial_position scalar 0-based index of the reference snapshot

Simple member functions:

Method Description
.n_snap() Returns the number of snapshots
.n_rx() Returns number of receive antennas; 0 if coefficients absent
.n_tx() Returns number of transmit antennas; 0 if coefficients absent
.n_path() Returns number of paths per snapshot as a vector
.empty() Returns true if the object contains no channel data
.is_valid() Returns empty string if valid, otherwise an error message

Complex member functions:


.add_paths - Append new propagation paths to an existing channel snapshot
Declaration:
void quadriga_lib::channel<dtype>::add_paths(
    arma::uword i_snap,
    const arma::Cube<dtype> *coeff_re_add = nullptr,
    const arma::Cube<dtype> *coeff_im_add = nullptr,
    const arma::Cube<dtype> *delay_add = nullptr,
    const arma::u32_vec *no_interact_add = nullptr,
    const arma::Mat<dtype> *interact_coord_add = nullptr,
    const arma::Col<dtype> *path_gain_add = nullptr,
    const arma::Col<dtype> *path_length_add = nullptr,
    const arma::Mat<dtype> *path_polarization_add = nullptr,
    const arma::Mat<dtype> *path_angles_add = nullptr,
    const arma::Mat<dtype> *path_fbs_pos_add = nullptr,
    const arma::Mat<dtype> *path_lbs_pos_add = nullptr);
Inputs:


.calc_effective_path_gain - Calculate the effective path gain per snapshot in linear scale
Declaration:
arma::Col<dtype> quadriga_lib::channel<dtype>::calc_effective_path_gain(bool assume_valid = false) const;
Inputs:
Returns:


.write_paths_to_obj_file - Export propagation paths to a Wavefront OBJ file for 3D visualization
Declaration:
void quadriga_lib::channel<dtype>::write_paths_to_obj_file(
    std::string fn,
    arma::uword max_no_paths = 0,
    dtype gain_max = -60.0,
    dtype gain_min = -140.0,
    std::string colormap = "jet",
    arma::uvec i_snap = {},
    dtype radius_max = 0.05,
    dtype radius_min = 0.01,
    arma::uword n_edges = 5) const;
Inputs:
See also:




Channel functions



any_type_id - Get type ID and raw access from a `std::any` object
Declaration:
int quadriga_lib::any_type_id(
    const std::any *data,
    unsigned long long *dims = nullptr,
    void **dataptr = nullptr);
Inputs:
Outputs:
Returns:
See also:


baseband_freq_response - Compute the baseband frequency response of a MIMO channel
Declaration:
void quadriga_lib::baseband_freq_response(
    const arma::Cube<dtype> *coeff_re,
    const arma::Cube<dtype> *coeff_im,
    const arma::Cube<dtype> *delay,
    const arma::Col<dtype> *pilot_grid,
    const double bandwidth,
    arma::Cube<dtype> *hmat_re,
    arma::Cube<dtype> *hmat_im,
    arma::Cube<std::complex<dtype>> *hmat = nullptr);
Inputs:
Outputs:
See also:


baseband_freq_response_multi - Compute the wideband frequency response of a MIMO channel with frequency-dependent coefficients
Declaration:
void quadriga_lib::baseband_freq_response_multi(
    const std::vector<arma::Cube<dtype>> &coeff_re,
    const std::vector<arma::Cube<dtype>> &coeff_im,
    const std::vector<arma::Cube<dtype>> &delay,
    const arma::Col<dtype> &freq_in,
    const arma::Col<dtype> &freq_out,
    arma::Cube<dtype> *hmat_re = nullptr,
    arma::Cube<dtype> *hmat_im = nullptr,
    arma::Cube<std::complex<dtype>> *hmat = nullptr,
    bool remove_delay_phase = true);
Inputs:
Outputs:
See also:


baseband_freq_response_vec - Compute the baseband frequency response of multiple MIMO channels
Declaration:
void quadriga_lib::baseband_freq_response_vec(
    const std::vector<arma::Cube<dtype>> *coeff_re,
    const std::vector<arma::Cube<dtype>> *coeff_im,
    const std::vector<arma::Cube<dtype>> *delay,
    const arma::Col<dtype> *pilot_grid,
    const double bandwidth,
    std::vector<arma::Cube<dtype>> *hmat_re,
    std::vector<arma::Cube<dtype>> *hmat_im,
    const arma::u32_vec *i_snap = nullptr);
Inputs:
Outputs:
See also:


get_HDF5_version - Return the HDF5 version string as defined by the compile-time header macros

Declaration:
std::string quadriga_lib::get_HDF5_version();
Returns:


hdf5_create - Create a new HDF5 channel file with a defined storage layout
Declaration:
void quadriga_lib::hdf5_create(
    std::string fn,
    unsigned nx = 65536,
    unsigned ny = 1,
    unsigned nz = 1,
    unsigned nw = 1);
Inputs:
See also:


hdf5_read_channel - Read a channel object from an HDF5 file at a specified 4D index
Declaration:
quadriga_lib::channel<dtype> quadriga_lib::hdf5_read_channel(
    std::string fn,
    unsigned ix = 0,
    unsigned iy = 0,
    unsigned iz = 0,
    unsigned iw = 0);
Inputs:
Returns:
See also:


hdf5_read_dset - Read an unstructured dataset from an HDF5 file at a specified 4D index
Declaration:
std::any quadriga_lib::hdf5_read_dset(
    std::string fn,
    std::string par_name,
    unsigned ix = 0,
    unsigned iy = 0,
    unsigned iz = 0,
    unsigned iw = 0,
    std::string prefix = "par_");
Inputs:
Returns:
See also:


hdf5_read_dset_names - Read names of unstructured datasets stored at a 4D slot in an HDF5 file
Declaration:
arma::uword quadriga_lib::hdf5_read_dset_names(
    std::string fn,
    std::vector<std::string> *par_names,
    unsigned ix = 0,
    unsigned iy = 0,
    unsigned iz = 0,
    unsigned iw = 0,
    std::string prefix = "par_");
Inputs:
Returns:
See also:


hdf5_read_layout - Read the storage layout of an HDF5 channel file
Declaration:
arma::u32_vec quadriga_lib::hdf5_read_layout(
    std::string fn,
    arma::u32_vec *channelID = nullptr);
Inputs:
Returns:
See also:


hdf5_reshape_layout - Reshape the 4D storage layout of an existing HDF5 channel file
Declaration:
void quadriga_lib::hdf5_reshape_layout(
    std::string fn,
    unsigned nx,
    unsigned ny = 1,
    unsigned nz = 1,
    unsigned nw = 1);
Inputs:
See also:


hdf5_write - Write a channel object to an HDF5 file at a specified 4D index
Declaration:
int quadriga_lib::hdf5_write(
    const quadriga_lib::channel<dtype> *ch,
    std::string fn,
    unsigned ix = 0,
    unsigned iy = 0,
    unsigned iz = 0,
    unsigned iw = 0,
    bool assume_valid = false);
Inputs:
Returns:
See also:


hdf5_write_dset - Write a single unstructured dataset to an HDF5 file at a specified 4D index
Declaration:
void quadriga_lib::hdf5_write_dset(
    std::string fn,
    std::string par_name,
    const std::any *par_data,
    unsigned ix = 0,
    unsigned iy = 0,
    unsigned iz = 0,
    unsigned iw = 0,
    std::string prefix = "par_");
Inputs:
See also:


qrt_file_parse - Read metadata from a QRT file
Declaration:
void quadriga_lib::qrt_file_parse(
    const std::string &fn,
    arma::uword *no_cir = nullptr,
    arma::uword *no_orig = nullptr,
    arma::uword *no_dest = nullptr,
    arma::uword *no_freq = nullptr,
    arma::uvec *cir_offset = nullptr,
    std::vector<std::string> *orig_names = nullptr,
    std::vector<std::string> *dest_names = nullptr,
    int *version = nullptr,
    arma::fvec *freq = nullptr,
    arma::fmat *cir_pos = nullptr,
    arma::fmat *cir_orientation = nullptr,
    arma::fmat *orig_pos = nullptr,
    arma::fmat *orig_orientation = nullptr,
    std::ifstream *file = nullptr);
Inputs:
Outputs:


qrt_file_read - Read ray-tracing CIR data from a QRT file
Declaration:
void quadriga_lib::qrt_file_read(
    const std::string &fn,
    arma::uword i_cir = 0,
    arma::uword i_orig = 0,
    bool downlink = true,
    arma::Col<dtype> *center_frequency = nullptr,
    arma::Col<dtype> *tx_pos = nullptr,
    arma::Col<dtype> *tx_orientation = nullptr,
    arma::Col<dtype> *rx_pos = nullptr,
    arma::Col<dtype> *rx_orientation = nullptr,
    arma::Mat<dtype> *fbs_pos = nullptr,
    arma::Mat<dtype> *lbs_pos = nullptr,
    arma::Mat<dtype> *path_gain = nullptr,
    arma::Col<dtype> *path_length = nullptr,
    arma::Cube<dtype> *M = nullptr,
    arma::Col<dtype> *aod = nullptr,
    arma::Col<dtype> *eod = nullptr,
    arma::Col<dtype> *aoa = nullptr,
    arma::Col<dtype> *eoa = nullptr,
    std::vector<arma::Mat<dtype>> *path_coord = nullptr,
    int normalize_M = 1,
    arma::u32_vec *no_int = nullptr,
    arma::fmat *coord = nullptr,
    std::ifstream *file = nullptr,
    const qrt_read_cache *cache = nullptr);
Inputs:
Outputs:
Example:
std::ifstream stream("scene.qrt", std::ios::in | std::ios::binary);
auto cache = quadriga_lib::qrt_read_cache_init("scene.qrt", &stream);
arma::vec center_freq, tx_pos, rx_pos, path_length;
arma::mat path_gain; arma::cube M;
for (arma::uword ic = 0; ic < cache.no_cir; ++ic)
    for (arma::uword io = 0; io < cache.no_orig; ++io)
        quadriga_lib::qrt_file_read<double>("", ic, io, true,
            &center_freq, &tx_pos, nullptr, &rx_pos, nullptr,
            nullptr, nullptr, &path_gain, &path_length, &M,
            nullptr, nullptr, nullptr, nullptr, nullptr, 1,
            nullptr, nullptr, &stream, &cache);
See also:


qrt_read_cache_init - Initialize a QRT read cache for fast repeated access
Declaration:
quadriga_lib::qrt_read_cache quadriga_lib::qrt_read_cache_init(
    const std::string &fn,
    std::ifstream *file = nullptr);
Inputs:
Returns:


quantize_delays - Map path delays to a fixed tap grid using two-tap power-weighted interpolation
Declaration:
void quadriga_lib::quantize_delays(
    const std::vector<arma::Cube<dtype>> *coeff_re,
    const std::vector<arma::Cube<dtype>> *coeff_im,
    const std::vector<arma::Cube<dtype>> *delay,
    std::vector<arma::Cube<dtype>> *coeff_re_quant,
    std::vector<arma::Cube<dtype>> *coeff_im_quant,
    std::vector<arma::Cube<dtype>> *delay_quant,
    dtype tap_spacing = (dtype)5.0e-9,
    arma::uword max_no_taps = 48,
    dtype power_exponent = (dtype)1.0,
    int fix_taps = 0);
Inputs:
Outputs:
Example:
std::vector<arma::Cube<double>> cre(2), cim(2), dl(2);
cre[0].set_size(1,1,3); cim[0].set_size(1,1,3); dl[0].set_size(1,1,3);
cre[1].set_size(1,1,2); cim[1].set_size(1,1,2); dl[1].set_size(1,1,2);
dl[0](0,0,1) = 12.5e-9; dl[0](0,0,2) = 33.4e-9; dl[1](0,0,1) = 10.0e-9;
std::vector<arma::Cube<double>> cre_q, cim_q, dl_q;
quadriga_lib::quantize_delays(&cre, &cim, &dl, &cre_q, &cim_q, &dl_q, 5.0e-9, 48, 1.0, 0);




Channel generation functions



get_channels_ieee_indoor - Generate indoor MIMO channel realizations for IEEE TGn/TGac/TGax/TGah models
Declaration:
std::vector<quadriga_lib::channel<double>> quadriga_lib::get_channels_ieee_indoor(
    const quadriga_lib::arrayant<double> &ap_array,
    const quadriga_lib::arrayant<double> &sta_array,
    std::string ChannelType,
    double CarrierFreq_Hz = 5.25e9,
    double tap_spacing_s = 10.0e-9,
    arma::uword n_users = 1,
    double observation_time = 0.0,
    double update_rate = 1.0e-3,
    double speed_station_kmh = 0.0,
    double speed_env_kmh = 1.2,
    arma::vec Dist_m = {4.99},
    arma::uvec n_floors = {0},
    bool uplink = false,
    arma::mat offset_angles = {},
    arma::uword n_subpath = 20,
    double Doppler_effect = 50.0,
    arma::sword seed = -1,
    double KF_linear = NAN,
    double XPR_NLOS_linear = NAN,
    double SF_std_dB_LOS = NAN,
    double SF_std_dB_NLOS = NAN,
    double dBP_m = NAN,
    arma::uvec n_walls = {0},
    double wall_loss = 5.0);
Inputs:
Returns:
See also:


get_channels_irs - Calculate MIMO channel coefficients for IRS-assisted communication
Declaration:
std::vector<bool> quadriga_lib::get_channels_irs(
    const quadriga_lib::arrayant<dtype> *tx_array,
    const quadriga_lib::arrayant<dtype> *rx_array,
    const quadriga_lib::arrayant<dtype> *irs_array,
    dtype Tx, dtype Ty, dtype Tz,
    dtype Tb, dtype Tt, dtype Th,
    dtype Rx, dtype Ry, dtype Rz,
    dtype Rb, dtype Rt, dtype Rh,
    dtype Ix, dtype Iy, dtype Iz,
    dtype Ib, dtype It, dtype Ih,
    const arma::Mat<dtype> *fbs_pos_1,
    const arma::Mat<dtype> *lbs_pos_1,
    const arma::Col<dtype> *path_gain_1,
    const arma::Col<dtype> *path_length_1,
    const arma::Mat<dtype> *M_1,
    const arma::Mat<dtype> *fbs_pos_2,
    const arma::Mat<dtype> *lbs_pos_2,
    const arma::Col<dtype> *path_gain_2,
    const arma::Col<dtype> *path_length_2,
    const arma::Mat<dtype> *M_2,
    arma::Cube<dtype> *coeff_re,
    arma::Cube<dtype> *coeff_im,
    arma::Cube<dtype> *delay,
    arma::uword i_irs = 0,
    dtype threshold_dB = -140.0,
    dtype center_frequency = 0.0,
    bool use_absolute_delays = false,
    arma::Cube<dtype> *aod = nullptr,
    arma::Cube<dtype> *eod = nullptr,
    arma::Cube<dtype> *aoa = nullptr,
    arma::Cube<dtype> *eoa = nullptr,
    const quadriga_lib::arrayant<dtype> *irs_array_2 = nullptr,
    const std::vector<bool> *active_path = nullptr);
Inputs:
Outputs:
Returns:
See also:


get_channels_multifreq - Compute channel coefficients for spherical waves across multiple frequencies
Declaration:
void quadriga_lib::get_channels_multifreq(
    const std::vector<arrayant<dtype>> &tx_array,
    const std::vector<arrayant<dtype>> &rx_array,
    dtype Tx, dtype Ty, dtype Tz,
    dtype Tb, dtype Tt, dtype Th,
    dtype Rx, dtype Ry, dtype Rz,
    dtype Rb, dtype Rt, dtype Rh,
    const arma::Mat<dtype> &fbs_pos,
    const arma::Mat<dtype> &lbs_pos,
    const arma::Mat<dtype> &path_gain,
    const arma::Col<dtype> &path_length,
    const arma::Cube<dtype> &M,
    const arma::Col<dtype> &freq_in,
    const arma::Col<dtype> &freq_out,
    std::vector<arma::Cube<dtype>> &coeff_re,
    std::vector<arma::Cube<dtype>> &coeff_im,
    std::vector<arma::Cube<dtype>> &delay,
    bool use_absolute_delays = false,
    bool add_fake_los_path = false,
    dtype propagation_speed = dtype(299792458.0));
Inputs:
Outputs:
See also:


get_channels_planar - Calculate MIMO channel coefficients for planar wave paths
Declaration:
void quadriga_lib::get_channels_planar(
    const quadriga_lib::arrayant<dtype> *tx_array,
    const quadriga_lib::arrayant<dtype> *rx_array,
    dtype Tx, dtype Ty, dtype Tz,
    dtype Tb, dtype Tt, dtype Th,
    dtype Rx, dtype Ry, dtype Rz,
    dtype Rb, dtype Rt, dtype Rh,
    const arma::Col<dtype> *aod,
    const arma::Col<dtype> *eod,
    const arma::Col<dtype> *aoa,
    const arma::Col<dtype> *eoa,
    const arma::Col<dtype> *path_gain,
    const arma::Col<dtype> *path_length,
    const arma::Mat<dtype> *M,
    arma::Cube<dtype> *coeff_re,
    arma::Cube<dtype> *coeff_im,
    arma::Cube<dtype> *delay,
    dtype center_frequency = dtype(0.0),
    bool use_absolute_delays = false,
    bool add_fake_los_path = false,
    arma::Col<dtype> *rx_Doppler = nullptr);
Inputs:
Outputs:
See also:


get_channels_spherical - Calculate MIMO channel coefficients and delays for spherical wave propagation
Declaration:
void quadriga_lib::get_channels_spherical(
    const quadriga_lib::arrayant<dtype> *tx_array,
    const quadriga_lib::arrayant<dtype> *rx_array,
    dtype Tx, dtype Ty, dtype Tz,
    dtype Tb, dtype Tt, dtype Th,
    dtype Rx, dtype Ry, dtype Rz,
    dtype Rb, dtype Rt, dtype Rh,
    const arma::Mat<dtype> *fbs_pos,
    const arma::Mat<dtype> *lbs_pos,
    const arma::Col<dtype> *path_gain,
    const arma::Col<dtype> *path_length,
    const arma::Mat<dtype> *M,
    arma::Cube<dtype> *coeff_re,
    arma::Cube<dtype> *coeff_im,
    arma::Cube<dtype> *delay,
    dtype center_frequency = dtype(0.0),
    bool use_absolute_delays = false,
    bool add_fake_los_path = false,
    arma::Cube<dtype> *aod = nullptr,
    arma::Cube<dtype> *eod = nullptr,
    arma::Cube<dtype> *aoa = nullptr,
    arma::Cube<dtype> *eoa = nullptr,
    bool use_avx2 = false);
Inputs:
Outputs:
See also:




Channel statistics



acdf - Calculate the empirical averaged cumulative distribution function (CDF)
Declaration:
void quadriga_lib::acdf(const arma::Mat<dtype> &data,
    arma::Col<dtype> *bins = nullptr,
    arma::Mat<dtype> *cdf_per_set = nullptr,
    arma::Col<dtype> *cdf_avg = nullptr,
    arma::Col<dtype> *mu = nullptr,
    arma::Col<dtype> *sig = nullptr,
    arma::uword n_bins = 201);
Inputs:
Outputs:


calc_angular_spreads_sphere - Calculate azimuth and elevation angular spreads with spherical wrapping
Declaration:
void quadriga_lib::calc_angular_spreads_sphere(
    const std::vector<arma::Col<dtype>> &az,
    const std::vector<arma::Col<dtype>> &el,
    const std::vector<arma::Col<dtype>> &powers,
    arma::Col<dtype> *azimuth_spread = nullptr,
    arma::Col<dtype> *elevation_spread = nullptr,
    arma::Mat<dtype> *orientation = nullptr,
    std::vector<arma::Col<dtype>> *phi = nullptr,
    std::vector<arma::Col<dtype>> *theta = nullptr,
    bool disable_wrapping = false,
    bool calc_bank_angle = true,
    dtype quantize = (dtype)0);
Inputs:
Outputs:


calc_cross_polarization_ratio - Calculate the cross-polarization ratio (XPR) for linear and circular polarization bases
Declaration:
void quadriga_lib::calc_cross_polarization_ratio(
    const std::vector<arma::Col<dtype>> &powers,
    const std::vector<arma::Mat<dtype>> &M,
    const std::vector<arma::Col<dtype>> &path_length,
    const arma::Mat<dtype> &tx_pos,
    const arma::Mat<dtype> &rx_pos,
    arma::Mat<dtype> *xpr = nullptr,
    arma::Col<dtype> *pg = nullptr,
    bool include_los = false,
    dtype window_size = 0.01);
Inputs:
Outputs:


calc_delay_spread - Calculates RMS delay spread from per-CIR delays and linear-scale powers
Declaration:
arma::Col<dtype> quadriga_lib::calc_delay_spread(
    const std::vector<arma::Col<dtype>> &delays,
    const std::vector<arma::Col<dtype>> &powers,
    dtype threshold = 100.0,
    dtype granularity = 0.0,
    arma::Col<dtype> *mean_delay = nullptr);
Inputs:
Outputs:
Returns:
See also:


calc_rician_k_factor - Calculate the Rician K-Factor from channel impulse response data
Declaration:
void quadriga_lib::calc_rician_k_factor(
    const std::vector<arma::Col<dtype>> &powers,
    const std::vector<arma::Col<dtype>> &path_length,
    const arma::Mat<dtype> &tx_pos,
    const arma::Mat<dtype> &rx_pos,
    arma::Col<dtype> *kf = nullptr,
    arma::Col<dtype> *pg = nullptr,
    dtype window_size = 0.01);
Inputs:
Outputs:




Math functions



calc_rotation_matrix - Calculate rotation matrices from Euler angles
Declaration:
arma::Cube<dtype> quadriga_lib::calc_rotation_matrix(
    const arma::Cube<dtype> &orientation,
    bool invert_y_axis = false, 
    bool transposeR = false);

arma::Mat<dtype> quadriga_lib::calc_rotation_matrix(
    const arma::Mat<dtype> &orientation,
    bool invert_y_axis = false, 
    bool transposeR = false);

arma::Col<dtype> quadriga_lib::calc_rotation_matrix(
    const arma::Col<dtype> &orientation,
    bool invert_y_axis = false, 
    bool transposeR = false);
Inputs:
Returns:


fast_acos - Compute elementwise approximate arc-cosine of a vector
Declaration:
void quadriga_lib::fast_acos(const arma::fvec &x, arma::fvec &c);
void quadriga_lib::fast_acos(const arma::vec &x,  arma::fvec &c);
Inputs:
Outputs:


fast_asin - Compute elementwise approximate arc-sine of a vector
Declaration:
void quadriga_lib::fast_asin(const arma::fvec &x, arma::fvec &s);
void quadriga_lib::fast_asin(const arma::vec &x,  arma::fvec &s);
Inputs:
Outputs:


fast_atan2 - Compute elementwise approximate two-argument arc-tangent of two vectors
Declaration:
void quadriga_lib::fast_atan2(const arma::fvec &y, const arma::fvec &x, arma::fvec &a);
void quadriga_lib::fast_atan2(const arma::vec &y,  const arma::vec &x,  arma::fvec &a);
Inputs:
Outputs:


fast_cart2geo - Convert elementwise Cartesian coordinates to azimuth/elevation angles and vector length
Declaration:
void quadriga_lib::fast_cart2geo(const arma::fvec &x, const arma::fvec &y, const arma::fvec &z,
                                 arma::fvec &az, arma::fvec &el, arma::fvec *len = nullptr, int use_kernel = 0);

void quadriga_lib::fast_cart2geo(const arma::vec &x, const arma::vec &y, const arma::vec &z,
                                 arma::vec &az, arma::vec &el, arma::vec *len = nullptr, int use_kernel = 0);
Inputs:
Outputs:
See also:


fast_geo2cart - Convert elementwise azimuth/elevation angles to Cartesian coordinates
Declaration:
void fast_geo2cart(
    const arma::Col<dtype> &az,
    const arma::Col<dtype> &el,
    arma::Col<dtype> &x,
    arma::Col<dtype> &y,
    arma::Col<dtype> &z,
    arma::Col<dtype> *sAZ = nullptr,
    arma::Col<dtype> *cAZ = nullptr,
    arma::Col<dtype> *sEL = nullptr,
    arma::Col<dtype> *cEL = nullptr,
    const arma::Col<dtype> *len = nullptr,
    int use_kernel = 0);
Inputs:
Outputs:
See also:


fast_sincos - Compute elementwise approximate sine and/or cosine of a vector
Declaration:
void quadriga_lib::fast_sincos(const arma::fvec &x, arma::fvec *s = nullptr, arma::fvec *c = nullptr);
void quadriga_lib::fast_sincos(const arma::vec &x,  arma::fvec *s = nullptr, arma::fvec *c = nullptr);
Inputs:
Outputs:


fast_slerp - Compute elementwise approximate SLERP interpolation between two complex-valued vectors
Declaration:
void quadriga_lib::fast_slerp(const arma::fvec &Ar, const arma::fvec &Ai,
                              const arma::fvec &Br, const arma::fvec &Bi,
                              const arma::fvec &w,
                              arma::fvec &Xr, arma::fvec &Xi);

void quadriga_lib::fast_slerp(const arma::vec &Ar, const arma::vec &Ai,
                              const arma::vec &Br, const arma::vec &Bi,
                              const arma::vec &w,
                              arma::fvec &Xr, arma::fvec &Xi);
Inputs:
Outputs:


interp_1D / interp_2D - Perform linear interpolation (1D or 2D) on single or multiple data sets.
Declarations:
void interp_2D(const arma::Cube<dtype> &input, const arma::Col<dtype> &xi, const arma::Col<dtype> &yi,
               const arma::Col<dtype> &xo, const arma::Col<dtype> &yo, arma::Cube<dtype> &output);

arma::Cube<dtype> interp_2D(const arma::Cube<dtype> &input, const arma::Col<dtype> &xi, const arma::Col<dtype> &yi,
                            const arma::Col<dtype> &xo, const arma::Col<dtype> &yo);

arma::Mat<dtype> interp_2D(const arma::Mat<dtype> &input, const arma::Col<dtype> &xi, const arma::Col<dtype> &yi,
                           const arma::Col<dtype> &xo, const arma::Col<dtype> &yo);

arma::Mat<dtype> interp_1D(const arma::Mat<dtype> &input, const arma::Col<dtype> &xi, const arma::Col<dtype> &xo);

arma::Col<dtype> interp_1D(const arma::Col<dtype> &input, const arma::Col<dtype> &xi, const arma::Col<dtype> &xo);
Arguments:
Input / Output size details:
Examples:
arma::cube input(5, 5, 2, arma::fill::randu); // example input data
arma::vec xi = arma::linspace(0, 4, 5);
arma::vec yi = arma::linspace(0, 4, 5);
arma::vec xo = arma::linspace(0, 4, 10);
arma::vec yo = arma::linspace(0, 4, 10);

arma::cube output;
quadriga_lib::interp_2D(input, xi, yi, xo, yo, output);
arma::vec input = arma::linspace(0, 1, 5);
arma::vec xi = arma::linspace(0, 4, 5);
arma::vec xo = arma::linspace(0, 4, 10);

auto output = quadriga_lib::interp_1D(input, xi, xo);




Site-specific simulation tools



calc_diffraction_gain - Calculate diffraction gain for multiple TX-RX pairs using a 3D triangular mesh
Declaration:
void calc_diffraction_gain(
    const arma::Mat<dtype> *orig,
    const arma::Mat<dtype> *dest,
    const arma::Mat<dtype> *mesh,
    const arma::Mat<dtype> *mtl_prop,
    dtype center_frequency,
    int lod = 2,
    arma::Col<dtype> *gain = nullptr,
    arma::Cube<dtype> *coord = nullptr,
    int verbose = 0,
    const arma::u32_vec *sub_mesh_index = nullptr,
    int use_kernel = 0,
    int gpu_id = 0);
Inputs:
Outputs:
See also:


colormap - Generate a colormap matrix with RGB values
Declaration:
arma::uchar_mat quadriga_lib::colormap(std::string map, bool high_res = false);
Inputs:
Returns:


combine_irs_coord - Combine path interaction coordinates for IRS-assisted TX → RX channels
Declaration:
void quadriga_lib::combine_irs_coord(
    dtype Ix, dtype Iy, dtype Iz,
    const arma::u32_vec *no_interact_1,
    const arma::Mat<dtype> *interact_coord_1,
    const arma::u32_vec *no_interact_2,
    const arma::Mat<dtype> *interact_coord_2,
    arma::u32_vec *no_interact,
    arma::Mat<dtype> *interact_coord,
    bool reverse_segment_1 = false,
    bool reverse_segment_2 = false,
    const std::vector<bool> *active_path = nullptr);
Inputs:
Outputs:
See also:


coord2path - Convert path interaction coordinates into FBS/LBS positions, path length, and angles
Declaration:
void quadriga_lib::coord2path(
    dtype Tx, dtype Ty, dtype Tz,
    dtype Rx, dtype Ry, dtype Rz,
    const arma::u32_vec *no_interact,
    const arma::Mat<dtype> *interact_coord,
    arma::Col<dtype> *path_length = nullptr,
    arma::Mat<dtype> *fbs_pos = nullptr,
    arma::Mat<dtype> *lbs_pos = nullptr,
    arma::Mat<dtype> *path_angles = nullptr,
    std::vector<arma::Mat<dtype>> *path_coord = nullptr,
    bool reverse_path = false);
Inputs:
Outputs:


generate_diffraction_paths - Generate elliptic propagation paths and weights for diffraction gain estimation
Declaration:
void generate_diffraction_paths(
    const arma::Mat<dtype> *orig,
    const arma::Mat<dtype> *dest,
    dtype center_frequency,
    int lod,
    arma::Cube<dtype> *ray_x,
    arma::Cube<dtype> *ray_y,
    arma::Cube<dtype> *ray_z,
    arma::Cube<dtype> *weight);
Inputs:
Outputs:
See also:


icosphere - Construct a geodesic polyhedron from recursive icosahedron subdivision
Declaration:
arma::uword quadriga_lib::icosphere(
    arma::uword n_div,
    dtype radius,
    arma::Mat<dtype> *center,
    arma::Col<dtype> *length = nullptr,
    arma::Mat<dtype> *vert = nullptr,
    arma::Mat<dtype> *direction = nullptr,
    bool direction_xyz = false);
Inputs:
Outputs:
Returns:


medium_gain - Linear gain of a ray traversing a homogeneous lossy medium
Declaration:
dtype quadriga_lib::medium_gain(
        const arma::Mat<dtype> &mtl_prop,
        arma::uword iM,
        dtype dist,
        dtype fGHz);
Inputs:
Returns:
See also:


mitsuba_xml_file_write - Write a triangular mesh to a Mitsuba 3 XML scene file
Declaration:
void quadriga_lib::mitsuba_xml_file_write(
    const std::string &fn,
    const arma::Mat<dtype> &vert_list,
    const arma::umat &face_ind,
    const arma::uvec &obj_ind,
    const arma::uvec &mtl_ind,
    const std::vector<std::string> &obj_names,
    const std::vector<std::string> &mtl_names,
    const arma::Mat<dtype> &bsdf = {},
    bool map_to_itu_materials = false);
Inputs:
See also:


obj_file_read - Read a Wavefront .obj file and extract geometry and material information
Declaration:
arma::uword quadriga_lib::obj_file_read(
    std::string fn,
    arma::Mat<dtype> *mesh = nullptr,
    arma::Mat<dtype> *mtl_prop = nullptr,
    arma::Mat<dtype> *vert_list = nullptr,
    arma::umat *face_ind = nullptr,
    arma::uvec *obj_ind = nullptr,
    arma::uvec *mtl_ind = nullptr,
    std::vector<std::string> *obj_names = nullptr,
    std::vector<std::string> *mtl_names = nullptr,
    arma::Mat<dtype> *bsdf = nullptr,
    const std::string &materials_csv = "");
Inputs:
Outputs:
Returns:
Default material table:
See also:


obj_overlap_test - Detect overlapping 3D objects in a triangular mesh
Declaration:
arma::uvec quadriga_lib::obj_overlap_test(
    const arma::Mat<dtype> *mesh,
    const arma::uvec *obj_ind,
    std::vector<std::string> *reason = nullptr,
    dtype tolerance = 0.0005);
Inputs:
Returns:
See also:


path_to_tube - Convert a 3D path into a tube surface mesh for visualization
Declaration:
void quadriga_lib::path_to_tube(
    const arma::Mat<dtype> *path_coord,
    arma::Mat<dtype> *vert,
    arma::umat *faces,
    dtype radius = 1.0,
    arma::uword n_edges = 5);
Inputs:
Outputs:


point_cloud_aabb - Compute the axis-aligned bounding boxes (AABB) of a 3D point cloud
Declaration:
arma::Mat<dtype> quadriga_lib::point_cloud_aabb(
    const arma::Mat<dtype> *points,
    const arma::u32_vec *sub_cloud_index = nullptr,
    arma::uword vec_size = 1);
Inputs:
Returns:
See also:


point_cloud_segmentation - Reorganize a point cloud into spatial sub-clouds for efficient processing
Declaration:
arma::uword quadriga_lib::point_cloud_segmentation(
    const arma::Mat<dtype> *points,
    arma::Mat<dtype> *pointsR,
    arma::u32_vec *sub_cloud_index,
    arma::uword target_size = 1024,
    arma::uword vec_size = 1,
    arma::u32_vec *forward_index = nullptr,
    arma::u32_vec *reverse_index = nullptr);
Inputs:
Outputs:
Returns:
See also:


point_cloud_split - Split a point cloud into two sub-clouds along a spatial axis
Declaration:
int quadriga_lib::point_cloud_split(
    const arma::Mat<dtype> *points,
    arma::Mat<dtype> *pointsA,
    arma::Mat<dtype> *pointsB,
    int axis = 0,
    arma::Col<int> *split_ind = nullptr);
Inputs:
Outputs:
Returns:
See also:


point_inside_mesh - Test whether 3D points are inside a triangle mesh using raycasting
Declaration:
arma::u32_vec quadriga_lib::point_inside_mesh(
    const arma::Mat<dtype> *points,
    const arma::Mat<dtype> *mesh,
    const arma::u32_vec *obj_ind = nullptr,
    dtype distance = 0.0);
Inputs:
Returns:
See also:


ray_mesh_interact - Calculates reflection, transmission, or refraction of EM/acoustic waves at mesh surfaces
Declaration:
void quadriga_lib::ray_mesh_interact(
    int interaction_type, dtype center_frequency,
    const arma::Mat<dtype> *orig, const arma::Mat<dtype> *dest,
    const arma::Mat<dtype> *fbs, const arma::Mat<dtype> *sbs,
    const arma::Mat<dtype> *mesh, const arma::Mat<dtype> *mtl_prop,
    const arma::u32_vec *fbs_ind, const arma::u32_vec *sbs_ind,
    const arma::Mat<dtype> *trivec = nullptr,
    const arma::Mat<dtype> *tridir = nullptr,
    const arma::Col<dtype> *orig_length = nullptr,
    arma::Mat<dtype> *origN = nullptr, arma::Mat<dtype> *destN = nullptr,
    arma::Col<dtype> *gainN = nullptr, arma::Mat<dtype> *xprmatN = nullptr,
    arma::Mat<dtype> *trivecN = nullptr, arma::Mat<dtype> *tridirN = nullptr,
    arma::Col<dtype> *orig_lengthN = nullptr,
    arma::Col<dtype> *fbs_angleN = nullptr,
    arma::Col<dtype> *thicknessN = nullptr,
    arma::Col<dtype> *edge_lengthN = nullptr,
    arma::Mat<dtype> *normal_vecN = nullptr,
    arma::s32_vec *out_typeN = nullptr);
Inputs:
Outputs:
See also:


ray_point_intersect - Calculate intersections of ray beams with points in 3D space
Declaration:
std::vector<arma::u32_vec> quadriga_lib::ray_point_intersect(
    const arma::Mat<dtype> *points,
    const arma::Mat<dtype> *orig,
    const arma::Mat<dtype> *trivec,
    const arma::Mat<dtype> *tridir,
    const arma::u32_vec *sub_cloud_index = nullptr,
    arma::u32_vec *hit_count = nullptr,
    int use_kernel = 0,
    int gpu_id = 0);
Inputs:
Optional output:
Returns:
See also:


ray_triangle_intersect - Compute ray-triangle intersections in 3D using the Möller–Trumbore algorithm
Declaration:
void quadriga_lib::ray_triangle_intersect(
    const arma::Mat<dtype> *orig,
    const arma::Mat<dtype> *dest,
    const arma::Mat<dtype> *mesh,
    arma::Mat<dtype> *fbs = nullptr,
    arma::Mat<dtype> *sbs = nullptr,
    arma::u32_vec *no_interact = nullptr,
    arma::u32_vec *fbs_ind = nullptr,
    arma::u32_vec *sbs_ind = nullptr,
    const arma::u32_vec *sub_mesh_index = nullptr,
    const arma::Mat<dtype> *aabb = nullptr,
    int use_kernel = 0,
    int gpu_id = 0);
Inputs:
Outputs:
See also:


subdivide_rays - Subdivide ray beams into four smaller sub-beams
Declaration:
arma::uword quadriga_lib::subdivide_rays(
    const arma::Mat<dtype> *orig,
    const arma::Mat<dtype> *trivec,
    const arma::Mat<dtype> *tridir,
    const arma::Mat<dtype> *dest = nullptr,
    arma::Mat<dtype> *origN = nullptr,
    arma::Mat<dtype> *trivecN = nullptr,
    arma::Mat<dtype> *tridirN = nullptr,
    arma::Mat<dtype> *destN = nullptr,
    const arma::u32_vec *index = nullptr,
    const double ray_offset = 0.0);
Inputs:
Outputs:
Returns:
See also:


subdivide_triangles - Subdivide triangles into smaller triangles
Declaration:
arma::uword quadriga_lib::subdivide_triangles(
    arma::uword n_div,
    const arma::Mat<dtype> *triangles_in,
    arma::Mat<dtype> *triangles_out,
    const arma::Mat<dtype> *mtl_prop = nullptr,
    arma::Mat<dtype> *mtl_prop_out = nullptr);
Inputs:
Outputs:
Returns:


triangle_mesh_aabb - Calculate the axis-aligned bounding box (AABB) of a triangle mesh and its sub-meshes
Declaration:
arma::Mat<dtype> quadriga_lib::triangle_mesh_aabb(
    const arma::Mat<dtype> *mesh,
    const arma::u32_vec *sub_mesh_index = nullptr,
    arma::uword vec_size = 1);
Inputs:
Returns:
See also:


triangle_mesh_segmentation - Reorganize a 3D triangular mesh into spatially clustered sub-meshes for faster processing
Declaration:
arma::uword triangle_mesh_segmentation(
    const arma::Mat<dtype> *mesh,
    arma::Mat<dtype> *meshR,
    arma::u32_vec *sub_mesh_index,
    arma::uword target_size = 1024,
    arma::uword vec_size = 1,
    const arma::Mat<dtype> *mtl_prop = nullptr,
    arma::Mat<dtype> *mtl_propR = nullptr,
    arma::u32_vec *mesh_index = nullptr);
Inputs:
Outputs:
Returns:
See also:


triangle_mesh_split - Split a 3D triangular mesh into two sub-meshes along a given axis
Declaration:
int triangle_mesh_split(
    const arma::Mat<dtype> *mesh,
    arma::Mat<dtype> *meshA,
    arma::Mat<dtype> *meshB,
    int axis = 0,
    arma::Col<int> *split_ind = nullptr);
Inputs:
Outputs:
Returns:
See also:


write_png - Write a data matrix to a color-coded PNG file
Declaration:
void quadriga_lib::write_png(
    const arma::Mat<dtype> &data,
    std::string fn,
    std::string colormap = "jet",
    dtype min_val = NAN,
    dtype max_val = NAN,
    bool log_transform = false);
Inputs: