MALAB / Octave API Documentation for Quadriga-Lib v0.10.1
Overview
Array antenna functions
Channel functions
Channel generation functions
Miscellaneous / Tools
| calc_rotation_matrix | Calculates a 3x3 rotation matrix from a 3-element orientation vector |
| cart2geo | Transform Cartesian (x,y,z) coordinates to Geographic (az, el, length) coordinates |
| fast_sincos | Fast, approximate sine/cosine for MATLAB numeric arrays |
| geo2cart | Transform Geographic (az, el, length) to Cartesian (x,y,z) coordinates coordinates |
| version | Returns the quadriga-lib version number |
| interp | 2D and 1D linear interpolation. |
| write_png | Write data to a PNG file |
Site-Specific Simulation Tools
Array antenna functions
arrayant_calc_directivity - Calculates the directivity (in dBi) of array antenna elements
-
Description:
Directivity is a parameter of an antenna or which measures the degree to which the radiation emitted
is concentrated in a single direction. It is the ratio of the radiation intensity in a given direction
from the antenna to the radiation intensity averaged over all directions. Therefore, the directivity
of a hypothetical isotropic radiator is 1, or 0 dBi.
-
Usage:
% Input as struct
directivity = quadriga_lib.arrayant_calc_directivity(arrayant);
directivity = quadriga_lib.arrayant_calc_directivity(arrayant, i_element);
% Separate inputs
directivity = quadriga_lib.arrayant_calc_directivity(e_theta_re, e_theta_im, e_phi_re, ...
e_phi_im, azimuth_grid, elevation_grid);
directivity = quadriga_lib.arrayant_calc_directivity(e_theta_re, e_theta_im, e_phi_re, ...
e_phi_im, azimuth_grid, elevation_grid, i_element);
-
Examples:
% Generate dipole antenna
ant = quadriga_lib.arrayant_generate('dipole');
% Calculate directivity
directivity = quadriga_lib.arrayant_calc_directivity(ant);
-
Input arguments for struct mode:
arrayant [1]
Struct containing a array antenna pattern with at least the following fields:
e_theta_re |
Real part of e-theta field component |
Size: [n_elevation, n_azimuth, n_elements] |
e_theta_im |
Imaginary part of e-theta field component |
Size: [n_elevation, n_azimuth, n_elements] |
e_phi_re |
Real part of e-phi field component |
Size: [n_elevation, n_azimuth, n_elements] |
e_phi_im |
Imaginary part of e-phi field component |
Size: [n_elevation, n_azimuth, n_elements] |
azimuth_grid |
Azimuth angles in [rad] -pi to pi, sorted |
Size: [n_azimuth] |
elevation_grid |
Elevation angles in [rad], -pi/2 to pi/2, sorted |
Size: [n_elevation] |
i_element [2] (optional)
Element index, 1-based. If not provided or empty, the directivity is calculated for all elements in the
array antenna.
Size: [n_out] or empty
-
Input arguments for separate inputs:
e_theta_re [1]
Real part of e-theta field component, Size: [n_elevation, n_azimuth, n_elements]
e_theta_im [2]
Imaginary part of e-theta field component, Size: [n_elevation, n_azimuth, n_elements]
e_phi_re [3]
Real part of e-phi field component, Size: [n_elevation, n_azimuth, n_elements]
e_phi_im [4]
Imaginary part of e-phi field component, Size: [n_elevation, n_azimuth, n_elements]
azimuth_grid [5]
Azimuth angles in [rad] -pi to pi, sorted, Size: [n_azimuth]
elevation_grid [6]
Elevation angles in [rad], -pi/2 to pi/2, sorted, Size: [n_elevation]
i_element [7] (optional)
Element index, 1-based. If not provided or empty, the directivity is calculated for all elements in the
array antenna.
Size: [n_out] or empty
-
Output Argument:
directivity
Directivity of the antenna pattern in dBi, double precision,
Size: [n_out] or [n_elements]
arrayant_combine_pattern - Calculate effective radiation patterns for array antennas
-
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. By integrating the element radiation patterns, their 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. The function arrayant_combine_pattern is designed to compute these effective
radiation patterns.
-
Usage:
% Minimal example (input/output = struct)
arrayant_out = quadriga_lib.arrayant_combine_pattern(arrayant_in);
% Optional inputs: freq, azimuth_grid, elevation_grid
arrayant_out = quadriga_lib.arrayant_combine_pattern(arrayant_in, freq, azimuth_grid, elevation_grid);
% Separate outputs, struct input
[e_theta_re, e_theta_im, e_phi_re, e_phi_im, azimuth_grid, elevation_grid, element_pos, ...
coupling_re, coupling_im, freq, name] = quadriga_lib.arrayant_combine_pattern(arrayant_in);
% Separate inputs
arrayant_out = quadriga_lib.arrayant_combine_pattern([], freq, azimuth_grid, elevation_grid, ...
e_theta_re, e_theta_im, e_phi_re, e_phi_im, azimuth_grid, elevation_grid, element_pos, ...
coupling_re, coupling_im, freq, name);
-
Examples:
The following example creates a unified linear array of 4 dipoles, spaced at half-wavelength. The
elements are then coupled with each other (i.e., they receive the same signal). The effective pattern
is calculated using arrayant_combine_pattern.
% Generate dipole pattern
ant = quadriga_lib.arrayant_generate('dipole');
% Duplicate 4 times
ant.e_theta_re = repmat(ant.e_theta_re, [1,1,4]);
ant.e_theta_im = repmat(ant.e_theta_im, [1,1,4]);
ant.e_phi_re = repmat(ant.e_phi_re, [1,1,4]);
ant.e_phi_im = repmat(ant.e_phi_im, [1,1,4]);
ant.element_pos = repmat(ant.element_pos, [1,4]);
% Set element positions and coupling matrix
ant.element_pos(2,:) = [ -0.75, -0.25, 0.25, 0.75]; % lambda, along y-axis
ant.coupling_re = [ 1 ; 1 ; 1 ; 1 ]/sqrt(4);
ant.coupling_im = [ 0 ; 0 ; 0 ; 0 ];
% Calculate effective pattern
ant_c = quadriga_lib.arrayant_combine_pattern( ant );
% Plot gain
plot( ant.azimuth_grid*180/pi, [ 10*log10( ant.e_theta_re(91,:,1).^2 ); 10*log10( ant_c.e_theta_re(91,:).^2 ) ]);
axis([-180 180 -20 15]); ylabel('Gain (dBi)'); xlabel('Azimth angle (deg)'); legend('Dipole','Array')
-
Input Arguments:
arrayant_in [1]
Struct containing the arrayant data with the following fields:
e_theta_re |
e-theta field component, real part |
Size: [n_elevation, n_azimuth, n_elements] |
e_theta_im |
e-theta field component, imaginary part |
Size: [n_elevation, n_azimuth, n_elements] |
e_phi_re |
e-phi field component, real part |
Size: [n_elevation, n_azimuth, n_elements] |
e_phi_im |
e-phi field component, imaginary part |
Size: [n_elevation, n_azimuth, n_elements] |
azimuth_grid |
Azimuth angles in [rad], -pi to pi, sorted |
Size: [n_azimuth] |
elevation_grid |
Elevation angles in [rad], -pi/2 to pi/2, sorted |
Size: [n_elevation] |
element_pos |
Antenna element (x,y,z) positions, optional |
Size: [3, n_elements] |
coupling_re |
Coupling matrix, real part, optional |
Size: [n_elements, n_ports] |
coupling_im |
Coupling matrix, imaginary part, optional |
Size: [n_elements, n_ports] |
center_freq |
Center frequency in [Hz], optional |
Scalar |
name |
Name of the array antenna object, optional |
String |
If an empty array is passed, array antenna data is provided as separate inputs (Inputs 5-15)
freq [2] (optional)
An alternative value for the center frequency. Overwrites the value given in arrayant_in. If
neither freq not arrayant_in["center_freq"] are given, an error is thrown.
azimuth_grid [3] (optional)
Alternative azimuth angles for the output in [rad], -pi to pi, sorted, Size: [n_azimuth_out],
If not given, arrayant_in.azimuth_grid is used instead.
elevation_grid [4] (optional)
Alternative elevation angles for the output in [rad], -pi/2 to pi/2, sorted, Size: [n_elevation_out],
If not given, arrayant_in.elevation_grid is used instead.
-
Output Arguments:
arrayant_out
Struct containing the arrayant data with the following fields:
e_theta_re |
e-theta field component, real part |
Size: [n_elevation, n_azimuth, n_elements] |
e_theta_im |
e-theta field component, imaginary part |
Size: [n_elevation, n_azimuth, n_elements] |
e_phi_re |
e-phi field component, real part |
Size: [n_elevation, n_azimuth, n_elements] |
e_phi_im |
e-phi field component, imaginary part |
Size: [n_elevation, n_azimuth, n_elements] |
azimuth_grid |
Azimuth angles in [rad] -pi to pi, sorted |
Size: [n_azimuth] |
elevation_grid |
Elevation angles in [rad], -pi/2 to pi/2, sorted |
Size: [n_elevation] |
element_pos |
Antenna element (x,y,z) positions |
Size: [3, n_elements] |
coupling_re |
Coupling matrix, real part |
Size: [n_elements, n_ports] |
coupling_im |
Coupling matrix, imaginary part |
Size: [n_elements, n_ports] |
center_freq |
Center frequency in [Hz], default = 0.3 GHz |
Scalar |
name |
Name of the array antenna object |
String |
Can be returned as separate outputs.
arrayant_copy_element - Create copies of array antenna elements
arrayant_export_obj_file - Creates a Wavefront OBJ file for visualizing the shape of the antenna pattern
-
Usage:
quadriga_lib.arrayant_export_obj_file( fn, arrayant, directivity_range, colormap, ...
object_radius, icosphere_n_div, i_element );
-
Input Arguments:
fn [1]
Filename of the OBJ file, string
arrayant [2]
Struct containing the arrayant data with the following fields:
e_theta_re |
e-theta field component, real part |
Size: [n_elevation, n_azimuth, n_elements] |
e_theta_im |
e-theta field component, imaginary part |
Size: [n_elevation, n_azimuth, n_elements] |
e_phi_re |
e-phi field component, real part |
Size: [n_elevation, n_azimuth, n_elements] |
e_phi_im |
e-phi field component, imaginary part |
Size: [n_elevation, n_azimuth, n_elements] |
azimuth_grid |
Azimuth angles in [rad], -pi to pi, sorted |
Size: [n_azimuth] |
elevation_grid |
Elevation angles in [rad], -pi/2 to pi/2, sorted |
Size: [n_elevation] |
element_pos |
Antenna element (x,y,z) positions, optional |
Size: [3, n_elements] |
name |
Name of the array antenna object, optional |
String |
directivity_range [3]
Directivity range of the antenna pattern visualization in dB
colormap [4]
Colormap for the visualization, string, supported are 'jet', 'parula', 'winter', 'hot', 'turbo',
'copper', 'spring', 'cool', 'gray', 'autumn', 'summer', Optional, default = 'jet'
object_radius [5]
Radius in meters of the exported object
icosphere_n_div [6]
Map pattern to an Icosphere with given number of subdivisions
element [7]
Antenna element indices, 1-based, empty = export all
arrayant_generate - Generates predefined array antenna models
-
Description:
This functions can be used to generate a variety of pre-defined array antenna models, including 3GPP
array antennas used for 5G-NR simulations. The first argument is the array type. The following input
arguments are then specific to this type.
-
Usage:
% Simple antenna models, output as struct
ant = quadriga_lib.arrayant_generate('omni', res); % Isotropic radiator, v-pol
ant = quadriga_lib.arrayant_generate('dipole', res); % Short dipole, v-pol
ant = quadriga_lib.arrayant_generate('half-wave-dipole', res); % Half-wave dipole, v-pol
ant = quadriga_lib.arrayant_generate('xpol', res); % Cross-polarized isotropic radiator
% An antenna with a custom 3dB beam with (in degree)
ant = quadriga_lib.arrayant_generate('custom', res, freq, az_3dB, el_3db, rear_gain_lin);
% A unified linear array with isotropic patterns
ant = quadriga_lib.arrayant_generate('ula', res, freq, [], [], [], 1, N, [], [], spacing);
% Antenna model for the 3GPP-NR channel model with 3GPP default pattern
ant = quadriga_lib.arrayant_generate('3GPP', res, freq, [], [], [],
M, N, pol, tilt, spacing, Mg, Ng, dgv, dgh);
% Antenna model for the 3GPP-NR channel model with a custom beam width
ant = quadriga_lib.arrayant_generate('3GPP', res, freq, az_3dB, el_3db, rear_gain_lin,
M, N, pol, tilt, spacing, Mg, Ng, dgv, dgh);
% Antenna model for the 3GPP-NR channel model with a custom antenna pattern
ant = quadriga_lib.arrayant_generate('3GPP', res, freq, [], [], [],
M, N, pol, tilt, spacing, Mg, Ng, dgv, dgh, pattern);
% Multi-beam antenna (single beam serving multiple directions)
ant = quadriga_lib.arrayant_generate('multibeam', res, freq, az_3dB, el_3db, rear_gain_lin,
M, N, pol, beam_angles, spacing);
% Multi-beam antenna (one beam per direction)
ant = quadriga_lib.arrayant_generate('multibeam_sep', res, freq, az_3dB, el_3db, rear_gain_lin,
M, N, pol, beam_angles, spacing);
% Optional for all types: output as separate variables, (must have exactly 11 outputs)
[e_theta_re, e_theta_im, e_phi_re, e_phi_im, azimuth_grid, elevation_grid, element_pos, ...
coupling_re, coupling_im, freq, name] = quadriga_lib.arrayant_generate( ... );
-
Input Arguments:
type [1]
Antenna model type, string
res [2]
Pattern resolution in [deg], scalar, default = 1 deg
freq [3]
The center frequency in [Hz], scalar, default = 299792458 Hz
-
Input arguments for type 'custom' and '3GPP' (custom beam width):
az_3dB [4]
3dB beam width in azimuth direction in [deg], scalar,
default for custom = 90 deg, default for 3gpp = 67 deg
el_3db [5]
3dB beam width in elevation direction in [deg], scalar,
default for custom = 90 deg, default for 3gpp = 67 deg
rear_gain_lin [6]
Isotropic gain (linear scale) at the back of the antenna, scalar, default = 0.0
-
Input arguments for type '3GPP':
M [7]
Number of vertically stacked elements, scalar, default = 1
N [8]
Number of horizontally stacked elements, scalar, default = 1
pol [9]
Polarization indicator to be applied for each of the M elements:
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 |
Polarization indicator is ignored when a custom pattern is provided.
tilt [10]
The electric downtilt angle in [deg], Only relevant for pol = 4/5/6, scalar, default = 0
spacing [11]
Element spacing in [λ], scalar, default = 0.5
Mg [12]
Number of nested panels in a column, scalar, default = 1
Ng [13]
Number of nested panels in a row, scalar, default = 1
dgv [14]
Panel spacing in vertical direction in [λ], scalar, default = 0.5
dgh [15]
Panel spacing in horizontal direction in [λ], scalar, default = 0.5
pattern [16]
Struct containing a custom pattern (default = empty) with at least the following fields:
e_theta_re_c |
Real part of e-theta field component |
Size: [n_elevation, n_azimuth, n_elements_c] |
e_theta_im_c |
Imaginary part of e-theta field component |
Size: [n_elevation, n_azimuth, n_elements_c] |
e_phi_re_c |
Real part of e-phi field component |
Size: [n_elevation, n_azimuth, n_elements_c] |
e_phi_im_c |
Imaginary part of e-phi field component |
Size: [n_elevation, n_azimuth, n_elements_c] |
azimuth_grid_c |
Azimuth angles in [rad] -pi to pi, sorted |
Size: [n_azimuth] |
elevation_grid_c |
Elevation angles in [rad], -pi/2 to pi/2, sorted |
Size: [n_elevation] |
If custom pattern data is not provided, the pattern is generated internally (either with a custom
beam width if az_3dB and el_3db are given or using the default 3GPP pattern).
-
Input arguments for type 'multibeam' and 'multibeam_sep':
beam_angles [10]
Matrix containing the beam steering angles in [deg] for the multi-beam antenna, size = [3, n_beams],
Rows are: [azimuth_deg, elevation_deg, weight].
-
Output Arguments:
ant
Struct containing the arrayant data with the following fields:
e_theta_re |
e-theta field component, real part |
Size: [n_elevation, n_azimuth, n_elements] |
e_theta_im |
e-theta field component, imaginary part |
Size: [n_elevation, n_azimuth, n_elements] |
e_phi_re |
e-phi field component, real part |
Size: [n_elevation, n_azimuth, n_elements] |
e_phi_im |
e-phi field component, imaginary part |
Size: [n_elevation, n_azimuth, n_elements] |
azimuth_grid |
Azimuth angles in [rad] -pi to pi, sorted |
Size: [n_azimuth] |
elevation_grid |
Elevation angles in [rad], -pi/2 to pi/2, sorted |
Size: [n_elevation] |
element_pos |
Antenna element (x,y,z) positions |
Size: [3, n_elements] |
coupling_re |
Coupling matrix, real part |
Size: [n_elements, n_ports] |
coupling_im |
Coupling matrix, imaginary part |
Size: [n_elements, n_ports] |
center_freq |
Center frequency in [Hz], optional, default = 0.3 GHz |
Scalar |
name |
Name of the array antenna object |
String |
arrayant_interpolate - Interpolate array antenna field patterns
-
Description:
This function interpolates polarimetric antenna field patterns for a given set of azimuth and
elevation angles.
-
Usage:
% Arrayant as struct
[V_re, V_im, H_re, H_im, dist, azimuth_loc, elevation_loc, gamma] = ...
quadriga_lib.arrayant_interpolate( arrayant, azimuth, elevation, element, orientation, element_pos );
% Arrayant as separate inputs
[V_re, V_im, H_re, H_im, dist, azimuth_loc, elevation_loc, gamma] = ...
quadriga_lib.arrayant_interpolate( [], azimuth, elevation, element, orientation, element_pos,
e_theta_re, e_theta_im, e_phi_re, e_phi_im, azimuth_grid, elevation_grid );
-
Input Arguments:
arrayant [1] (optional)
Struct containing the arrayant data with the following fields:
e_theta_re |
e-theta field component, real part |
Size: [n_elevation, n_azimuth, n_elements] |
e_theta_im |
e-theta field component, imaginary part |
Size: [n_elevation, n_azimuth, n_elements] |
e_phi_re |
e-phi field component, real part |
Size: [n_elevation, n_azimuth, n_elements] |
e_phi_im |
e-phi field component, imaginary part |
Size: [n_elevation, n_azimuth, n_elements] |
azimuth_grid |
Azimuth angles in [rad], -pi to pi, sorted |
Size: [n_azimuth] |
elevation_grid |
Elevation angles in [rad], -pi/2 to pi/2, sorted |
Size: [n_elevation] |
element_pos |
Antenna element (x,y,z) positions, optional |
Size: [3, n_elements] |
If an empty array is passed, array antenna data is provided as separate inputs (Inputs 7-12)
azimuth [2] (required)
Azimuth angles in [rad] for which the field pattern should be interpolated. Values must be
between -pi and pi, single or double precision.
| 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] |
elevation [3] (required)
Elevation angles in [rad] for which the field pattern should be interpolated. Values must be
between -pi/2 and pi/2, single or double precision.
| 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] |
element [4] (optional)
The element indices for which the interpolation should be done. Optional parameter. 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 [1:n_elements]. In this case, n_out = n_elements. Allowed types: uint32 or double.
Size: [1, n_out] or [n_out, 1] or empty []
orientation [5] (optional)
This (optional) 3-element vector describes the orientation of the array antenna or of individual
array elements. The The first value describes the ”bank angle”, the second value describes the
”tilt angle”, (positive values point upwards), the third value describes the bearing or ”heading
angle”, in mathematic sense. Values must be given in [rad]. East corresponds to 0, and the
angles increase counter-clockwise, so north is pi/2, south is -pi/2, and west is equal to pi. By
default, the orientation is [0,0,0]', i.e. the broadside of the antenna points at the horizon
towards the East. Single or double precision
Size: [3, 1] or [3, n_out] or [3, 1, n_ang] or [3, n_out, n_ang] or empty []
element_pos [6] (optional)
Alternative positions of the array antenna elements in local cartesian coordinates (using units of [m]).
If this parameter is not given, element positions arrayant are used. If the arrayant has no
positions, they are initialzed to [0,0,0]. For example, when duplicating the fist element by setting
element = [1,1], different element positions can be set for the two elements in the output.
Size: [3, n_out] or empty []
-
Derived inputs:
n_azimuth |
Number of azimuth angles in the filed pattern |
n_elevation |
Number of elevation angles in the filed pattern |
n_elements |
Number of antenna elements filed pattern of the array antenna |
n_ang |
Number of interpolation angles |
n_out |
Number of antenna elements in the generated output (may differ from n_elements) |
-
Output Arguments:
vr
Real part of the interpolated e-theta (vertical) field component. Size [n_out, n_ang]
vi
Imaginary part of the interpolated e-theta (vertical) field component. Size [n_out, n_ang]
hr
Real part of the interpolated e-phi (horizontal) field component. Size [n_out, n_ang]
hi
Imaginary part of the interpolated e-phi (horizontal) field component. Size [n_out, n_ang]
dist (optional)
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 [n_out, n_ang]
azimuth_loc (optional)
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 [n_out, n_ang]
elevation_loc (optional)
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 [n_out, n_ang]
gamma (optional)
Polarization rotation angles in [rad], Size [n_out, n_ang]
arrayant_qdant_read - Reads array antenna data from QDANT files
-
Description:
The QuaDRiGa array antenna exchange format (QDANT) is a file format used to store antenna pattern
data in XML. This function reads pattern data from the specified file.
-
Usage:
% Read as struct
[ ant, layout ] = quadriga_lib.arrayant_qdant_read( fn, id );
% Read as separate fields
[e_theta_re, e_theta_im, e_phi_re, e_phi_im, azimuth_grid, elevation_grid, element_pos, coupling_re,
coupling_im, center_freq, name, layout ] = quadriga_lib.arrayant_qdant_read( fn, id );
-
Input Arguments:
fn
Filename of the QDANT file, string
id (optional)
ID of the antenna to be read from the file, optional, Default: Read first
-
Output Arguments:
ant
Struct containing the arrayant data with the following fields:
e_theta_re |
e-theta field component, real part |
Size: [n_elevation, n_azimuth, n_elements] |
e_theta_im |
e-theta field component, imaginary part |
Size: [n_elevation, n_azimuth, n_elements] |
e_phi_re |
e-phi field component, real part |
Size: [n_elevation, n_azimuth, n_elements] |
e_phi_im |
e-phi field component, imaginary part |
Size: [n_elevation, n_azimuth, n_elements] |
azimuth_grid |
Azimuth angles in [rad] -pi to pi, sorted |
Size: [n_azimuth] |
elevation_grid |
Elevation angles in [rad], -pi/2 to pi/2, sorted |
Size: [n_elevation] |
element_pos |
Antenna element (x,y,z) positions |
Size: [3, n_elements] |
coupling_re |
Coupling matrix, real part |
Size: [n_elements, n_ports] |
coupling_im |
Coupling matrix, imaginary part |
Size: [n_elements, n_ports] |
center_freq |
Center frequency in [Hz], optional, default = 0.3 GHz |
Scalar |
name |
Name of the array antenna object |
String |
layout
Layout of multiple array antennas. Contain element ids that are present in the file
-
See also:
arrayant_qdant_write - Writes array antenna data to QDANT files
-
Description:
The QuaDRiGa array antenna exchange format (QDANT) is a file format used to store antenna pattern
data in XML. This function writes pattern data to the specified file.
-
Usage:
% Arrayant as struct
id_in_file = quadriga_lib.arrayant_qdant_write( fn, arrayant, id, layout);
% Arrayant as separate inputs
id_in_file = quadriga_lib.arrayant_qdant_write( fn, [], id, layout, e_theta_re, e_theta_im, e_phi_re,
e_phi_im, azimuth_grid, elevation_grid, element_pos, coupling_re, coupling_im, center_freq, name);
-
Caveat:
- 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.
-
Input Arguments:
fn [1]
Filename of the QDANT file, string
arrayant [2] (optional)
Struct containing the arrayant data with the following fields:
e_theta_re |
e-theta field component, real part |
Size: [n_elevation, n_azimuth, n_elements] |
e_theta_im |
e-theta field component, imaginary part |
Size: [n_elevation, n_azimuth, n_elements] |
e_phi_re |
e-phi field component, real part |
Size: [n_elevation, n_azimuth, n_elements] |
e_phi_im |
e-phi field component, imaginary part |
Size: [n_elevation, n_azimuth, n_elements] |
azimuth_grid |
Azimuth angles in [rad], -pi to pi, sorted |
Size: [n_azimuth] |
elevation_grid |
Elevation angles in [rad], -pi/2 to pi/2, sorted |
Size: [n_elevation] |
element_pos |
Antenna element (x,y,z) positions, optional |
Size: [3, n_elements] |
coupling_re |
Coupling matrix, real part, optional |
Size: [n_elements, n_ports] |
coupling_im |
Coupling matrix, imaginary part, optional |
Size: [n_elements, n_ports] |
center_freq |
Center frequency in [Hz], optional, default = 0.3 GHz |
Scalar |
name |
Name of the array antenna object, optional |
String |
If an empty array is passed, array antenna data is provided as separate inputs (Inputs 5-15)
id [3] (optional)
ID of the antenna to be written to the file, optional, Default: Max-ID in existing file + 1
layout [4] (optional)
Layout of multiple array antennas. Must only contain element ids that are present in the file. optional
-
Output Argument:
id_in_file
ID of the antenna in the file after writing
-
See also:
arrayant_rotate_pattern - Rotates antenna patterns
-
Description:
This MATLAB function transforms the radiation patterns of array antenna elements, allowing for
precise rotations around the three principal axes (x, y, z) of the local Cartesian coordinate system.
This is essential in antenna design and optimization, enabling engineers to tailor the radiation
pattern for enhanced performance. The function also adjusts the sampling grid for non-uniformly
sampled antennas, such as parabolic antennas with small apertures, ensuring accurate and efficient
computations. The 3 rotations are applies in the order: 1. rotation around the x-axis (bank angle);
2. rotation around the y-axis (tilt angle), 3. rotation around the z-axis (heading angle)
-
Usage:
% Minimal example (input/output = struct)
arrayant_out = quadriga_lib.arrayant_rotate_pattern(arrayant_in, bank, tilt, head, usage, element);
% Separate outputs, struct input
[e_theta_re, e_theta_im, e_phi_re, e_phi_im, azimuth_grid, elevation_grid, element_pos, ...
coupling_re, coupling_im, freq, name] = quadriga_lib.arrayant_rotate_pattern(arrayant_in, ...
bank, tilt, head, usage, element);
% Separate inputs
arrayant_out = quadriga_lib.arrayant_rotate_pattern([], bank, tilt, head, usage, element, ...
e_theta_re, e_theta_im, e_phi_re, e_phi_im, azimuth_grid, elevation_grid, element_pos, ...
coupling_re, coupling_im, freq, name);
-
Input Arguments:
arrayant_in [1]
Struct containing the arrayant data with the following fields:
e_theta_re |
e-theta field component, real part |
Size: [n_elevation, n_azimuth, n_elements] |
e_theta_im |
e-theta field component, imaginary part |
Size: [n_elevation, n_azimuth, n_elements] |
e_phi_re |
e-phi field component, real part |
Size: [n_elevation, n_azimuth, n_elements] |
e_phi_im |
e-phi field component, imaginary part |
Size: [n_elevation, n_azimuth, n_elements] |
azimuth_grid |
Azimuth angles in [rad], -pi to pi, sorted |
Size: [n_azimuth] |
elevation_grid |
Elevation angles in [rad], -pi/2 to pi/2, sorted |
Size: [n_elevation] |
element_pos |
Antenna element (x,y,z) positions, optional |
Size: [3, n_elements] |
coupling_re |
Coupling matrix, real part, optional |
Size: [n_elements, n_ports] |
coupling_im |
Coupling matrix, imaginary part, optional |
Size: [n_elements, n_ports] |
center_freq |
Center frequency in [Hz], optional |
Scalar |
name |
Name of the array antenna object, optional |
String |
If an empty array is passed, array antenna data is provided as separate inputs (Inputs 7-17)
x_deg [2] (optional)
The rotation angle around x-axis (bank angle) in [degrees]
y_deg [3] (optional)
The rotation angle around y-axis (tilt angle) in [degrees]
z_deg [4] (optional)
The rotation angle around z-axis (heading angle) in [degrees]
usage [5] (optional)
The optional parameter 'usage' can limit the rotation procedure either to the pattern or polarization.
usage = 0 |
Rotate both, pattern and polarization, adjusts sampling grid (default) |
usage = 1 |
Rotate only pattern, adjusts sampling grid |
usage = 2 |
Rotate only polarization |
usage = 3 |
Rotate both, but do not adjust the sampling grid |
element [6] (optional)
The element indices for which the pattern should be transformed. Optional parameter. Values must
be between 1 and n_elements. If this parameter is not provided (or an empty array is passed),
all elements will be rotated by the same angles. Size: [1, n_elements] or empty []
-
Output Arguments:
arrayant_out
Struct containing the arrayant data with the following fields:
e_theta_re |
e-theta field component, real part |
Size: [n_elevation, n_azimuth, n_elements] |
e_theta_im |
e-theta field component, imaginary part |
Size: [n_elevation, n_azimuth, n_elements] |
e_phi_re |
e-phi field component, real part |
Size: [n_elevation, n_azimuth, n_elements] |
e_phi_im |
e-phi field component, imaginary part |
Size: [n_elevation, n_azimuth, n_elements] |
azimuth_grid |
Azimuth angles in [rad] -pi to pi, sorted |
Size: [n_azimuth] |
elevation_grid |
Elevation angles in [rad], -pi/2 to pi/2, sorted |
Size: [n_elevation] |
element_pos |
Antenna element (x,y,z) positions |
Size: [3, n_elements] |
coupling_re |
Coupling matrix, real part |
Size: [n_elements, n_ports] |
coupling_im |
Coupling matrix, imaginary part |
Size: [n_elements, n_ports] |
center_freq |
Center frequency in [Hz], default = 0.3 GHz |
Scalar |
name |
Name of the array antenna object |
String |
Can be returned as separate outputs.
Channel functions
baseband_freq_response - Transforms the channel into frequency domain and returns the frequency response
-
Usage:
[ hmat_re, hmat_im ] = quadriga_lib.baseband_freq_response( coeff_re, coeff_im, delay, pilot_grid, bandwidth, i_snap );
-
Input Arguments:
coeff_re
Channel coefficients, real part, Size: [ n_rx, n_tx, n_path, n_snap ]
coeff_im
Channel coefficients, imaginary part, Size: [ n_rx, n_tx, n_path, n_snap ]
delays
Propagation delay in seconds, Size: [ n_rx, n_tx, n_path, n_snap ] or [ 1, 1, n_path, n_snap ]
or [ n_path, n_snap ]
pilot_grid
Sub-carrier positions relative to the bandwidth. The carrier positions are given relative to the
bandwidth where '0' is the begin of the spectrum (i.e., the center frequency f0) and '1' is
equal to f0+bandwidth. To obtain the channel frequency response centered around f0, the
input variable 'pilot_grid' must be set to '(-N/2:N/2)/N', where N is the number of sub-
carriers. Vector of length: [ n_carriers ]
bandwidth
The baseband bandwidth in [Hz], scalar
i_snap (optional)
Snapshot indices for which the frequency response should be generated (1-based index). If this
variable is not given, all snapshots are processed. Length: [ n_out ]
-
Output Argument:
hmat_re
Freq. domain channel matrices (H), real part, Size [ n_rx, n_tx, n_carriers, n_out ]
hmat_im
Freq. domain channel matrices (H), imaginary part, Size [ n_rx, n_tx, n_carriers, n_out ]
channel_export_obj_file - Export path data to a Wavefront OBJ file for visualization in Blender
-
Description:
This function exports path data to a Wavefront OBJ file, which can be used for visualization in 3D
software such as Blender. It supports various colormaps for color-coding the paths based on their
gain values. In addition, the function allows you to control the maximum number of paths displayed,
set gain thresholds for color-coding and selection.
-
Usage:
quadriga_lib.channel_export_obj_file( fn, max_no_paths, gain_max, gain_min, colormap, radius_max,
radius_min, n_edges, rx_position, tx_position, no_interact, interact_coord, center_freq,
coeff_re, coeff_im, i_snap )
-
Input Arguments:
fn
Filename of the OBJ file, string, required
max_no_paths (optional)
Maximum number of paths to be shown, optional, default: 0 = export all above gain_min
gain_max (optional)
Maximum path gain in dB (only for color-coding), optional, default = -60.0
gain_min (optional)
Minimum path gain in dB (for color-coding and path selection), optional, default = -140.0
colormap (optional)
Colormap for the visualization, string, supported are 'jet', 'parula', 'winter', 'hot', 'turbo',
'copper', 'spring', 'cool', 'gray', 'autumn', 'summer', optional, default = 'jet'
radius_max (optional)
Maximum tube radius in meters, optional, default = 0.05
radius_min (optional)
Minimum tube radius in meters, optional, default = 0.01
n_edges (optional)
Number of vertices in the circle building the tube, must be >= 3, optional, default = 5
rx_position
Receiver positions, required, size [3, n_snap] or [3, 1]
tx_position
Transmitter positions, required, size [3, n_snap] or [3, 1]
no_interact
Number interaction points of paths with the environment, required, uint32, Size [n_path, n_snap]
interact_coord
Interaction coordinates, required, Size [3, max(sum(no_interact)), n_snap]
center_freq
Center frequency in [Hz], required, Size [n_snap, 1] or scalar
coeff_re
Channel coefficients, real part, Size: [ n_rx, n_tx, n_path, n_snap ]
coeff_im
Channel coefficients, imaginary part, Size: [ n_rx, n_tx, n_path, n_snap ]
i_snap
(optional)
Snapshot indices, optional, 1-based, range [1 ... n_snap]
-
Output Argument:
This function does not return a value. It writes the OBJ file directly to disk.
hdf5_create_file - Create a new HDF5 channel file with a custom storage layout
hdf5_read_channel - Reads channel data from HDF5 files
-
Description:
Quadriga-Lib provides an HDF5-based solution for storing and organizing channel data. This data
comprises various well-defined sets, including channel coefficients, positions of transmitters and
receivers, as well as path data that reflects the interaction of radio waves with the environment.
Typically, these datasets are multi-dimensional, encompassing data for n_rx receive antennas,
n_tx transmit antennas, n_path propagation paths, and n_snap snapshots. Snapshots are
particularly useful for recording data across different locations (such as along a trajectory) or
various frequencies. It is important to note that not all datasets include all these dimensions.
The library also supports the addition of extra datasets of any type or shape, which can be useful
for incorporating descriptive data or analysis results. To facilitate data access, the function
quadriga_lib.hdf5_read_channel is designed to read both structured and unstructured data from the
file.
-
Usage:
[ par, rx_position, tx_position, coeff_re, coeff_im, delay, center_freq, name, initial_pos, ...
path_gain, path_length, path_polarization, path_angles, path_fbs_pos, path_lbs_pos, no_interact, ...
interact_coord, rx_orientation, tx_orientation ] = quadriga_lib.hdf5_read_channel( fn, location, snap );
-
Input Arguments:
fn
Filename of the HDF5 file, string
location (optional)
Storage location inside the file; 1-based; vector with 1-4 elements, i.e. [ix], [ix, iy],
[ix,iy,iz] or [ix,iy,iz,iw]; Default: ix = iy = iz = iw = 1
snap (optional)
Snapshot range; optional; vector, default = read all
-
Output Arguments:
par
Unstructured data as struct, may be empty if no unstructured data is present
- Structured data: (outputs 2-19, single precision)
rx_position |
Receiver positions |
[3, n_snap] or [3, 1] |
tx_position |
Transmitter positions |
[3, n_snap] or [3, 1] |
coeff_re |
Channel coefficients, real part |
[n_rx, n_tx, n_path, n_snap] |
coeff_im |
Channel coefficients, imaginary part |
[n_rx, n_tx, n_path, n_snap] |
delay |
Propagation delays in seconds |
[n_rx, n_tx, n_path, n_snap] or [1, 1, n_path, n_snap] |
center_freq |
Center frequency in [Hz] |
[n_snap, 1] or scalar |
name |
Name of the channel |
String |
initial_pos |
Index of reference position, 1-based |
uint32, scalar |
path_gain |
Path gain before antenna, linear scale |
[n_path, n_snap] |
path_length |
Path length from TX to RX phase center in m |
[n_path, n_snap] |
polarization |
Polarization transfer function, interleaved complex |
[8, n_path, n_snap] |
path_angles |
Departure and arrival angles {AOD, EOD, AOA, EOA} in rad |
[n_path, 4, n_snap] |
path_fbs_pos |
First-bounce scatterer positions |
[3, n_path, n_snap] |
path_lbs_pos |
Last-bounce scatterer positions |
[3, n_path, n_snap] |
no_interact |
Number interaction points of paths with the environment |
uint32, [n_path, n_snap] |
interact_coord |
Interaction coordinates |
[3, max(sum(no_interact)), n_snap] |
rx_orientation |
Transmitter orientation |
[3, n_snap] or [3, 1] |
tx_orientation |
Receiver orientation |
[3, n_snap] or [3, 1] |
-
Caveat:
- Empty outputs are returned if data set does not exist in the file
- All structured data is stored in single precision. Hence, outputs are also in single precision.
- Unstructured datatypes are returned as stored in the HDF file (same type, dimensions and storage order)
- Typically,
n_path may vary for each snapshot. In such cases, n_path is set to the maximum value found
within the range of snapshots, and any missing paths are padded with zeroes.
hdf5_read_dset - Read a single unstructured dataset from an HDF5 file
hdf5_read_dset_names - Read the names of unstructured data fields from an HDF5 file
hdf5_read_layout - Read the storage layout of channel data inside an HDF5 file
hdf5_reshape_layout - Reshapes the storage layout inside an existing HDF5 file
hdf5_write_channel - Writes channel data to HDF5 files
-
Description:
Quadriga-Lib provides an HDF5-based solution for storing and organizing channel data. This function
can be used to write structured and unstructured data to an HDF5 file.
-
Usage:
storage_dims = quadriga_lib.hdf5_write_channel( fn, location, rx_position, tx_position, ...
coeff_re, coeff_im, delay, center_freq, name, initial_pos, path_gain, path_length, ...
path_polarization, path_angles, path_fbs_pos, path_lbs_pos, no_interact, interact_coord, ...
rx_orientation, tx_orientation )
-
Input Arguments:
fn
Filename of the HDF5 file, string
location (optional)
Storage location inside the file; 1-based; vector with 1-4 elements, i.e. [ix], [ix, iy],
[ix,iy,iz] or [ix,iy,iz,iw]; Default: ix = iy = iz = iw = 1
par
Unstructured data as struct, can be empty if no unstructured data should be written
- Structured data: (inputs 4-21, single or double precision)
rx_position |
Receiver positions |
[3, n_snap] or [3, 1] |
tx_position |
Transmitter positions |
[3, n_snap] or [3, 1] |
coeff_re |
Channel coefficients, real part |
[n_rx, n_tx, n_path, n_snap] |
coeff_im |
Channel coefficients, imaginary part |
[n_rx, n_tx, n_path, n_snap] |
delay |
Propagation delays in seconds |
[n_rx, n_tx, n_path, n_snap] or [1, 1, n_path, n_snap] |
center_freq |
Center frequency in [Hz] |
[n_snap, 1] or scalar |
name |
Name of the channel |
String |
initial_pos |
Index of reference position, 1-based |
uint32, scalar |
path_gain |
Path gain before antenna, linear scale |
[n_path, n_snap] |
path_length |
Path length from TX to RX phase center in m |
[n_path, n_snap] |
polarization |
Polarization transfer function, interleaved complex |
[8, n_path, n_snap] |
path_angles |
Departure and arrival angles {AOD, EOD, AOA, EOA} in rad |
[n_path, 4, n_snap] |
path_fbs_pos |
First-bounce scatterer positions |
[3, n_path, n_snap] |
path_lbs_pos |
Last-bounce scatterer positions |
[3, n_path, n_snap] |
no_interact |
Number interaction points of paths with the environment |
uint32, [n_path, n_snap] |
interact_coord |
Interaction coordinates |
[3, max(sum(no_interact)), n_snap] |
rx_orientation |
Transmitter orientation |
[3, n_snap] or [3, 1] |
tx_orientation |
Receiver orientation |
[3, n_snap] or [3, 1] |
-
Output Arguments:
storage_dims
Size of the dimensions of the storage space, vector with 4 elements, i.e. [nx,ny,nz,nw].
-
Caveat:
- If the file exists already, the new data is added to the exisiting file
- If a new file is created, a storage layout is created to store the location of datasets in the file
- For
location = [ix] storage layout is [65536,1,1,1] or [ix,1,1,1] if (ix > 65536)
- For
location = [ix,iy] storage layout is [1024,64,1,1]
- For
location = [ix,iy,iz] storage layout is [256,16,16,1]
- For
location = [ix,iy,iz,iw] storage layout is [128,8,8,8]
- You can create a custom storage layout by creating the file first using "
hdf5_create_file"
- You can reshape the storage layout by using "
hdf5_reshape_storage", but the total number of elements must not change
- Inputs can be empty or missing.
- All structured data is written in single precision (but can can be provided as single 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_write_dset - Writes unstructured data to a HDF5 file
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.
- For 3D antenna models (default models from arrayant_generate), only the azimuth cut at
elevation_grid = 0 is used
- 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:
chan = quadriga_lib.get_channels_ieee_indoor(ap_array, sta_array, ChannelType, CarrierFreq_Hz, ...
tap_spacing_s, n_users, observation_time, update_rate, speed_station_kmh, speed_env_kmh, ...
Dist_m, n_floors, uplink, offset_angles, n_subpath, Doppler_effect, seed, ...
KF_linear, XPR_NLOS_linear, SF_std_dB_LOS, SF_std_dB_NLOS, dBP_m);
-
ap_array:
ap_array [1]
Struct containing the access point array antenna with n_tx elements (= ports after element coupling)
e_theta_re |
Real part of e-theta field component |
Size: [n_elevation_ap, n_azimuth_ap, n_elements_ap] |
e_theta_im |
Imaginary part of e-theta field component |
Size: [n_elevation_ap, n_azimuth_ap, n_elements_ap] |
e_phi_re |
Real part of e-phi field component |
Size: [n_elevation_ap, n_azimuth_ap, n_elements_ap] |
e_phi_im |
Imaginary part of e-phi field component |
Size: [n_elevation_ap, n_azimuth_ap, n_elements_ap] |
azimuth_grid |
Azimuth angles in [rad] -pi to pi, sorted |
Size: [n_azimuth_ap] |
elevation_grid |
Elevation angles in [rad], -pi/2 to pi/2, sorted |
Size: [n_elevation_ap] |
sta_array [2]
Struct containing the mobile station array antenna with n_rx elements (= ports after element coupling)
e_theta_re |
Real part of e-theta field component |
Size: [n_elevation_sta, n_azimuth_sta, n_elements_sta] |
e_theta_im |
Imaginary part of e-theta field component |
Size: [n_elevation_sta, n_azimuth_sta, n_elements_sta] |
e_phi_re |
Real part of e-phi field component |
Size: [n_elevation_sta, n_azimuth_sta, n_elements_sta] |
e_phi_im |
Imaginary part of e-phi field component |
Size: [n_elevation_sta, n_azimuth_sta, n_elements_sta] |
azimuth_grid |
Azimuth angles in [rad] -pi to pi, sorted |
Size: [n_azimuth_sta] |
elevation_grid |
Elevation angles in [rad], -pi/2 to pi/2, sorted |
Size: [n_elevation_sta] |
ChannelType [3]
Channel model type as defined by TGn. String. Supported: A, B, C, D, E, F.
CarrierFreq_Hz = 5.25e9 [4] (optional)
Carrier frequency in Hz.
tap_spacing_s = 10e-9 [5] (optional)
Tap spacing in seconds. Must be equal to 10 ns / 2^k (TGn default = 10e-9).
n_users = 1 [6] (optional)
Number of users (only for TGac, TGah). Output struct array length equals n_users.
observation_time = 0 [7] (optional)
Channel observation time in seconds. 0 creates a static channel.
update_rate = 1e-3 [8] (optional)
Channel update interval in seconds (only relevant when observation_time > 0).
speed_station_kmh = 0 [9] (optional)
Station movement speed in km/h. Movement direction is AoA_offset. Only relevant when observation_time > 0.
speed_env_kmh = 1.2 [10] (optional)
Environment movement speed in km/h. Default 1.2 for TGn, use 0.089 for TGac. Only relevant when observation_time > 0.
vector Dist_m = [4.99] [11] (optional)
TX-to-RX distance(s) in meters. Length n_users or length 1 (same distance for all users).
vector n_floors = [0] [12] (optional)
Number of floors for TGah model (per user), up to 4 floors. Length n_users or length 1.
uplink = false [13] (optional)
Channel direction flag. Default is downlink; set to true to generate reverse (uplink) direction.
offset_angles = [] [14] (optional)
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.
n_subpath = 20 [15] (optional)
Number of sub-paths per path/cluster used for Laplacian angular spread mapping.
Doppler_effect = 50 [16] (optional)
Special Doppler effects: models D, E (fluorescent lights, value = mains freq.) and F (moving vehicle speed in km/h).
Use 0 to disable.
seed = -1 [17] (optional)
Numeric seed for repeatability. -1 disables the fixed seed and uses the system random device.
KF_linear = [] [18] (optional)
Overwrites the model-specific KF-value. If this parameter is empty (default), NAN 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.
XPR_NLOS_linear = [] [19] (optional)
Overwrites the model-specific Cross-polarization ratio. If this parameter is empty (default), NAN or negative,
the model default of 2 (3 dB) is used. XPR is applied to all NLOS taps.
SF_std_dB_LOS = [] [20] (optional)
Overwrites the model-specific shadow fading for LOS channels. If this parameter is empty (default) or NAN,
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.
SF_std_dB_NLOS = [] [21] (optional)
Overwrites the model-specific shadow fading for LOS channels. If this parameter is empty (default) or NAN,
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.
dBP_m = [] [22] (optional)
Overwrites the model-specific breakpoint distance. If this parameter is empty (default), NAN or negative,
the model defaults are A/B/C: 5 m, D: 10 m, E: 20 m, F: 30 m.
-
Returns:
chan
Struct array of length n_users containing the channel data with the following fields.
name |
Channel name |
String |
tx_position |
Transmitter positions (AP for downlink, STA for uplink) |
Size: [3, 1] or [3, n_snap] |
rx_position |
Receiver positions (STA for downlink, AP for uplink) |
Size: [3, 1] or [3, n_snap] |
tx_orientation |
Transmitter orientation, Euler angles (AP for downlink, STA for uplink) |
Size: [3, 1] or [3, n_snap] |
rx_orientation |
Receiver orientation, Euler angles (STA for downlink, AP for uplink) |
Size: [3, 1] or [3, n_snap] |
coeff_re |
Channel coefficients, real part |
Size: [n_rx, n_tx, n_path, n_snap] |
coeff_im |
Channel coefficients, imaginary part |
Size: [n_rx, n_tx, n_path, n_snap] |
delay |
Propagation delays in seconds |
Size: [n_rx, n_tx, n_path, n_snap] |
path_gain |
Path gain before antenna, linear scale |
Size: [n_path, n_snap] |
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.
-
Usage:
[ coeff_re, coeff_im, delays, active_path, aod, eod, aoa, eoa ] = quadriga_lib.get_channels_irs( ...
ant_tx, ant_rx, ant_irs, ...
fbs_pos_1, lbs_pos_1, path_gain_1, path_length_1, M_1, ...
fbs_pos_2, lbs_pos_2, path_gain_2, path_length_2, M_2, ...
tx_pos, tx_orientation, rx_pos, rx_orientation, irs_pos, irs_orientation, ...
i_irs, threshold_dB, center_freq, use_absolute_delays, active_path, ant_irs_2 );
-
Input Arguments:
ant_tx [1] (required)
Struct containing the transmit (TX) arrayant data with the following fields:
e_theta_re |
e-theta field component, real part |
Size: [n_elevation_tx, n_azimuth_tx, n_elements_tx] |
e_theta_im |
e-theta field component, imaginary part |
Size: [n_elevation_tx, n_azimuth_tx, n_elements_tx] |
e_phi_re |
e-phi field component, real part |
Size: [n_elevation_tx, n_azimuth_tx, n_elements_tx] |
e_phi_im |
e-phi field component, imaginary part |
Size: [n_elevation_tx, n_azimuth_tx, n_elements_tx] |
azimuth_grid |
Azimuth angles in [rad], -pi to pi, sorted |
Size: [n_azimuth_tx] |
elevation_grid |
Elevation angles in [rad], -pi/2 to pi/2, sorted |
Size: [n_elevation_tx] |
element_pos |
Antenna element (x,y,z) positions, optional |
Size: [3, n_elements_tx] |
coupling_re |
Coupling matrix, real part, optional |
Size: [n_elements_tx, n_ports_tx] |
coupling_im |
Coupling matrix, imaginary part, optional |
Size: [n_elements_tx, n_ports_tx] |
ant_rx [2] (required)
Struct containing the receive (RX) arrayant data with the following fields:
e_theta_re |
e-theta field component, real part |
Size: [n_elevation_rx, n_azimuth_rx, n_elements_rx] |
e_theta_im |
e-theta field component, imaginary part |
Size: [n_elevation_rx, n_azimuth_rx, n_elements_rx] |
e_phi_re |
e-phi field component, real part |
Size: [n_elevation_rx, n_azimuth_rx, n_elements_rx] |
e_phi_im |
e-phi field component, imaginary part |
Size: [n_elevation_rx, n_azimuth_rx, n_elements_rx] |
azimuth_grid |
Azimuth angles in [rad], -pi to pi, sorted |
Size: [n_azimuth_rx] |
elevation_grid |
Elevation angles in [rad], -pi/2 to pi/2, sorted |
Size: [n_elevation_rx] |
element_pos |
Antenna element (x,y,z) positions, optional |
Size: [3, n_elements_rx] |
coupling_re |
Coupling matrix, real part, optional |
Size: [n_elements_rx, n_ports_rx] |
coupling_im |
Coupling matrix, imaginary part, optional |
Size: [n_elements_rx, n_ports_rx] |
ant_irs [3] (required)
Struct containing the intelligent reflective surface (IRS) model:
e_theta_re |
e-theta field component, real part |
Size: [n_elevation_irs, n_azimuth_irs, n_elements_irs] |
e_theta_im |
e-theta field component, imaginary part |
Size: [n_elevation_irs, n_azimuth_irs, n_elements_irs] |
e_phi_re |
e-phi field component, real part |
Size: [n_elevation_irs, n_azimuth_irs, n_elements_irs] |
e_phi_im |
e-phi field component, imaginary part |
Size: [n_elevation_irs, n_azimuth_irs, n_elements_irs] |
azimuth_grid |
Azimuth angles in [rad], -pi to pi, sorted |
Size: [n_azimuth_irs] |
elevation_grid |
Elevation angles in [rad], -pi/2 to pi/2, sorted |
Size: [n_elevation_irs] |
element_pos |
Antenna element (x,y,z) positions, optional |
Size: [3, n_elements_irs] |
coupling_re |
Coupling matrix, real part, optional |
Size: [n_elements_irs, n_ports_irs] |
coupling_im |
Coupling matrix, imaginary part, optional |
Size: [n_elements_irs, n_ports_irs] |
fbs_pos_1 [4] (required)
First-bounce scatterer positions of TX → IRS paths, Size: [ 3, n_path_1 ].
lbs_pos_1 [5] (required)
Last-bounce scatterer positions of TX → IRS paths, Size [3, n_path_1].
path_gain_1 [6] (required)
Path gains (linear) for TX → IRS paths, Length n_path_1.
path_length_1 [7] (required)
Path lengths for TX → IRS paths, Length n_path_1.
M_1 [8] (required)
Polarization transfer matrix for TX → IRS paths, Size [8, n_path_1].
fbs_pos_2 [9] (required)
First-bounce scatterer positions of IRS → RX paths, Size: [ 3, n_path_2 ]
lbs_pos_2 [10] (required)
Last-bounce scatterer positions of IRS → RX paths, Size [3, n_path_2]
path_gain_2 [11] (required)
Path gains (linear) for IRS → RX paths, Length n_path_2.
path_length_2 [12] (required)
Path lengths for IRS → RX paths, Length n_path_2.
M_2 [13] (required)
Polarization transfer matrix for IRS → RX paths, Size [8, n_path_2].
tx_pos [14] (required)
Transmitter position in 3D Cartesian coordinates, Size: [3,1] or [1,3]
tx_orientation [15] (required)
3-element vector describing the orientation of the transmit antenna in Euler angles (bank, tilt, heading),
Size: [3,1] or [1,3]
rx_pos [16] (required)
Receiver position in 3D Cartesian coordinates, Size: [3,1] or [1,3]
rx_orientation [17] (required)
3-element vector describing the orientation of the receive antenna, Size: [3,1] or [1,3]
irs_pos [18] (required)
IRS position in 3D Cartesian coordinates, Size: [3,1] or [1,3]
irs_orientation [19] (required)
3-element (Euler) vector in Radians describing the orientation of the IRS, Size: [3,1] or [1,3]
i_irs [20] (optional)
Index of IRS codebook entry (port number), Scalar, Default: 0.
threshold_dB [21] (optional)
Threshold (in dB) below which paths are discarded, Scalar, Default: -140.0.
center_freq [22] (optional)
Center frequency in [Hz]; optional; If the value is not provided or set to 0, phase calculation
in coefficients is disabled, i.e. that path length has not influence on the results. This can be
used to calculate the antenna response for a specific angle and polarization. Scalar value
use_absolute_delays [23] (optional)
If true, the LOS delay is included for all paths; Default is false, i.e. delays are normalized
to the LOS delay.
active_path [24] (optional)
Optional bitmask for selecting active TX-IRS and IRS-RX path pairs. Ignores threshold_dB when provided.
ant_irs_2 [25] (optional)
Optional second IRS array (TX side for IRS → RX paths) for asymmetric IRS behavior. Same structure as for ant_irs
-
Output Arguments:
coeff_re
Channel coefficients, real part, Size: [ n_ports_rx, n_ports_tx, n_path ]
coeff_im
Channel coefficients, imaginary part, Size: [ n_ports_rx, n_ports_tx, n_path ]
delays
Propagation delay in seconds, Size: [ n_ports_rx, n_ports_tx, n_path ]
active_path (optional)
Boolean mask of length n_path_1 * n_path_2, indicating which path combinations were used.
aod (optional)
Azimuth of Departure angles in [rad], Size: [ n_ports_rx, n_ports_tx, n_path ]
eod (optional)
Elevation of Departure angles in [rad], Size: [ n_ports_rx, n_ports_tx, n_path ]
aoa (optional)
Azimuth of Arrival angles in [rad], Size: [ n_ports_rx, n_ports_tx, n_path ]
eoa (optional)
Elevation of Arrival angles in [rad], Size: [ n_ports_rx, n_ports_tx, n_path ]
get_channels_planar - Calculate channel coefficients for planar waves
-
Description:
In this function, the wireless propagation channel between a transmitter and a receiver is calculated,
based on a single transmit and receive position. Additionally, interaction points with the environment,
which are derived from either Ray Tracing or Geometric Stochastic Models such as QuaDRiGa, are
considered. The calculation is performed under the assumption of planar wave propagation. For accurate
execution of this process, several pieces of input data are required:
- The 3D Cartesian (local) coordinates of both the transmitter and the receiver.
- The azimuth/elevation departure and arrval angles.
- The polarization transfer matrix for each propagation path.
- Antenna models for both the transmitter and the receiver.
- The orientations of the antennas.
-
Usage:
[ coeff_re, coeff_im, delays, rx_Doppler ] = quadriga_lib.get_channels_planar( ant_tx, ant_rx, ...
aod, eod, aoa, eoa, path_gain, path_length, M, tx_pos, tx_orientation, rx_pos, rx_orientation, ...
center_freq, use_absolute_delays, add_fake_los_path );
-
Input Arguments:
ant_tx [1] (required)
Struct containing the transmit (TX) arrayant data with the following fields:
e_theta_re |
e-theta field component, real part |
Size: [n_elevation_tx, n_azimuth_tx, n_elements_tx] |
e_theta_im |
e-theta field component, imaginary part |
Size: [n_elevation_tx, n_azimuth_tx, n_elements_tx] |
e_phi_re |
e-phi field component, real part |
Size: [n_elevation_tx, n_azimuth_tx, n_elements_tx] |
e_phi_im |
e-phi field component, imaginary part |
Size: [n_elevation_tx, n_azimuth_tx, n_elements_tx] |
azimuth_grid |
Azimuth angles in [rad], -pi to pi, sorted |
Size: [n_azimuth_tx] |
elevation_grid |
Elevation angles in [rad], -pi/2 to pi/2, sorted |
Size: [n_elevation_tx] |
element_pos |
Antenna element (x,y,z) positions, optional |
Size: [3, n_elements_tx] |
coupling_re |
Coupling matrix, real part, optional |
Size: [n_elements_tx, n_ports_tx] |
coupling_im |
Coupling matrix, imaginary part, optional |
Size: [n_elements_tx, n_ports_tx] |
ant_rx [2] (required)
Struct containing the receive (RX) arrayant data with the following fields:
e_theta_re |
e-theta field component, real part |
Size: [n_elevation_rx, n_azimuth_rx, n_elements_rx] |
e_theta_im |
e-theta field component, imaginary part |
Size: [n_elevation_rx, n_azimuth_rx, n_elements_rx] |
e_phi_re |
e-phi field component, real part |
Size: [n_elevation_rx, n_azimuth_rx, n_elements_rx] |
e_phi_im |
e-phi field component, imaginary part |
Size: [n_elevation_rx, n_azimuth_rx, n_elements_rx] |
azimuth_grid |
Azimuth angles in [rad], -pi to pi, sorted |
Size: [n_azimuth_rx] |
elevation_grid |
Elevation angles in [rad], -pi/2 to pi/2, sorted |
Size: [n_elevation_rx] |
element_pos |
Antenna element (x,y,z) positions, optional |
Size: [3, n_elements_rx] |
coupling_re |
Coupling matrix, real part, optional |
Size: [n_elements_rx, n_ports_rx] |
coupling_im |
Coupling matrix, imaginary part, optional |
Size: [n_elements_rx, n_ports_rx] |
aod [3] (required)
Departure azimuth angles in [rad], Size: [ 1, n_path ] or [ n_path, 1 ]
eod [4] (required)
Departure elevation angles in [rad], Size: [ 1, n_path ] or [ n_path, 1 ]
aoa [5] (required)
Arrival azimuth angles in [rad], Size: [ 1, n_path ] or [ n_path, 1 ]
eoa [6] (required)
Arrival elevation angles in [rad], Size: [ 1, n_path ] or [ n_path, 1 ]
path_gain [7] (required)
Path gain (linear scale), Size: [ 1, n_path ] or [ n_path, 1 ]
path_length [8] (required)
Total path length in meters, Size: [ 1, n_path ] or [ n_path, 1 ]
M [9] (required)
Polarization transfer matrix, interleaved complex values (ReVV, ImVV, ReVH, ImVH, ReHV, ImHV, ReHH, ImHH),
Size: [ 8, n_path ]
tx_pos [10] (required)
Transmitter position in 3D Cartesian coordinates; Size: [3,1] or [1,3]
tx_orientation [11] (required)
3-element vector describing the orientation of the transmit antenna in Euler angles (bank, tilt, heading),
Size: [3,1] or [1,3]
rx_pos [12] (required)
Receiver position in 3D Cartesian coordinates, Size: [3,1] or [1,3]
rx_orientation [13 (required)]
3-element vector describing the orientation of the receive antenna in Euler angles,
Size: [3,1] or [1,3]
center_freq [14] (optional)
Center frequency in [Hz]; optional; If the value is not provided or set to 0, phase calculation
in coefficients is disabled, i.e. that path length has not influence on the results. This can be
used to calculate the antenna response for a specific angle and polarization. Scalar value
use_absolute_delays [15] (optional)
If true, the LOS delay is included for all paths; Default is false, i.e. delays are normalized
to the LOS delay.
add_fake_los_path [16] (optional)
If true, adds a zero-power LOS path as the first path in case where no LOS path was present.
Default: false
-
Derived inputs:
n_azimuth_tx |
Number of azimuth angles in the TX antenna pattern |
n_elevation_tx |
Number of elevation angles in the TX antenna pattern |
n_elements_tx |
Number of physical antenna elements in the TX array antenna |
n_ports_tx |
Number of ports (after coupling) in the TX array antenna |
n_azimuth_rx |
Number of azimuth angles in the RX antenna pattern |
n_elevation_rx |
Number of elevation angles in the RX antenna pattern |
n_elements_rx |
Number of physical antenna elements in the RX array antenna |
n_ports_rx |
Number of ports (after coupling) in the RX array antenna |
n_path |
Number of propagation paths |
-
Output Arguments:
coeff_re
Channel coefficients, real part, Size: [ n_ports_tx, n_ports_rx, n_path ]
coeff_im
Channel coefficients, imaginary part, Size: [ n_ports_tx, n_ports_rx, n_path ]
delays
Propagation delay in seconds, Size: [ n_ports_tx, n_ports_rx, n_path ]
rx_Doppler
Doppler weights for moving RX, Size: [ 1, n_path ]
-
Caveat:
- Input data is directly accessed from MATLAB / Octave memory, without copying if it is provided in
double precision.
- Other formats (e.g. single precision inputs) will be converted to double automatically, causing
additional computation steps.
- To improve performance of repeated computations (e.g. in loops), consider preparing the data
accordingly to avoid unecessary computations.
get_channels_spherical - Calculate channel coefficients from path data and antenna patterns
-
Description:
In this function, the wireless propagation channel between a transmitter and a receiver is calculated,
based on a single transmit and receive position. Additionally, interaction points with the environment,
which are derived from either Ray Tracing or Geometric Stochastic Models such as QuaDRiGa, are
considered. The calculation is performed under the assumption of spherical wave propagation. For accurate
execution of this process, several pieces of input data are required:
- The 3D Cartesian (local) coordinates of both the transmitter and the receiver.
- The specific interaction positions of the propagation paths within the environment.
- The polarization transfer matrix for each propagation path.
- Antenna models for both the transmitter and the receiver.
- The orientations of the antennas.
-
Usage:
[ coeff_re, coeff_im, delays, aod, eod, aoa, eoa ] = quadriga_lib.get_channels_spherical( ant_tx, ant_rx, ...
fbs_pos, lbs_pos, path_gain, path_length, M, tx_pos, tx_orientation, rx_pos, rx_orientation, ...
center_freq, use_absolute_delays, add_fake_los_path );
-
Input Arguments:
ant_tx [1] (required)
Struct containing the transmit (TX) arrayant data with the following fields:
e_theta_re |
e-theta field component, real part |
Size: [n_elevation_tx, n_azimuth_tx, n_elements_tx] |
e_theta_im |
e-theta field component, imaginary part |
Size: [n_elevation_tx, n_azimuth_tx, n_elements_tx] |
e_phi_re |
e-phi field component, real part |
Size: [n_elevation_tx, n_azimuth_tx, n_elements_tx] |
e_phi_im |
e-phi field component, imaginary part |
Size: [n_elevation_tx, n_azimuth_tx, n_elements_tx] |
azimuth_grid |
Azimuth angles in [rad], -pi to pi, sorted |
Size: [n_azimuth_tx] |
elevation_grid |
Elevation angles in [rad], -pi/2 to pi/2, sorted |
Size: [n_elevation_tx] |
element_pos |
Antenna element (x,y,z) positions, optional |
Size: [3, n_elements_tx] |
coupling_re |
Coupling matrix, real part, optional |
Size: [n_elements_tx, n_ports_tx] |
coupling_im |
Coupling matrix, imaginary part, optional |
Size: [n_elements_tx, n_ports_tx] |
ant_rx [2] (required)
Struct containing the receive (RX) arrayant data with the following fields:
e_theta_re |
e-theta field component, real part |
Size: [n_elevation_rx, n_azimuth_rx, n_elements_rx] |
e_theta_im |
e-theta field component, imaginary part |
Size: [n_elevation_rx, n_azimuth_rx, n_elements_rx] |
e_phi_re |
e-phi field component, real part |
Size: [n_elevation_rx, n_azimuth_rx, n_elements_rx] |
e_phi_im |
e-phi field component, imaginary part |
Size: [n_elevation_rx, n_azimuth_rx, n_elements_rx] |
azimuth_grid |
Azimuth angles in [rad], -pi to pi, sorted |
Size: [n_azimuth_rx] |
elevation_grid |
Elevation angles in [rad], -pi/2 to pi/2, sorted |
Size: [n_elevation_rx] |
element_pos |
Antenna element (x,y,z) positions, optional |
Size: [3, n_elements_rx] |
coupling_re |
Coupling matrix, real part, optional |
Size: [n_elements_rx, n_ports_rx] |
coupling_im |
Coupling matrix, imaginary part, optional |
Size: [n_elements_rx, n_ports_rx] |
fbs_pos [3] (required)
First interaction point of the rays and the environment, Size: [ 3, n_path ]
lbs_pos [4] (required)
Last interaction point of the rays and the environment; For single-bounce models, this must be
identical to fbs_pos, Size: [ 3, n_path ]
path_gain [5] (required)
Path gain (linear scale), Size: [ 1, n_path ] or [ n_path, 1 ]
path_length [6] (required)
Total path length in meters; If path_length is shorter than the shortest possible path from TX to
FBS to LBS to RX, it is replaced by the shortest path length, Size: [ 1, n_path ] or [ n_path, 1 ]
M [7] (required)
Polarization transfer matrix; interleaved complex values (ReVV, ImVV, ReVH, ImVH, ReHV, ImHV, ReHH, ImHH),
Size: [ 8, n_path ]
tx_pos [8] (required)
Transmitter position in 3D Cartesian coordinates, Size: [3,1] or [1,3]
tx_orientation [9] (required)
3-element vector describing the orientation of the transmit antenna in Euler angles (bank, tilt, heading),
Size: [3,1] or [1,3]
rx_pos [10] (required)
Receiver position in 3D Cartesian coordinates, Size: [3,1] or [1,3]
rx_orientation [11] (required)
3-element vector describing the orientation of the receive antenna, Size: [3,1] or [1,3]
center_freq [12] (optional)
Center frequency in [Hz]; optional; If the value is not provided or set to 0, phase calculation
in coefficients is disabled, i.e. that path length has not influence on the results. This can be
used to calculate the antenna response for a specific angle and polarization. Scalar value
use_absolute_delays [13] (optional)
If true, the LOS delay is included for all paths; Default is false, i.e. delays are normalized
to the LOS delay.
add_fake_los_path [14] (optional)
If true, adds a zero-power LOS path as the first path in case where no LOS path was present.
Default: false
-
Derived inputs:
n_azimuth_tx |
Number of azimuth angles in the TX antenna pattern |
n_elevation_tx |
Number of elevation angles in the TX antenna pattern |
n_elements_tx |
Number of physical antenna elements in the TX array antenna |
n_ports_tx |
Number of ports (after coupling) in the TX array antenna |
n_azimuth_rx |
Number of azimuth angles in the RX antenna pattern |
n_elevation_rx |
Number of elevation angles in the RX antenna pattern |
n_elements_rx |
Number of physical antenna elements in the RX array antenna |
n_ports_rx |
Number of ports (after coupling) in the RX array antenna |
n_path |
Number of propagation paths |
-
Output Arguments:
coeff_re
Channel coefficients, real part, Size: [ n_ports_rx, n_ports_tx, n_path ]
coeff_im
Channel coefficients, imaginary part, Size: [ n_ports_rx, n_ports_tx, n_path ]
delays
Propagation delay in seconds, Size: [ n_ports_rx, n_ports_tx, n_path ]
aod (optional)
Azimuth of Departure angles in [rad], Size: [ n_ports_rx, n_ports_tx, n_path ]
eod (optional)
Elevation of Departure angles in [rad], Size: [ n_ports_rx, n_ports_tx, n_path ]
aoa (optional)
Azimuth of Arrival angles in [rad], Size: [ n_ports_rx, n_ports_tx, n_path ]
eoa (optional)
Elevation of Arrival angles in [rad], Size: [ n_ports_rx, n_ports_tx, n_path ]
-
Caveat:
- Input data is directly accessed from MATLAB / Octave memory, without copying if it is provided in
double precision.
- Other formats (e.g. single precision inputs) will be converted to double automatically, causing
additional computation steps.
- To improve performance of repeated computations (e.g. in loops), consider preparing the data
accordingly to avoid unecessary computations.
Miscellaneous / Tools
calc_rotation_matrix - Calculates a 3x3 rotation matrix from a 3-element orientation vector
-
Description:
In linear algebra, a rotation matrix is a transformation matrix that is used to perform a rotation
in Euclidean space. The rotation of a rigid body (or three-dimensional coordinate system with a
fixed origin) is described by a single rotation about some axis. Such a rotation may be uniquely
described by three real-valued parameters. The idea behind Euler rotations is to split the complete
rotation of the coordinate system into three simpler constitutive rotation. This function calculates
the 3x3 rotation matrix R from the intrinsic Euler angles.
-
Usage:
rotation = quadriga_lib.calc_rotation_matrix( orientation, invert_y_axis, transpose )
-
Example:
The following example obtains the 3x3 matrix R for a 45 degree rotation around the z-axis:
bank = 0;
tilt = 0;
heading = 45 * pi/180;
orientation = [ bank; tilt; heading ];
rotation = quadriga_lib.calc_rotation_matrix( orientation );
R = reshape( rotation, 3, 3 );
-
Input Arguments:
orientation
Euler angles (bank, tilt, head) in [rad]; Single or double precision, Size: [3, n_row, n_col]
invert_y_axis
Optional parameter. If set to 1, the rotation around the y-axis is inverted.
transpose
Optional parameter. If set to 1, the output is transposed.
-
Output Argument:
rotation
The rotation matrix, i.e. a transformation matrix that is used to perform a rotation in 3D
Euclidean space; Size: [9, n_row, n_col]
cart2geo - Transform Cartesian (x,y,z) coordinates to Geographic (az, el, length) coordinates
fast_sincos - Fast, approximate sine/cosine for MATLAB numeric arrays
-
Description:
Computes elementwise sine and/or cosine for input angles in radians.
- Works on vectors, matrices, and 3-D arrays
- Accepts any numeric input class; best performance with single precision
- Outputs are always single precision
- Results are approximate and may differ from MATLAB
sin / cos
- 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)
- Request one or two outputs to control which results are returned
- With one output, set the optional
cosineOnly flag to true to return cosine instead of sine
-
Usage:
[s, c] = arrayant_lib.fast_sincos(x);
s = arrayant_lib.fast_sincos(x);
c = arrayant_lib.fast_sincos(x, true);
-
Input Arguments:
x (input)
Numeric array of angles in radians. Any size/shape.
cosineOnly = false (optional input)
Logical scalar. When requesting a single output, set to true to return cos(x); otherwise returns
sin(x).
-
Output Arguments:
s
Single-precision sin(x). Same size as x.
c
Single-precision cos(x). Same size as x.
-
Examples:
% Input as single for best performance
x = single(linspace(0, 2*pi, 1000));
% Compute sine and cosine
[s, c] = arrayant_lib.fast_sincos(x);
% Compute only sine (single output)
s = arrayant_lib.fast_sincos(x);
% Compute only cosine (single output with flag)
c = arrayant_lib.fast_sincos(x, true);
% Double input is accepted; outputs remain single
xd = linspace(0, 2*pi, 8);
s_only = arrayant_lib.fast_sincos(xd); % class(s_only) == 'single'
c_only = arrayant_lib.fast_sincos(xd, true); % class(c_only) == 'single'
geo2cart - Transform Geographic (az, el, length) to Cartesian (x,y,z) coordinates coordinates
version - Returns the quadriga-lib version number
interp - 2D and 1D linear interpolation.
-
Description:
This function implements 2D and 1D linear interpolation.
-
Usage:
dataI = quadriga_lib.interp( x, y, data, xI, yI ); % 2D case
dataI = quadriga_lib.interp( x, [], data, xI ); % 1D case
-
Input Arguments:
x
Vector of sample points in x direction for which data is provided; single or double; Length: [nx]
y
Vector of sample points in y direction for which data is provided; single or double; Length: [ny]
Must be an empty array [] in case of 1D interpolation.
data
The input data tensor; single or double; Size: [ny, nx, ne] or [1, nx, ne] for 1D case
The 3rd dimension enables interpolation for mutiple datasets simultaneously.
xI
Vector of sample points in x direction for which data should be interpolated; single or double;
Length: [nxI]
yI
Vector of sample points in y direction for which data should be interpolated; single or double;
Length: [nyI]
-
Output Arguments:
dataI
The interpolated dat; single or double (same as data);
Size: [nyI, nxI, ne] or [1, nxI, ne] for 1D case
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
- Uses the LodePNG library for PNG writing
-
Declaration:
quadriga_lib.write_png( fn, data, colormap, min_val, max_val, log_transform )
-
Arguments:
fn
Filename of the PNG file, string, required
data
Data matrix, required, size [N, M]
colormap (optional)
Colormap for the visualization, string, supported are 'jet', 'parula', 'winter', 'hot', 'turbo',
'copper', 'spring', 'cool', 'gray', 'autumn', 'summer', optional, default = 'jet'
min_val (optional)
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.
max_val (optional)
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.
log_transform (optional)
If enabled, the data values are transformed to the log-domain (10*log10(data)) before processing.
Default: false (disabled)
Site-Specific Simulation Tools
generate_diffraction_paths - Generate propagation paths for estimating the diffraction gain
-
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 generates the elliptic propagation paths and corresponding weights
necessary for this calculation.
-
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
-
Usage:
[ rays, weights ] = quadriga_lib.generate_diffraction_paths( orig, dest, center_frequency, lod );
-
Input Arguments:
orig
Origin point of the propagation ellipsoid (e.g. transmitter positions). Size: [ n_pos, 3 ]
dest
Destination point of the propagation ellipsoid (e.g. receiver positions). Size: [ n_pos, 3 ]
center_freq
The center frequency in [Hz], scalar, default = 299792458 Hz
lod
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) |
-
Output Arguments:
rays
Coordinates of the generated rays; Size: [ n_pos, n_path, n_seg-1, 3 ]
weights
Weights; Size: [ n_pos, n_path, n_seg ]
icosphere - Construct a geodesic polyhedron (icosphere), a convex polyhedron made from triangles
-
Description:
An icosphere is constructed by subdividing faces of an icosahedron, a polyhedron with 20 faces,
12 vertices and 30 edges, and then projecting the new vertices onto the surface of a sphere. The
resulting mesh has 6 triangles at each vertex, except for 12 vertices which have 5 triangles.
The approximate equilateral triangles have roughly the same edge length and surface area.
-
Usage:
[ center, length, vert, direction ] = quadriga_lib.icosphere( no_div, radius, direction_xyz );
-
Input Arguments:
no_div
Number of divisions per edge of the generating icosahedron. The resulting number of faces is
equal to no_face = 20 · no_div^2
radius
Radius of the sphere in meters
direction_xyz
Direction format indicator: 0 = Spherical (default), 1 = Cartesian
-
Output Arguments:
center
Position of the center point of each triangle; Size: [ no_face, 3 ]
length
Length of the vector pointing from the origin to the center point. This number is smaller than
1 since the triangles are located inside the unit sphere;
Size: [ no_face ]
vert
The 3 vectors pointing from the center point to the vertices of each triangle; the values are
in the order [ v1x, v1y, v1z, v2x, v2y, v2z, v3x, v3y, v3z ];
Size: [ no_face, 9 ]
direction
The directions of the vertex-rays. If the format indicator direction_xyz is set to 0, the
output is in geographic coordinates (azimuth and elevation angle in rad); the values are in the
order [ v1az, v1el, v2az, v2el, v3az, v3el ];Size: [ no_face, 6 ] If the format indicator
direction_xyz is set to 1, the output is in Cartesian coordinates and the values are in the
order [ v1x, v1y, v1z, v2x, v2y, v2z, v3x, v3y, v3z ]; Size: [ no_face, 9 ]
obj_file_read - Reads a triangulated 3D polygon mesh from a Wavefront OBJ file
-
Description:
The function imports a polygon mesh from an OBJ file. The OBJ file format is a straightforward data
format, representing 3D geometry. It details the position of each vertex and defines polygons as lists
of these vertices. By default, vertices are arranged in a counter-clockwise order, eliminating the
need for explicit declaration of face normals. When exporting the mesh from software like Blender,
it's essential to triangulate the mesh and include material definitions. If the material name
exists in the material database, the function loads the corresponding properties.
Otherwise, it defaults to using standard properties.
-
Usage:
[ mesh, mtl_prop, vert_list, face_ind, obj_ind, mtl_ind, obj_names, mtl_names, bsdf ] = ...
quadriga_lib.obj_file_read( fn );
-
Input Arguments:
fn
Filename of the OBJ file, string
-
Output Arguments:
mesh
Vertices of the triangular mesh in global Cartesian coordinates. Each face is described by 3 points
in 3D-space. Hence, a face has 9 values in the order [ v1x, v1y, v1z, v2x, v2y, v2z, v3x, v3y, v3z ];
Size: [ no_mesh, 9 ]
mtl_prop
Material properties of each mesh element; If no material is defined for an object, the properties
for vacuum are used. Size: [ no_mesh, 5 ]
vert_list
List of vertices found in the OBJ file; Size: [ no_vert, 3 ]
face_ind
Triangular faces are defined by three vertices. Vertex indices match the corresponding vertex elements
of the previously defined vert_list (1-based).
uint32; Size: [ no_mesh, 3 ]
obj_id
Mesh elements in the OBJ file can be grouped into objects (e.g. 12 triangles define the walls of a
cube). Each object is identified by a unique ID (1-based).
uint32; Size: [ no_mesh, 1 ]
mtl_id
Each mesh element gets assigned a material and each unique material gets assigned an ID. Different
faces of an object can have different materials. If no material is defined in the OBJ file, the
id is set to 0 and no entry is made in mtl_names.
uint32; Size: [ no_mesh, 1 ]
obj_names
Names of the objects in the OBJ file; Cell array of strings
mtl_names
Names of the materials in the OBJ file; Cell array of strings
bsdf
Principled BSDF (Bidirectional Scattering Distribution Function) values extracted from the
.MTL file. Size [mtl_names.size(), 17]. Columns are:
| 1 |
Base Color Red |
Range 0-1 |
Default = 0.8 |
| 2 |
Base Color Green |
Range 0-1 |
Default = 0.8 |
| 3 |
Base Color Blue |
Range 0-1 |
Default = 0.8 |
| 4 |
Transparency (alpha) |
Range 0-1 |
Default = 1.0 (fully opaque) |
| 5 |
Roughness |
Range 0-1 |
Default = 0.5 |
| 6 |
Metallic |
Range 0-1 |
Default = 0.0 |
| 7 |
Index of refraction (IOR) |
Range 0-4 |
Default = 1.45 |
| 8 |
Specular Adjustment to the IOR |
Range 0-1 |
Default = 0.5 (no adjustment) |
| 9 |
Emission Color Red |
Range 0-1 |
Default = 0.0 |
| 10 |
Emission Color Green |
Range 0-1 |
Default = 0.0 |
| 11 |
Emission Color Blue |
Range 0-1 |
Default = 0.0 |
| 12 |
Sheen |
Range 0-1 |
Default = 0.0 |
| 13 |
Clearcoat |
Range 0-1 |
Default = 0.0 |
| 14 |
Clearcoat roughness |
Range 0-1 |
Default = 0.0 |
| 15 |
Anisotropic |
Range 0-1 |
Default = 0.0 |
| 16 |
Anisotropic rotation |
Range 0-1 |
Default = 0.0 |
| 17 |
Transmission |
Range 0-1 |
Default = 0.0 |
-
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.
In addition, custom properties can be set by assigning adding the 5 properties after the material
name, separated by :, e.g.:
usemtl custom::2.1:0.1:0.1:0.5:20
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 |
|
point_cloud_aabb - Calculates the axis-aligned bounding box (AABB) for a 3D point cloud or a set of sub-clouds
point_cloud_segmentation - Rearranges elements of a point cloud into smaller sub-clouds
-
Description:
This function processes a large 3D point cloud by clustering closely spaced points and recursively
partitioning it into smaller sub-clouds, each below a specified size threshold. It minimizes the
axis-aligned bounding box of each sub-cloud while striving to maintain a target number of points
per cluster.
Sub-clouds are aligned to a specified SIMD vector size (e.g., for AVX or CUDA), with padding applied
as needed. The function outputs a reorganized version of the input points (pointsR), where points
are grouped by sub-cloud, and provides forward and reverse index maps to track the reordering. This
organization is particularly useful for optimizing spatial processing tasks such as bounding volume
hierarchies or GPU batch execution.
-
Usage:
[ points_out, sub_cloud_index, forward_index, reverse_index ] = ...
quadriga_lib.point_cloud_segmentation( points_in, target_size, vec_size );
-
Input Arguments:
points_in
Points in 3D-Cartesian space; Size: [ n_points_in, 3 ]
target_size (optional)
The target number of elements of each sub-cloud. Default value = 1024. For best performance, the
value should be around 10 * sgrt( n_points_in )
vec_size (optional)
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-cloud 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-cloud.
-
Output Arguments:
points_out
Points in 3D-Cartesian space; singe or double precision; Size: [ n_points_out, 3 ]
sub_cloud_index
Start indices of the sub-clouds in 0-based notation. Type: uint32; Vector of length [ n_sub_cloud ]
forward_index
Indices for mapping elements of "points_in" to "points_out"; 1-based;
Length: [ n_points_out ]; For vec_size > 1, the added elements not contained in the input
are indicated by zeros.
reverse_index
Indices for mapping elements of "points_out" to "points_in"; 1-based; Length: [ n_points_in ]
point_inside_mesh - Test whether 3D points are inside a triangle mesh using raycasting
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.
- 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.
-
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 (σ). 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
-
Usage:
[ origN, destN, gainN, xprmatN, trivecN, tridirN, orig_lengthN, fbs_angleN, thicknessN, edge_lengthN, normal_vecN] = ...
quadriga_lib.ray_mesh_interact( interaction_type, center_freq, orig, dest, fbs, sbs, mesh, mtl_prop,
fbs_ind, sbs_ind, trivec, tridir, orig_length );
-
Input Arguments:
interaction_type
Interaction type: (0) Reflection, (1) Transmission, (2) Refraction
center_freq
Center frequency in [Hz]; Scalar value
orig
Ray origins in 3D Cartesian coordinates; Size: [ no_ray, 3 ]
dest
Ray destinations in 3D Cartesian coordinates; Size: [ no_ray, 3 ]
fbs
First interaction point between the rays and the triangular mesh. Size: [ no_ray, 3 ]
sbs
Second interaction point between the rays and the triangular mesh. Size: [ no_ray, 3 ]
mesh
Vertices of the triangular mesh in global Cartesian coordinates. Each face is described by 3 points
in 3D-space. Hence, a face has 9 values in the order [ v1x, v1y, v1z, v2x, v2y, v2z, v3x, v3y, v3z ];
Size: [ no_mesh, 9 ]
mtl_prop
Material properties of each mesh element (see above); Size: [ no_mesh, 5 ]
fbs_ind
Index of the triangle that was hit by the ray at the FBS location; 1-based; Length: [ no_ray ]
sbs_ind
Index of the triangle that was hit by the ray at the SBS location; 1-based; Length: [ no_ray ]
trivec (optional)
The 3 vectors pointing from the center point of the ray at the origin to the vertices of a triangular
propagation tube, the values are in the order [ v1x, v1y, v1z, v2x, v2y, v2z, v3x, v3y, v3z ];
Size: [ no_ray, 9 ]
tridir (optional)
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 ]
orig_length (optional)
Path length at origin point, default is 0, Size: [ n_ray ]
-
Output Arguments:
origN
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 ]
destN
New ray destinaion after the interaction with the medium, taking the change of direction into account;
Size: [ no_rayN, 3 ]
gainN
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. Linear scale.
Size: [ no_rayN ]
xprmatN
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 ]
trivecN
The 3 vectors pointing from the center point of the ray at origN to the vertices of a triangular
propagation tube. Size: [ no_rayN, 9 ]
tridirN
The directions of the vertex-rays after interaction with the medium;
Size: [ no_rayN, 6 ] or [ no_rayN, 9 ] depending on input size
orig_lengthN
Length of the ray from orig to origN. If orig_length is given as input, its value is added.
Size: [ no_rayN ]
fbs_angleN
Angle between incoming ray and FBS in [rad], Size [ n_rayN ]
thicknessN
Material thickness in meters calculated from the difference between FBS and SBS, Size [ n_rayN ]
edge_lengthN
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 ]
normal_vecN
Normal vector of FBS and SBS, Size [ n_rayN, 6 ]
out_typeN
A numeric indicator describing the type of the interaction. The total refection indicator is only
set in refraction mode.
| 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 |
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.
Ray beams are determined by an origin point, three vectors pointing from the origin to the three
vertices of a triangle that defines the shape of the tube and the three direction of the rays at
the vertices.
-
Usage:
[ hit_count, ray_ind ] = quadriga_lib.ray_point_intersect( orig, trivec, tridir, points, ...
max_no_hit, sub_cloud_index, target_size );
-
Input Arguments:
orig
Ray origins in 3D Cartesian coordinates; Size: [ no_ray, 3 ]
trivec
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 ]
tridir
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.
points
Points in 3D-Cartesian space; Size: [ n_points_in, 3 ]
max_no_hit (optional)
Max. number of hits in the output ray_ind, default = 32
sub_cloud_index (optional)
Start indices of the sub-clouds in 0-based notation. Type: uint32; Vector of length [ n_sub_cloud ]
If this optional input is not given, the sub-could index is calculated automatically. Passing a
value of uint32(0) will disable the sub-cloud calculation.
target_size (optional)
Target value for the sub-cloud size, only evaluated if 'sub_cloud_index' is not given.
-
Output Arguments:
hit_count
Number of rays that hit a point, unit32, Length: [ n_points ]
ray_ind
Ray indices that hit the points, 1-based, 0 = no hit, Size [ n_points, max_no_hit ]
ray_triangle_intersect - Calculates the intersection of rays and triangles in three dimensions
-
Description:
- This function implements the Möller–Trumbore ray-triangle intersection algorithm, known for its
efficiency in calculating the intersection of a ray and a triangle in three-dimensional space.
This method achieves its speed by eliminating the need for precomputed plane equations of the plane
containing the triangle.
- The algorithm defines the ray using two points: an origin and a destination. Similarly, the triangle
is specified by its three vertices.
- To enhance performance, this implementation leverages AVX2 intrinsic functions and OpenMP, when
available, to speed up the computational process.
-
Usage:
[ fbs, sbs, no_interact, fbs_ind, sbs_ind ] = quadriga_lib.ray_triangle_intersect( orig, dest, mesh, sub_mesh_index );
-
Input Arguments:
orig
Ray origins in 3D Cartesian coordinates; Size: [ no_ray, 3 ]
dest
Ray destinations in 3D Cartesian coordinates; Size: [ no_ray, 3 ]
mesh
Vertices of the triangular mesh in global Cartesian coordinates. Each face is described by 3 points
in 3D-space. Hence, a face has 9 values in the order [ v1x, v1y, v1z, v2x, v2y, v2z, v3x, v3y, v3z ];
Size: [ no_mesh, 9 ]
sub_mesh_index (optional)
Start indices of the sub-meshes in 0-based notation. If this parameter is not given, intersections
are calculated for each mesh element, leading to poor performance for large meshes.
Type: uint32; Vector of length [ n_sub_mesh ]
-
Output Arguments:
fbs
First interaction point between the rays and the triangular mesh. If no interaction was found, the
FBS location is equal to dest. Size: [ no_ray, 3 ]
sbs
Second interaction point between the rays and the triangular mesh. If no interaction was found, the
SBS location is equal to dest. Size: [ no_ray, 3 ]
no_interact
Total number of interactions between the origin point and the destination; uint32; Length: [ no_ray ]
fbs_ind
Index of the triangle that was hit by the ray at the FBS location; 1-based; uint32; Length: [ no_ray ]
sbs_ind
Index of the triangle that was hit by the ray at the SBS location; 1-based; uint32; Length: [ no_ray ]
-
Caveat:
orig, dest, and mesh can be provided in single or double precision; fbs and lbs will have
the same type.
- All internal computation are done in single precision to achieve an additional 2x improvement in
speed compared to double precision when using AVX2 intrinsic instructions
subdivide_triangles - Subdivide the faces of a triangle mesh into smaller faces
-
Description:
This function splits the triangles of a mesh into smaller triangles by subdividing the edges
into n_div sub-edges. This creates n_div^2 sub-faces per face.
-
Usage:
triangles_out = quadriga_lib.subdivide_triangles( triangles_in, no_div );
[ triangles_out, mtl_prop_out ] = quadriga_lib.subdivide_triangles( triangles_in, no_div, mtl_prop_in );
-
Input Arguments:
triangles_in
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 ]; single or double precision;
Size: [ n_triangles_in, 9 ]
no_div
Number of divisions per edge of the input mesh. The resulting number of faces is equal to
n_triangles_out = n_triangles_in * n_div^2
mtl_prop_in (optional)
Material properties of each mesh element; Size: [ n_triangles_in, 5 ]
-
Output Arguments:
triangles_out
Vertices of the sub-divided mesh in global Cartesian coordinates; singe or double precision;
Size: [ n_triangles_out, 9 ]
mtl_prop_out
Material properties for the sub-divided triangle mesh elements. The values for the new faces are
copied from mtl_prop_in; Size: [ n_triangles_out, 5 ]
triangle_mesh_aabb - Calculate the axis-aligned bounding box (AABB) of a triangle mesh and its sub-meshes
triangle_mesh_segmentation - Rearranges elements of a triangle mesh into smaller sub-meshes
-
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.
This approach is particularly useful in computer graphics and simulation applications where
managing computational resources efficiently is crucial. By organizing the mesh elements into
compact clusters, the function enhances rendering performance and accelerates computational
tasks, such as collision detection and physics simulations. It allows for quicker processing
and reduced memory usage, making it an essential technique in both real-time graphics rendering
and complex simulation environments.
-
Usage:
[ triangles_out, sub_mesh_index, mesh_index ] = quadriga_lib.triangle_mesh_segmentation( ...
triangles_in, target_size, vec_size );
[ triangles_out, sub_mesh_index, mesh_index, mtl_prop_out ] = ...
quadriga_lib.triangle_mesh_segmentation( triangles_in, target_size, vec_size, mtl_prop_in );
-
Input Arguments:
triangles_in
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 ]; single or double precision;
Size: [ n_triangles_in, 9 ]
target_size (optional)
The target number of elements of each sub-mesh. Default value = 1024. For best performance, the
value should be around sgrt( n_triangles_in )
vec_size (optional)
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.
mtl_prop_in (optional)
Material properties of each mesh element; Size: [ n_triangles_in, 5 ]
-
Output Arguments:
triangles_out
Vertices of the clustered mesh in global Cartesian coordinates; singe or double precision;
Size: [ n_triangles_out, 9 ]
sub_mesh_index
Start indices of the sub-meshes in 0-based notation. Type: uint32; Vector of length [ n_sub_mesh ]
mesh_index
Indices for mapping elements of "triangles_in" to "triangles_out"; 1-based;
Length: [ n_triangles_out ]; For vec_size > 1, the added elements not contained in the input
are indicated by zeros.
mtl_prop_out
Material properties for the sub-divided triangle mesh elements. The values for the new faces are
copied from mtl_prop_in; Size: [ n_triangles_out, 5 ]; For vec_size > 1, the added elements
will contain the vacuum / air material.