virtualscanner.server.simulation.bloch

Numerical MR signal simulation based on the Bloch equation

caller_script_blochsim.py

run_blochsim(seqinfo, phtinfo, pat_id)

Caller function that runs Bloch simulation

This function parses parameters for the requested simulation and then executes the main simulation script, pulseq_bloch_simulator.py, passing on the parameters

Parameters
  • seqinfo (dict) – Acquire page payload

  • phtinfo (dict) – Register page payload

  • pat_id (str) – Patient ID of current session

Returns

num

Return type

int=1

phantom.py

Numerical phantom generation and access

class BrainwebPhantom(filename, dsf=1, make2d=False, loc=0, dir='z', dBmap=0)

Bases: virtualscanner.server.simulation.bloch.phantom.Phantom

This phantom is in development.

class DTTPhantom(type_map, type_params, vsize, dBmap=0, loc=(0, 0, 0))

Bases: virtualscanner.server.simulation.bloch.phantom.Phantom

Discrete tissue type phantom

Phantom constructed from a finite set of tissue types and their parameters

Parameters
  • type_map (numpy.ndarray) – Matrix of integers that map to tissue types

  • type_params (dict) – Dictionary that maps tissue type number to tissue type parameters (PD,T1,T2)

  • vsize (float) – Voxel size in meters (isotropic)

  • dBmap (numpy.ndarray, optional) – Matrix of B0 magnetic field variation across phantom The default is 0 and means no variation

  • loc (tuple, optional) – Overall location of phantom; default is (0,0,0)

class Phantom(T1map, T2map, PDmap, vsize, dBmap=0, loc=(0, 0, 0))

Bases: object

Generic numerical phantom for MRI simulations

The phantom is mainly defined by three matrices of T1, T2, and PD values, respectively. At the moment, each index in the matrix corresponds to a single spin group. The overall physical size is determined by vsize; phantom voxels must be isotropic.

Parameters
  • T1map (numpy.ndarray) – Matrix of T1 values in seconds

  • T2map (numpy.ndarray) – Matrix of T2 values in seconds

  • PDmap (numpy.ndarray) – Matrix PD values between 0 and 1

  • vsize (float) – Voxel size in meters (isotropic)

  • dBmap (numpy.ndarray, optional) – Matrix of B0 magnetic field variation across phantom The default is 0 and means no variation

  • loc (tuple, optional) – Overall location of phantom in meters; default is (0,0,0)

fov

[fov_x, fov_y, fov_z] 1 x 3 array of fields-of-view in x, y, and z directions

Type

numpy.ndarray

Xs

1D array of all x locations in phantom

Type

numpy.ndarray

Ys

1D array of all y locations in phantom

Type

numpy.ndarray

Zs

1D array of all z locations in phantom

Type

numpy.ndarray

get_list_inds()

Returns a flattened 1D array of all indices in phantom [(u1,v1,w1),…,(uk,vk,wk)]

Returns

list_inds

Return type

list

get_list_locs()

Returns a flattened 1D array of all location vectors [(x1,y1,z1),…,(xk,yk,zk)]

Returns

list_locs

Return type

list

get_location(indx)

Returns (x,y,z) physical location in meters at given indices

Parameters

indx (tuple or array_like) – (ind1, ind2, ind3) Index for querying

Returns

x, y, z – physical location corresponding to index

Return type

float

get_params(indx)

Returns PD, T1, and T2 at given indices

Parameters

indx (tuple) – Index for querying

Returns

PD, T1, T2 – Tissue parameters corresponding to the queried index

Return type

float

get_shape()

Returns the phantom’s matrix size

Returns

shape – The matrix size in three dimensions

Return type

tuple

class SpheresArrayPlanarPhantom(centers, radii, type_params, fov, n, dir='z', R=0, loc=(0, 0, 0))

Bases: virtualscanner.server.simulation.bloch.phantom.DTTPhantom

2D phantom extracted from a cylinder containing spheres

Regardless of dir, this will be an axial slice of a cylinder That is, the plane is constructed as a z-slice and then re-indexed to lie in the x or y plane The centers of spheres will correspond to locations before re-indexing

Parameters
  • centers (list or array_like) – List of 3D physical locations of the spheres’ centers

  • radii (list or array_like) – List of radii for the spheres

  • type_params (dict) – Dictionary that maps tissue type number to tissue type parameters (PD,T1,T2)

  • fov (float) – Field of view (isotropic)

  • n (int) – Matrix size

  • dir (str, optional {'z','x','y'}) – Orientation of plane; default is z

  • R (float, optional) – Cylinder’s cross-section radius; default is half of fov

  • loc (tuple, optional) – Overall location (x,y,z) of phantom from isocenter in meters Default is (0,0,0)

makeCylindricalPhantom(dim=2, n=16, dir='z', loc=0, fov=0.24)

Makes a cylindrical phantom with fixed geometry and T1, T2, PD but variable resolution and overall size

The cylinder’s diameter is the same as its height; three layers of spheres represent T1, T2, and PD variation.

Parameters
  • dim (int, optional {2,3}) – Dimension of phantom created

  • n (int) – Number of spin groups in each dimension; default is 16

  • dir (str, optional {'z', 'x', 'y'}) – Direction (norm) of plane in the case of 2D phantom

  • loc (float, optional) – Location of plane relative to isocenter; default is 0

  • fov (float, optional) – Physical length for both diameter and height of cylinder

Returns

phantom

Return type

DTTPhantom

makePlanarPhantom(n, fov, T1s, T2s, PDs, radii, dir='z', loc=(0, 0, 0))

Make a circular 2D phantom with concentric layers

Parameters
  • n (int) – Matrix size of phantom (isotropic)

  • fov (float) – Field of view of phantom (isotropic)

  • T1s (numpy.ndarray or list) – List of T1s in seconds for the layers, going outward

  • T2s (numpy.ndarray or list) – List of T2s in seconds for the layers, going outward

  • PDs (numpy.ndarray or list) – List of PDs between 0 and 1 for the layers, going outward

  • radii (numpy.ndarray) – List of radii that define the layers Note that the radii are expected to go from smallest to largest If not, they will be sorted first without sorting the parameters

  • dir (str, optional {'z','x','y'}) – Orientation of the plane; default is z, axial

  • loc (tuple, optional) – Overall (x,y,z) location of phantom; default is (0,0,0)

Returns

phantom

Return type

DTTPhantom

makeSphericalPhantom(n, fov, T1s, T2s, PDs, radii, loc=(0, 0, 0))

Make a spherical phantom with concentric layers

Parameters
  • n (int) – Matrix size of phantom (isotropic)

  • fov (float) – Field of view of phantom (isotropic)

  • T1s (numpy.ndarray or list) – List of T1s in seconds for the layers, going outward

  • T2s (numpy.ndarray or list) – List of T2s in seconds for the layers, going outward

  • PDs (numpy.ndarray or list) – List of PDs between 0 and 1 for the layers, going outward

  • radii (numpy.ndarray) – List of radii that define the layers Note that the radii are expected to go from smallest to largest If not, they will be sorted first without sorting the parameters

  • loc (tuple, optional) – Overall (x,y,z) location of phantom; default is (0,0,0)

Returns

phantom

Return type

DTTPhantom

pulseq_bloch_simulator.py

Main Bloch simulation script

This script is called by run_blochsim() with command line arguments to execute a particular simulation request, including the construction of a phantom, a pulseq sequence, parallel simulation, reconstruction, and data & image output

Notes

See beginning of the main section for explanations of the command line arguments (line 41-71)

pulseq_blochsim_methods.py

Methods to help Bloch simulation from pulseq objects

apply_pulseq_commands(isc, seq_info)

Imposes sequence commands on a single spin group

This is the key simulation function that goes through the commands and applies each to the spin group

Parameters
  • isc (SpinGroup) – The affected spin group

  • seq_info (dict) – Commands generated by store_pulseq_commands() from a pulseq object

apply_pulseq_old(isc, seq)

Deprecated function for applying a seq on a spin group and retrieving the signal

combine_gradient_areas(blk)

Helper function that combines gradient areas in a pulseq block

Parameters

blk (dict) – Pulseq block obtained from seq.get_block()

Returns

grad_areas – [Gx_area, Gy_area, Gz_area] Gradient areas converted into units of seconds*Tesla/meter

Return type

numpy.ndarray

combine_gradients(blk, dt=0, timing=(), delay=0)

Helper function that merges multiple gradients into a format for simulation

Interpolate x, y, and z gradients starting from time 0 at dt intervals, for as long as the longest gradient lasts and combine them into a 3 x N array

Parameters
  • blk (dict) – Pulseq block obtained from seq.get_block()

  • dt (float, optional) – Raster time used in interpolating gradients, in seconds Default is 0 - in this case, timing is supposed to be inputted

  • timing (numpy.ndarray, optional) – Time points at which gradients are interpolated, in seconds Default is () - in this case, dt is supposed to be inputted

  • delay (float, optional) – Adds an additional time interval in seconds at the beginning of the interpolation Default is 0; when nonzero it is only used in ADC sampling to realize ADC delay

Returns

  • grad (numpy.ndarray) – Gradient shape in Tesla/meter

  • grad_timing (numpy.ndarray) – Gradient timing in seconds

  • duration (float) – Duration of input block in seconds

Notes

Only input one argument between dt and timing and not the other

find_precessing_time(blk, dt)

Helper function that finds and returns longest duration among Gx, Gy, and Gz for use in SpinGroup.fpwg()

Parameters
  • blk (dict) – Pulseq Block obtained from seq.get_block()

  • dt (float) – Gradient raster time for calculating duration of only arbitrary gradients (‘grad’ instead of ‘trap’)

Returns

max_time – Maximum gradient time, in seconds, among the three gradients Gx, Gy, and Gz

Return type

float

get_dB0_map(maptype=0)

Returns a predefined B0 map for simulating effects of B0 inhomogeneity

Parameters

maptype (int) – Index for retrieving a map. 1 - linear map (center out) 2 - quadratic map (center out) Others - uniform map

Returns

dB0_map – This function takes location (x,y,z) as a single parameter and returns delta B0 in Tesla

Return type

function

sim_single_spingroup(loc_ind, freq_offset, phantom, seq_info)

Function for applying a seq on a spin group and retrieving the signal

Parameters
  • loc_ind (tuple) – Index in phantom of the specific spin group

  • freq_offset (float) – Off-resonance in Hertz

  • phantom (Phantom) – Phantom where spin group is located

  • seq_info (dict) – Commands generated by store_pulseq_commands() from a pulseq object

Returns

signal – Complex signal consisting of all readouts stored in the SpinGroup object

Return type

numpy.ndarray

sim_single_spingroup_old(loc_ind, freq_offset, phantom, seq)

Deprecated function for applying a seq on a spin group and retrieving the signal

store_pulseq_commands(seq)

Converts seq file into set of commands for more efficient simulation

Parameters

seq (Sequence) – Pulseq object to parse from

Returns

seq_info – Pulseq commands used by apply_pulseq_commands()

Return type

dict

pulseq_library.py

Library for generating pulseq sequences: GRE, SE, IRSE, EPI

combine_trap_grad_xyz(gradients, system, dur)

Helper function that merges multiple gradients

A list of gradients are combined into one set of 3 oblique gradients (Gx, Gy, Gz) with equivalent areas Note that the waveforms are not preserved : the outputs will always be trapezoidal gradients

Parameters
  • gradients (list) – List of gradients to be combined; there can be any number of x, y, or z gradients

  • system (Opts) – Pulseq object that indicates system constraints for gradient parameters

  • dur (float) – Duration of the output oblique gradients

Returns

gtx, gty, gtz – Oblique pulseq gradients with equivalent areas to all input gradients combined

Return type

Gradient

make_oblique_gradients(gradient, unit_grad)

Helper function to make oblique gradients

(Gx, Gy, Gz) are generated from a single orthogonal gradient and a direction indicated by unit vector

Parameters
  • gradient (Gradient) – Pulseq gradient object

  • unit_grad (array_like) – Length-3 unit vector indicating direction of resulting oblique gradient

Returns

ngx, ngy, ngz – Oblique gradients in x, y, and z directions

Return type

Gradient

make_pulseq_epi_oblique(fov, n, thk, fa, tr, te, enc='xyz', slice_locs=None, echo_type='se', n_shots=1, seg_type='blocked', write=False)

Makes an Echo Planar Imaging (EPI) sequence in any plane

2D oblique multi-slice EPI pulse sequence with Cartesian encoding Oblique means that each of slice-selection, phase encoding, and frequency encoding can point in any specified direction

Parameters
  • fov (array_like) – Isotropic field-of-view, or length-2 list [fov_readout, fov_phase], in meters

  • n (array_like) – Isotropic matrix size, or length-2 list [n_readout, n_phase]

  • thk (float) – Slice thickness in meters

  • fa (float) – Flip angle in degrees

  • tr (float) – Repetition time in seconds

  • te (float) – Echo time in seconds

  • enc (str or array_like, optional) – Spatial encoding directions 1st - readout; 2nd - phase encoding; 3rd - slice select - Use str with any permutation of x, y, and z to obtain orthogonal slices e.g. The default ‘xyz’ means axial(z) slice with readout in x and phase encoding in y - Use list to indicate vectors in the encoding directions for oblique slices They should be perpendicular to each other, but not necessarily unit vectors e.g. [(2,1,0),(-1,2,0),(0,0,1)] rotates the two in-plane encoding directions for an axial slice

  • slice_locs (array_like, optional) – Slice locations from isocenter in meters Default is None which means a single slice at the center

  • echo_type (str, optional {'se','gre'}) – Type of echo generated se (default) - spin echo (an 180 deg pulse is used) gre - gradient echo

  • n_shots (int, optional) – Number of shots used to encode each slicel; default is 1

  • seg_type (str, optional {'blocked','interleaved'}) – Method to divide up k-space in the case of n_shots > 1; default is ‘blocked’ ‘blocked’ - each shot covers a rectangle, with no overlap between shots ‘interleaved’ - each shot samples the full k-space but with wider phase steps

  • write (bool, optional) – Whether to write seq into file; default is False

Returns

  • seq (Sequence) – Pulse sequence as a Pulseq object

  • ro_dirs (numpy.ndarray) – List of 0s and 1s indicating direction of readout 0 - left to right 1 - right to left (needs to be reversed at recon)

  • ro_order (numpy.ndarray) – Order in which to re-arrange the readout lines It is [] for blocked acquisition (retain original order)

make_pulseq_gre(fov, n, thk, fa, tr, te, enc='xyz', slice_locs=None, write=False)

Makes a gradient-echo sequence

2D orthogonal multi-slice gradient-echo pulse sequence with Cartesian encoding Orthogonal means that each of slice-selection, phase encoding, and frequency encoding aligns with the x, y, or z directions

Parameters
  • fov (float) – Field-of-view in meters (isotropic)

  • n (int) – Matrix size (isotropic)

  • thk (float) – Slice thickness in meters

  • fa (float) – Flip angle in degrees

  • tr (float) – Repetition time in seconds

  • te (float) – Echo time in seconds

  • enc (str, optional) – Spatial encoding directions 1st - readout; 2nd - phase encoding; 3rd - slice select Default ‘xyz’ means axial(z) slice with readout in x and phase encoding in y

  • slice_locs (array_like, optional) – Slice locations from isocenter in meters Default is None which means a single slice at the center

  • write (bool, optional) – Whether to write seq into file; default is False

Returns

seq – Pulse sequence as a Pulseq object

Return type

Sequence

make_pulseq_gre_oblique(fov, n, thk, fa, tr, te, enc='xyz', slice_locs=None, write=False)

Makes a gradient-echo sequence in any plane

2D oblique multi-slice gradient-echo pulse sequence with Cartesian encoding Oblique means that each of slice-selection, phase encoding, and frequency encoding can point in any specified direction

Parameters
  • fov (array_like) – Isotropic field-of-view, or length-2 list [fov_readout, fov_phase], in meters

  • n (array_like) – Isotropic matrix size, or length-2 list [n_readout, n_phase]

  • thk (float) – Slice thickness in meters

  • fa (float) – Flip angle in degrees

  • tr (float) – Repetition time in seconds

  • te (float) – Echo time in seconds

  • enc (str or array_like, optional) – Spatial encoding directions 1st - readout; 2nd - phase encoding; 3rd - slice select - Use str with any permutation of x, y, and z to obtain orthogonal slices e.g. The default ‘xyz’ means axial(z) slice with readout in x and phase encoding in y - Use list to indicate vectors in the encoding directions for oblique slices They should be perpendicular to each other, but not necessarily unit vectors e.g. [(2,1,0),(-1,2,0),(0,0,1)] rotates the two in-plane encoding directions for an axial slice

  • slice_locs (array_like, optional) – Slice locations from isocenter in meters Default is None which means a single slice at the center

  • write (bool, optional) – Whether to write seq into file; default is False

Returns

seq – Pulse sequence as a Pulseq object

Return type

Sequence

make_pulseq_irse(fov, n, thk, fa, tr, te, ti, enc='xyz', slice_locs=None, write=False)

Makes an Inversion Recovery Spin Echo (IRSE) sequence

2D orthogonal multi-slice IRSE pulse sequence with Cartesian encoding Orthogonal means that each of slice-selection, phase encoding, and frequency encoding aligns with the x, y, or z directions

Parameters
  • fov (float) – Field-of-view in meters (isotropic)

  • n (int) – Matrix size (isotropic)

  • thk (float) – Slice thickness in meters

  • fa (float) – Flip angle in degrees

  • tr (float) – Repetition time in seconds

  • te (float) – Echo time in seconds

  • ti (float) – Inversion time in seconds

  • enc (str, optional) – Spatial encoding directions 1st - readout; 2nd - phase encoding; 3rd - slice select Default ‘xyz’ means axial(z) slice with readout in x and phase encoding in y

  • slice_locs (array_like, optional) – Slice locations from isocenter in meters Default is None which means a single slice at the center

  • write (bool, optional) – Whether to write seq into file; default is False

Returns

seq – Pulse sequence as a Pulseq object

Return type

Sequence

make_pulseq_irse_oblique(fov, n, thk, fa, tr, te, ti, enc='xyz', slice_locs=None, write=False)

Makes an Inversion Recovery Spin Echo (IRSE) sequence in any plane

2D oblique multi-slice IRSE pulse sequence with Cartesian encoding Oblique means that each of slice-selection, phase encoding, and frequency encoding can point in any specified direction

Parameters
  • fov (array_like) – Isotropic field-of-view, or length-2 list [fov_readout, fov_phase], in meters

  • n (array_like) – Isotropic matrix size, or length-2 list [n_readout, n_phase]

  • thk (float) – Slice thickness in meters

  • fa (float) – Flip angle in degrees

  • tr (float) – Repetition time in seconds

  • te (float) – Echo time in seconds

  • ti (float) – Inversion time in seconds

  • enc (str or array_like, optional) – Spatial encoding directions 1st - readout; 2nd - phase encoding; 3rd - slice select - Use str with any permutation of x, y, and z to obtain orthogonal slices e.g. The default ‘xyz’ means axial(z) slice with readout in x and phase encoding in y - Use list to indicate vectors in the encoding directions for oblique slices They should be perpendicular to each other, but not necessarily unit vectors e.g. [(2,1,0),(-1,2,0),(0,0,1)] rotates the two in-plane encoding directions for an axial slice

  • slice_locs (array_like, optional) – Slice locations from isocenter in meters Default is None which means a single slice at the center

  • write (bool, optional) – Whether to write seq into file; default is False

Returns

seq – Pulse sequence as a Pulseq object

Return type

Sequence

make_pulseq_se(fov, n, thk, fa, tr, te, enc='xyz', slice_locs=None, write=False)

Makes a Spin Echo (SE) sequence

2D orthogonal multi-slice Spin-Echo pulse sequence with Cartesian encoding Orthogonal means that each of slice-selection, phase encoding, and frequency encoding aligns with the x, y, or z directions

Parameters
  • fov (float) – Field-of-view in meters (isotropic)

  • n (int) – Matrix size (isotropic)

  • thk (float) – Slice thickness in meters

  • fa (float) – Flip angle in degrees

  • tr (float) – Repetition time in seconds

  • te (float) – Echo time in seconds

  • enc (str, optional) – Spatial encoding directions 1st - readout; 2nd - phase encoding; 3rd - slice select Default ‘xyz’ means axial(z) slice with readout in x and phase encoding in y

  • slice_locs (array_like, optional) – Slice locations from isocenter in meters Default is None which means a single slice at the center

  • write (bool, optional) – Whether to write seq into file; default is False

Returns

seq – Pulse sequence as a Pulseq object

Return type

Sequence

make_pulseq_se_oblique(fov, n, thk, fa, tr, te, enc='xyz', slice_locs=None, write=False)

Makes a Spin Echo (SE) sequence in any plane

2D oblique multi-slice Spin-Echo pulse sequence with Cartesian encoding Oblique means that each of slice-selection, phase encoding, and frequency encoding can point in any specified direction

Parameters
  • fov (array_like) – Isotropic field-of-view, or length-2 list [fov_readout, fov_phase], in meters

  • n (array_like) – Isotropic matrix size, or length-2 list [n_readout, n_phase]

  • thk (float) – Slice thickness in meters

  • fa (float) – Flip angle in degrees

  • tr (float) – Repetition time in seconds

  • te (float) – Echo time in seconds

  • enc (str or array_like, optional) – Spatial encoding directions 1st - readout; 2nd - phase encoding; 3rd - slice select - Use str with any permutation of x, y, and z to obtain orthogonal slices e.g. The default ‘xyz’ means axial(z) slice with readout in x and phase encoding in y - Use list to indicate vectors in the encoding directions for oblique slices They should be perpendicular to each other, but not necessarily unit vectors e.g. [(2,1,0),(-1,2,0),(0,0,1)] rotates the two in-plane encoding directions for an axial slice

  • slice_locs (array_like, optional) – Slice locations from isocenter in meters Default is None which means a single slice at the center

  • write (bool, optional) – Whether to write seq into file; default is False

Returns

seq – Pulse sequence as a Pulseq object

Return type

Sequence

modify_gradient(gradient, scale, channel=None)

Helper function to modify the strength and channel of an existing gradient

Parameters
  • gradient (Gradient) – Pulseq gradient object to be modified

  • scale (float) – Scalar to multiply the gradient strength by

  • channel (str, optional {None, 'x','y','z'}) – Channel to switch gradient into Default is None which keeps the original channel

parse_enc(enc)

Helper function for decoding enc parameter

Parameters

enc (str or array_like) – Inputted encoding scheme to parse

Returns

  • ug_fe (numpy.ndarray) – Length-3 vector of readout direction

  • ug_pe (numpy.ndarray) – Length-3 vector of phase encoding direction

  • ug_ss (numpy.ndarray) – Length-3 vector of slice selecting direction

spingroup_ps.py

class SpinGroup(loc=(0, 0, 0), pdt1t2=(1, 0, 0), df=0)

Bases: object

Basic magnetization unit for Bloch simulation

A SpinGroup has a defined location, relaxation times, proton density, and off-resonance.

Parameters
  • loc (tuple,optional) – (x,y,z) Location from isocenter in meters; default is (0,0,0)

  • pdt1t2 (tuple, optional) – (PD,T1,T2) Proton density between 0 and 1 T1, T2 in seconds; if zero it signifies no relaxation Default is PD = 1 and T1, T2 = 0

  • df (float, optional) – Off-resonance in Hertz; default is 0

m

[[Mx],[My],[Mz]] 3 x 1 array of magnetization

Type

numpy.ndarray

PD

Proton density between 0 and 1 that scales the signal

Type

float

T1

Longitudinal relaxation time in seconds

Type

float

T2

Transverse relaxation time in seconds

Type

float

loc

(x,y,z) Location of spin group from isocenter in meters

Type

tuple

df

Off-resonance in Hertz

Type

float

signal

Complex signal produced from this spin group; only generated by self.readout()

Type

array

apply_rf(pulse_shape, grads_shape, dt)

Applies an RF pulse

Euler’s method numerical integration of Bloch equation with both B1(RF) field and gradient field

Parameters
  • pulse_shape – 1 x n complex array (B1)[tesla]

  • grads_shape – 3 x n real array [tesla/meter]

  • dt – raster time for both shapes [seconds]

delay(t)

Applies a time passage to the spin group

The spin group free-precesses with only T1, T2 effects (no RF, no gradients). This function changes self.m and has no output.

Parameters

t (float) – Delay interval in seconds

fpwg(grad_area, t)

Apply only gradients to the spin group

The spin group precesses with gradients applied. This function changes self.m and has no output.

Parameters
  • grad_area (numpy.ndarray) – [Gx_area, Gy_area, Gz_area] Total area under Gx, Gy, and Gz in seconds*Tesla/meter

  • t (float) – Total time of precession in seconds

get_m_signal()

Gets spin group’s transverse magnetization

Returns

m_signal – Single complex number representing transverse magnetization Real part = Mx Imaginary part = My

Return type

numpy.ndarray

readout(dwell, n, delay, grad, timing)

ADC sampling for single spin group

Samples spin group’s magnetization while playing an arbitrary gradient This data is then stored in self.signal

Parameters
  • dwell (float) – Constant sampling interval in seconds

  • n (int) – Number of samples

  • delay (float) – Delay of the first point sampled relative to beginning of gradient waveform

  • grad (numpy.ndarray) – 2D array with shape 3 x m (i.e. m samples of the 3D gradient (Gx, Gy, Gz)) Arbitrary gradient waveform in Tesla/meter

  • timing (numpy.ndarray) – 1D array with length m Timing of gradient waveform

anyrot(v)

Helper method that generates rotational matrix from Rodrigues’s formula

Rodrigue’s formula gives a 3 x 3 matrix for arbitrary rotational axis and angle

Parameters

v (tuple) – (vx, vy, vz) 3D angular rotation vector in Cartesian space whose direction is the rotational axis and whose length the angle rotated in radians

Returns

R – 3 x 3 rotational matrix

Return type

numpy.ndarray

utest_bloch_pulseq.py

Unit test for bloch simulation on GRE, SE, and IRSE sequences (2D Cartesian, line-by-line, rectangular trajectory) Run the script to generated a simulated image. Modify the code directly to set the phantom and acquisition parameters.

utest_bloch_pulseq_epi.py

Unit test for bloch simulation on pulseq EPI sequence (2D Cartesian, rectangular trajectory) Run the script to generated a simulated image. Modify the code directly to set the phantom and acquisition parameters.