C++ API Documentation for Quadriga-Lib v0.10.1
Overview
Array antenna class
Array antenna functions
Channel class
Channel functions
Channel generation functions
Miscellaneous tools
| calc_rotation_matrix | Calculate rotation matrices from Euler angles |
| cart2geo | Convert Cartesian coordinates to geographic coordinates (azimuth, elevation, distance) |
| colormap | Generate colormap |
| fast_sincos | Fast, approximate sine/cosine |
| geo2cart | Transform geographic (azimuth, elevation, length) to Cartesian coordinates |
| interp_1d / interp_2d | Perform linear interpolation (1D or 2D) on single or multiple data sets. |
| write_png | Write data to a PNG file |
Site-Specific Simulation Tools
Array antenna class
arrayant - Class for storing and manipulating array antenna models
-
Description:
- An array antenna consists of multiple individual elements.
- Each element occupies a specific position relative to the array's phase-center, its local origin.
- Elements can also be inter-coupled, represented by a coupling matrix.
-
Attributes:
arma::Cube<dtype> e_theta_re |
Vertical component of the electric field, real part |
arma::Cube<dtype> e_theta_im |
Vertical component of the electric field, imaginary part |
arma::Cube<dtype> e_phi_re |
Horizontal component of the electric field, real part |
arma::Cube<dtype> e_phi_im |
Horizontal component of the electric field, imaginary part |
arma::Col<dtype> azimuth_grid |
Azimuth angles in pattern (theta) in [rad], between -pi and pi, sorted |
arma::Col<dtype> elevation_grid |
Elevation angles in pattern (phi) in [rad], between -pi/2 and pi/2, sorted |
arma::Mat<dtype> element_pos |
Element positions (optional), Size: Empty or [3, n_elements] |
arma::Mat<dtype> coupling_re |
Coupling matrix, real part (optional), Size: [n_elements, n_ports] |
arma::Mat<dtype> coupling_im |
Coupling matrix, imaginary part (optional), Size: [n_elements, n_ports] |
dtype center_frequency |
Center frequency in [Hz] |
- Allowed datatypes (
dtype): float and double
e_theta_re, e_theta_im, e_phi_re, e_phi_im must have size [n_elevation, n_azimuth, n_elements]
-
Example:
float pi = arma::datum::pi;
quadriga_lib::arrayant<float> ant;
ant.azimuth_grid = {-0.75f * pi, 0.0f, 0.75f * pi, pi};
ant.elevation_grid = {-0.45f * pi, 0.0f, 0.45f * pi};
arma::mat A = arma::linspace(1.0, 12.0, 12);
A.reshape(3, 4);
arma::fcube B;
B.zeros(3, 4, 1);
B.slice(0) = arma::conv_to<arma::fmat>::from(A);
ant.e_theta_re = B * 0.5f;
ant.e_theta_im = B * 0.002f;
ant.e_phi_re = -B;
ant.e_phi_im = -B * 0.001f;
arma::fmat C = {1.0f, 2.0f, 4.0f};
ant.element_pos = C.t();
ant.coupling_re = {1.0f};
ant.coupling_im = {0.1f};
ant.center_frequency = 2.0e9f;
ant.name = "name";
-
Simple member functions:
.n_elevation() |
Returns number of elevation angles as 64bit integer |
.n_azimuth() |
Returns number of azimuth angles as 64bit integer |
.n_elements() |
Returns number of antenna elements as 64bit integer |
.n_ports() |
Returns number of ports (after coupling) as 64bit integer |
.copy() |
Creates a copy of the array antenna object |
.reset() |
Reset the size to zero (the arrayant object will contain no data) |
.is_valid() |
Returns an empty string if arrayant object is valid or an error message otherwise |
-
Complex member fuctions:
.append - Append elements of another antenna array
-
Description:
- Combines elements of another antenna array (
new_arrayant) with the current antenna array object.
- Returns a new
arrayant object containing elements from both antenna arrays.
- Throws an error if the sampling grids of the two antenna arrays do not match.
- Member function of arrayant
- Allowed datatypes (
dtype): float or double
-
Declaration:
quadriga_lib::arrayant<dtype> quadriga_lib::arrayant<dtype>::append(const arrayant<dtype> *new_arrayant) const;
-
Arguments:
const arrayant<dtype> *new_arrayant (input)
Pointer to an antenna array object whose elements will be added to the current object. Sampling grids must match exactly.
-
Returns:
quadriga_lib::arrayant<dtype>
A new antenna array object combining the current and new antenna elements.
-
Example:
quadriga_lib::arrayant<double> ant1 = quadriga_lib::generate_arrayant_custom<double>(90.0, 90.0, 0.0);
quadriga_lib::arrayant<double> ant2 = quadriga_lib::generate_arrayant_custom<double>(120.0, 60.0, 0.0);
quadriga_lib::arrayant<double> combined_ant = ant1.append(&ant2);
.calc_directivity_dbi - Calculate the directivity (in dBi) of array antenna elements
.combine_pattern - Calculate effective radiation patterns for array antennas
-
Description:
- Member function of arrayant
- By integrating element radiation patterns, element positions, and the coupling weights, one can
determine an effective radiation pattern observable by a receiver in the antenna's far field.
- Leveraging these effective patterns is especially beneficial in antenna design, beamforming
applications such as in 5G systems, and in planning wireless communication networks in complex
environments like urban areas. This streamlined approach offers a significant boost in computation
speed when calculating MIMO channel coefficients, as it reduces the number of necessary operations.
- Allowed datatypes (
dtype): float and double
-
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;
-
Arguments:
arma::Col<dtype> *azimuth_grid_new (optional)
Azimuth angle grid of the output array antenna in [rad], between -pi and pi, sorted
arma::Col<dtype> *elevation_grid_new (optional)
Elevation angle grid of the output array antenna in [rad], between -pi/2 and pi/2, sorted
-
Example:
auto ant = quadriga_lib::generate_arrayant_omni<double>(); // Generate omni antenna
ant.copy_element(0, 1); // Duplicate the first element
ant.element_pos.row(1) = {-0.25, 0.25}; // Set element positions (in lambda)
ant.coupling_re.ones(2, 1); // Set coupling matrix (real part)
ant.coupling_im.reset(); // Remove imaginary part
ant = ant.combine_pattern(); // Calculate the combined pattern
.copy_element - Creates a copy of a single array antenna element
-
Description:
- Member function of arrayant
- Allowed datatypes (
dtype): float and double
-
Declaration:
void quadriga_lib::arrayant<dtype>::copy_element(arma::uword source, arma::uvec destination);
void quadriga_lib::arrayant<dtype>::copy_element(arma::uword source, arma::uword destination);
-
Arguments:
arma::uword source (optional)
Index of the source element (0-based)
arma::uvec destination or arma::uword destination
Index of the destinations element (0-based), either as a vector or as a scalar.
-
Example:
auto ant = quadriga_lib::generate_arrayant_omni<double>(); // Generate omni antenna
ant.copy_element(0, 1); // Duplicate the first element
ant.copy_element(1, {2,3}); // Duplicate multiple times
.export_obj_file - Export antenna pattern geometry to Wavefront OBJ file
-
Description:
- This function exports the antenna pattern geometry to a Wavefront OBJ file, useful for visualization in 3D software such as Blender.
- Member function of arrayant
- Allowed datatypes (
dtype): float or double
-
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;
-
Arguments:
std::string fn (input)
Filename of the OBJ file to which the antenna pattern will be exported. Cannot be empty.
dtype directivity_range = 30.0 (optional input)
Directivity range in decibels (dB) for visualizing the antenna pattern. This value defines the
dynamic range of the visualized directivity pattern. Default: 30.0
std::string colormap = "jet" (optional input)
Colormap used for visualizing the antenna directivity. Supported colormaps are: jet,
parula, winter, hot, turbo, copper, spring, cool, gray, autumn, summer.
Default: "jet"
dtype object_radius = 1.0 (optional input)
Radius of the exported antenna pattern geometry object, specified in meters. Default: 1.0
arma::uword icosphere_n_div = 4 (optional input)
Number of subdivisions used to map the antenna pattern onto an icosphere. Higher values yield finer
mesh resolution. Default: 4
arma::uvec i_element = {} (optional input)
Antenna element indices for which the pattern geometry is exported. Indices are 0-based. Providing
an empty vector {} (default) exports the geometry for all elements of the antenna array.
-
Example:
auto ant = quadriga_lib::generate_arrayant_custom<double>(90.0, 90.0, 0.0);
ant.export_obj_file("antenna_pattern.obj", 40.0, "turbo", 1.5, 5);
.interpolate - Interpolate array antenna field patterns
-
Description:
- This function interpolates polarimetric antenna field patterns for a given set of azimuth and
elevation angles.
- Member function of arrayant
- Allowed datatypes (
dtype): float or double
-
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,
const arma::Mat<dtype> *element_pos_i,
arma::Mat<dtype> *dist,
arma::Mat<dtype> *azimuth_loc, arma::Mat<dtype> *elevation_loc,
arma::Mat<dtype> *gamma) const;
-
Arguments:
const arma::Mat<dtype> *azimuth (input)
Azimuth angles in [rad] for which the field pattern should be interpolated. Values must be
between -pi and pi, cannot be NULL
| Option 1: |
Use the same angles for all antenna elements (planar wave approximation) |
|
Size: [1, n_ang] |
| Option 2: |
Provide different angles for each array element (e.g. for spherical waves) |
|
Size: [n_out, n_ang] |
const arma::Mat<dtype> *elevation (input)
Elevation angles in [rad] for which the field pattern should be interpolated. Values must be
between -pi/2 and pi/2, cannot be NULL
| Option 1: |
Use the same angles for all antenna elements (planar wave approximation) |
|
Size: [1, n_ang] |
| Option 2: |
Provide different angles for each array element (e.g. for spherical waves) |
|
Size: [n_out, n_ang] |
arma::Mat<dtype> *V_re (output)
Real part of the interpolated e-theta (vertical) field component, Size [n_out, n_ang],
will be resized if it does not match the required size (invalidates data pointers), cannot be NULL
arma::Mat<dtype> *V_im (output)
Imaginary part of the interpolated e-theta (vertical) field component, Size [n_out, n_ang]
will be resized if it does not match the required size (invalidates data pointers), cannot be NULL
arma::Mat<dtype> *H_re (output)
Real part of the interpolated e-phi (horizontal) field component, Size [n_out, n_ang]
will be resized if it does not match the required size (invalidates data pointers), cannot be NULL
arma::Mat<dtype> *H_im (output)
Imaginary part of the interpolated e-phi (horizontal) field component, Size [n_out, n_ang]
will be resized if it does not match the required size (invalidates data pointers), cannot be NULL
arma::uvec i_element = {} (optional input)
The element indices for which the interpolation should be done, optional argument,
values must be between 1 and n_elements. It is possible to duplicate elements, i.e. by passing
{1,1,2}. If this parameter is not provided (or an empty array is passed), i_element is initialized
to include all elements of the array antenna. In this case, n_out = n_elements,
Length: n_out or empty {}
const arma::Cube<dtype> *orientation = nullptr (optional input)
This (optional) 3-element vector allows for setting orientation of the array antenna or
of individual elements using Euler angles (bank, tilt, heading); values must be given in [rad];
By default, the orientation is {0,0,0}, i.e. the broadside of the antenna points at the horizon
towards the East. Sizes: nullptr (use default), [3, 1] (set orientation for entire array),
[3, n_out] (set orientation for individual elements), or [3, 1, n_ang] (set orientation for
individual angles) or [3, n_out, n_ang] (set orientation for individual elements and angles)
const arma::Mat<dtype> *element_pos_i = nullptr (optional input)
Positions of the array antenna elements in local cartesian coordinates (using units
of [m]). If this parameter is not given, the element positions from the arrayant object are used.
Sizes: nullptr (use arrayant.element_pos), [3, n_out] (set alternative positions)
arma::Mat<dtype> *dist = nullptr (optional output)
The effective distances between the antenna elements when seen from the direction
of the incident path. The distance is calculated by an projection of the array positions on the normal
plane of the incident path. This is needed for calculating the phase of the antenna response.
Size: nullptr (do not calculate this) or [n_out, n_ang] (argument be resized if it does not already
match this size)
arma::Mat<dtype> *azimuth_loc = nullptr (optional output)
The azimuth angles in [rad] for the local antenna coordinate system, i.e., after
applying the orientation. If no orientation vector is given, these angles are identical to the input
azimuth angles. Size: nullptr or [n_out, n_ang]
arma::Mat<dtype> *elevation_loc = nullptr (optional output)
The elevation angles in [rad] for the local antenna coordinate system, i.e., after
applying the orientation. If no orientation vector is given, these angles are identical to the input
elevation angles. Size: nullptr or [n_out, n_ang]
arma::Mat<dtype> *gamma = nullptr (optional output)
Polarization rotation angles in [rad]. Size: nullptr or [n_out, n_ang]
-
Example:
double pi = arma::datum::pi;
// Directional antenna, pointing east
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}; // Azimuth angles: East, North, South, West
arma::mat elevation(1, azimuth.n_elem); // Initialize to 0
arma::mat V_re, V_im, H_re, H_im; // Output variables (uninitialized)
ant.interpolate(&azimuth, &elevation, &V_re, &V_im, &H_re, &H_im);
V_re.print();
.qdant_write - Write array antenna object and layout to QDANT file
-
Description:
- This function writes array antenna patterns and their layout into the QuaDRiGa array antenna exchange
format (QDANT), an XML-based file format
- Multiple array antennas can be stored in the same file using the
id parameter.
- If writing to an exisiting file without specifying an
id, the data gests appended at the end.
The output id_in_file identifies the location inside the file.
- An optional storage
layout can be provided to organize data inside the file.
- Member function of arrayant
- Allowed datatypes (
dtype): float or double
-
Declaration:
unsigned quadriga_lib::arrayant<dtype>::qdant_write(
std::string fn,
unsigned id = 0,
arma::u32_mat layout = {}) const;
-
Arguments:
std::string fn (input)
Filename of the QDANT file to write the antenna pattern data. Cannot be empty.
unsigned id = 0 (optional input)
ID of the antenna to write into the file. If not provided or set to 0, the antenna pattern is appended with a new ID equal to the maximum existing ID in the file plus one.
arma::u32_mat layout = {} (optional input)
Layout specifying the organization of multiple antenna elements inside the file. This matrix must only contain element IDs present within the file. Default: empty matrix {}.
-
Returns:
unsigned
Returns the ID assigned to the antenna pattern within the file after writing.
-
Example:
quadriga_lib::arrayant<double> ant = quadriga_lib::generate_arrayant_custom<double>(90.0, 90.0, 0.0);
unsigned ant_id = ant.qdant_write("antenna_data.qdant");
-
See also:
.remove_zeros - Remove zeros from antenna pattern data
-
Description:
- This function removes zeros from the antenna pattern data, altering its size accordingly.
- If called without an argument, the function modifies the antenna array properties in place.
- Member function of arrayant
- Allowed datatypes (
dtype): float or double
-
Declaration:
void quadriga_lib::arrayant<dtype>::remove_zeros(arrayant<dtype> *output = nullptr);
-
Arguments:
arrayant<dtype> *output = nullptr (optional output)
Pointer to an antenna array object where the modified pattern data is should be written to. If set
to nullptr (default), the modifications are applied directly to the calling antenna object.
-
Example:
auto ant = quadriga_lib::generate_arrayant_custom<double>(90.0, 90.0, 0.0);
ant.remove_zeros(); // Modifies ant in-place
.rotate_pattern - Adjust orientation of antenna patterns
-
Description:
* Adjusts the orientation of antenna radiation patterns by performing precise rotations around the
three principal axes (x, y, z) of the local Cartesian coordinate system (Euler rotations)
* Transforms both uniformly and non-uniformly sampled antenna patterns, useful for precise adjustments
in antennas like parabolic antennas with small apertures.
* Member function of arrayant
* Allowed datatypes (dtype): float or double
-
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);
-
Arguments:
dtype x_deg = 0.0 (optional input)
Rotation angle around the x-axis (bank angle), specified in degrees. Default: 0.0
dtype y_deg = 0.0 (optional input)
Rotation angle around the y-axis (tilt angle), specified in degrees. Default: 0.0
dtype z_deg = 0.0 (optional input)
Rotation angle around the z-axis (heading angle), specified in degrees. Default: 0.0
unsigned usage = 0 (optional input)
Rotation usage model, specifying which components to rotate: (0): Rotate both pattern and polarization,
(1): Rotate only pattern, (2): Rotate only polarization, (3): Rotate both pattern and polarization without adjusting the grid
unsigned element = -1 (optional input)
Index of the antenna element (0-based) to rotate. Default (-1) applies rotation to all elements.
arrayant<dtype> output = nullptr (optional output)
Pointer to an antenna array object to store the modified pattern data. If nullptr (default), modifications are applied directly to the calling antenna object.
-
Example:
auto ant = quadriga_lib::generate_arrayant_custom<double>(90.0, 90.0, 0.0);
ant.rotate_pattern(0.0, 0.0, 45.0);
.set_size - Change size of antenna array object
-
Description:
- Changes the size of an antenna array (
arrayant) without explicitly preserving existing data.
- Resets
element_pos to zero and sets coupling_re and coupling_im to identity matrices.
- Other properties may contain undefined or garbage data after resizing
- Size update is performed only if the existing size differs from the specified new size
- Function returns an error if the antenna object is marked as read-only
- Member function of arrayant
- Allowed datatypes (
dtype): float or double
-
Declaration:
void quadriga_lib::arrayant<dtype>::set_size(
arma::uword n_elevation,
arma::uword n_azimuth,
arma::uword n_elements,
arma::uword n_ports);
-
Arguments:
arma::uword n_elevation (input)
Number of elevation angles to resize to.
arma::uword n_azimuth (input)
Number of azimuth angles to resize to.
arma::uword n_elements (input)
Number of antenna elements in the array after resizing.
arma::uword n_ports (input)
Number of ports (after coupling of elements) in the resized antenna array.
-
Example:
quadriga_lib::arrayant<double> ant;
ant.set_size(180, 360, 4, 2);
.is_valid - Validate integrity of antenna array object
-
Description:
- Checks the integrity of an antenna array (
arrayant) object.
- Returns an empty string if the antenna object is valid.
- Provides an error message describing any issue if the object is invalid.
- A quick integrity check can be performed for efficiency.
- Member function of arrayant
- Allowed datatypes (
dtype): float or double
-
Declaration:
std::string quadriga_lib::arrayant<dtype>::is_valid(bool quick_check = true) const;
-
Arguments:
bool quick_check = true (optional input)
If set to true (default), performs a quick validation check. Setting it to false performs a more thorough validation.
-
Returns:
std::string
Returns an empty string ("") if the antenna object passes the integrity check; otherwise, returns an error message detailing the issue.
-
Example:
quadriga_lib::arrayant<double> ant;
std::string result = ant.is_valid();
if(result.empty()) {
std::cout << "Antenna array is valid." << std::endl;
} else {
std::cout << "Error: " << result << std::endl;
}
Array antenna functions
qdant_read - Reads array antenna data from QDANT files
-
Description:
- Reads antenna pattern data and layout from a QuaDRiGa array antenna exchange format (QDANT) file, which stores antenna pattern data in XML format.
- Returns an antenna array (
arrayant) object constructed from the specified data in the QDANT file.
- Allowed datatypes (
dtype): float or double
-
Declaration:
arrayant<dtype> quadriga_lib::qdant_read(std::string fn, unsigned id = 1, arma::u32_mat *layout = nullptr);
-
Arguments:
std::string fn (input)
Filename of the QDANT file from which antenna pattern data will be read. Cannot be empty.
unsigned id = 1 (optional input)
ID of the antenna within the QDANT file to read. Default is 1.
arma::u32_mat layout = nullptr (optional output)
Pointer to a matrix that will store the layout information of multiple antenna elements from the file. The layout contains element IDs present in the QDANT file.
-
Returns:
arrayant<dtype>
Antenna array object containing data from the specified antenna ID within the QDANT file.
-
Example:
arma::u32_mat layout;
auto ant = quadriga_lib::qdant_read<double>("antenna_data.qdant", 2, &layout);
-
See also:
generate_arrayant_omni - Generate isotropic radiator with vertical polarization
generate_arrayant_xpol - Generate cross-polarized isotropic radiator
generate_arrayant_dipole - Generate short dipole with vertical polarization
generate_arrayant_half_wave_dipole - Generate half-wave dipole with vertical polarization
generate_arrayant_custom - Generate antenna with custom 3dB beamwidth
-
Description:
- Generates an antenna array object with a custom-defined 3dB beamwidth (FWHM) in azimuth and elevation.
- Allows control over the rear-side attenuation using a front-to-back gain ratio.
- Supports specification of the pattern sampling resolution.
- Allowed datatypes (
dtype): float or double
-
Declaration:
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)
-
Arguments:
dtype az_3dB = 90.0 (optional input)
3dB beamwidth in azimuth, specified in degrees. Default: 90.0
dtype el_3db = 90.0 (optional input)
3dB beamwidth in elevation, specified in degrees. Default: 90.0
dtype rear_gain_lin = 0.0 (optional input)
Front-to-back gain ratio as a linear value. Default: 0.0 (no rear gain)
dtype res = 1.0 (optional input)
Resolution of the antenna pattern sampling grid, specified in degrees. Default: 1.0
-
Returns:
arrayant<dtype>
Antenna array object with the specified 3dB beamwidth and rear gain.
-
Example:
auto ant = quadriga_lib::generate_arrayant_custom<double>(60.0, 45.0, 0.01, 1.0);
generate_arrayant_ula - Generate an unified linear array
-
Description:
- Supports horizontal stacking of elements
- Default pattern: Isotropic omnidirectional radiator, v-polarized
- Optionally, a custom per-element pattern can be provided (element positions, coupling, and center frequency are ignored).
- Allowed datatypes (
dtype): float or double
-
Declaration:
arrayant<dtype> quadriga_lib::generate_arrayant_ula(
arma::uword N = 1, dtype center_freq = 299792458.0, dtype spacing = 0.5,
const arrayant<dtype> *pattern = nullptr, dtype res = 1.0)
-
Arguments:
arma::uword N = 1 (optional input)
Number of horizontal elements in the array. Default: 1
dtype center_freq = 299792458.0 (optional input)
Center frequency of the antenna array in Hz. Default: 299792458.0
dtype spacing = 0.5 (optional input)
Spacing between elements in wavelengths (λ). Default: 0.5
const arrayant<dtype> pattern = nullptr (optional input)
Optional pointer to a custom per-element antenna pattern. If provided, it overrides default generation. Only the shape of the pattern is used; element positions, coupling, and frequency are ignored.
dtype res = 1.0 (optional input)
Pattern resolution in degrees. Ignored if a custom pattern is provided. Default: 1.0
-
Returns:
arrayant<dtype>
Generated antenna array object according to 3GPP-NR specifications.
generate_arrayant_3gpp - Generate 3GPP-NR compliant antenna model
-
Description:
- Generates an antenna array based on the 3GPP-NR channel model specification.
- Supports vertical and horizontal stacking of elements and panels with a range of polarization configurations.
- Optionally, a custom per-element pattern can be provided (element positions, coupling, and center frequency are ignored).
- Allowed datatypes (
dtype): float or double
-
Declaration:
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 arrayant<dtype> *pattern = nullptr, dtype res = 1.0)
-
Arguments:
arma::uword M = 1 (optional input)
Number of vertical elements in the array. Default: 1
arma::uword N = 1 (optional input)
Number of horizontal elements in the array. Default: 1
dtype center_freq = 299792458.0 (optional input)
Center frequency of the antenna array in Hz. Default: 299792458.0
unsigned pol = 1 (optional input)
Polarization configuration:
pol = 1 |
vertical polarization (default value) |
pol = 2 |
H/V polarized elements, results in 2NM elements |
pol = 3 |
+/-45° polarized elements, results in 2NM elements |
pol = 4 |
vertical polarization, combines elements in vertical direction, results in N elements |
pol = 5 |
H/V polarization, combines elements in vertical direction, results in 2N elements |
pol = 6 |
+/-45° polarization, combines elements in vertical direction, results in 2N elements |
dtype tilt = 0.0 (optional input)
Electrical downtilt angle in degrees. Used for pol values 4, 5, and 6. Default: 0.0
dtype spacing = 0.5 (optional input)
Spacing between elements in wavelengths (λ). Default: 0.5
arma::uword Mg = 1 (optional input)
Number of vertically stacked panels (columns). Default: 1
arma::uword Ng = 1 (optional input)
Number of horizontally stacked panels (rows). Default: 1
dtype dgv = 0.5 (optional input)
Panel spacing in vertical direction in wavelengths (λ). Default: 0.5
dtype dgh = 0.5 (optional input)
Panel spacing in horizontal direction in wavelengths (λ). Default: 0.5
const arrayant<dtype> pattern = nullptr (optional input)
Optional pointer to a custom per-element antenna pattern. If provided, it overrides default generation. Only the shape of the pattern is used; element positions, coupling, and frequency are ignored.
dtype res = 1.0 (optional input)
Pattern resolution in degrees. Ignored if a custom pattern is provided. Default: 1.0
-
Returns:
arrayant<dtype>
Generated antenna array object according to 3GPP-NR specifications.
-
Example:
// Generate 3GPP antenna array with 4x4 elements, H/V polarization
auto ant = quadriga_lib::generate_arrayant_3GPP<double>(4, 4, 3e9, 2);
generate_arrayant_multibeam - Generate a planar multi-element antenna with support for multiple beam directions.
-
Description:
This function generates a planar array antenna with M rows and N columns.
It allows customization of the per-element radiation pattern, polarization, and spacing.
Multiple beam directions can be specified via azimuth and elevation angles.
Beamforming uses maximum-ratio transmission (MRT), which is optimal for a single
beam and approximate when multiple beams are specified.
Supported data types for dtype: float or double.
-
Declaration:
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 );
-
Arguments:
arma::uword M = 1 (optional input)
Number of vertical (rows) elements in the array. Default: 1
arma::uword N = 1 (optional input)
Number of horizontal (columns) elements in the array. Default: 1
arma::Col<dtype> az = {0.0} (optional input)
Azimuth beam angles (degrees). Vector of length n_beams. Default: {0.0}
arma::Col<dtype> el = {0.0} (optional input)
Elevation beam angles (degrees). Vector of length n_beams. Default: {0.0}
arma::Col<dtype> weight = {1.0} (optional input)
Scaling factors for each beam. The vector must have the same length as az and el.
Values are normalized so that their sum equals 1. Can be used to prioritize beams.
Default: {1.0}
dtype center_freq = 299792458.0 (optional input)
Center frequency of the antenna array in Hz. Default: 299792458.0
unsigned pol = 1 (optional input)
Polarization configuration:
pol = 1 |
vertical polarization (default value) |
pol = 2 |
H/V polarized elements, results in 2NM elements |
pol = 3 |
+/-45° polarized elements, results in 2NM elements |
dtype spacing = 0.5 (optional input)
Spacing between elements in wavelengths (λ). Default: 0.5
dtype az_3dB = 120.0 (optional input)
3dB beamwidth in azimuth, specified in degrees. Default: 120.0
dtype el_3db = 120.0 (optional input)
3dB beamwidth in elevation, specified in degrees. Default: 120.0
dtype rear_gain_lin = 0.0 (optional input)
Front-to-back gain ratio as a linear value. Default: 0.0 (no rear gain)
dtype res = 1.0 (optional input)
Resolution of the antenna pattern sampling grid, specified in degrees. Default: 1.0
bool separate_beams = false (optional input)
If set to true, create a separate beam for each angle pair (ignores weights)
bool apply_weights = false (optional input)
Switch to apply the beam-forming weights
-
Returns:
arrayant<dtype>
Array antenna object containing the generated multibeam pattern.
-
Example:
// Generate a pattern with 2 weighted beams
double freq = 3.75e9;
arma::vec az = {20.0, 0.0};
arma::vec el = {-7.0, 30.0};
arma::vec weight = {2.0, 1.0};
auto ant = quadriga_lib::generate_arrayant_multibeam<double>(6, 6, az, el, weight, freq);
Channel class
channel - Class for storing and managing MIMO channel data and associated metadata
-
Description:
- A channel object represents MIMO path-level channel data between antenna arrays over multiple time snapshots.
- Each snapshot may have a different number of propagation paths.
- Contains structured fields for positions, delays, gains, angles, coefficients, and more.
- Supports optional metadata via
par_names and par_data.
- Allowed datatypes (
dtype): float and double
-
Attributes:
std::string name |
Name of the channel object |
arma::Col<dtype> center_frequency |
Center frequency in [Hz], size [1] or [n_snap] or [] |
arma::Mat<dtype> tx_pos |
Transmitter positions, size [3, n_snap] or [3, 1] |
arma::Mat<dtype> rx_pos |
Receiver positions, size [3, n_snap] or [3, 1] |
arma::Mat<dtype> tx_orientation |
Transmitter orientation (Euler), size [3, n_snap], [3, 1], or [] |
arma::Mat<dtype> rx_orientation |
Receiver orientation (Euler), size [3, n_snap], [3, 1], or [] |
std::vector<arma::Cube<dtype>> coeff_re |
Channel coefficients, real part, size [n_rx, n_tx, n_path] per snapshot |
std::vector<arma::Cube<dtype>> coeff_im |
Channel coefficients, imaginary part, same size as coeff_re |
std::vector<arma::Cube<dtype>> delay |
Path delays [s], size [n_rx, n_tx, n_path] or [1,1,n_path] per snapshot |
std::vector<arma::Col<dtype>> path_gain |
Path gains (pre-pattern), length [n_path] per snapshot |
std::vector<arma::Col<dtype>> path_length |
Path lengths TX→RX [m], length [n_path] per snapshot |
std::vector<arma::Mat<dtype>> path_polarization |
Polarization matrix, size [8, n_path] per snapshot |
std::vector<arma::Mat<dtype>> path_angles |
Departure/arrival angles, size [n_path, 4], columns: AOD, EOD, AOA, EOA |
std::vector<arma::Mat<dtype>> path_fbs_pos |
First-bounce scatterer positions, size [3, n_path] |
std::vector<arma::Mat<dtype>> path_lbs_pos |
Last-bounce scatterer positions, size [3, n_path] |
std::vector<arma::Col<unsigned>> no_interact |
Number of interactions per path, length [n_path] per snapshot |
std::vector<arma::Mat<dtype>> interact_coord |
Coordinates of all interactions, size [3, sum(no_interact)] per snapshot |
std::vector<std::string> par_names |
Names of unstructured parameters |
std::vector<std::any> par_data |
Unstructured metadata fields (e.g., string, scalar, matrix) |
int initial_position |
Index of reference snapshot (0-based) |
-
Simple member functions:
.n_snap() |
Returns the number of snapshots |
.n_rx() |
Returns the number of receive antennas (0 if coeffs not present) |
.n_tx() |
Returns the number of transmit antennas (0 if coeffs not present) |
.n_path() |
Returns the number of paths per snapshot as vector |
.empty() |
Returns true if the object has no channel data |
.is_valid() |
Returns an empty string if object is valid, else an error message |
-
Complex member functions:
.add_paths - Append new propagation paths to an existing channel snapshot
-
Description:
- Adds path-level channel data to a specific snapshot (
i_snap) in a channel object.
- All fields provided must be consistent in length (
n_path_add) and structure.
- The number of antennas must match existing entries for the snapshot.
- Existing fields in the channel object must also be provided to this method if relevant.
- Member function of channel
- Allowed datatypes (
dtype): float or double
-
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);
-
Arguments:
arma::uword i_snap (input)
Index of the snapshot to which the new paths should be added (0-based).
const arma::Cube<dtype> *coeff_re_add (optional input)
Real part of channel coefficients. Size: [n_rx, n_tx, n_path_add].
const arma::Cube<dtype> *coeff_im_add (optional input)
Imaginary part of channel coefficients. Size: [n_rx, n_tx, n_path_add].
const arma::Cube<dtype> *delay_add (optional input)
Propagation delay in seconds. Size: [n_rx, n_tx, n_path_add] or [1, 1, n_path_add].
const arma::u32_vec *no_interact_add (optional input)
Number of interaction points per path. Length: [n_path_add].
const arma::Mat<dtype> *interact_coord_add (optional input)
Coordinates of interaction points. Size: [3, sum(no_interact)].
const arma::Col<dtype> *path_gain_add (optional input)
Path gains before antenna effects. Length: [n_path_add].
const arma::Col<dtype> *path_length_add (optional input)
Path lengths from TX to RX phase center. Length: [n_path_add].
const arma::Mat<dtype> *path_polarization_add (optional input)
Polarization transfer matrices (interleaved). Size: [8, n_path_add].
const arma::Mat<dtype> *path_angles_add (optional input)
Departure and arrival angles. Size: [n_path_add, 4], format: {AOD, EOD, AOA, EOA}.
const arma::Mat<dtype> *path_fbs_pos_add (optional input)
First-bounce scatterer positions. Size: [3, n_path_add].
const arma::Mat<dtype> *path_lbs_pos_add (optional input)
Last-bounce scatterer positions. Size: [3, n_path_add].
-
Notes:
- Any provided input must match the snapshot structure and existing fields of the
channel object.
- This method does not update
tx_pos, rx_pos, or orientation fields.
.calc_effective_path_gain - Calculate the effective path gain for each snapshot (in linear scale)
.write_paths_to_obj_file - Export propagation paths to a Wavefront OBJ file
-
Description:
- Writes ray-traced propagation paths to a
.obj file for 3D visualization (e.g., in Blender).
- Each path is represented as a tube, optionally colored by path gain using a selected colormap.
- The function supports filtering by path gain, maximum number of paths, and snapshot index.
- Tube radius and detail can be customized.
- Member function of channel
- Allowed datatypes (
dtype): float or double
-
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;
-
Arguments:
std::string fn (input)
Path to the output .obj file.
arma::uword max_no_paths = 0 (optional input)
Maximum number of paths to be visualized. A value of 0 includes all available paths above gain_min.
Default: 0.
dtype gain_max = -60.0 (optional input)
Maximum path gain (in dB) used for color-coding. Paths with higher gain are clipped. Default: -60.0.
dtype gain_min = -140.0 (optional input)
Minimum path gain (in dB) for color-coding and optional path filtering. Default: -140.0.
std::string colormap = "jet" (optional input)
Name of the colormap to be used for path coloring.
Supported maps: jet, parula, winter, hot, turbo, copper, spring, cool, gray,
autumn, summer. Default: "jet".
arma::uvec i_snap = {} (optional input)
Indices of the snapshots to be included (0-based). Empty vector exports all snapshots. Default: {}.
dtype radius_max = 0.05 (optional input)
Maximum radius (in meters) of the visualized tube geometry. Default: 0.05.
dtype radius_min = 0.01 (optional input)
Minimum radius (in meters) of the visualized tube geometry. Default: 0.01.
arma::uword n_edges = 5 (optional input)
Number of vertices used to create each tube cross-section. Must be ≥ 3. Default: 5.
-
See also:
Channel functions
any_type_id - Get type ID and raw access from a 'std::any' object
baseband_freq_response - Compute the baseband frequency response of a MIMO channel
-
Description:
- Computes the frequency-domain response of a time-domain MIMO channel using a discrete Fourier transform (DFT).
- Input consists of real and imaginary channel coefficients and corresponding delays for each MIMO sub-link.
- Outputs the complex channel response matrix
H at given sub-carrier frequency positions.
- Internally uses AVX2 instructions for fast parallel computation of 8 carriers at once.
- Can be efficiently called in a loop (e.g., over snapshots) and parallelized with OpenMP.
- Internal arithmetic is performed in single precision for speed.
- Allowed datatypes (
dtype): float or double
-
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);
-
Arguments:
const arma::Cube<dtype> *coeff_re (input)
Real part of channel coefficients in time domain, Size [n_rx, n_tx, n_path].
const arma::Cube<dtype> *coeff_im (input)
Imaginary part of channel coefficients in time domain, Size [n_rx, n_tx, n_path].
const arma::Cube<dtype> *delay (input)
Path delays in seconds. Size [n_rx, n_tx, n_path] or broadcastable shape [1, 1, n_path].
const arma::Col<dtype> *pilot_grid (input)
Normalized sub-carrier positions relative to bandwidth. Range: 0.0 (center freq) to 1.0 (center + bandwidth). Length: n_carriers.
const double bandwidth (input)
Total baseband bandwidth in Hz (defines absolute frequency spacing of the pilot grid).
arma::Cube<dtype> *hmat_re (output)
Output: Real part of the frequency-domain channel matrix, Size [n_rx, n_tx, n_carriers].
arma::Cube<dtype> *hmat_im (output)
Output: Imaginary part of the frequency-domain channel matrix, Size [n_rx, n_tx, n_carriers].
baseband_freq_response_vec - Compute the baseband frequency response of multiple MIMO channels
-
Description:
- Computes the frequency-domain response of a batch of time-domain MIMO channels using a discrete Fourier transform (DFT).
- This function wraps
quadriga_lib::baseband_freq_response and applies it across multiple snapshots in parallel using OpenMP.
- Input consists of vectors of real/imaginary coefficients and delay Cubes for each snapshot.
- Output is a vector of frequency-domain channel matrices
H (one per snapshot).
- Can optionally compute a selected subset of snapshots using
i_snap.
- Internal arithmetic is performed in single precision for performance.
- Allowed datatypes (
dtype): float or double
-
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);
-
Arguments:
const std::vector<arma::Cube<dtype>> *coeff_re (input)
Real part of channel coefficients, vector of length n_snap. Each cube has shape [n_rx, n_tx, n_path].
const std::vector<arma::Cube<dtype>> *coeff_im (input)
Imaginary part of channel coefficients, same structure as coeff_re.
const std::vector<arma::Cube<dtype>> *delay (input)
Path delays in seconds, same structure as coeff_re, shape can be broadcasted [1, 1, n_path].
const arma::Col<dtype> *pilot_grid (input)
Normalized sub-carrier positions relative to bandwidth. Range: 0.0 (center freq) to 1.0 (center + bandwidth). Length: n_carriers.
const double bandwidth (input)
Total baseband bandwidth in Hz, used to compute sub-carrier frequencies.
std::vector<arma::Cube<dtype>> *hmat_re (output)
Output: Real part of the frequency-domain channel matrices. Vector of length n_out. Each cube is [n_rx, n_tx, n_carriers].
std::vector<arma::Cube<dtype>> *hmat_im (output)
Output: Imaginary part of the frequency-domain channel matrices. Same structure as hmat_re.
const arma::u32_vec *i_snap = nullptr (optional input)
Optional subset of snapshot indices to process. If omitted, all n_snap snapshots are processed. Length: n_out.
-
See also:
get_hdf5_version - Get the version of the linked HDF5 library
hdf5_create - Create a new HDF5 channel file with a defined storage layout
-
Description:
- Quadriga-Lib offers an HDF5-based method for storing and managing channel data. A key feature of this
library is its ability to organize multiple channels within a single HDF5 file while enabling access
to individual data sets without the need to read the entire file.
- This function initializes a new HDF5 file for storing wireless channel data.
- Defines a multi-dimensional array layout for organizing channels (up to 4 dimensions).
- Typical usage: map base stations (BS), user equipment (UE), and frequencies to dimensions.
- The layout can later be reshaped using
hdf5_reshape_layout, provided the total number of entries remains constant.
- Each combination of indices corresponds to a storage slot that can hold a channel object.
-
Declaration:
void quadriga_lib::hdf5_create(
std::string fn,
unsigned nx = 65536,
unsigned ny = 1,
unsigned nz = 1,
unsigned nw = 1);
-
Arguments:
std::string fn (input)
Filename (including path) of the HDF5 file to be created. If the file exists, it will be overwritten.
unsigned nx = 65536 (input)
Number of entries in the x-dimension, e.g., for base stations (BSs). Default: 65536.
unsigned ny = 1 (input)
Number of entries in the y-dimension, e.g., for user equipment (UEs). Default: 1.
unsigned nz = 1 (input)
Number of entries in the z-dimension, e.g., for frequency points. Default: 1.
unsigned nw = 1 (input)
Number of entries in the w-dimension, e.g., for repetitions, scenarios, or configurations. Default: 1.
-
Example:
// Create a file with layout: [10 BSs, 4 UEs, 2 frequencies]
quadriga_lib::hdf5_create("channels.hdf5", 10, 4, 2);
hdf5_read_layout - Read the HDF5 channel storage layout
hdf5_write - Write channel data to HDF5 file
-
Description:
- Quadriga-Lib provides an HDF5-based solution for storing and organizing channel data.
- This function rites a
channel object to a specified HDF5 file at the given 4D index location.
- If the file does not exist, a new file is created with default layout
(65535 × 1 × 1 × 1).
-
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);
-
Arguments:
const channel<dtype> *ch (input)
Pointer to the channel object to be stored.
std::string fn (input)
Path to the HDF5 file.
unsigned ix = 0 (input)
Index in the x-dimension (e.g., Base Station ID). Default: 0.
unsigned iy = 0 (input)
Index in the y-dimension (e.g., User Equipment ID). Default: 0.
unsigned iz = 0 (input)
Index in the z-dimension (e.g., Frequency point). Default: 0.
unsigned iw = 0 (input)
Index in the w-dimension (e.g., Repetition/Scenario). Default: 0.
bool assume_valid = false (input)
If true, skips channel integrity validation before writing. Default: false.
-
Returns:
0 if a new dataset was created.
1 if an existing dataset was overwritten or extended.
-
Caveat:
- If the file exists already, the new data is added to the exisiting file
- If the index already contains data, it will be overwritten
- Use
assume_valid = true to skip internal validation (faster, but unsafe if data may be corrupted).
- Throws an error if the index was not reserved during
hdf5_create.
- All structured data is written in single precision (but can can be provided as float or double)
- Unstructured datatypes are maintained in the HDF file
- Supported unstructured types: string, double, float, (u)int32, (u)int64
- Supported unstructured size: up to 3 dimensions
- Storage order of the unstructured data is maintained
hdf5_read_channel - Read a channel object from an HDF5 file
-
Description:
- Loads a
quadriga_lib::channel<dtype> object from the specified index in a previously created HDF5 file.
- If the selected index does not contain a valid channel, an empty channel object is returned (with
no_snapshots == 0).
- Allowed datatypes (
dtype): float or double
- All structured data (e.g., channel coefficients, delays, positions) is stored in single precision in the
file but will be converted to the appropriate precision (
float or double) depending on the provided template parameter.
- Unstructured or user-defined fields stored in
std::any containers are not converted and retain their original type
-
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);
-
Arguments:
std::string fn (input)
Filename of the HDF5 file containing the channel data.
unsigned ix = 0 (input)
x-Index in the file's 4D storage layout. Default: 0.
unsigned iy = 0 (input)
y-Index in the file's 4D storage layout. Default: 0.
unsigned iz = 0 (input)
z-Index in the file's 4D storage layout. Default: 0.
unsigned iw = 0 (input)
w-Index in the file's 4D storage layout. Default: 0.
-
Returns:
quadriga_lib::channel<dtype>
A channel object containing the channel data at the specified index. If no data is found, an empty channel object is returned.
hdf5_reshape_layout - Reshape the storage layout of an HDF5 channel file
-
Description:
- Updates the multi-dimensional storage layout of an existing HDF5 file that contains channel data.
- The layout consists of up to four dimensions:
{nx, ny, nz, nw}.
- This is useful for reorganizing data after initial creation (e.g., grouping entries by BS/UE/frequency).
- The total number of entries must remain unchanged, i.e.,
nx × ny × nz × nw must equal the original layout.
- Throws an error if the reshaped layout violates this constraint.
-
Declaration:
void quadriga_lib::hdf5_reshape_layout(
std::string fn,
unsigned nx,
unsigned ny = 1,
unsigned nz = 1,
unsigned nw = 1);
-
Arguments:
std::string fn (input)
Filename of the HDF5 file to reshape.
unsigned nx (input)
Number of entries in the first dimension (e.g., base stations).
unsigned ny = 1 (input)
Number of entries in the second dimension (e.g., user equipment). Default: 1.
unsigned nz = 1 (input)
Number of entries in the third dimension (e.g., carrier frequency). Default: 1.
unsigned nw = 1 (input)
Number of entries in the fourth dimension. Default: 1.
hdf5_read_dset - Read an unstructured dataset from an HDF5 file
-
Description:
- Reads a user-defined, unstructured dataset stored in an HDF5 file at the specified index.
- Unstructured datasets are typically used to store additional parameters or metadata and are stored under a name prefix (e.g.
"par_") followed by the dataset name.
- Returns the dataset as a
std::any object, which can hold any supported type (e.g., scalar values, vectors, matrices, cubes).
- Use
quadriga_lib::any_type_id to determine the contained type and obtain a raw pointer for direct access.
- If the dataset does not exist at the specified index or name, an empty
std::any object is returned.
-
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_");
-
Arguments:
std::string fn (input)
Filename of the HDF5 file containing the dataset.
std::string par_name (input)
Name of the dataset, without the prefix (e.g., "carrier_frequency").
unsigned ix = 0 (input)
x-Index in the HDF5 file’s 4D storage layout. Default: 0.
unsigned iy = 0 (input)
y-Index in the HDF5 file’s 4D storage layout. Default: 0.
unsigned iz = 0 (input)
z-Index in the HDF5 file’s 4D storage layout. Default: 0.
unsigned iw = 0 (input)
w-Index in the HDF5 file’s 4D storage layout. Default: 0.
std::string prefix = "par_" (input)
Optional dataset name prefix. Default is "par_".
-
Returns:
- A
std::any object containing the dataset. If the dataset is not present, the return value is an empty std::any.
hdf5_read_dset_names - Read names of unstructured datasets from an HDF5 file
-
Description:
- Retrieves the names of all unstructured datasets stored at a specified 4D index in an HDF5 file.
- Dataset names are identified by a common prefix (default:
"par_") and the actual parameter name follows the prefix.
- The returned names in
par_names exclude the prefix for convenience.
- Returns the number of datasets found at the given index.
-
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_");
-
Arguments:
std::string fn (input)
Filename of the HDF5 file containing the datasets.
std::vector<std::string> *par_names (output)
Pointer to a vector that will be filled with all dataset names found at the specified index (excluding the prefix).
unsigned ix = 0 (input)
x-Index in the HDF5 file’s 4D layout. Default: 0.
unsigned iy = 0 (input)
y-Index in the HDF5 file’s 4D layout. Default: 0.
unsigned iz = 0 (input)
z-Index in the HDF5 file’s 4D layout. Default: 0.
unsigned iw = 0 (input)
w-Index in the HDF5 file’s 4D layout. Default: 0.
std::string prefix = "par_" (input)
Optional prefix string used to identify unstructured datasets. Default: "par_".
-
Returns:
- The number of unstructured datasets found at the specified location with the given prefix.
hdf5_write_dset - Write a single unstructured dataset to an HDF5 file
-
Description:
- Writes a single unstructured data field to an HDF5 file at the specified 4D index.
- Supported scalar types:
std::string, unsigned, int, long long, unsigned long long, float, double
- Supported Armadillo types:
arma::Col, arma::Row, arma::Mat, and arma::Cube with float, double, int, unsigned, sword, uword, unsigned long long
arma::Row vectors are converted to column vectors (arma::Col) before writing.
- The dataset name is prefixed with
"par_" (default) unless another prefix is specified.
- Throws an error for unsupported types.
- Dataset names may only include alphanumeric characters and underscores.
-
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_");
-
Arguments:
std::string fn (input)
Filename of the HDF5 file to which the dataset will be written.
std::string par_name (input)
Name of the parameter to store (without prefix). Must contain only letters, digits, and underscores.
const std::any *par_data (input)
Pointer to the data to be written. Type must be supported (see above).
unsigned ix = 0 (input)
x-Index in the HDF5 file’s 4D layout. Default: 0.
unsigned iy = 0 (input)
y-Index in the HDF5 file’s 4D layout. Default: 0.
unsigned iz = 0 (input)
z-Index in the HDF5 file’s 4D layout. Default: 0.
unsigned iw = 0 (input)
w-Index in the HDF5 file’s 4D layout. Default: 0.
std::string prefix = "par_" (input)
Optional prefix for the dataset name. Default: "par_".
Channel generation functions
get_channels_ieee_indoor - Generate indoor MIMO channel realizations for IEEE TGn/TGac/TGax/TGah models
-
Description:
- Generates one or multiple indoor channel realizations based on IEEE TGn/TGac/TGax/TGah model definitions.
- 2D model: no elevation angles are used; azimuth angles and planar motion are considered.
- Supports channel model types
A, B, C, D, E, F (as defined by TGn) via ChannelType.
- Can generate MU-MIMO channels (
n_users > 1) with per-user distances/floors and optional angle offsets according to TGac
- Optional time evolution via
observation_time, update_rate, and mobility parameters.
-
Declaration:
std::vector<quadriga_lib::channel<double>> quadriga_lib::get_channels_ieee_indoor(
const arrayant<double> &ap_array,
const 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 );
-
Arguments:
const arrayant<double> ap_array (input)
Access point array antenna with n_tx elements (= ports after element coupling).
const arrayant<double> sta_array (input)
Mobile station array antenna with n_rx elements (= ports after element coupling).
std::string ChannelType (input)
Channel model type as defined by TGn. Supported: A, B, C, D, E, F.
double CarrierFreq_Hz = 5.25e9 (optional input)
Carrier frequency in Hz.
double tap_spacing_s = 10.0e-9 (optional input)
Tap spacing in seconds. Must be equal to 10 ns / 2^k (TGn default = 10e-9).
arma::uword n_users = 1 (optional input)
Number of users (only for TGac, TGah). Output vector length equals n_users.
double observation_time = 0.0 (optional input)
Channel observation time in seconds. 0.0 creates a static channel.
double update_rate = 1.0e-3 (optional input)
Channel update interval in seconds (only relevant when observation_time > 0).
double speed_station_kmh = 0.0 (optional input)
Station movement speed in km/h. Movement direction is AoA_offset. Only relevant when observation_time > 0.
double speed_env_kmh = 1.2 (optional input)
Environment movement speed in km/h. Default 1.2 for TGn, use 0.089 for TGac. Only relevant when
observation_time > 0.
arma::vec Dist_m = {4.99} (optional input)
TX-to-RX distance(s) in meters. Length n_users or length 1 (same distance for all users). Size
[n_users] or [1].
arma::uvec n_floors = {0} (optional input)
Number of floors for TGah model (per user), up to 4 floors. Length n_users or length 1. Size
[n_users] or [1].
bool uplink = false (optional input)
Channel direction flag. Default is downlink; set to true to generate reverse (uplink) direction.
arma::mat offset_angles = {} (optional input)
Offset angles in degree for MU-MIMO channels. Empty uses model defaults (TGac auto for n_users > 1).
Size [4, n_users] with rows: AoD LOS, AoD NLOS, AoA LOS, AoA NLOS.
arma::uword n_subpath = 20 (optional input)
Number of sub-paths per path/cluster used for Laplacian angular spread mapping.
double Doppler_effect = 50.0 (optional input)
Special Doppler effects: models D, E (fluorescent lights, value = mains freq.) and F (moving
vehicle speed in km/h). Use 0.0 to disable.
arma::sword seed = -1 (optional input)
Numeric seed for repeatability. -1 disables the fixed seed and uses the system random device.
double KF_linear = NAN (optional input)
Overwrites the model-specific KF-value. If this parameter is NAN (default) or negative, model defaults are used:
A/B/C (KF = 1 for d < dBP, 0 otherwise); D (KF = 2 for d < dBP, 0 otherwise); E/F (KF = 4 for d < dBP, 0 otherwise).
KF is applied to the first tap only. Breakpoint distance is ignored for KF_linear >= 0.
double XPR_NLOS_linear = NAN (optional input)
Overwrites the model-specific Cross-polarization ratio. If this parameter is NAN (default) or negative,
the model default of 2 (3 dB) is used. XPR is applied to all NLOS taps.
double SF_std_dB_LOS = NAN (optional input)
Overwrites the model-specific shadow fading for LOS channels. If this parameter is NAN (default),
the model default of 3 dB is used. SF_std_dB_LOS is applied to all LOS channels, where the
AP-STA distance d < dBP.
double SF_std_dB_NLOS = NAN (optional input)
Overwrites the model-specific shadow fading for LOS channels. If this parameter is NAN (default),
the model defaults are A/B: 4 dB, C/D: 5 dB, E/F: 6 dB. SF_std_dB_NLOS is applied to all NLOS channels,
where the AP-STA distance d >= dBP.
double dBP_m = NAN (optional input)
Overwrites the model-specific breakpoint distance. If this parameter is NAN (default) or negative,
the model defaults are A/B/C: 5 m, D: 10 m, E: 20 m, F: 30 m.
-
Returns:
std::vector<quadriga_lib::channel<double>> (output)
Vector of channel objects with length n_users. Each entry contains the generated indoor channel
realization for one user (including direction determined by uplink).
get_channels_irs - Calculate channel coefficients for intelligent reflective surfaces (IRS)
-
Description:
- Calculates MIMO channel coefficients and delays for IRS-assisted communication using two channel segments:
1. TX → IRS; 2. IRS → RX
- The IRS is modeled as a passive antenna array with phase shifts defined via its coupling matrix.
- IRS codebook entries can be selected via a port index (
i_irs).
- Supports combining paths from both segments to form
n_path_irs valid output paths, subject to a gain threshold.
- Optional second IRS array allows different antenna behavior for TX-IRS and IRS-RX directions.
- Returns a boolean vector indicating which path combinations are included in the output.
- Allowed datatypes (
dtype): float or double
-
Declaration:
std::vector<bool> quadriga_lib::get_channels_irs(
const arrayant<dtype> *tx_array,
const arrayant<dtype> *rx_array,
const 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 arrayant<dtype> *irs_array_2 = nullptr,
const std::vector<bool> *active_path = nullptr);
-
Arguments:
const arrayant<dtype> *tx_array (input)
Pointer to the transmit antenna array object (with n_tx elements).
const arrayant<dtype> *rx_array (input)
Pointer to the receive antenna array object (with n_rx elements).
const arrayant<dtype> *irs_array (input)
Pointer to the IRS array antenna (with n_irs elements).
dtype Tx, Ty, Tz (input)
Transmitter position in Cartesian coordinates [m].
dtype Tb, Tt, Th (input)
Transmitter orientation (Euler angles) in radians.
dtype Rx, Ry, Rz (input)
Receiver position in Cartesian coordinates [m].
dtype Rb, Rt, Rh (input)
Receiver orientation (Euler angles) in radians.
dtype Ix, Iy, Iz (input)
IRS position in Cartesian coordinates [m].
dtype Ib, It, Ih (input)
IRS orientation (Euler angles) in radians.
const arma::Mat<dtype> *fbs_pos_1 (input)
First-bounce scatterer positions of TX → IRS paths, Size [3, n_path_1].
const arma::Mat<dtype> *lbs_pos_1 (input)
Last-bounce scatterer positions of TX → IRS paths, Size [3, n_path_1].
const arma::Col<dtype> *path_gain_1 (input)
Path gains (linear) for TX → IRS paths, Length n_path_1.
const arma::Col<dtype> *path_length_1 (input)
Path lengths for TX → IRS paths, Length n_path_1.
const arma::Mat<dtype> *M_1 (input)
Polarization transfer matrix for TX → IRS paths, Size [8, n_path_1].
const arma::Mat<dtype> *fbs_pos_2 (input)
First-bounce scatterer positions of IRS → RX paths, Size [3, n_path_2].
const arma::Mat<dtype> *lbs_pos_2 (input)
Last-bounce scatterer positions of IRS → RX paths, Size [3, n_path_2].
const arma::Col<dtype> *path_gain_2 (input)
Path gains (linear) for IRS → RX paths, Length n_path_2.
const arma::Col<dtype> *path_length_2 (input)
Path lengths for IRS → RX paths, Length n_path_2.
const arma::Mat<dtype> *M_2 (input)
Polarization transfer matrix for IRS → RX paths, Size [8, n_path_2].
arma::Cube<dtype> *coeff_re (output)
Real part of resulting IRS-assisted channel coefficients, Size [n_rx, n_tx, n_path_irs].
arma::Cube<dtype> *coeff_im (output)
Imaginary part of channel coefficients, Size [n_rx, n_tx, n_path_irs].
arma::Cube<dtype> *delay (output)
Propagation delays in seconds, Size [n_rx, n_tx, n_path_irs].
arma::uword i_irs = 0 (optional input)
Index of IRS codebook entry (port number), Default: 0.
dtype threshold_dB = -140.0 (optional input)
Threshold (in dB) below which paths are discarded. Default: -140.0.
dtype center_frequency = 0.0 (optional input)
Center frequency in Hz; 0.0 disables phase computation. Default: 0.0.
bool use_absolute_delays = false (optional input)
If true, includes LOS delay in all paths. Default: false.
arma::Cube<dtype> *aod = nullptr (optional output)
Azimuth of Departure angles [rad], Size [n_rx, n_tx, n_path_irs].
arma::Cube<dtype> *eod = nullptr (optional output)
Elevation of Departure angles [rad], Size [n_rx, n_tx, n_path_irs].
arma::Cube<dtype> *aoa = nullptr (optional output)
Azimuth of Arrival angles [rad], Size [n_rx, n_tx, n_path_irs].
arma::Cube<dtype> *eoa = nullptr (optional output)
Elevation of Arrival angles [rad], Size [n_rx, n_tx, n_path_irs].
const arrayant<dtype> *irs_array_2 = nullptr (optional input)
Optional second IRS array (TX side) for asymmetric IRS behavior.
const std::vector<bool> *active_path = nullptr (optional input)
Optional bitmask for selecting active TX-IRS and IRS-RX path pairs. Ignores threshold_dB when provided.
-
Returns:
std::vector<bool>
Boolean mask of length n_path_1 * n_path_2, indicating which path combinations were used.
-
See also:
get_channels_planar - Calculate channel coefficients for planar waves
-
Description:
- Calculates MIMO channel coefficients and delays for a set of planar wave paths between two antenna arrays.
- Interpolates antenna patterns (including orientation and polarization) for both transmitter and receiver arrays.
- Supports LOS path identification based on distance (angles are ignored).
- Polarization transfer matrix models polarization coupling and must be normalized.
- Doppler weights can optionally be calculated from receiver motion relative to path direction.
- Element positions and antenna orientation are fully considered for delay and phase.
- Allowed datatypes (
dtype): float or double
-
Declaration:
void quadriga_lib::get_channels_planar(
const arrayant<dtype> *tx_array,
const 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);
-
Arguments:
const arrayant<dtype> *tx_array (input)
Pointer to the transmit antenna array object (with n_tx elements).
const arrayant<dtype> *rx_array (input)
Pointer to the receive antenna array object (with n_rx elements).
dtype Tx, Ty, Tz (input)
Transmitter position in Cartesian coordinates [m].
dtype Tb, Tt, Th (input)
Transmitter orientation (Euler) angles (bank, tilt, head) in [rad].
dtype Rx, Ry, Rz (input)
Receiver position in Cartesian coordinates [m].
dtype Rb, Rt, Rh (input)
Receiver orientation (Euler) angles (bank, tilt, head) in [rad].
const arma::Col<dtype> *aod (input)
Departure azimuth angles in radians, Length n_path.
const arma::Col<dtype> *eod (input)
Departure elevation angles in radians, Length n_path.
const arma::Col<dtype> *aoa (input)
Arrival azimuth angles in radians, Length n_path.
const arma::Col<dtype> *eoa (input)
Arrival elevation angles in radians, Length n_path.
const arma::Col<dtype> *path_gain (input)
Path gains in linear scale, Length n_path.
const arma::Col<dtype> *path_length (input)
Path lengths from TX to RX phase center, Length n_path.
const arma::Mat<dtype> *M (input)
Polarization transfer matrix of size [8, n_path], interleaved: (ReVV, ImVV, ReVH, ImVH, ReHV, ImHV, ReHH, ImHH).
arma::Cube<dtype> *coeff_re (output)
Real part of channel coefficients, Size [n_rx, n_tx, n_path(+1)].
arma::Cube<dtype> *coeff_im (output)
Imaginary part of channel coefficients, Size [n_rx, n_tx, n_path(+1)].
arma::Cube<dtype> *delay (output)
Propagation delays in seconds, Size [n_rx, n_tx, n_path(+1)].
dtype center_frequency = 0.0 (optional input)
Center frequency in Hz; set to 0 to disable phase calculation. Default: 0.0
bool use_absolute_delays = false (optional input)
If true, includes LOS delay in all paths. Default: false
bool add_fake_los_path = false (optional input)
Adds a zero-power LOS path if no LOS is present. Default: false
arma::Col<dtype> *rx_Doppler = nullptr (optional output)
Doppler weights for moving RX, Length n_path(+1). Positive = towards path, Negative = away.
get_channels_spherical - Calculate channel coefficients for spherical waves
-
Description:
- Calculates MIMO channel coefficients and delays for a set of spherical wave paths between two antenna arrays.
- Interpolates antenna patterns (including orientation and polarization) for both transmitter and receiver arrays.
- Accurately models path-based propagation using provided scatterer positions.
- Supports LOS path identification and handles complex polarization coupling.
- Element positions and antenna orientation are fully considered for delay and phase.
- Allowed datatypes (
dtype): float or double
-
Declaration:
void quadriga_lib::get_channels_spherical(
const arrayant<dtype> *tx_array,
const 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)
-
Arguments:
const arrayant<dtype> *tx_array (input)
Pointer to the transmit antenna array object (with n_tx elements).
const arrayant<dtype> *rx_array (input)
Pointer to the receive antenna array object (with n_rx elements).
dtype Tx, Ty, Tz (input)
Transmitter position in Cartesian coordinates [m].
dtype Tb, Tt, Th (input)
Transmitter orientation (Euler) angles (bank, tilt, head) in [rad].
dtype Rx, Ry, Rz (input)
Receiver position in Cartesian coordinates [m].
dtype Rb, Rt, Rh (input)
Receiver orientation (Euler) angles (bank, tilt, head) in [rad].
const arma::Mat<dtype> *fbs_pos (input)
First-bounce scatterer positions, Size: [3, n_path].
const arma::Mat<dtype> *lbs_pos (input)
Last-bounce scatterer positions, Size: [3, n_path].
const arma::Col<dtype> *path_gain (input)
Path gains in linear scale, Length n_path.
const arma::Col<dtype> *path_length (input)
Path lengths from TX to RX phase center Length n_path.
const arma::Mat<dtype> *M (input)
Polarization transfer matrix of size [8, n_path], interleaved: (ReVV, ImVV, ReVH, ImVH, ReHV, ImHV, ReHH, ImHH).
arma::Cube<dtype> *coeff_re (output)
Real part of channel coefficients, Size [n_rx, n_tx, n_path].
arma::Cube<dtype> *coeff_im (output)
Imaginary part of channel coefficients, Size [n_rx, n_tx, n_path].
arma::Cube<dtype> *delay (output)
Propagation delays in seconds, Size [n_rx, n_tx, n_path].
dtype center_frequency = 0.0 (optional input)
Center frequency in Hz; set to 0 to disable phase calculation. Default: 0.0
bool use_absolute_delays = false (optional input)
If true, includes LOS delay in all paths. Default: false
bool add_fake_los_path = false (optional input)
Adds a zero-power LOS path if no LOS is present. Default: false
arma::Cube<dtype> *aod = nullptr (optional output)
Azimuth of Departure angles in radians, Size [n_rx, n_tx, n_path].
arma::Cube<dtype> *eod = nullptr (optional output)
Elevation of Departure angles in radians, Size [n_rx, n_tx, n_path].
arma::Cube<dtype> *aoa = nullptr (optional output)
Azimuth of Arrival angles in radians, Size [n_rx, n_tx, n_path].
arma::Cube<dtype> *eoa = nullptr (optional output)
Elevation of Arrival angles in radians, Size [n_rx, n_tx, n_path].
Miscellaneous tools
calc_rotation_matrix - Calculate rotation matrices from Euler angles
-
Description:
- Computes 3D rotation matrices from input Euler angles (bank, tilt, head).
- The result is returned in column-major order as a 3×3 matrix per input orientation vector.
- Calculations are internally performed in double precision for improved numerical accuracy, even if
dtype is float.
- Supports optional inversion of the y-axis and optional transposition of the output matrix.
- Allowed datatypes (
dtype): float or double
-
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);
-
Arguments:
const arma::Cube<dtype> &orientation or const arma::Mat<dtype> &orientation or const arma::Col<dtype> &orientation (input)
Input Euler angles (bank, tilt, head) in radians, Size [3, n_row, n_col] or [3, n_mat] or Size [3].
bool invert_y_axis = false (optional input)
If true, the y-axis of the rotation is inverted. Default: false.
bool transposeR = false (optional input)
If true, the transpose of the rotation matrix is returned. Default: false.
-
Returns:
arma::Cube<dtype> or arma::Mat<dtype> or arma::Col<dtype>
Rotation matrices in column-major ordering. Size [9, n_row, n_col] or [9, n_mat] or [9].
-
Example:
arma::cube ori(3, 1, 1);
ori(0, 0, 0) = 0.0; // bank
ori(1, 0, 0) = 0.0; // tilt
ori(2, 0, 0) = 1.5708; // head
auto R = quadriga_lib::calc_rotation_matrix(ori);
cart2geo - Convert Cartesian coordinates to geographic coordinates (azimuth, elevation, distance)
-
Description:
- Transforms 3D Cartesian coordinates
(x, y, z) into geographic coordinates:
- Azimuth angle [rad]
- Elevation angle [rad]
- Distance (vector norm)
- Azimuth is measured in the x-y plane from the x-axis; elevation is from the x-y plane upward.
- Allowed datatypes (
dtype): float or double
-
Declaration:
arma::Cube<dtype> quadriga_lib::cart2geo(const arma::Cube<dtype> &cart);
arma::Mat<dtype> quadriga_lib::cart2geo(const arma::Mat<dtype> &cart);
arma::Col<dtype> quadriga_lib::cart2geo(const arma::Col<dtype> &cart);
-
Arguments:
const arma::Cube<dtype> cart or const arma::Mat<dtype> cart or const arma::Col<dtype> cart (input)
Cartesian coordinate vectors (x, y, z), Size [3, n_row, n_col] or [3, n_row] or [3].
-
Returns:
arma::Cube<dtype> or arma::Mat<dtype> or arma::Col<dtype>
Geographic coordinate vectors (azimuth, elevation, distance), Size [n_row, n_col, 3] or [n_row, 3] or [3].
-
Example:
arma::vec cart = {1.0, 1.0, 1.0};
auto geo = quadriga_lib::cart2geo(cart);
colormap - Generate colormap
-
Description:
- Returns a 64x3 or 256x3 colormap matrix with RGB values in unsigned char format.
- Each row corresponds to an RGB color entry of the selected colormap.
- Useful for visualization purposes (e.g., heatmaps or 3D rendering).
- Available color maps include:
jet, parula, winter, hot, turbo, copper, spring, cool, gray, autumn, summer.
-
Declaration:
arma::uchar_mat quadriga_lib::colormap(std::string map, bool high_res = false)
-
Arguments:
std::string map (input)
Name of the desired colormap. Must be one of:
"jet", "parula", "winter", "hot", "turbo", "copper", "spring", "cool", "gray", "autumn", "summer".
bool high_res (input)
Enables 256 color steps
-
Returns:
arma::uchar_mat
A matrix of size [64 x 3] or [256 x 3] containing RGB color values as unsigned chars in the range [0, 255].
-
Example:
arma::uchar_mat cm = quadriga_lib::colormap("turbo");
fast_sincos - Fast, approximate sine/cosine
-
Description:
Computes elementwise sine and/or cosine for an Armadillo vector. Designed for high throughput on modern CPUs.
- Operates on input angles in radians
- AVX2-optimized (8 floats per lane); scalar fallback without AVX2 or on non-AVX2 CPUs
- Parallelizes across cores with OpenMP when enabled
- Results are approximate and may differ from
std::sinf / std::cosf
- For x in [-pi, pi], the maximum absolute error is 2^(-22.1), and larger otherwise
- For x in [-500, 500], the maximum absolute error is 2^(-16.0)
- Either output (
s or c) may be nullptr to skip its computation
- Allowed input datatype:
float (Armadillo fvec) or double (Armadillo vec)
-
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);
-
Arguments:
const arma::fvec x or const arma::vec x (input)
Input angles in radians. Size [n].
arma::fvec *s = nullptr (optional output)
If non-null, set to sin(x). Resized to length n if needed. Size [n] or nullptr.
arma::fvec *c = nullptr (optional output)
If non-null, set to cos(x). Resized to length n if needed. Size [n] or nullptr.
-
Returns:
void (output)
No return value. Results written via output pointers.
-
Example:
arma::fvec x = arma::linspace[arma::fvec](arma::fvec)(0.0f, 6.2831853f, 1000);
arma::fvec s, c;
quadriga_lib::fast_sincos(x, &s, &c); // compute both
quadriga_lib::fast_sincos(x, &s, nullptr); // compute sine only
quadriga_lib::fast_sincos(x, nullptr, &c); // compute cosine only
geo2cart - Transform geographic (azimuth, elevation, length) to Cartesian coordinates
-
Description:
- Converts azimuth and elevation angles (in radians) into 3D Cartesian coordinates.
- Optional radial distance (
length) can be provided; otherwise, unit vectors are returned.
- Useful for converting spherical direction data into vector representations.
- Allowed datatypes (
dtype): float or double
-
Declaration:
arma::Cube<dtype> quadriga_lib::geo2cart(
const arma::Mat<dtype> &azimuth,
const arma::Mat<dtype> &elevation,
const arma::Mat<dtype> &length = {})
-
Arguments:
const arma::Mat<dtype> azimuth (input)
Azimuth angles in radians. Size [n_row, n_col].
const arma::Mat<dtype> elevation (input)
Elevation angles in radians. Size [n_row, n_col].
const arma::Mat<dtype> length = {} (optional input)
Radial distance (length). Same size as azimuth/elevation or empty for unit vectors. Size [n_row, n_col] or [0, 0].
-
Returns:
arma::Cube<dtype> (output)
Cartesian coordinates with dimensions [3, n_row, n_col], representing (x, y, z) for each input direction.
-
Example:
arma::mat az(2, 2), el(2, 2), len(2, 2, arma::fill::ones);
auto cart = quadriga_lib::geo2cart(az, el, len);
interp_1d / interp_2d - Perform linear interpolation (1D or 2D) on single or multiple data sets.
-
Description:
- Interpolates given input data at specified output points.
- Supports single and multiple data sets.
- Returns interpolated results either directly or through reference argument.
- Data types (
dtype): float or double
-
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: Input data array/matrix (size details below)
xi: Input x-axis sampling points, vector of length nx
yi: Input y-axis sampling points (for 2D only), vector of length ny
xo: Output x-axis sampling points, vector of length mx
yo: Output y-axis sampling points (for 2D only), vector of length my
output: Interpolated data cube (modified in-place for one variant)
-
Input / Output size details:
- 2D interpolation of multiple datasets (
arma::Cube):
Input size: [ny, nx, ne], Output size: [my, mx, ne]
- 2D interpolation of single dataset (
arma::Mat):
Input size: [ny, nx], Output size: [my, mx]
- 1D interpolation of multiple datasets (
arma::Mat):
Input size: [nx, ne], Output size: [mx, ne]
- 1D interpolation of single dataset (
arma::Col):
Input length: [nx], Output length: [mx]
-
Examples:
- 2D interpolation example:
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);
- 1D interpolation example:
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);
write_png - Write data to a PNG file
-
Description:
- Converts input data into a color-coded PNG file for visualization
- Support optional selection of a colormap, as well a minimum and maximum value limits
- Allowed datatypes (
dtype): float or double
- Uses the LodePNG library for PNG writing
-
Declaration:
void 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);
-
Arguments:
const arma::Mat<dtype> &data
Data matrix
std::string fn
Path to the .png file to be written
std::string colormap
Name of the desired colormap. Must be one of:
"jet", "parula", "winter", "hot", "turbo", "copper", "spring", "cool", "gray", "autumn", "summer".
dtype min_val
Minimum value. Values below this value will have be encoded with the color of the smallest value.
If NAN is provided (default), the lowest values is determined from the data.
dtype max_val
Maximum value. Values above this value will have be encoded with the color of the largest value.
If NAN is provided (default), the largest values is determined from the data.
bool log_transform
If enabled, the data values are transformed to the log-domain (10*log10(data)) before processing.
Default: false (disabled)
-
See also:
Site-Specific Simulation Tools
calc_diffraction_gain - Calculate diffraction gain for multiple transmit and receive positions
-
Description:
Diffraction refers to the phenomenon where waves bend or interfere around the edges of an obstacle,
extending into the region that would otherwise be in the obstacle's geometrical shadow. The object
causing the diffraction acts as a secondary source for the wave's propagation. A specific example of
this is the knife-edge effect, or knife-edge diffraction, where a sharp, well-defined obstacle—like
a mountain range or a building wall—partially truncates the incident radiation.
To estimate the diffraction gain in a three-dimensional space, one can assess the extent to which the
Fresnel ellipsoid is obstructed by objects, and then evaluate the impact of this obstruction on the
received power. This method presupposes that diffracted waves travel along slightly varied paths
before arriving at a receiver. These waves may reach the receiver out of phase with the primary wave
due to their different travel lengths, leading to either constructive or destructive interference.
The process of estimating the gain involves dividing the wave propagation from a transmitter to a
receiver into n_path paths. These paths are represented by elliptic arcs, which are further
approximated using n_seg line segments. Each segment can be individually blocked or attenuated
by environmental objects. To determine the overall diffraction gain, a weighted sum of these
individual path contributions is calculated. The weighting is adjusted to align with the uniform
theory of diffraction (UTD) coefficients in two dimensions, but the methodology is adapted for
any 3D object shape.
- This function computes the diffraction gain between multiple transmit and receive positions using a 3D triangular mesh.
- Supports accelerated computation via mesh segmentation: large meshes can be split into smaller sub-meshes to skip irrelevant geometry.
- Optionally returns the gain per link and the coordinates along the diffracted path.
- Allowed datatypes (
dtype): float or double.
-
Declaration:
void quadriga_lib::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);
-
Arguments:
const arma::Mat<dtype> *orig (input)
Matrix of origin (transmitter) points, size [n_pos, 3].
const arma::Mat<dtype> *dest (input)
Matrix of destination (receiver) points, size [n_pos, 3].
const arma::Mat<dtype> *mesh (input)
Triangular mesh geometry. Each row contains 3 vertices, flattened as [X1, Y1, Z1, X2, Y2, Z2, X3, Y3, Z3], size [no_mesh, 9].
const arma::Mat<dtype> *mtl_prop (input)
Material properties per triangle, size [no_mesh, 5], See obj_file_read for details
dtype center_frequency (input)
Center frequency in Hz used for diffraction calculations.
arma::Col<dtype> *gain (output, optional)
Calculated diffraction gain in linear scale, vector of size [n_pos].
arma::Cube<dtype> *coord (output, optional)
Coordinates of the diffracted path, size [3, n_seg-1, n_pos]. Includes one point per segment excluding endpoints.
int verbose = 0 (optional input)
Verbosity level for debug and diagnostic output. Higher values enable more log messages. Default: 0.
const arma::u32_vec *sub_mesh_index (input, optional)
Optional sub-mesh indexing for acceleration. Vector mapping triangles to sub-meshes, length [no_mesh].
-
Technical Notes:
- If
sub_mesh_index is provided, spatial bounding boxes for each sub-mesh are automatically computed.
Any TX-RX path not intersecting a bounding box is excluded from detailed edge evaluation.
- The diffraction algorithm supports multiple edges and considers the geometric configuration and material properties of each triangle.
- The LOD parameter controls both accuracy and computational cost. Lower values use simplified heuristics;
higher values may result in significantly more computation time.
- Each ellipsoid consists of
n_path diffraction paths. The number of paths is determined by the
level of detail (lod). See generate_diffraction_paths for details
-
See also:
combine_irs_coord - Combine path interaction coordinates for channels with intelligent reflective surfaces (IRS)
-
Description:
- Merges two propagation segments — (1) TX → IRS and (2) IRS → RX — into complete TX → RX paths via an IRS.
- Interaction coordinates of both segments are stored in a compressed format where
no_interact is a vector of
length n_path storing the number of interactions per path and interact_coord stores all
interaction coordinates in path order.
- Each combined path includes interaction points from both segments, optionally reversed for each segment.
- The number of output paths
n_path_irs is at most n_path_1 × n_path_2, but can be reduced by specifying active_path.
- Output includes the number and coordinates of interaction points per IRS-reflected path.
- Allowed datatypes (
dtype): float or double.
-
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);
-
Arguments:
dtype Ix, Iy, Iz (input)
IRS position in Cartesian coordinates [m].
const arma::u32_vec *no_interact_1 (input)
Number of interaction points in segment 1 (TX → IRS), vector of length [n_path_1]
const arma::Mat<dtype> *interact_coord_1 (input)
Interaction coordinates for segment 1, matrix of size [3, sum(no_interact_1)].
const arma::u32_vec *no_interact_2 (input)
Number of interaction points in segment 2 (IRS → RX), vector of length [n_path_2]
const arma::Mat<dtype> *interact_coord_2 (input)
Interaction coordinates for segment 2, matrix of size [3, sum(no_interact_2)].
arma::u32_vec *no_interact (output)
Combined number of interaction points per IRS-based path, vector of length [n_path_irs].
arma::Mat<dtype> *interact_coord (output)
Combined interaction coordinates, matrix of size [3, sum(no_interact)].
bool reverse_segment_1 = false (optional input)
If true, reverses the interaction coordinates of segment 1. TX and IRS positions are not swapped. Default: false.
bool reverse_segment_2 = false (optional input)
If true, reverses the interaction coordinates of segment 2. IRS and RX positions are not swapped. Default: false.
const std::vector<bool> *active_path = nullptr (optional input)
Optional mask vector of length [n_path_1 × n_path_2]. If provided, only paths with true are combined.
This is generated as the output of get_channels_irs
-
Technical Notes:
- Paths are created by appending interaction coordinates from both segments. Reversing only affects the order of these coordinates.
- Output matrix
interact_coord is built sequentially and should be preallocated only if the total number of interaction points is known in advance.
- The function supports sparsely activated paths via
active_path, this is generated as a return value of get_channels_irs
and filters out paths that have very little power after being reflected by the IRS.
- This function is typically used as a complementary step to delay or coefficient calculations involving IRS-based channels.
It calculated the required data for visualizing paths, e.g. in Blender or other visualization tools.
-
See also:
- get_channels_irs (for computing IRS channels)
- coord2path (for processing coordinates, calculating departure and arrival angels, etc.)
coord2path - Convert path interaction coordinates into FBS/LBS positions, path length and angles
-
Description:
- Interaction coordinates can be stored in a compressed format where
no_interact is a vector of
length n_path storing the number of interactions per path and interact_coord stores all
interaction coordinates in path order.
- This function processes the interaction coordinates of to extract relevant propagation metrics.
- Calculates absolute path lengths and first/last bounce scatterer positions and angles.
- LOS paths are assigned a virtual FBS/LBS position at the midpoint between TX and RX.
- Automatically resizes output arguments if necessary.
- Optionally reverses TX and RX to simulate the reverse link.
- Allowed datatypes (
dtype): float or double
-
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,
arma::Mat<dtype> *fbs_pos,
arma::Mat<dtype> *lbs_pos,
arma::Mat<dtype> *path_angles,
std::vector<arma::Mat<dtype>> *path_coord,
bool reverse_path);
-
Arguments:
dtype Tx, Ty, Tz (input)
Transmitter position in Cartesian coordinates in [m]
dtype Rx, Ry, Rz (input)
Receiver position in Cartesian coordinates in [m]
const arma::u32_vec *no_interact (input)
Vector of length n_path indicating the number of interactions per path (0 = LOS)
const arma::Mat<dtype> *interact_coord (input)
Matrix of size [3, sum(no_interact)] containing interaction coordinates in path order
arma::Col<dtype> *path_length (output, optional)
Absolute path lengths from TX to RX, Length [n_path]
arma::Mat<dtype> *fbs_pos (output, optional)
First-bounce scatterer positions, Size [3, n_path]; For LOS, set to midpoint between TX and RX
arma::Mat<dtype> *lbs_pos (output, optional)
Last-bounce scatterer positions, Size [3, n_path]. For LOS, set to midpoint between TX and RX
arma::Mat<dtype> *path_angles (output, optional)
Departure and arrival angles: columns {AOD, EOD, AOA, EOA}, size [n_path, 4]
std::vector<arma::Mat<dtype>> *path_coord (output, optional)
Full path coordinates (including TX and RX), vector of size n_path with each entry of size [3, n_interact+2]
bool reverse_path = false (optional input)
If true, swaps TX and RX and reverses all interaction sequences.
-
Example:
// Suppose we have two paths: one LOS, one single-bounce.
// Path 0: LOS (no_interact[0] = 0)
// Path 1: single bounce at {5,0,2}.
arma::u32_vec no_int = {0, 1};
arma::mat interact(3, 1);
interact.col(0) = {5.0, 0.0, 2.0};
arma::vec lengths;
arma::mat fbs, lbs;
arma::mat angles;
std::vector<arma::mat> coords;
quadriga_lib::coord2path<double>(
0.0, 0.0, 0.0, // TX at origin
10.0, 0.0, 0.0, // RX at x = 10 m
&no_int, &interact, &lengths, &fbs, &lbs, &angles, &coords);
// After the call:
// lengths: [10.0, 10.77]
// fbs: Path 1 [5,0,0], Path 2 [5,0,2]
// lbs: Path 1 [5,0,0], Path 2 [5,0,2]
// angles: Path 1 [0, 0, pi, p], Path 2 [0, 22°, pi, 22°]
// coords[0]: [ [0,0,0], [10,0,0] ]
// coords[1]: [ [0,0,0], [5,0,2], [10,0,0] ]
generate_diffraction_paths - Generate propagation paths for estimating the diffraction gain
-
Description:
This function generates the elliptic propagation paths and corresponding weights necessary for the
calculation of the diffraction gain in calc_diffraction_gain.
-
Caveat:
- Each ellipsoid consists of
n_path diffraction paths. The number of paths is determined by the
level of detail (lod).
- All diffraction paths of an ellipsoid originate at
orig and arrive at dest
- Each diffraction path has
n_seg segments
- Points
orig and dest lay on the semi-major axis of the ellipsoid
- The generated rays sample the volume of the ellipsoid
- Weights are calculated from the Knife-edge diffraction model when parts of the ellipsoid are shadowed
- Initial weights are normalized such that
sum(prod(weights,3),2) = 1
- Inputs
orig and dest may be provided as double or single precision
- Supported datatypes
dtype are float or double
-
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);
-
Arguments:
const arma::Mat<dtype> *orig (input)
Pointer to Armadillo matrix containing the origin points of the propagation ellipsoid (e.g.
transmitter positions). Size: [ n_pos, 3 ]
const arma::Mat<dtype> *dest (input)
Pointer to Armadillo matrix containing the destination point of the propagation ellipsoid (e.g.
receiver positions). Size: [ n_pos, 3 ]
dtype center_frequency (input)
The center frequency in [Hz], scalar, default = 299792458 Hz
int lod (input)
Level of detail, scalar value
lod = 1 |
results in n_path = 7 and n_seg = 3 |
lod = 2 |
results in n_path = 19 and n_seg = 3 |
lod = 3 |
results in n_path = 37 and n_seg = 4 |
lod = 4 |
results in n_path = 61 and n_seg = 5 |
lod = 5 |
results in n_path = 1 and n_seg = 2 (for debugging) |
lod = 6 |
results in n_path = 2 and n_seg = 2 (for debugging) |
arma::Cube<dtype> *ray_x (output)
Pointer to an Armadillo cube for the x-coordinates of the generated rays; Size: [ n_pos, n_path, n_seg-1 ]
Size will be adjusted if not set correctly.
arma::Cube<dtype> *ray_y (output)
Pointer to an Armadillo cube for the y-coordinates of the generated rays; Size: [ n_pos, n_path, n_seg-1 ]
Size will be adjusted if not set correctly.
arma::Cube<dtype> *ray_z (output)
Pointer to an Armadillo cube for the z-coordinates of the generated rays; Size: [ n_pos, n_path, n_seg-1 ]
Size will be adjusted if not set correctly.
arma::Cube<dtype> *weight (output)
Pointer to an Armadillo cube for the weights; Size: [ n_pos, n_path, n_seg ]
Size will be adjusted if not set correctly.
-
See also:
icosphere - Construct a geodesic polyhedron (icosphere) from triangles
-
Description:
- Generates a convex polyhedral surface (icosphere) made entirely of triangles, based on recursive subdivision of an icosahedron.
- Useful for sampling directions uniformly on a sphere, for applications like ray tracing, antenna pattern evaluation, or spatial grids.
- Each triangular face points outward from the center and has an associated normal.
- Optionally returns vertex directions and vector lengths in either Cartesian or spherical coordinates.
- Allowed datatypes (
dtype): float or double.
-
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);
-
Arguments:
arma::uword n_div (input)
Number of subdivisions per triangle edge. The total number of faces will be n_faces = 20 × n_div².
dtype radius (input)
Radius of the icosphere in meters. All triangle vertices lie on this sphere.
arma::Mat<dtype> *center (output)
Unit vectors pointing from the origin to the center of each triangle face. Usually a bit shorter than radius, Size [n_faces, 3].
arma::Col<dtype> *length = nullptr (optional output)
Vector magnitudes of each center vector (usually slightly less than radius). Vector size [n_faces].
arma::Mat<dtype> *vert = nullptr (optional output)
Vertex vectors from each triangle’s center to its three vertices, flattened as [x1, y1, z1, x2, y2, z2, x3, y3, z3]. Size [n_faces, 9].
arma::Mat<dtype> *direction = nullptr (optional output)
Direction vectors of the three triangle edges. Format depends on direction_xyz: If false
(spherical): [v1az, v1el, v2az, v2el, v3az, v3el], If true (Cartesian): [v1x, v1y, v1z, v2x, v2y, v2z, v3x, v3y, v3z],
Size [n_faces, 6] or [n_faces, 9].
bool direction_xyz = false (optional input)
If true, output directions in Cartesian coordinates. If false, output in spherical azimuth/elevation. Default: false.
-
Returns:
arma::uword
The number of triangular faces generated: n_faces = 20 × n_div².
-
Technical Notes:
- The generated mesh is well-suited for uniform angular sampling on a sphere.
- Triangle vertices are calculated relative to the center to ensure they lie on the desired sphere.
- The radius parameter scales the final structure without changing angular spacing.
-
Example:
arma::fmat center, vert, direction;
arma::fvec length;
// 4 subdivisions → 320 faces, map to unit sphere, output Cartesian directions
auto n = quadriga_lib::icosphere<float>(4, 1.0, ¢er, &length, &vert, &direction, true);
mitsuba_xml_file_write - Write geometry and material data to a Mitsuba 3 XML scene file.
-
Description:
This routine converts a triangular surface mesh stored in quadriga-lib data structures into the
XML format understood by Mitsuba 3 www.mitsuba-renderer.org.
The generated file can be loaded directly by NVIDIA Sionna RT for differentiable radio-propagation
simulations.
- Converts a 3D geometry mesh into Mitsuba 3 XML format for use with rendering tools.
- Enables exporting models from
quadriga-lib to be used with Mitsuba 3 or Sionna RT:
- Mitsuba 3: Research-oriented retargetable rendering system.
- NVIDIA Sionna: Hardware-accelerated differentiable ray tracer for wireless propagation, built on Mitsuba 3.
- Supports grouping faces into named objects and assigning materials by name.
- Optionally maps materials to ITU default presets used by Sionna RT.
- Allowed datatypes (
dtype): float or double
-
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);
-
Arguments:
const std::string fn (input)
Output file name (including path and .xml extension).
const arma::Mat<dtype> vert_list (input)
Vertex list, size [n_vert, 3], each row is a vertex (x, y, z) in Cartesian coordinates [m].
const arma::umat face_ind (input)
Face indices (0-based), size [n_mesh, 3], each row defines a triangle via vertex indices.
const arma::uvec obj_ind (input)
Object indices (1-based), size [n_mesh]. Assigns each triangle to an object.
const arma::uvec mtl_ind (input)
Material indices (1-based), size [n_mesh]. Assigns each triangle to a material.
const std::vector<std::string> obj_names (input)
Names of objects. Length must be equal to max(obj_ind).
const std::vector<std::string> mtl_names (input)
Names of materials. Length must be equal to max(mtl_ind).
const arma::Mat<dtype> bsdf = {} (optional input)
Material reflectivity data (BSDF parameters), size [mtl_names.size(), 17]. If omitted, the null BSDF is used.
Note that Sionna RT ignores all BSDF parameters. They are only used by the Mitsuma rendering system.
See obj_file_read for a definition of the data fields.
bool map_to_itu_materials = false (optional input)
If true, maps material names to ITU-defined presets used by Sionna RT. Default: false
-
See also:
obj_file_read - Read Wavefront `.obj` file and extract geometry and material information
-
Description:
- Parses a Wavefront
.obj file containing triangularized 3D geometry.
- Extracts triangle face data, material properties, vertex indices, and optional metadata such as object/material names.
- Multiple triangles belonging to the same object are grouped together by
obj_ind.
- Supports default and custom ITU-compliant materials encoded via the
usemtl tag.
- Automatically resizes output matrices/vectors as needed to match the file content.
- Returns the number of triangular mesh elements found in the file.
- Allowed datatypes (
dtype): float or double.
-
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);
-
Arguments:
std::string fn (input)
Path to the .obj file to be read.
arma::Mat<dtype> *mesh = nullptr (optional output)
Flattened triangle mesh data. Each row holds [X1, Y1, Z1, X2, Y2, Z2, X3, Y3, Z3]. Size: [n_mesh, 9].
arma::Mat<dtype> *mtl_prop = nullptr (optional output)
Material properties for each triangle. Size: [n_mesh, 5].
arma::Mat<dtype> *vert_list = nullptr (optional output)
List of all vertex positions found in the .obj file. Size: [n_vert, 3].
arma::umat *face_ind = nullptr (optional output)
Indices into vert_list for each triangle (0-based). Size: [n_mesh, 3].
arma::uvec *obj_ind = nullptr (optional output)
Object index (1-based) for each triangle. Size: [n_mesh].
arma::uvec *mtl_ind = nullptr (optional output)
Material index (1-based) for each triangle. Size: [n_mesh].
std::vector<std::string> *obj_names = nullptr (optional output)
Names of objects found in the file. Length: max(obj_ind).
std::vector<std::string> *mtl_names = nullptr (optional output)
Names of materials found in the file. Length: max(mtl_ind).
arma::Mat<dtype> *bsdf = nullptr (optional output)
Principled BSDF (Bidirectional Scattering Distribution Function) values extracted from the
.MTL file. Size [mtl_names.size(), 17]. Values are:
| 0 |
Base Color Red |
Range 0-1 |
Default = 0.8 |
| 1 |
Base Color Green |
Range 0-1 |
Default = 0.8 |
| 2 |
Base Color Blue |
Range 0-1 |
Default = 0.8 |
| 3 |
Transparency (alpha) |
Range 0-1 |
Default = 1.0 (fully opaque) |
| 4 |
Roughness |
Range 0-1 |
Default = 0.5 |
| 5 |
Metallic |
Range 0-1 |
Default = 0.0 |
| 6 |
Index of refraction (IOR) |
Range 0-4 |
Default = 1.45 |
| 7 |
Specular Adjustment to the IOR |
Range 0-1 |
Default = 0.5 (no adjustment) |
| 8 |
Emission Color Red |
Range 0-1 |
Default = 0.0 |
| 9 |
Emission Color Green |
Range 0-1 |
Default = 0.0 |
| 10 |
Emission Color Blue |
Range 0-1 |
Default = 0.0 |
| 11 |
Sheen |
Range 0-1 |
Default = 0.0 |
| 12 |
Clearcoat |
Range 0-1 |
Default = 0.0 |
| 13 |
Clearcoat roughness |
Range 0-1 |
Default = 0.0 |
| 14 |
Anisotropic |
Range 0-1 |
Default = 0.0 |
| 15 |
Anisotropic rotation |
Range 0-1 |
Default = 0.0 |
| 16 |
Transmission |
Range 0-1 |
Default = 0.0 |
-
Returns:
arma::uword
Number of mesh triangles found in the file (n_mesh).
-
Technical Notes:
- Unknown or missing materials default to
"vacuum" (ε_r = 1, σ = 0).
- Materials are applied per triangle via the
usemtl tag in the .obj file.
- Input geometry must be fully triangulated—quads and n-gons are not supported.
- File parsing is case-sensitive for material names.
-
Material Tag Format:
- Default materials (ITU-R P.2040-3 Table 3):
"usemtl itu_concrete", "itu_brick", "itu_wood", "itu_water", etc.
- Frequency range: 1–40 GHz (limited to 1–10 GHz for ground materials)
- Custom materials syntax:
"usemtl Name::A:B:C:D:att" with A, B: Real permittivity ε_r = A * fGHz^B,
C, D: Conductivity σ = C * fGHz^D, att: Penetration loss in dB (fixed, per interaction)
-
Material properties:
Each material is defined by its electrical properties. Radio waves that interact with a building will
produce losses that depend on the electrical properties of the building materials, the material
structure and the frequency of the radio wave. The fundamental quantities of interest are the electrical
permittivity (ϵ) and the conductivity (σ). A simple regression model for the frequency dependence is
obtained by fitting measured values of the permittivity and the conductivity at a number of frequencies.
The five parameters returned in mtl_prop then are:
- Real part of relative permittivity at f = 1 GHz (a)
- Frequency dependence of rel. permittivity (b) such that ϵ = a · f^b
- Conductivity at f = 1 GHz (c)
- Frequency dependence of conductivity (d) such that σ = c· f^d
- Fixed attenuation in dB applied to each transition
A more detailed explanation together with a derivation can be found in ITU-R P.2040. The following
list of material is currently supported and the material can be selected by using the usemtl tag
in the OBJ file. When using Blender, the simply assign a material with that name to an object or face.
The following materials are defined by default:
| Name |
a |
b |
c |
d |
Att |
max fGHz |
|
| vacuum / air |
1.0 |
0.0 |
0.0 |
0.0 |
0.0 |
100 |
|
| textiles |
1.5 |
0.0 |
5e-5 |
0.62 |
0.0 |
100 |
|
| plastic |
2.44 |
0.0 |
2.33e-5 |
1.0 |
0.0 |
100 |
|
| ceramic |
6.5 |
0.0 |
0.0023 |
1.32 |
0.0 |
100 |
|
| sea_water |
80.0 |
-0.25 |
4.0 |
0.58 |
0.0 |
100 |
|
| sea_ice |
3.2 |
-0.022 |
1.1 |
1.5 |
0.0 |
100 |
|
| water |
80.0 |
-0.18 |
0.6 |
1.52 |
0.0 |
20 |
|
| water_ice |
3.17 |
-0.005 |
5.6e-5 |
1.7 |
0.0 |
20 |
|
| itu_concrete |
5.24 |
0.0 |
0.0462 |
0.7822 |
0.0 |
100 |
|
| itu_brick |
3.91 |
0.0 |
0.0238 |
0.16 |
0.0 |
40 |
|
| itu_plasterboard |
2.73 |
0.0 |
0.0085 |
0.9395 |
0.0 |
100 |
|
| itu_wood |
1.99 |
0.0 |
0.0047 |
1.0718 |
0.0 |
100 |
|
| itu_glass |
6.31 |
0.0 |
0.0036 |
1.3394 |
0.0 |
100 |
|
| itu_ceiling_board |
1.48 |
0.0 |
0.0011 |
1.075 |
0.0 |
100 |
|
| itu_chipboard |
2.58 |
0.0 |
0.0217 |
0.78 |
0.0 |
100 |
|
| itu_plywood |
2.71 |
0.0 |
0.33 |
0.0 |
0.0 |
40 |
|
| itu_marble |
7.074 |
0.0 |
0.0055 |
0.9262 |
0.0 |
60 |
|
| itu_floorboard |
3.66 |
0.0 |
0.0044 |
1.3515 |
0.0 |
100 |
|
| itu_metal |
1.0 |
0.0 |
1.0e7 |
0.0 |
0.0 |
100 |
|
| itu_very_dry_ground |
3.0 |
0.0 |
0.00015 |
2.52 |
0.0 |
10 |
|
| itu_medium_dry_ground |
15.0 |
-0.1 |
0.035 |
1.63 |
0.0 |
10 |
|
| itu_wet_ground |
30.0 |
-0.4 |
0.15 |
1.3 |
0.0 |
10 |
|
| itu_vegetation |
1.0 |
0.0 |
1.0e-4 |
1.1 |
0.0 |
100 |
|
| irr_glass |
6.27 |
0.0 |
0.0043 |
1.1925 |
23.0 |
100 |
|
-
Example:
arma::mat mesh, mtl_prop, vert_list;
arma::umat face_ind;
arma::uvec obj_ind, mtl_ind;
std::vector<std::string> obj_names, mtl_names;
quadriga_lib::obj_file_read<double>("cube.obj", &mesh, &mtl_prop, &vert_list, &face_ind, &obj_ind, &mtl_ind, &obj_names, &mtl_names);
obj_overlap_test - Detect overlapping 3D objects in a triangular mesh
-
Description:
- Tests whether any objects in a triangular mesh overlap by checking for shared volume or intersection.
- Touching faces or edges are not considered overlapping
- Returns the indices (1-based) of all objects that intersect with at least one other object.
- Can optionally output a list of overlap reasons for diagnostic purposes.
- Uses a configurable geometric tolerance to account for numerical precision.
- Allowed datatypes (
dtype): float or double.
-
Declaration:
arma::uvec quadriga_lib::obj_overlap_test(
const arma::Mat<dtype> *mesh,
const arma::uvec *obj_ind,
std::vector[std::string](std::string) *reason = nullptr,
dtype tolerance = 0.0005);
-
Arguments:
const arma::Mat<dtype> *mesh (input)
Triangular mesh geometry. Each row contains 3 vertices flattened as [X1, Y1, Z1, X2, Y2, Z2, X3, Y3, Z3]. Size: [n_mesh, 9].
const arma::uvec *obj_ind (input)
Object indices (1-based) that map multiple triangles in mesh to objects; Size: [n_mesh];
This is an output generated by obj_file_read.
std::vector<std::string> *reason = nullptr (optional output)
Human-readable list of overlap reasons corresponding to each overlapping object. Length: [n_overlap].
dtype tolerance = 0.0005 (optional input)
Geometric tolerance (in meters) used to determine intersections. Default: 0.0005 (0.5 mm).
-
Returns:
arma::uvec: Vector of unique object indices (1-based) that were found to overlap, size [n_overlap].
-
Technical Notes:
- Overlap detection includes checks for: Intersecting triangle faces (shared volume), Vertices or edges penetrating another object’s bounding volume.
- The
tolerance accounts for modeling inaccuracies and numerical instability—small overlaps below this threshold are ignored.
- This function does not modify the mesh or attempt to repair overlapping geometry — it only reports it.
-
See also:
path_to_tube - Convert a 3D path into a tube surface for visualization
-
Description:
- Converts a sequence of 3D points (path) into a tubular surface using a ring of vertices around each path segment.
- Produces a quad-based mesh suitable for rendering in 3D tools such as Blender or MeshLab.
- Uses circular cross-sections with configurable radius and edge count.
- Internal calculations are performed in double precision to ensure numerical stability along complex paths.
- Allowed datatypes (
dtype): float or double.
-
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);
-
Arguments:
const arma::Mat<dtype> *path_coord (input)
Ordered list of 3D coordinates along the path, matrix of size [3, n_coord].
arma::Mat<dtype> *vert (output)
Generated tube vertex positions. Size: [3, n_coord × n_edges].
arma::umat *faces (output)
Quad face indices into vert. Each row contains 4 indices forming a quad. Size: [4, (n_coord - 1) × n_edges].
dtype radius = 1.0 (optional input)
Radius of the tube cross-section (in meters). Default: 1.0.
arma::uword n_edges = 5 (optional input)
Number of vertices used to approximate each circular cross-section. Must be ≥ 3. Default: 5.
-
Technical Notes:
- The generated tube is centered around the path with perpendicular circular cross-sections.
- Orientation between path segments is handled with continuous frame alignment to reduce twisting.
- Quad faces are generated by connecting adjacent rings along the path.
- Output
faces can be directly exported to formats like .obj or .ply.
point_cloud_aabb - Compute the Axis-Aligned Bounding Boxes (AABB) of a 3D point cloud
point_cloud_segmentation - Reorganize a point cloud into spatial sub-clouds for efficient processing
-
Description:
- Recursively partitions a 3D point cloud into smaller sub-clouds, each below a given size threshold.
- Sub-clouds are aligned to a specified SIMD vector size (e.g., for AVX or CUDA), with padding if necessary.
- Outputs (
pointsR) a reorganized version of the input points that groups points by sub-cloud.
- Also produces forward and reverse index maps to track the reordering of points.
- Useful for optimizing spatial processing tasks such as bounding volume hierarchies or GPU batch execution.
- Allowed datatypes (
dtype): float or double
-
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);
-
Arguments:
const arma::Mat<dtype> *points (input)
Original 3D point cloud to be segmented. Size: [n_points, 3].
arma::Mat<dtype> *pointsR (output)
Reorganized point cloud with points grouped by sub-cloud. Size: [n_pointsR, 3].
arma::u32_vec *sub_cloud_index (output)
Vector of starting indices (0-based) for each sub-cloud within pointsR. Length: [n_sub].
arma::uword target_size = 1024 (optional input)
Maximum number of elements allowed per sub-cloud (before padding). Default: 1024.
arma::uword vec_size = 1 (optional input)
Vector alignment size for SIMD or CUDA. The number of points in each sub-cloud is padded to a multiple of this value. Default: 1.
arma::u32_vec *forward_index = nullptr (optional output)
Index map from original points to reorganized pointsR (1-based). Size: [n_pointsR]. Padding indices are 0.
arma::u32_vec *reverse_index = nullptr (optional output)
Index map from pointsR back to original points (0-based). Size: [n_points].
-
Returns:
arma::uword
Number of generated sub-clouds, n_sub.
-
Technical Notes:
- Sub-clouds are formed using recursive spatial splitting (e.g., median-split along bounding box axes).
- Padding points are placed at the AABB center of the corresponding sub-cloud and can be ignored in processing.
- This function is typically used as a preprocessing step for GPU acceleration or bounding volume hierarchy (BVH) generation.
- If
vec_size = 1, no padding is applied and all output maps contain valid indices only.
-
See also:
point_cloud_split - Split a point cloud into two sub-clouds along a spatial axis
-
Description:
- Divides a 3D point cloud into two sub-clouds along the specified axis.
- Attempts to split the data at the median value to balance the number of points in each half.
- Returns the axis used for the split, or a negative value if a valid split was not possible (e.g., all points fall on one side).
- Output point clouds are written into
pointsA and pointsB, and their size is adjusted accordingly.
- An optional indicator vector identifies the target sub-cloud (A or B) for each input point.
- Used in recursive spatial partitioning such as building BVH structures.
- Allowed datatypes (
dtype): float or double
-
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);
-
Arguments:
const arma::Mat<dtype> *points (input)
Input point cloud. Size: [n_points, 3].
arma::Mat<dtype> *pointsA (output)
First sub-cloud after split. Size: [n_pointsA, 3].
arma::Mat<dtype> *pointsB (output)
Second sub-cloud after split. Size: [n_pointsB, 3].
int axis = 0 (optional input)
Axis to split along: 0 = longest extent (default), 1 = x-axis, 2 = y-axis, 3 = z-axis.
arma::Col<int> *split_ind = nullptr (optional output)
Vector of length [n_points], where each element is: 1 if the point goes to pointsA, 2 if it goes to pointsB, 0 if error.
-
Returns:
int
Axis used for splitting: 1 = x, 2 = y, 3 = z, or -1, -2, -3 if the split failed (no points assigned to one of the outputs).
-
Notes:
- The function does not modify
pointsA or pointsB if the split fails.
- The selected axis is based on the bounding box if
axis == 0
- This function is a building block for spatial acceleration structures (e.g., BVH, KD-trees), see point_cloud_segmentation
-
See also:
point_inside_mesh - Test whether 3D points are inside a triangle mesh using raycasting
-
Description:
- Uses raycasting to determine whether each 3D point lies inside a triangle mesh.
- Requires that the mesh is watertight and all normals are pointing outwards.
- For each point, multiple rays are cast in various directions.
- If any ray intersects a mesh element with a negative incidence angle, the point is classified as inside.
- Output can be binary (0 = outside, 1 = inside) or labeled with object indices.
- Allowed datatypes (
dtype): float or double.
-
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);
-
Arguments:
const arma::Mat<dtype> *points (input)
3D point coordinates to test, size [n_points, 3].
const arma::Mat<dtype> *mesh (input)
Triangular mesh faces. Each row represents a triangle using 3 vertices in row-major format (x1,y1,z1,x2,y2,z2,x3,y3,z3), size [n_mesh, 9].
const arma::u32_vec *obj_ind = nullptr (optional input)
Optional object index for each mesh element (1-based), size [n_mesh]. If provided, the return vector will contain the index of the enclosing object instead of binary values.
dtype distance = 0.0 (optional input)
Optional distance in meters from objects that should be considered as inside the object.
Possible range: 0 - 20 m. Using this parameter significantly increases computation time.
-
Returns:
arma::u32_vec, size [n_points]
For each point: Returns 0 if the point is outside the mesh (or all objects), 1 if inside (or close to) any mesh object
(if obj_ind not given), or returns the 1-based object index if obj_ind is provided.
ray_mesh_interact - Calculates interactions (reflection, transmission, refraction) of radio waves with objects
-
Description:
- Radio waves that interact with a building, or other objects in the environment, will produce losses
that depend on the electrical properties of the materials and material structure.
- This function calculates the interactions of radio-waves with objects in a 3D environment.
- It considers a plane wave incident upon a planar interface between two homogeneous and isotropic
media of differing electric properties. The media extend sufficiently far from the interface such
that the effect of any other interface is negligible.
- Supports beam-based modeling using triangular ray tubes defined by
trivec and tridir.
- Air to media transition is assumed if front side of a face is hit and FBS != SBS
- Media to air transition is assumed if back side of a face is hit and FBS != SBS
- Media to media transition is assumed if FBS = SBS with opposing face orientations
- Order of the vertices determines side (front or back) of a mesh element
- Overlapping geometry in the triangle mesh must be avoided, since materials are transparent to radio
waves.
- Implementation is done according to ITU-R P.2040-1.
- Rays that do not interact with the environment (i.e. for which
fbs_ind = 0) are omitted from the output.
- Maintains consistency of input and output ray formats, including direction encoding (spherical or Cartesian).
- Assigns interaction type codes to output rays, allowing classification (e.g., single hit, total
reflection, material-to-material).
- Allowed datatypes (
dtype): float or double.
-
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);
-
Arguments:
int interaction_type (input)
Type of interaction: 0 = reflection, 1 = transmission, 2 = refraction.
dtype center_frequency (input)
Center frequency in Hz.
const arma::Mat<dtype> orig, dest (input)
Ray origin and destination points in global Cartesian space. Size: [n_ray, 3].
const arma::Mat<dtype> fbs, sbs (input)
First and second interaction points. Size: [n_ray, 3].
const arma::Mat<dtype> *mesh (input)
Triangle mesh faces. Size: [n_mesh, 9], See obj_file_read for details
const arma::Mat<dtype> *mtl_prop (input)
Material properties per face. Size: [n_mesh, 5], See obj_file_read for details
const arma::u32_vec fbs_ind, sbs_ind (input)
Mesh indices of interaction points (1-based). Size: [n_ray].
const arma::Mat<dtype> *trivec (optional input)
The 3 vectors pointing from the center point of the ray at the origin to the vertices of a triangular
propagation tube (=beam wavefront), the values are in the order [ v1x, v1y, v1z, v2x, v2y, v2z, v3x, v3y, v3z ];
Size: [no_ray, 9]
const arma::Mat<dtype> *tridir (optional input)
The directions of the vertex-rays. Size: [ n_ray, 6 ] or [ n_ray, 9 ],
For 6 columns, values are in geographic coordinates (azimuth and elevation angle in rad); the
values are in the order [ v1az, v1el, v2az, v2el, v3az, v3el ];
For 9 columns, the input is in Cartesian coordinates and the values are in the order
[ v1x, v1y, v1z, v2x, v2y, v2z, v3x, v3y, v3z ]
const arma::Col<dtype> *orig_length (optional input)
Path length at origin. Size: [n_ray]. Default is 0
arma::Mat<dtype> *origN (output)
New ray origins after the interaction with the medium, usually placed close to the FBS location.
A small offset of 0.001 m in the direction of travel after the interaction with the medium is added
to avoid getting stuct inside a mesh element. Size: [ no_rayN, 3 ]
arma::Mat<dtype> *destN (output)
New ray destinaion after the interaction with the medium, taking the change of direction into account;
Size: [no_rayN, 3]
arma::Col<dtype> *gainN (output)
Gain (negative loss) caused by the interaction with the medium, averaged over both polarization
directions. This value includes the in-medium attenuation, but does not account for FSPL.
Values are in linear scale (not dB); Size: [no_rayN]
arma::Mat<dtype> *xprmatN (output)
Polarization transfer matrix; Interleaved complex values (ReVV, ImVV, ReVH, ImVH, ReHV, ImHV, ReHH, ImHH);
The values account for the following effects: (1) gain caused by the interaction with the medium,
(2) different reflection/transmission coefficients for transverse electric (TE) and transverse
magnetic (TM) polarisation, (3) orientation of the incidence plane, (4) in-medium attenuation.
FSPL is excluded, Size: [ no_rayN, 8 ]
arma::Mat<dtype> trivecN, tridirN (output)
New ray tube geometry and direction. Format matches input.
arma::Col<dtype> *orig_lengthN (output)
Length of the ray from orig to origN. If orig_length is given as input, its value is added.
Size: [ no_rayN ]
arma::Col<dtype> *fbs_angleN (output)
Angle between incoming ray and FBS in [rad], Size [ n_rayN ]
arma::Col<dtype> *thicknessN (output)
Material thickness in meters calculated from the difference between FBS and SBS, Size [ n_rayN ]
arma::Col<dtype> *edge_lengthN (output)
Max. edge length of the ray tube triangle at the new origin. A value of infinity indicates that only
a part of the ray tube hits the object. Size [ n_rayN, 3 ]
arma::Mat<dtype> *normal_vecN (output)
Normal vector of FBS and SBS, Size [ n_rayN, 6 ]
arma::s32_vec *out_typeN (output)
A numeric indicator describing the type of the interaction. The total refection indicator is only
set in refraction mode. Size [ n_rayN ]
| No |
θF<0 |
θS<0 |
dFS=0 |
TotRef |
iSBS=0 |
NF=-NS |
NF=NS |
startIn |
endIn |
Meaning |
| 0 |
|
|
|
|
|
|
|
|
|
Undefined |
| 1 |
no |
N/A |
no |
N/A |
yes |
N/A |
N/A |
no |
yes |
Single Hit o-i |
| 2 |
yes |
N/A |
no |
no |
yes |
N/A |
N/A |
yes |
no |
Single Hit i-o |
| 3 |
yes |
N/A |
no |
yes |
yes |
N/A |
N/A |
yes |
no |
Single Hit i-o, TR |
| --- |
----- |
------ |
------- |
-------- |
-------- |
-------- |
------- |
--------- |
------- |
---------------------------- |
| 4 |
no |
yes |
yes |
no |
no |
yes |
no |
yes |
yes |
M2M, M2 hit first |
| 5 |
yes |
no |
yes |
no |
no |
yes |
no |
yes |
yes |
M2M, M1 hit first |
| 6 |
yes |
no |
yes |
yes |
no |
yes |
no |
yes |
yes |
M2M, M1 hit first, TR |
| --- |
----- |
------ |
------- |
-------- |
-------- |
-------- |
------- |
--------- |
------- |
---------------------------- |
| 7 |
no |
no |
yes |
N/A |
no |
no |
yes |
no |
yes |
Overlapping Faces, o-i |
| 8 |
yes |
yes |
yes |
no |
no |
no |
yes |
yes |
no |
Overlapping Faces, i-o |
| 9 |
yes |
yes |
yes |
yes |
no |
no |
yes |
yes |
no |
Overlapping Faces, i-o, TR |
| --- |
----- |
------ |
------- |
-------- |
-------- |
-------- |
------- |
--------- |
------- |
---------------------------- |
| 10 |
no |
yes |
yes |
N/A |
no |
no |
no |
no |
no |
Edge Hit, o-i-o |
| 11 |
yes |
no |
yes |
no |
no |
no |
no |
yes |
yes |
Edge Hit, i-o-i |
| 12 |
yes |
no |
yes |
yes |
no |
no |
no |
yes |
yes |
Edge Hit, i-o-i, TR |
| 13 |
no |
no |
yes |
N/A |
no |
no |
no |
no |
yes |
Edge Hit, o-i |
| 14 |
yes |
yes |
yes |
no |
no |
no |
no |
yes |
no |
Edge Hit, i-o |
| 15 |
yes |
yes |
yes |
yes |
no |
no |
no |
yes |
no |
Edge Hit, i-o, TR |
-
See also:
ray_point_intersect - Calculates the intersection of ray beams with points in three dimensions
-
Description:
Unlike traditional ray tracing, where rays do not have a physical size, beam tracing models rays as
beams with volume. Beams are defined by triangles whose vertices diverge as the beam extends. This
approach is used to simulate a kind of divergence or spread in the beam, reminiscent of how radio
waves spreads as they travel from a point source. The volumetric nature of the beams allows for more
realistic modeling of energy distribution. As beams widen, the energy they carry can be distributed
across their cross-sectional area, affecting the intensity of the interaction with surfaces.
Unlike traditional ray tracing where intersections are line-to-geometry tests, beam tracing requires
volumetric intersection tests.
- This function computes whether points in 3D Cartesian space are intersected by ray beams.
- A ray beam is defined by a ray origin and a triangular wavefront formed by three directional vectors.
- Returns a list of ray indices (0-based) that intersect with each point in the input point cloud.
- Optional support for pre-segmented point clouds (e.g., using point_cloud_segmentation)
to reduce computational cost.
- All internal computations use single precision for speed.
- Uses Advanced Vector Extensions (AVX2) for efficient computation, if supported by the CPU
- Recommended to use with small tube radius and well-distributed points for optimal accuracy.
-
Declaration:
template <typename dtype>
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);
-
Arguments:
const arma::Mat<dtype> *points (input)
3D coordinates of the point cloud. Size: [n_points, 3].
const arma::Mat<dtype> *orig (input)
Ray origin positions in global coordinate system. Size: [n_ray, 3].
const arma::Mat<dtype> *trivec (input)
The 3 vectors pointing from the center point of the ray at the ray origin to the vertices of
a triangular propagation tube (the beam), the values are in the order
[ v1x, v1y, v1z, v2x, v2y, v2z, v3x, v3y, v3z ]; Size: [ no_ray, 9 ]
const arma::Mat<dtype> *tridir (input)
The directions of the vertex-rays. Size: [ n_ray, 9 ], Values must be given in Cartesian
coordinates in the order [ d1x, d1y, d1z, d2x, d2y, d2z, d3x, d3y, d3z ]; The vector does
not need to be normalized.
const arma::u32_vec *sub_cloud_index = nullptr (optional input)
Index vector to mark boundaries between point cloud segments (e.g. for SIMD optimization). Length: [n_sub].
arma::u32_vec *hit_count = nullptr (optional output)
Output array with number of rays that intersected each point. Length: [n_points].
-
Returns:
std::vector<arma::u32_vec>
List of ray indices that intersected each point (0-based). Each entry in the returned vector corresponds to one point.
-
See also:
ray_triangle_intersect - Calculates the intersection of rays and triangles in three dimensions
-
Description:
- Implements the Möller–Trumbore algorithm to compute intersections between rays and triangles in 3D.
- Efficient SIMD acceleration using AVX2 intrinsics to process 8 mesh triangles in parallel.
- Can detect first and second intersections (FBS/SBS), number of intersections, and intersection indices.
- Supports a compact input format where origin and destination coordinates are stored together.
- Allowed datatypes (
dtype): float or double
- Internal computations are carried out in single precision regardless of
dtype.
-
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,
bool transpose_inputs = false);
-
Arguments:
const arma::Mat<dtype> *orig (input)
Ray origins in global coordinate system (GCS). Size: [n_ray, 3], or [n_ray, 6] if dest == nullptr.
In the latter case, columns must be ordered as {xo, yo, zo, xd, yd, zd}.
const arma::Mat<dtype> *dest (input)
Ray destinations in GCS. Size: [n_ray, 3]. Set to nullptr if orig has 6 columns and contains
both origin and destination.
const arma::Mat<dtype> *mesh (input)
Triangular surface mesh. Size: [n_mesh, 9], where each row contains the 3 vertices
{x1 y1 z1 x2 y2 z2 x3 y3 z3}.
arma::Mat<dtype> *fbs (optional output)
First-bounce surface intersection points (FBS). Size: [n_ray, 3].
arma::Mat<dtype> *sbs (optional output)
Second-bounce surface intersection points (SBS). Size: [n_ray, 3].
arma::u32_vec *no_interact (optional output)
Number of intersections per ray (0, 1, or 2). Size: [n_ray].
arma::u32_vec *fbs_ind (optional output)
1-based index of the first intersected mesh element, 0 = no intersection. Size: [n_ray].
arma::u32_vec *sbs_ind (optional output)
1-based index of the second intersected mesh element, 0 = no second intersection. Size: [n_ray].
const arma::u32_vec *sub_mesh_index (optional input)
Indexes indicating start of sub-meshes in mesh. Size: [n_sub]. Enables faster processing via segmentation.
bool transpose_inputs = false (optional input)
If true, treats orig/dest as [3, n_ray] and mesh as [9, n_mesh].
-
See also:
subdivide_rays - Subdivide ray beams into four smaller sub-beams
-
Description:
- Subdivides each ray beam (defined by a triangular wavefront) into four new beams with adjusted origin, shape, and direction.
- Supports input in Spherical or Cartesian direction format.
- When
dest is not provided, the corresponding output destN is omitted.
- Useful for hierarchical ray tracing or angular resolution refinement.
- Allowed datatypes (
dtype): float or double.
-
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);
-
Arguments:
const arma::Mat<dtype> *orig (input)
Ray origin points in global coordinate system (GCS).
Size: [n_ray, 3].
const arma::Mat<dtype> *trivec (input)
Vectors pointing from the ray origin to the three triangle vertices.
Size: [n_ray, 9], order: [x1 y1 z1 x2 y2 z2 x3 y3 z3].
const arma::Mat<dtype> *tridir (input)
Directions of the three vertex-rays.
Format can be Spherical [n_ray, 6] as [v1az v1el v2az v2el v3az v3el],
or Cartesian [n_ray, 9] as [v1x v1y v1z v2x v2y v2z v3x v3y v3z].
const arma::Mat<dtype> *dest = nullptr (input)
Ray destination points. If nullptr, the output destN will remain empty.
Size: [n_ray, 3].
arma::Mat<dtype> *origN (output)
New ray origins after subdivision.
Size: [n_rayN, 3].
arma::Mat<dtype> *trivecN (output)
Updated vectors for each subdivided triangle beam.
Size: [n_rayN, 9].
arma::Mat<dtype> *tridirN (output)
New directions of the subdivided vertex-rays, in the same format as input.
Size: [n_rayN, 6] (spherical) or [n_rayN, 9] (Cartesian).
arma::Mat<dtype> *destN (output)
Updated destination points.
Size: [n_rayN, 3], empty if input dest was nullptr.
const arma::u32_vec *index (optional input)
List of ray indices to be subdivided (0-based). Only the specified rays are subdivided.
Size: [n_ind].
const double ray_offset = 0.0 (optional input)
Offset (in meters) applied to the origin of each subdivided ray along its propagation direction.
Default: 0.0.
-
Returns:
arma::uword n_rayN
Number of output rays, typically 4 × n_ray or 4 × n_ind if index is provided.
-
See also:
subdivide_triangles - Subdivide triangles into smaller triangles
-
Description:
- Uniformly subdivides each input triangle into
n_div × n_div smaller triangles.
- Increases spatial resolution for mesh-based processing (e.g., ray tracing or visualization).
- Optional input/output material properties are duplicated across subdivided triangles.
- Allowed datatypes (
dtype): float or double.
-
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);
-
Arguments:
arma::uword n_div (input)
Number of subdivisions per triangle edge;
total output triangles: n_triangles_out = n_triangles_in × n_div × n_div.
const arma::Mat<dtype> *triangles_in (input)
Vertices of the triangular mesh in global Cartesian coordinates; each face is described by 3
points in 3D-space: [ v1x, v1y, v1z, v2x, v2y, v2z, v3x, v3y, v3z ]; Size: [n_triangles_in, 9]
arma::Mat<dtype> *triangles_out (output)
Vertices of the sub-divided mesh in global Cartesian coordinates; Size: [n_triangles_out, 9]
const arma::Mat<dtype> *mtl_prop = nullptr (optional input)
Material properties associated for the input triangles; Size: [n_triangles_in, 5].
arma::Mat<dtype> *mtl_prop_out = nullptr (optional output)
Material properties for the subdivided triangles, copied from the parent triangle,
Size: [n_triangles_out, 5].
-
Returns:
arma::uword n_triangles_out
Number of generated triangles (equals n_triangles_in × n_div × n_div).
triangle_mesh_aabb - Calculate the axis-aligned bounding box (AABB) of a triangle mesh and its sub-meshes
triangle_mesh_segmentation - Reorganize a 3D mesh into smaller sub-meshes for faster processing
-
Description:
This function processes the elements of a large triangle mesh by clustering those that are
closely spaced. The resulting mesh retains the same elements but rearranges their order.
The function aims to minimize the size of the axis-aligned bounding box around each cluster,
referred to as a sub-mesh, while striving to maintain a specific number of elements within
each cluster.
- Subdivision is recursive and based on bounding box partitioning until each sub-mesh contains no more than
target_size triangles.
- Sub-meshes are aligned to
vec_size for SIMD or GPU optimization; padded with dummy triangles at the center of each sub-mesh if needed.
- If material properties are provided, these are also reorganized and padded accordingly.
- The function returns the number of created sub-meshes, and reorders the triangles and materials.
- Allowed datatypes (
dtype): float or double.
-
Declaration:
arma::uword quadriga_lib::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);
-
Arguments:
const arma::Mat<dtype> *mesh (input)
Vertices of the triangular mesh in global Cartesian coordinates. Each face is described by 3
points in 3D-space: [ v1x, v1y, v1z, v2x, v2y, v2z, v3x, v3y, v3z ]; Size: [n_mesh, 9]
arma::Mat<dtype> *meshR (output)
Vertices of the clustered mesh in global Cartesian coordinates; Size: [n_meshR, 9]
arma::u32_vec *sub_mesh_index (output)
Start indices of the sub-meshes in 0-based notation; Vector of length [n_sub]
arma::uword target_size = 1024 (input)
The target number of elements of each sub-mesh. Default value = 1024. For best performance, the
value should be around sgrt(n_mesh)
arma::uword vec_size = 1 (input)
Vector size for SIMD processing (e.g. 8 for AVX2, 32 for CUDA). Default value = 1.
For values > 1,the number of rows for each sub-mesh in the output is increased to a multiple
of vec_size. For padding, zero-sized triangles are placed at the center of the AABB of
the corresponding sub-mesh.
const arma::Mat<dtype> *mtl_prop = nullptr (optional input)
Material properties corresponding to the original mesh; Size: [n_mesh, 5]; See obj_file_read
arma::Mat<dtype> *mtl_propR = nullptr (optional output)
Reorganized material properties, aligned and padded if necessary, Size: [n_meshR, 5]
arma::u32_vec *mesh_index = nullptr (optional output)
1-based mapping from the original mesh to the reorganized mesh; Size: [n_meshR] (0 = padding)
-
Returns:
arma::uword n_sub
The number of created sub-meshes. Output matrices are resized accordingly.
triangle_mesh_split - Split a 3D mesh into two sub-meshes along a given axis
-
Description:
- Divides a triangular mesh into two sub-meshes along a selected axis (or automatically the longest).
- The function chooses a split point based on the bounding box center of the selected axis.
- Returns the axis used for the split:
1 = x, 2 = y, 3 = z; or negative values if the split failed.
- An optional indicator vector identifies the target sub-mesh (A or B) for each input point.
- On failure (i.e., all triangles fall into one side), outputs
meshA and meshB remain unchanged.
- Allowed datatypes (
dtype): float or double
-
Declaration:
int quadriga_lib::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);
-
Arguments:
const arma::Mat<dtype> *mesh (input)
Triangle mesh input; each row contains one triangle as [x1 y1 z1 x2 y2 z2 x3 y3 z3]
Size: [n_mesh, 9]
arma::Mat<dtype> *meshA (output)
First resulting sub-mesh; triangles with centroid below split threshold. Size: [n_meshA, 9]
arma::Mat<dtype> *meshB (output)
Second resulting sub-mesh; triangles with centroid above split threshold. Size: [n_meshB, 9]
int axis = 0 (optional input)
Axis to split along: 0 = longest extent (default), 1 = x-axis, 2 = y-axis, 3 = z-axis.
arma::Col<int> *split_ind = nullptr (optional output)
Output vector indicating assignment of each triangle: 1 = meshA, 2 = meshB, 0 = not assigned (on failure)
Length: [n_mesh]
-
Returns:
int
The axis used for the split (1, 2, or 3), or negative value on failure (-1, -2, or -3).