Reference

The Parameter Class

class precon.Parameter(file: Path)

This class contains all parameters and labels of the raw data.

Depending on the input file, the parameters are read from different sources:

1) .raw file:

1a): if the .raw file was acquired with the ReconFrame patch, the parameters are read from it

1b): if it was acquired without patch a .sin file is searched in the same folder and if found the parameters are read from it.

2) .sin file: the parameters are read from the .sin file

3) .lab file: a .raw file is searched in the same folder and then continued as in 1)

4) .rc file: the parameters are read from the .rc file (.rc and .sin contain the same content)

5) .cpx file: a .rc file is searched in the same folder and then continued as in 4)

Variables:
  • rawfile (Path) – The path to the current raw file

  • labfile (Path) – The path to the current lab file

  • sinfile (Path) – The path to the current sin file

  • rcfile (Path) – The path to the current rc file

  • cpxfile (Path) – The path to the current cpx file

  • goal (SinPars) – A GoalPars class which holds all the Philips Goal parameters

  • sin – A SinPars class which holds all the Philips sin parameters

  • scanner_labels (List[LabelStructStandardC]) – The unprocessed scanner labels as stored in the labfile

  • labels (List[LabelStructC]) – The processed and formatted scanner labels which are used in the reconstruction

  • coil_info (CoilInfo) – Holds information about the connected coil arrays

Examples

import precon as pr
pars = pr.Parameter(Path('/data/my_rawfile.raw'))
get_concom_factors()

Returns the concomitant field correction factors.

Returns:

A list of concomitant field correction factors. If the list is empty or has more than 383 elements, it is truncated to 384 elements.

Return type:

list

Examples

concom_factors = pars.get_concom_factors()
get_coordinate_system(system: int = 1)

Returns the coordinate system of the MPS system expressed in RL, AP, FH axes.

Parameters:

system (int, optional) – The source coordinate system (1=MPS, 2=RL-AP-FH, 3=XYZ). Only MPS is implemented at the moment.

Returns:

The coordinate system in the specified format.

Return type:

list

Raises:

NotImplementedError – If the source coordinate system is not Enums.MPS.

Examples

coord_system = pars.get_coordinate_system()
get_geo_corr_pars()

Returns the geometric correction parameters.

Returns:

A tuple containing the following parameters:
  • r (float or None): The reference radius.

  • gys (ndarray or None): The field S coefficients for the y-axis.

  • gxc (ndarray or None): The field C coefficients for the x-axis.

  • gz (ndarray or None): The field coefficients for the z-axis.

Return type:

tuple

Examples

r, gys, gxc, gz = pars.get_geo_corr_pars()
get_in_plane_transformation(mix=0, stack=0)

Get the in-plane transformation for a given mix and stack. The in-plane transformation is used to rotate the imaged into the radiological convention.

Parameters:
  • mix (int) – The mix index. Defaults to 0.

  • stack (int) – The stack index. Defaults to 0.

Returns:

The in-plane transformation for the given mix and stack, or None if the transformation is not available.

Examples

in_plane_trafo = pars.get_in_plane_transformation()
get_loca_per_stack(stack=0)

Get the location indices for a given stack.

Parameters:

stack (int) – The index of the stack. Default is 0.

Returns:

A list of locations numbers

Return type:

list or None

Examples

locas = pars.get_loca_per_stack(0)
get_location_spectrum_offset(loca=0, stack=0)

Get the spectrum offset for a given location and stack.

Parameters:
  • enc (int) – The encoding direction (0=readout, 1=phase encoding, 2=slice encoding). Default is 0.

  • loca (int) – The location index. Default is 0.

  • stack (int) – The stack index. Default is 0.

Returns:

The location spectrum offset if it exists, otherwise None.

Return type:

int or None

Examples

spectrum_offset = pars.get_location_spectrum_offset(loca=0)
get_nr_phases()

Returns the number of cardiac phases.

Returns:

The number of cardiac phases.

Return type:

int

Examples

nr_phases = pars.get_nr_phases()
get_nus_enc_nrs()

Get the non-uniform sampling coordinates for an EPI scan. Every k-space profile must be interpolated on these coordinates to make sure the samples lay on a rectangular grid.

Returns:

A list of NUS sampling coordinates.

Return type:

list

Examples

nus_enc_nrs = pars.get_nus_enc_nrs()
get_oversampling(enc=0, mix=0, echo=0)

Get the oversampling factor for a specific encoding direction, mix and echo.

Parameters:
  • enc (int) – The encoding direction (0=readout, 1=phase encoding, 2=slice encoding). Default is 0.

  • mix (int) – The mix index. Default is 0.

  • echo (int) – The echo index. Default is 0.

Returns:

The oversampling factor.

Return type:

float

Examples

xovs = pars.get_oversampling(enc=0)
yovs = pars.get_oversampling(enc=1)
zovs = pars.get_oversampling(enc=2)
get_range(space=0, enc=0, mix=0, echo=0, stack=0, ovs=True, symmetric=False)

Get the sampling range of the current scan in k-space or image-space for a given encoding direction, mix, echo and stack.

Parameters:
  • space (str, optional) – The space in which the range is calculated (0=k-space, 1=image-space). Defaults to k-space.

  • enc (int, optional) – The encoding direction (0=readout, 1=phase encoding, 2=slice encoding). Default is 0.

  • mix (int, optional) – The mix index. Defaults to 0.

  • echo (int, optional) – The echo index. Defaults to 0.

  • stack (int, optional) – The stack index. Defaults to 0.

  • ovs (bool, optional) – Whether to include oversampling in the returned range. Defaults to True.

Returns:

The range of encoding numbers as a tuple (min_enc_number, max_enc_number) if the conditions are met, otherwise None.

Return type:

tuple or None

Examples

kx_range = pars.get_range(enc=0, ovs=False)
ky_range = pars.get_range(enc=1)
kz_range = pars.get_range(enc=2)
get_recon_resolution(mix=0, xovs=False, yovs=False, zovs=False, folded=False)

Returns the reconstruction resolutions (the resolution of the final reconstructed image).

Parameters:
  • mix (int) – The mix index. Default is 0.

  • xovs (bool) – If True returns the resolution including oversampling in readout direction. Default is False.

  • yovs (bool) – If True returns the resolution including oversampling in phase encoding direction. Default is False.. Default is False.

  • zovs (bool) – If True returns the resolution including oversampling in slice encoding direction. Default is False.. Default is False.

  • folded (bool) – If True returns the resolution of the folded images (the recon resolutions are divided by the sense factor). Default is False.

Returns:

A list of three integers representing the reconstruction resolutions in the x, y, and z directions.

Return type:

list

Examples

cur_recon_resolution = pars.get_recon_resolution(mix=mix, xovs=False, yovs=True, zovs=True)
get_sampled_size(enc=0, mix=0, echo=0, stack=0, ovs=None)

Returns the size of the sampled k-space without zero-filling. If half-scan or partial echo was enabled then the “symmetric” size is returned which includes the non-sampled part of k-space.

Parameters:
  • enc (int) – The encoding direction (0=readout, 1=phase encoding, 2=slice encoding). Default is 0.

  • mix (int) – The mixing value (default: 0).

  • echo (int) – The echo value (default: 0).

  • stack (int) – The stack value (default: 0).

  • ovs (bool, optional) – Whether to include oversampling in the returned size. Defaults to True.

Note

This is the size of the k-space returned by the “sort” function when no “output_size” is given.

Returns:

The sampled size.

Return type:

int

Examples

k_space_size = (pars.get_sampled_size(enc=0, ovs=False), pars.get_sampled_size(enc=1), pars.get_sampled_size(enc=2))
get_shift(enc=0, mix=0, stack=0, echo=0)

Get the image-space shift which needs to be applied after the fourier transform for a given encoding direction, mix, echo and stack.

Parameters:
  • enc (int, optional) – The encoding direction (0=readout, 1=phase encoding, 2=slice encoding). Default is 0.

  • mix (int, optional) – The mixing value. Defaults to 0.

  • stack (int, optional) – The stack value. Defaults to 0.

  • echo (int, optional) – The echo value. Defaults to 0.

Note

For EPI or non-cartesian scan the in-plane shifts need to be done in reconstruction. Also often the raw data is shifted by half a FOV. This needs to be corrected for the image to be correctly placed.

Returns:

The shift value, or None if the required parameters are missing.

Return type:

int or None

Examples

yshift = pars.get_shift(enc=1)
zshift = pars.get_shift(enc=2)
get_transformation_matrix(loca: List | int = 0, mix: int = 0, source: int = 1, target: int = 2) _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes]

Get the transformation matrix to transform between the different coordinate systems of the scanner (MPS, RL-AP-FH, XYZ)

Parameters:
  • loca (Union[List, int], optional) – The location(s) to get the transformation matrix for. Defaults to 0.

  • mix (int, optional) – The mix index. Defaults to 0.

  • source (int, optional) – The source system (1=MPS, 2=RL-AP-FH, 3=XYZ). Defaults to Enums.MPS.

  • target (int, optional) – The target system (1=MPS, 2=RL-AP-FH, 3=XYZ). Defaults to Enums.RAF.

Returns:

The transformation matrix.

Return type:

npt.ArrayLike

Raises:

RuntimeError – If the source system is not Enums.MPS or the target system is not Enums.RAF or Enums.XYZ.

Note

  • Only “Enums.MPS” is implemented as the source system right now.

  • Only “Enums.RAF” and “Enums.XYZ” are possible options for the target system.

Examples

locations = pr.utils.get_unique(labels, 'loca')
MPS_to_XYZ = pars.get_transformation_matrix(loca=locations, target=pr.Enums.XYZ)
get_value(name: str, default=None, at: int | None = None)

Get the value of a goal parameter.

Parameters:
  • name (str) – The name of the parameter.

  • default (Any, optional) – The default value to return if the parameter is not found. Defaults to None.

  • at (int, optional) – returns the value at the given index. Defaults to None.

Returns:

The value of the parameter.

Return type:

Any

Raises:

NotImplementedError – If no parameter dictionary is found.

Note

get_value always returns a list of values. If the parameter is not found then the specified default value is returned. If no default value is given None is returned.

Examples

flip_angle = pars.get_value('EX_ACQ_flip_angle', default=[90.0])
traj_type = goal.get_value('RC_k_space_traj_type', default=0, at=0)
get_venc(system: int = 1)

Get the VENC (Velocity Encoding) values for the given system.

Parameters:

system (int, optional) – The system to get the VENC values for (1=MPS, 2=RL-AP-FH). Defaults to Enums.MPS.

Returns:

A list of VENC values [vencRL, vencAP, vencFH].

Return type:

list

Raises:

RuntimeError – If the system is neither Enums.MPS nor Enums.RAF.

Examples

venc = pars.get_venc()
get_voxel_sizes(mix=0)

Get the voxel sizes of the current scan in mm.

Parameters:

mix (int) – The index of the mix. Defaults to 0.

Returns:

A list containing the voxel sizes [vx, vy, vz].

Return type:

list or None

Examples

vs = pars.get_voxel_sizes(mix=0)
is_3d(mix: int = 0)

Check if the current dataset is a 3D scan.

Parameters:

mix (int, optional) – The index of the mix. Defaults to 0.

Returns:

True if 3D, False otherwise.

Return type:

bool

Examples

is3D = pars.is_3d()
is_epi()

Checks if the current dataset is an EPI scan.

Returns:

True if it is an EPI scan, False otherwise.

Return type:

bool

Examples

isEPI = pars.is_epi()
is_hadamard_encoding()

Check if Hadamard encoding is used in the current flow scan.

Returns:

True if the encoding type is Hadamard, False otherwise.

Return type:

bool

Examples

is_hadamard = pars.is_hadamard_encoding
static make_symmetric(krange: List[int])

Make the given input range symmetric.

Parameters:

krange (List[int]) – A range.

Returns:

The symmetric range.

Return type:

List[int]

Examples

symmetric_range = pars.make_symmetric([-6, 11])      # will return [-12, 11]

The Parameter2Read Class

class precon.Parameter2Read(labels: List[LabelStructC] | None = None, cpx_headers: List[CpxHeader] | None = None)

This Parameter2Read class defines what data to read from the rawfile. It is an input argument in the “read” function

Variables:
  • typ (int) – The data type (1=imaging data, 2=rejected data, 3=echo phase correction data, 4=junk data, 5=noise data, 6=frequency correction data, 7=trajectory (for spiral))

  • mix (int) – The mix number

  • dyn (int) – The dynamic scan number

  • card (int) – The cardiac phase number

  • echo (int) – The echo number

  • loca (int) – The location number

  • stack (int) – The stack

  • extr1 (int) – The extra attribute (flow segment, b-value etc)

  • extr2 (int) – The second extra attribute (diffusion direction etc)

  • ky (int) – The phase encoding numbers

  • kz (int) – The slice encoding numbers

  • aver (int) – The average numbers

  • rtop (int) – The trigger delay time

Examples

import precon as pr
pars = pr.Parameter(Path('/data/my_rawfile.raw'))

# define what to read
parameter2read = pr.Parameter2Read(pars.labels)
parameter2read.dyn = [15, 32, 37]
parameter2read.stack = 0

# read data
with open(args.rawfile, 'rb') as raw:
    data, labels = pr.read(raw, parameter2read, pars.labels, pars.coil_info)

The Sensitivity Class

class precon.Sensitivity

Holds the sensitivity maps (and calibration data) which are calculated from a SENSE reference scan

Variables:
  • sensitivity (npt.ArrayLike) – The sensitivity maps

  • bodycoil (npt.ArrayLike) – The reformatted bodycoil images

  • surfacecoil (npt.ArrayLike) – The reformatted surfacecoil images

  • target_scan (Path) – The path to the target scan

  • refscan (Path) – The path to the sense reference scan

  • psi (npt.ArrayLike) – The noise covariance matrix

Examples

import precon as pr
sens = pr.calculate_sensitivities(Path('/data/refscan.raw'), Path('/data/target_scan.raw'), stack=0, mix=0)

The Enums Class

class precon.Enums

This class defines constants representing various enums used in the reconstruction.

DATA_DIMS

The total number of data dimensions (12).

Type:

int

KSPACE

Enumeration for k-space.

Type:

int

IMAGESPACE

Enumeration for image space.

Type:

int

IPTRANSF_PLUS_A_PLUS_B

Enumeration for in-plane transformation: +A +B.

Type:

int

IPTRANSF_PLUS_A_MIN_B

Enumeration for in-plane transformation: +A -B.

Type:

int

IPTRANSF_MIN_A_PLUS_B

Enumeration for in-plane transformation: -A +B.

Type:

int

IPTRANSF_MIN_A_MIN_B

Enumeration for in-plane transformation: -A -B.

Type:

int

IPTRANSF_PLUS_B_PLUS_A

Enumeration for in-plane transformation: +B +A.

Type:

int

IPTRANSF_PLUS_B_MIN_A

Enumeration for in-plane transformation: +B -A.

Type:

int

IPTRANSF_MIN_B_PLUS_A

Enumeration for in-plane transformation: -B +A.

Type:

int

IPTRANSF_MIN_B_MIN_A

Enumeration for in-plane transformation: -B -A.

Type:

int

X_DIM

Enumeration for X dimension.

Type:

int

Y_DIM

Enumeration for Y dimension.

Type:

int

Z_DIM

Enumeration for Z dimension.

Type:

int

COIL_DIM

Enumeration for coil dimension (also CHANNEL_DIM).

Type:

int

DYN_DIM

Enumeration for dynamic dimension.

Type:

int

HP_DIM

Enumeration for high-pass dimension (also PHASE_DIM).

Type:

int

ECHO_DIM

Enumeration for echo dimension.

Type:

int

LOCA_DIM

Enumeration for localization dimension.

Type:

int

MIX_DIM

Enumeration for mix dimension.

Type:

int

EXTR1_DIM

Enumeration for extra dimension 1.

Type:

int

FLOW_SEGMENT_DIM

Enumeration for flow segment dimension.

Type:

int

EXTR2_DIM

Enumeration for extra dimension 2.

Type:

int

AVER_DIM

Enumeration for average dimension.

Type:

int

NORMAL_DATA

Enumeration for normal data.

Type:

int

REJECTED_DATA

Enumeration for rejected data.

Type:

int

ECHO_PHASE_DATA

Enumeration for echo phase data.

Type:

int

JUNK_DATA

Enumeration for junk data.

Type:

int

NOISE_DATA

Enumeration for noise data.

Type:

int

FRC_DATA

Enumeration for FRC data.

Type:

int

TRAJECTORY

Enumeration for trajectory data.

Type:

int

MPS

Enumeration for MPS coordinate system.

Type:

int

RAF

Enumeration for RL-AP-FH coordinate system.

Type:

int

XYZ

Enumeration for XYZ coordinate system.

Type:

int

RL

Enumeration for right-left direction in patient coordinate system.

Type:

int

AP

Enumeration for anterior-posterior direction in patient coordinate system.

Type:

int

FH

Enumeration for foot-head direction in patient coordinate system.

Type:

int

REC_IMAGE_TYPE_M

Enumeration for Magnitude image type.

Type:

int

REC_IMAGE_TYPE_R

Enumeration for Real image type.

Type:

int

REC_IMAGE_TYPE_I

Enumeration for Imaginary image type.

Type:

int

REC_IMAGE_TYPE_P

Enumeration for Phase image type.

Type:

int

TRAJ_TYPE_CARTESIAN

Enumeration for cartesian trajectory type.

Type:

int

TRAJ_TYPE_RADIAL

Enumeration for radial trajectory type.

Type:

int

TRAJ_TYPE_SPIRAL

Enumeration for spiral trajectory type.

Type:

int

TRAJ_TYPE_PRPPELLER

Enumeration for propeller trajectory type.

Type:

int

Reconstruction Functions

calculate_sensitivities

precon.calculate_sensitivities(refscan: Path, target_scan: Path, stack: int = 0, mix: int = 0, match_target_size: bool = False, qbc_stack: int | None = None, fov: List[float] | None = None, output_size: List[int] | None = None) Sensitivity

Calculates the sensitivity maps from the SENSE reference scan in the geometry of the target scan.

Parameters:
  • refscan (Path) – The path to the rawfile of the SENSE reference scan

  • target_scan (Path) – The path to the rawfile of the target scan

  • stack (int, optional) – The stack index of the target scan which is used to determine the geometry (only relevant for multistack data). Defaults to 0.

  • mix (int, optional) – The mix index of the target scan which is used to determine the geometry (only relevant for mixed data). Defaults to 0.

  • match_target_size (bool, optional) – If True the size of the sensitivities is equal to the matrix size of the target scan. Defaults to False.

  • fov (List[float], optional) – The field of view of the reformatted data. Defaults to None. Only change if a different FOV than the target scan is needed (e.g. in simulations)

  • output_size (List[int], optional) – The matrix size of the reformatted data. Defaults to None. If not given then it is determined automatically. If match_target_size is True then this argumant is ignored

Returns:

The sensitivity object containing the reformatted SENSE reference scan.

Return type:

Sensitivity

Raises:
  • RuntimeError – If the “fov” does not contain exactly 3 elements.

  • RuntimeError – If the “output_size” does not contain exactly 3 elements.

Examples

sens = pr.calculate_sensitivities(Path(refscan), Path(target_scan), match_target_size=True)
# calculate the sensitivities  with a user defined FOV and matrix size
sens = pr.calculate_sensitivities(Path(refscan), Path(target_scan), fov=[320, 230, 50], output_size=[512, 320, 20])

compress_data

precon.compress_data(data: np.ArrayLike, A: np.ArrayLike) np.ArrayLike

Compresses the channels of the given data with the provided compression matrix A. The matrix can be calculated with the “get_array_compression_matrix” function. A more detailed description of the compression algorithm can be found in “Buehrer M, Pruessmann KP, Boesiger P, Kozerke S. Array compression for MRI with large coil arrays. Magn Reson Med 2007;57:1131–1139”

Parameters:
  • data (npt.ArrayLike) – The data to be compressed. The compression will be applied along the channel dimension (Enums.CHANNEL_DIM)

  • A (npt.ArrayLike) – The matrix used for compression. It must have the size: NR_VIRTUAL_CHANNELS x NR_CHANNELS

Returns:

The compressed data. It has NR_VIRTUAL_CHANNELS in the channel dimension

Return type:

npt.ArrayLike

compress_sensitivity

precon.compress_sensitivity(sensitivity: Sensitivity, A: np.ArrayLike) Sensitivity

Compresses the given sensitivity object with the provided compression matrix A. It will compress the member arrays: “sensitivity.sensitivity” and “sensitivity.surfacecoil”. This function should be called when using array compression together with SENSE A more detailed description of the compression algorithm can be found in “Buehrer M, Pruessmann KP, Boesiger P, Kozerke S. Array compression for MRI with large coil arrays. Magn Reson Med 2007;57:1131–1139”

Parameters:
  • sensitivity (Sensitivity) – The sensitivity object to compress.

  • A (npt.ArrayLike) – The matrix used for compression. It must have the size: NR_VIRTUAL_CHANNELS x NR_CHANNELS

Returns:

The compressed sensitivity object.

Return type:

Sensitivity

Examples

# calculate the sensitivities
sens = pr.calculate_sensitivities(Path(refscan), Path(target_scan))

# get the array compression matrix
A = pr.get_array_compression_matrix(sens.surfacecoil)

# crop the compression matrix to the number of virtual coils
virtual_coils = 8
A = A[0:args.virtual_coils, :]

# compress the sensitivities (for the SENSE recon)
sens = pr.compress_sensitivity(sens, A)

concomitant_field_correction

precon.concomitant_field_correction(data: np.ArrayLike, mps_to_xyz_transformation: list | np.ArrayLike, concom_factors: List[float], voxel_sizes: tuple, segments: List[float] | None = None) np.ArrayLike

Applies the concomitant field correction to flow data.

Parameters:
  • data (npt.ArrayLike) – The input data to be corrected.

  • mps_to_xyz_transformation (list) – The transformation matrix from the MPS to the XYZ system. It can be obtained with the “get_transformation_matrix” function which is a member of the Parameter class.

  • concom_factors (List[float]) – The concomitant field correction factors for each flow segment. It can be gotten with the “get_concom_factors” function which is a member of the Parameter class.

  • voxel_sizes (tuple) – The voxel sizes for each dimension.

  • segments (List[float], optional) – The flow segments of the current data. Defaults to None.

Returns:

The corrected data.

Return type:

npt.ArrayLike

Raises:

RuntimeError – If the number of flow segments in the data does not match the length of the segments list.

Examples

# flow segments are stored in the extr1 label
segments = parameter2read.extr1

# get the transformation matrices (MPS to XYZ) for every location.
locations = pr.utils.get_unique(labels, 'loca')
MPS_to_XYZ = pars.get_transformation_matrix(loca=locations, mix=mix, target=pr.Enums.XYZ)
voxel_sizes = pars.get_voxel_sizes(mix=mix)

# concommitant field correction (process every location separately)
concom_factors = pars.get_concom_factors()
data = pr.concomitant_field_correction(data, MPS_to_XYZ, concom_factors, voxel_sizes, segments)

crop

precon.crop(data: np.ArrayLike, axis: int | tuple = 0, factor: int | tuple = 1, n: int | tuple | None = None, where: str | int = 0) np.ArrayLike

Crop the input data along the specified axes.

Parameters:
  • data (npt.ArrayLike) – The input data to be cropped.

  • axis (Union[int, tuple], optional) – The axis or axes along which to crop the data. Defaults to 0.

  • factor (Union[int, tuple], optional) – The factor or factors by which to crop the data. Defaults to 1. If both “factor” and “n” are passed as arguments then “n” will be used.

  • n (Union[int, tuple], optional) – The number or numbers of elements to keep after cropping. Defaults to None. If both “factor” and “n” are passed as arguments then “n” will be used.

  • where (Union[str, int], optional) – Specifies where to crop. Options are: ‘symmetric’: the data will be cropped symmetrically on both sides; ‘left’: the data will be cropped on the left; ‘right’: the data will be cropped on the right

Returns:

The cropped data.

Return type:

npt.ArrayLike

Raises:
  • RuntimeError – If the ‘where’ option is unknown.

  • RuntimeError – If ‘axis’ and ‘factor’ do not have the same number of elements.

  • RuntimeError – If ‘axis’ and ‘n’ do not have the same number of elements.

Examples

data = pr.crop(data, axis=(1, 2), factor=(1.5, 2), where='symmetric')
data = pr.crop(data, axis=0, n=256, where='left')

divide_flow_segments

precon.divide_flow_segments(data: np.ArrayLike, hadamard=False) np.ArrayLike

Divides the flow segments of the input data (in dimension Enums.FLOW_SEGMENT_DIM)

Parameters:
  • data (npt.ArrayLike) – The input data.

  • hadamard (bool, optional) – Whether to use the Hadamard encoding scheme. Defaults to False.

Returns:

The resulting flow data.

Return type:

npt.ArrayLike

Examples

data = pr.divide_flow_segments(data, pars.is_hadamard_encoding())

epi_corr

precon.epi_corr(data: np.ArrayLike, labels: np.ArrayLike, slopes: List[float], offsets: List[float]) np.ArrayLike

Performs a linear EPI correction on the input data using the provided slopes and offsets for each channel. The slopes and offsets can be calculated with the “get_epi_corr_data” function. The input data must be fourier transformed in readout direction but not in phase-encoding direction

Parameters:
  • data (npt.ArrayLike) – The input data.

  • labels (npt.ArrayLike) – The current labels associated with the input data.

  • slopes (List[float]) – The slopes for each channel.

  • offsets (List[float]) – The offsets for each channel.

Returns:

The corrected data.

Return type:

npt.ArrayLike

Raises:
  • RuntimeError – If the lengths of offsets and slopes are not the same.

  • RuntimeError – If the length of offsets and slopes is not the same as the number of channels.

Examples

# read epi correction data
with open(args.rawfile, 'rb') as raw:
    parameter2read.typ = Label.TYPE_ECHO_PHASE
    epi_corr_data, epi_corr_labels = pr.read(raw, parameter2read, pars.labels, pars.coil_info, oversampling_removal=False)

# sort the epi correction data (since ky is always 0 set the grad label as ky)
cur_recon_resolution = pars.get_recon_resolution(mix=mix, xovs=True, yovs=True, zovs=True)
epi_corr_data, epi_corr_labels = pr.sort(epi_corr_data, epi_corr_labels, output_size=[cur_recon_resolution[0]], zeropad=(True, False, False), immediate_averaging=False, ky='grad')

# epi correction
slopes, offsets = pr.get_epi_corr_data(epi_corr_data, epi_corr_labels)
data = pr.epi_corr(data, labels, slopes, offsets)

fit_flow_phase

precon.fit_flow_phase(data: np.ArrayLike, order: int = 3, weights: np.ArrayLike | None = None) np.ArrayLike

Removes a static phase offset in the flow data (due to eddy-currents) by fitting the background phase with a polynom of a given order (default = 3)

Parameters:
  • data (npt.ArrayLike) – The input data (complex).

  • order (int, optional) – The order of the polynomial fit. Defaults to 3.

  • weights (npt.ArrayLike, optional) – The weights for the fit. Defaults to None. Must be a 2-dimensional array with the same matrix size as the data

Returns:

The data with the offset phase removed.

Return type:

npt.ArrayLike

Examples

data = pr.fit_flow_phase(data, order=3)

format

precon.format(data: np.ArrayLike, in_plane_transformation: int) np.ArrayLike

Transforms the images into the radiological convention.

Parameters:
  • data (npt.ArrayLike) – The images to be transformed.

  • in_plane_transformation (int) – The in-plane transformation to be applied.

Returns:

The images in radiological convention.

Return type:

npt.ArrayLike

Examples

data = pr.format(data, pars.get_in_plane_transformation(mix=0, stack=0))

format_flow

precon.format_flow(data: np.ArrayLike, coordinate_system: List[int], venc: List[float] | None = None, hadamard: bool = False) np.ArrayLike

Makes sure the flow encoding is along the RF-AP-FH axis.

Parameters:
  • data (npt.ArrayLike) – The flow data to be formatted.

  • coordinate_system (List[int]) – The MPS coordinate system of the current scan expressed in the RL-AP-FH system. Can be obtained with the “get_coordinate_system” function which is a member of the Parameter class.

  • venc (List[float], optional) – The velocity encoding values for each flow segment. Defaults to None.

  • hadamard (bool, optional) – Flag indicating whether Hadamard encoding is used. Defaults to False.

Returns:

The formatted flow data.

Return type:

npt.ArrayLike

Raises:
  • RuntimeError – If the number of flow segments in the data is greater than 3.

  • RuntimeError – If the number of flow segments in the data is less than 3 and venc is not provided.

  • RuntimeError – If the length of venc is not equal to 3.

Examples

if get_data_size(data, pr.Enums.FLOW_SEGMENT_DIM) <= 3:
    data = pr.format_flow(data, pars.get_coordinate_system(), pars.get_venc(), pars.is_hadamard_encoding())

get_array_compression_matrix

precon.get_array_compression_matrix(ref_data: np.ArrayLike, max_nr_voxels=None) np.ArrayLike

Calculate the array compression matrix for a fully sampled reference data.

Parameters:
  • ref_data (npt.ArrayLike) – The fully sampled reference data (e.g. the sense reference scan or sensitivities).

  • max_nr_voxels (int, optional) – The maximum number of voxels to include in the calculation. Defaults to None. Reduce it if you see a performance impact in the calculations

Returns:

The array compression matrix.

Return type:

np.ndarray

Raises:

RuntimeError – If the reference data has less than 2 coils.

Examples

# reconstruct refscan
ref_pars = pr.Parameter(Path(args.refscan))
qbc, coil = pr.reconstruct_refscan(ref_pars)

# calculate the sensitivities
sens = pr.reformat_refscan(qbc, coil, ref_pars, pars, stack=stack, mix=mix, match_target_size=True)

# calculate the array compression matrix
A = pr.get_array_compression_matrix(sens.surfacecoil)

get_epi_corr_data

precon.get_epi_corr_data(data: npt.ArrayLike, labels: npt.ArrayLike)

Calculate the slopes and offsets for the EPI correction. The resulting offsets and slopes can be passed to “epi_corr”

Parameters:
  • data (npt.ArrayLike) – The EPI data.

  • labels (npt.ArrayLike) – The current labels for the data.

Returns:

A tuple containing the slopes and offsets which can be passed to the “epi_corr” function.

Return type:

tuple

Example

# read epi correction data
with open(args.rawfile, 'rb') as raw:
    parameter2read.typ = Label.TYPE_ECHO_PHASE
    epi_corr_data, epi_corr_labels = pr.read(raw, parameter2read, pars.labels, pars.coil_info, oversampling_removal=False)

# sort the epi correction data (since ky is always 0 set the grad label as ky)
cur_recon_resolution = pars.get_recon_resolution(mix=mix, xovs=True, yovs=True, zovs=True)
epi_corr_data, epi_corr_labels = pr.sort(epi_corr_data, epi_corr_labels, output_size=[cur_recon_resolution[0]], zeropad=(True, False, False), immediate_averaging=False, ky='grad')

# epi correction
slopes, offsets = pr.get_epi_corr_data(epi_corr_data, epi_corr_labels)

geo_corr

precon.geo_corr(data: np.ArrayLike, mps_to_xyz_transformation: list | np.ArrayLike, r: float, gys: np.ArrayLike, gxc: np.ArrayLike, gz: np.ArrayLike, voxel_sizes: tuple) np.ArrayLike

Corrects image distortion due gradient non-linearities for large FOV’s

Parameters:
  • data (npt.ArrayLike) – The input dataset.

  • mps_to_xyz_transformation (Union[list, npt.ArrayLike]) – The transformation matrix from the MPS to the XYZ system. Can be obtained with the “get_transformation_matrix” function which is a member of the Parameter class.

  • r (float) – The threshold radius the correction sphere. Can be obtained with the “get_concom_factors” function which is a member of the Parameter class.

  • gys (npt.ArrayLike) – The x-correction parameters. Can be obtained with the “get_concom_factors” function which is a member of the Parameter class.

  • gxc (npt.ArrayLike) – The y-correction parameters. Can be obtained with the “get_concom_factors” function which is a member of the Parameter class.

  • gz (npt.ArrayLike) – The z-correction parameters. Can be obtained with the “get_concom_factors” function which is a member of the Parameter class.

  • voxel_sizes (tuple) – The voxel sizes in the XYZ coordinates.

Returns:

The geometry corrected images.

Return type:

npt.ArrayLike

Raises:

RuntimeError – If the number of transformation matrices does not match the number of locations in the data.

Example

# get geometry correction coefficients
r, gys, gxc, gz = pars.get_geo_corr_pars()

# get the transformation from MPS to XYZ
locations = pr.utils.get_unique(labels, 'loca')
MPS_to_XYZ = pars.get_transformation_matrix(loca=locations, mix=mix, target=pr.Enums.XYZ)

# get the voxel sizes
voxel_sizes = pars.get_voxel_sizes(mix=mix)

# perform geometry correction
data = pr.geo_corr(data, MPS_to_XYZ, r, gys, gxc, gz, voxel_sizes=voxel_sizes)

hamming_filter

precon.hamming_filter(data: np.ArrayLike, strength: int | tuple, axis: int | tuple = 0, sampled_size: int | tuple | None = None) np.ArrayLike

Applies a Hamming filter to the input data along the specified axis.

Parameters:
  • data (npt.ArrayLike) – The input data to be filtered.

  • strength (Union[int, tuple]) – The strength of the filter between (0,1). Can be a single value or a tuple of values, one for each axis.

  • axis (Union[int, tuple], optional) – The axis or axes along which to apply the filter. Defaults to 0.

  • sampled_size (Union[int, tuple], optional) – The size of the sampled data. Can be a single value or a tuple of values, one for each axis. Defaults to None.

Raises:
  • RuntimeError – If the number of elements in “axis” and “strength” are not the same.

  • RuntimeError – If the number of elements in “axis” and “sampled_size” are not the same.

Returns:

The filtered data.

Examples

sampled_size = (pars.get_sampled_size(enc=0, stack=0, ovs=False), pars.get_sampled_size(enc=1, stack=0), pars.get_sampled_size(enc=2, stack=0))
data = pr.hamming_filter(data, (0.25, 0.25, 0.25), axis=(0, 1, 2), sampled_size=sampled_size)

homodyne

precon.homodyne(data: np.ArrayLike, kx_range, ky_range, kz_range, is_kspace=False, use_torch=True, py_impl=True) np.ArrayLike

Performs a homodyne reconstruction on data with partially sampled k-space data. In partial fourier imaging only a fraction of k-space is acquired for improved scan- or echo-time. If the k-space is partly sampled in readout direction we speak of partial echo and halfscan for a partly sampled space in phase encoding direction. In the first case the echo time is reduced while in the second one the scantime is shortened. Performing a standard imaging reconstruction on the partly sampled k-space however leads to image bluring and reduced SNR. Therefore a partial fourier reconstruction is performed which utilizes the k-space symmetry to reconstruct images of better quality. In MRecon a so called homodyne reconstruction is implemented: Noll, D.C.; Nishimura, D.G.; Macovski, A.; “Homodyne detection in magnetic resonance imaging” Medical Imaging, IEEE Transactions on, vol.10, no.2, pp.154-163, Jun 1991.

Parameters:
  • data (npt.ArrayLike) – The input data.

  • kx_range – The range of the k-space in readout direction.

  • ky_range – The range of the k-space in phase encoding direction.

  • kz_range – The range of the k-space in slice direction.

  • is_kspace (bool, optional) – Whether the data is in k-space. Defaults to False.

  • use_torch (bool, optional) – Whether to use the pytorch implementation. Defaults to True.

  • py_impl (bool, optional) – Whether to use the Python implementation. Defaults to True.

Returns:

The corrected data.

Return type:

npt.ArrayLike

Example

# get the k-space sampling ranges
kx_range = pars.get_range(enc=0, mix=mix, stack=stack, ovs=False)
ky_range = pars.get_range(enc=1, mix=mix, stack=stack)
kz_range = pars.get_range(enc=2, mix=mix, stack=stack)

# perform homodyne reconstruction if k-space is partially sampled
if pr.is_partial_fourier(kx_range) or pr.is_partial_fourier(ky_range) or pr.is_partial_fourier(kz_range):
    data = pr.homodyne(data, kx_range, ky_range, kz_range)

is_partial_fourier

precon.is_partial_fourier(krange) bool

Check if the given k-space range is partially sampled.

Parameters:

krange (list) – A k-space range.

Returns:

True if the krange is a partially sampled, False otherwise.

Return type:

bool

Example

# get the k-space sampling range in readout direction
kx_range = pars.get_range(enc=0, mix=mix, stack=stack, ovs=False)

# check if partially sampled
is_partial_echo = pr.is_partial_fourier(kx_range)

i2k

precon.i2k(data: np.ArrayLike, axis: int | tuple = 0, use_torch=True) np.ArrayLike

Calculates the fourier transform from image-space to k-space of the given data.

Parameters:
  • data (npt.ArrayLike) – The image-space data.

  • axis (Union[int, tuple], optional) – The axis or axes along which to calculate the fft. Defaults to 0.

  • use_torch (bool, optional) – Whether to use pytorch implementation of the fft. Defaults to True. If False then numpy is used

Returns:

The k-space data

Raises:

TypeError – If the data is not array-like.

Examples

k_data = pr.i2k(im_data, axis=(0, 1, 2))

k2i

precon.k2i(data: np.ArrayLike, axis: int | tuple = 0, use_torch=True) np.ArrayLike

Calculates the fourier transform from k-space to image-space of the given data.

Parameters:
  • data (npt.ArrayLike) – The k-space data.

  • axis (Union[int, tuple], optional) – The axis or axes along which to calculate the fft. Defaults to 0.

  • use_torch (bool, optional) – Whether to use pytorch implementation of the fft. Defaults to True. If False then numpy is used

Returns:

The image-space data

Raises:

TypeError – If the data is not array-like.

Examples

im_data = pr.k2i(k_data, axis=(0, 1, 2))

kshift

precon.kshift(data: np.ArrayLike, shift: int | tuple, axis: int | tuple | None = None) np.ArrayLike

Applies a shift to the images by adding a linear phase in k-space (fourier-shift-theorem).

Parameters:
  • data (npt.ArrayLike) – The input array (k-space data).

  • shift (Union[int, tuple]) – The shift in pixels. If a tuple is provided, each element corresponds to the shift along the corresponding axis.

  • axis (Union[int, tuple], optional) – The axis or axes along which to shift the elements. If None, the default axis is 0. If a tuple is provided, each element corresponds to the axis along which to shift the elements.

Raises:
  • RuntimeError – If the number of elements in “axis” and “shift” is not the same.

  • - RuntimeError – If the number of elements in “axis” and “shift” is not the same.

Returns:

The k-space with an added linear phase.

Examples

datak = pr.kshift(datak, (10, 20), (0, 1))

read

precon.read(raw: BinaryIO, parameter2read: Parameter2Read, labels: List[LabelStructC], coil_info: dict, normalization=True, nonlin_correction=True, pda_correction=True, random_phase_correction=True, meas_phase_correction=True, dc_offset_correction=True, sign_reversal=True, oversampling_removal=True, array_compression: np.ArrayLike | None = None) np.ArrayLike
Reads data from a Philips .raw file.

The read function returns the data together with the corresponding scanner labels. The data contains the acquired profiles in temporal order and is always a 3-dimensional array in the form:

samples x coils x profiles

By default, the function already performs all basic-corrections and removes the oversampling along the readout direction. However the default behavior can be modified with optional input arguments (see below)

Parameters:
  • raw (BinaryIO) – The binary file to read from as File-like object.

  • parameter2read (Parameter2Read) – A Parameter2Read class specifying what to read.

  • labels (List[LabelStructC]) – The scanner labels.

  • coil_info (dict) – The information about the connected receive coils.

  • normalization (bool, optional) – Whether to apply normalization. Defaults to True. The scanner provides a normalization factor for each acquired profile which is usually applied after reading.

  • nonlin_correction (bool, optional) – Whether to apply non-linear correction. Defaults to True.

  • pda_correction (bool, optional) – Whether to apply PDA correction. Defaults to True. In 3D-scans each profile might have different receiver gain values to account for the largely different signal levels between the center and outer k-space lines. This is corrected in the PDA-correction.

  • random_phase_correction (bool, optional) – Whether to apply random phase correction. Defaults to True. The scanner adds a random phase to each acquired profile to correct for systematic errors. This phase is subtracted in the random phase correction.

  • meas_phase_correction (bool, optional) – Whether to apply measurement phase correction. Defaults to True. Different echoes might have a phase offset which is subtracted in this correction.

  • dc_offset_correction (bool, optional) – Whether to apply DC offset correction. Defaults to True. The acquired data might have constant offset, which is calculated and subtracted in this correction.

  • sign_reversal (bool, optional) – Whether to apply sign reversal. Defaults to True. In EPI sequences k-space lines are acquired in different directions. This correction reverses all the odd lines.

  • oversampling_removal (bool, optional) – Whether to remove oversampling. Defaults to True. Usually the oversampling is removed to save memory (normally a factor of 2.

  • array_compression (npt.ArrayLike, optional) – The array compression matrix. Defaults to None. If given then the data is compressed to the given number of virtual channels. The number of virtual channels after compression is given by the number of rows of the compression matrix. See: Buehrer, M., Pruessmann, K. P., Boesiger, P. and Kozerke, S. (2007), Array compression for MRI with large coil arrays. Magnetic Resonance in Medicine, 57: 1131–1139. doi: 10.1002/mrm.21237

Returns:

The read data in the form samples x coils x profiles.

Return type:

npt.ArrayLike

Note

The read function only reads data of the same size and with the same geometry. Therefore, you should always loop over the number of mixes and stacks. If PARAMETER2READ specifies data of different stacks, mixes or sizes then an error will be raised.

Examples

# read parameter
pars = pr.Parameter(Path(rawfile))

# define what to read
parameter2read = pr.Parameter2Read(pars.labels)
parameter2read.stack = 0
parameter2read.mix = 0

# read data
with open(rawfile, 'rb') as raw:
    data, labels = pr.read(raw, parameter2read, pars.labels, pars.coil_info)

reconstruct_refscan

precon.reconstruct_refscan(pars: Parameter, qbc_stack: int = None)

Reconstructs the SENSE reference scan.

Parameters:
  • pars (Parameter) – The parameter object of the refscan.

  • qbc_stack (int, optional) – The stack of the bodycoil. If not given it is determined automatically. Only change if the automatic determination fails.

Returns:

A tuple containing the reconstructed bodycoil and surfacecoil data.

Return type:

tuple

Example

ref_pars = pr.Parameter(Path(refscan))
qbc, coil = reconstruct_refscan(ref_pars)

reformat_refscan

precon.reformat_refscan(qbc: np.ArrayLike, coil: np.ArrayLike, ref_pars: Parameter, target_pars: Parameter, stack: int = 0, mix: int = 0, match_target_size: bool = False, fov: List[float] | None = None, output_size: List[int] | None = None) Sensitivity

Reformats the SENSE reference scan to the geometry of the target scan.

Parameters:
  • qbc (npt.ArrayLike) – The bodycoil data of the reference scan.

  • coil (npt.ArrayLike) – The surfacecoil data of the reference scan.

  • ref_pars (Parameter) – The parameters of the sense reference scan.

  • target_pars (Parameter) – The parameters of the target scan.

  • stack (int, optional) – The stack index of the target scan which is used to determine the geometry (only relevant for multistack data). Defaults to 0.

  • mix (int, optional) – The mix index of the target scan which is used to determine the geometry (only relevant for mixed data). Defaults to 0.

  • match_target_size (bool, optional) – If True the size of the sensitivities is equal to the matrix size of the target scan. Defaults to False.

  • fov (List[float], optional) – The field of view of the reformatted data. Defaults to None. Only change if a different FOV than the target scan is needed (e.g. in simulations)

  • output_size (List[int], optional) – The matrix size of the reformatted data. Defaults to None. If not given then it is determined automatically. If match_target_size is True then this argumant is ignored

Returns:

The sensitivity object containing the reformatted SENSE reference scan.

Return type:

Sensitivity

Raises:
  • RuntimeError – If the “fov” does not contain exactly 3 elements.

  • RuntimeError – If the “output_size” does not contain exactly 3 elements.

Examples

# reconstruct refscan
ref_pars = pr.Parameter(Path(refscan))
qbc, coil = pr.reconstruct_refscan(ref_pars)

# calculate the sensitivities
sens = pr.reformat_refscan(qbc, coil, ref_pars, pars, match_target_size=True)

# calculate the sensitivities  with a user defined FOV and matrix size
sens = pr.reformat_refscan(qbc, coil, ref_pars, pars, fov=[320, 230, 50], output_size=[512, 320, 20])

retro_binning

precon.retro_binning(labels: List[LabelStructC], nr_phases: int)

Calculates the cardiac phase label from the rtop and rr label in a retrospectively triggered cardiac scan. The cardiac phase is calculated as:

label.card = int(label.rtop / label.rr * nr_phases)

Parameters: - labels (List[LabelStructC]): The scanner labels - nr_phases (int): The number of cardiac phases after binning. If not provided, the default value is 20.

Returns: - List[LabelStructC]: The new labels with the cardiac phase set.

Raises: - RuntimeError: If the labels parameter is not a list. - RuntimeError: If all rtop labels are 0. - RuntimeError: If all rr labels are 0.

Note: Make sure to call this function before sorting the data.

Examples

with open(rawfile, 'rb') as raw:
    data, labels = pr.read(raw, parameter2read, pars.labels, pars.coil_info)

# retrospective cardiac binning
labels = pr.retro_binning(labels, 30)

retro_fill_holes

precon.retro_fill_holes(data: np.ArrayLike)

sense_unfold

precon.sense_unfold(data: np.ArrayLike, sensitivity: Sensitivity, output_size: List[int], regularization_factor=2, slice_shifts: List[int] | None = None, psi_ignore_off_diagonals: bool = True, use_torch: bool = True) np.ArrayLike

Performs a SENSE unfolding on the given data.

Parameters:
  • data (npt.ArrayLike) – The sensitivity-encoded k-space data.

  • sensitivity (Sensitivity) – The sensitivity object containing the sensitivity map, noise covariance matrix and calibration data.

  • output_size (List[int]) – The desired output size of the unfolded image. The ratio between output_size and the size of the input data determines the SENSE factor

  • regularization_factor (int, optional) – The regularization factor. Defaults to 2.

  • slice_shifts (List[int], optional) – The slice shifts. Defaults to None. Slice shifts should only be given as input for a multiband reconstruction.

  • psi_ignore_off_diagonals (bool, optional) – Whether to ignore off-diagonal elements in the PSI matrix. Defaults to True.

  • use_torch (bool, optional) – Whether to use the SENSE implementation in pytorch. Defaults to True.

Returns:

The unfolded image.

Return type:

npt.ArrayLike

Example

# calculate the sensitivities
sens = pr.calculate_sensitivities(Path(refscan), Path(target_scan))

#  get the regularization factor from the scanner parameters
regularization_factor = pars.get_value(pars.SENSE_REGULARIZATION_FACTOR, at=0, default=2)

# get the unfolded size
output_size = pars.get_recon_resolution(mix=mix, xovs=False, yovs=True, zovs=True, folded=False)

# perform sense reconstruction
data = pr.sense_unfold(data, sensitivity, output_size, regularization_factor=regularization_factor)

sort

precon.sort(data: np.ArrayLike, labels: List[LabelStructC], output_size: List[int] | None = None, zeropad=True, immediate_averaging=True, **kwargs) np.ArrayLike

Sorts the input data according to its scanner labels and creates a proper k-space.

After reading the data, the acquired profiles are stored in temporal order. The sort function sorts the profiles by their image attributes (labels) and creates properly organized k-space data which can be Fourier transformed to obtain the images.

Parameters:
  • data (npt.ArrayLike) – The input data to be sorted.

  • labels (List[LabelStructC]) – The labels corresponding to the input data (output of the “read” function).

  • output_size (List[int]) – The desired size of the output data. If the output_size is larger than the sampled range, the resulting k-space is zero-filled. If it is smaller, the resulting k-space will be cropped, effectively lowering the resolution. If the output_size is not given, the resulting k-space will have the same size as the sampled range.

  • zeropad (tuple) – A tuple specifying whether to zero-pad the output data along the x, y, or z dimension. If zeropad is False, then the output_size along this direction is ignored.

  • immediate_averaging (bool, optional) – Whether to perform immediate averaging. Defaults to True.

Returns:

A tuple containing the sorted data and its corresponding labels.

Return type:

tuple

The output data has a fixed size of 12 dimensions, where each dimension corresponds to a predefined image attribute: (kx, ky, kz, channels, dynamics, cardiac phases, echoes, locations, mixes, extra1, extra2, averages)

The corresponding labels have the same size as the sorted data, except for the following dimensions:

  • The first dimension (kx or readout) is always 1 in the labels because the label information is the same for all samples along the readout direction.

  • The coil dimension (channels) is also always 1 in the labels, as the label information is the same for all coils.

Note that the labels for the zero-filled entries in the k-space will be None.

Examples

# read data
with open(args.rawfile, 'rb') as raw:
    data, labels = pr.read(raw, parameter2read, pars.labels, pars.coil_info)

# sort and zero-fill data (create k-space)
cur_recon_resolution = pars.get_recon_resolution(mix=mix, xovs=False, yovs=True, zovs=True)
data, labels = pr.sort(data, labels, output_size=cur_recon_resolution)

# Access the label for the central k-space line
ky, kz = 64, 0  # Example indices for ky and kz (assuming a 128x128 matrix)
dynamic, cardiac_phase, echo, location, mix, extra1, extra2, average = 0, 0, 0, 0, 0, 0, 0, 0

label = labels[0, ky, kz, 0, dynamic, cardiac_phase, echo, location, mix, extra1, extra2, average]

sos

precon.sos(data: np.ArrayLike, axis: int = 3) np.ArrayLike

Calculate the sum of squares along a specified axis.

Parameters:
  • data (npt.ArrayLike) – The input data array.

  • axis (int, optional) – The axis along which to calculate the sum of squares. Defaults to 3.

Returns:

The square root of the sum of squares along the specified axis.

Return type:

npt.ArrayLike

Examples

data = pr.sos(data, axis=3)

zeropad

precon.zeropad(data: np.ArrayLike, size: int | tuple, axis: int | tuple = 0) np.ArrayLike

Zero-pad an array along the specified axis or axes. The zeros will be padded symmetrically to both sides of the axis

Parameters:
  • data (npt.ArrayLike) – The input array to be zero-padded.

  • size (Union[int, tuple]) – The desired size of the output array along the specified axis or axes. If a single integer is provided, the same size will be used for all axes. If a tuple is provided, each element corresponds to the size along the corresponding axis.

  • axis (Union[int, tuple], optional) – The axis or axes along which the array should be zero-padded. If a single integer is provided, the same axis will be used for all sizes. If a tuple is provided, each element corresponds to the axis along which the corresponding size should be applied. Defaults to 0.

Raises:

RuntimeError – If the number of elements in “axis” and “size” is not the same.

Returns:

The zero-padded array.

Return type:

np.ndarray

Note

  • If the size of the input array along an axis is smaller than the desired size, zeros will be added to both ends of the axis.

  • If the size of the input array along an axis is larger than the desired size, the array will be cropped along that axis.

Example

data_padded = pr.zeropad(data, (512, 320), axis=(0, 1))