sequence_fixed.py 68 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147114811491150115111521153115411551156115711581159116011611162116311641165116611671168116911701171117211731174117511761177117811791180118111821183118411851186118711881189119011911192119311941195119611971198119912001201120212031204120512061207120812091210121112121213121412151216121712181219122012211222122312241225122612271228122912301231123212331234123512361237123812391240124112421243124412451246124712481249125012511252125312541255125612571258125912601261126212631264126512661267126812691270127112721273127412751276127712781279128012811282128312841285128612871288128912901291129212931294129512961297129812991300130113021303130413051306130713081309131013111312131313141315131613171318131913201321132213231324132513261327132813291330133113321333133413351336133713381339134013411342134313441345134613471348134913501351135213531354135513561357135813591360136113621363136413651366136713681369137013711372137313741375137613771378137913801381138213831384138513861387138813891390139113921393139413951396139713981399140014011402140314041405140614071408140914101411141214131414141514161417141814191420142114221423142414251426142714281429143014311432143314341435143614371438143914401441144214431444144514461447144814491450145114521453145414551456145714581459146014611462146314641465146614671468146914701471147214731474147514761477147814791480148114821483148414851486148714881489149014911492149314941495149614971498149915001501150215031504150515061507150815091510151115121513151415151516151715181519152015211522152315241525152615271528152915301531153215331534153515361537153815391540154115421543154415451546154715481549155015511552155315541555155615571558155915601561156215631564156515661567156815691570157115721573157415751576157715781579158015811582158315841585158615871588158915901591159215931594159515961597159815991600160116021603160416051606160716081609161016111612161316141615161616171618161916201621162216231624162516261627162816291630163116321633163416351636163716381639164016411642164316441645164616471648164916501651165216531654165516561657165816591660166116621663166416651666166716681669167016711672167316741675167616771678167916801681168216831684168516861687168816891690169116921693169416951696169716981699170017011702170317041705170617071708
  1. import itertools
  2. import math
  3. from collections import OrderedDict
  4. from copy import copy
  5. from itertools import chain
  6. from types import SimpleNamespace
  7. from typing import Tuple, List, Iterable
  8. from typing import Union
  9. from warnings import warn
  10. import matplotlib as mpl
  11. import numpy as np
  12. from matplotlib import pyplot as plt
  13. from pypulseq import eps
  14. from pypulseq.Sequence import block, parula
  15. from pypulseq.Sequence.ext_test_report import ext_test_report
  16. from pypulseq.Sequence.read_seq import read
  17. from pypulseq.Sequence.write_seq import write as write_seq
  18. from pypulseq.calc_rf_center import calc_rf_center
  19. from pypulseq.check_timing import check_timing as ext_check_timing
  20. from pypulseq.decompress_shape import decompress_shape
  21. from pypulseq.event_lib import EventLibrary
  22. from pypulseq.opts import Opts
  23. from pypulseq.supported_labels_rf_use import get_supported_labels
  24. from version import major, minor, revision
  25. class Sequence:
  26. """
  27. Generate sequences and read/write sequence files. This class defines properties and methods to define a complete MR
  28. sequence including RF pulses, gradients, ADC events, etc. The class provides an implementation of the open MR
  29. sequence format defined by the Pulseq project. See http://pulseq.github.io/.
  30. See also `demo_read.py`, `demo_write.py`.
  31. """
  32. read_seq = None
  33. version_major = major
  34. version_minor = minor
  35. version_revision = revision
  36. def __init__(self, system=Opts()):
  37. # =========
  38. # EVENT LIBRARIES
  39. # =========
  40. self.adc_library = EventLibrary() # Library of ADC events
  41. self.delay_library = EventLibrary() # Library of delay events
  42. # Library of extension events. Extension events form single-linked zero-terminated lists
  43. self.extensions_library = EventLibrary()
  44. self.grad_library = EventLibrary() # Library of gradient events
  45. self.label_inc_library = (
  46. EventLibrary()
  47. ) # Library of Label(inc) events (reference from the extensions library)
  48. self.label_set_library = (
  49. EventLibrary()
  50. ) # Library of Label(set) events (reference from the extensions library)
  51. self.rf_library = EventLibrary() # Library of RF events
  52. self.shape_library = EventLibrary() # Library of compressed shapes
  53. self.trigger_library = EventLibrary() # Library of trigger events
  54. # =========
  55. # OTHER
  56. # =========
  57. self.system = system
  58. self.block_events = OrderedDict() # Event table
  59. self.definitions = dict() # Optional sequence definitions
  60. self.rf_raster_time = (
  61. self.system.rf_raster_time
  62. ) # RF raster time (system dependent)
  63. self.grad_raster_time = (
  64. self.system.grad_raster_time
  65. ) # Gradient raster time (system dependent)
  66. self.adc_raster_time = (
  67. self.system.adc_raster_time
  68. ) # ADC raster time (system dependent)
  69. self.block_duration_raster = self.system.block_duration_raster
  70. self.set_definition("AdcRasterTime", self.adc_raster_time)
  71. self.set_definition("BlockDurationRaster", self.block_duration_raster)
  72. self.set_definition("GradientRasterTime", self.grad_raster_time)
  73. self.set_definition("RadiofrequencyRasterTime", self.rf_raster_time)
  74. self.signature_type = ""
  75. self.signature_file = ""
  76. self.signature_value = ""
  77. self.block_durations = [] # Cache of block durations
  78. self.extension_numeric_idx = [] # numeric IDs of the used extensions
  79. self.extension_string_idx = [] # string IDs of the used extensions
  80. def __str__(self) -> str:
  81. s = "Sequence:"
  82. s += "\nshape_library: " + str(self.shape_library)
  83. s += "\nrf_library: " + str(self.rf_library)
  84. s += "\ngrad_library: " + str(self.grad_library)
  85. s += "\nadc_library: " + str(self.adc_library)
  86. s += "\ndelay_library: " + str(self.delay_library)
  87. s += "\nextensions_library: " + str(
  88. self.extensions_library
  89. ) # inserted for trigger support by mveldmann
  90. s += "\nrf_raster_time: " + str(self.rf_raster_time)
  91. s += "\ngrad_raster_time: " + str(self.grad_raster_time)
  92. s += "\nblock_events: " + str(len(self.block_events))
  93. return s
  94. def add_block(self, *args: SimpleNamespace) -> None:
  95. """
  96. Add a new block/multiple events to the sequence. Adds a sequence block with provided as a block structure
  97. See also:
  98. - `pypulseq.Sequence.sequence.Sequence.set_block()`
  99. - `pypulseq.make_adc.make_adc()`
  100. - `pypulseq.make_trapezoid.make_trapezoid()`
  101. - `pypulseq.make_sinc_pulse.make_sinc_pulse()`
  102. Parameters
  103. ----------
  104. args : SimpleNamespace
  105. Block structure or events to be added as a block to `Sequence`.
  106. """
  107. block.set_block(self, len(self.block_events) + 1, *args)
  108. def calculate_kspace(
  109. self, trajectory_delay: int = 0
  110. ) -> Tuple[np.array, np.array, np.array, np.array, np.array]:
  111. """
  112. Calculates the k-space trajectory of the entire pulse sequence.
  113. Parameters
  114. ----------
  115. trajectory_delay : int, default=0
  116. Compensation factor in seconds (s) to align ADC and gradients in the reconstruction.
  117. Returns
  118. -------
  119. k_traj_adc : numpy.array
  120. K-space trajectory sampled at `t_adc` timepoints.
  121. k_traj : numpy.array
  122. K-space trajectory of the entire pulse sequence.
  123. t_excitation : numpy.array
  124. Excitation timepoints.
  125. t_refocusing : numpy.array
  126. Refocusing timepoints.
  127. t_adc : numpy.array
  128. Sampling timepoints.
  129. """
  130. if np.any(np.abs(trajectory_delay) > 100e-6):
  131. raise Warning(
  132. f"Trajectory delay of {trajectory_delay * 1e6} us is suspiciously high"
  133. )
  134. # Initialise the counters and accumulator objects
  135. count_excitation = 0
  136. count_refocusing = 0
  137. count_adc_samples = 0
  138. # Loop through the blocks to prepare preallocations
  139. for block_counter in range(len(self.block_events)):
  140. block = self.get_block(block_counter + 1)
  141. if block.rf is not None:
  142. if (
  143. not hasattr(block.rf, "use")
  144. or block.rf.use == "excitation"
  145. or block.rf.use == "undefined"
  146. ):
  147. count_excitation += 1
  148. elif block.rf.use == "refocusing":
  149. count_refocusing += 1
  150. if block.adc is not None:
  151. count_adc_samples += int(block.adc.num_samples)
  152. t_excitation = np.zeros(count_excitation)
  153. t_refocusing = np.zeros(count_refocusing)
  154. k_time = np.zeros(count_adc_samples)
  155. current_duration = 0
  156. count_excitation = 0
  157. count_refocusing = 0
  158. kc_outer = 0
  159. traj_recon_delay = trajectory_delay
  160. # Go through the blocks and collect RF and ADC timing data
  161. for block_counter in range(len(self.block_events)):
  162. block = self.get_block(block_counter + 1)
  163. if block.rf is not None:
  164. rf = block.rf
  165. rf_center, _ = calc_rf_center(rf)
  166. t = rf.delay + rf_center
  167. if (
  168. not hasattr(block.rf, "use")
  169. or block.rf.use == "excitation"
  170. or block.rf.use == "undefined"
  171. ):
  172. t_excitation[count_excitation] = current_duration + t
  173. count_excitation += 1
  174. elif block.rf.use == "refocusing":
  175. t_refocusing[count_refocusing] = current_duration + t
  176. count_refocusing += 1
  177. if block.adc is not None:
  178. _k_time = np.arange(block.adc.num_samples) + 0.5
  179. _k_time = (
  180. _k_time * block.adc.dwell
  181. + block.adc.delay
  182. + current_duration
  183. + traj_recon_delay
  184. )
  185. k_time[kc_outer : kc_outer + block.adc.num_samples] = _k_time
  186. kc_outer += block.adc.num_samples
  187. current_duration += self.block_durations[block_counter]
  188. # Now calculate the actual k-space trajectory based on the gradient waveforms
  189. gw = self.gradient_waveforms()
  190. i_excitation = np.round(t_excitation / self.grad_raster_time)
  191. i_refocusing = np.round(t_refocusing / self.grad_raster_time)
  192. i_periods = np.sort(
  193. [1, *(i_excitation + 1), *(i_refocusing + 1), gw.shape[1] + 1]
  194. ).astype(np.int32)
  195. # i_periods -= 1 # Python is 0-indexed
  196. ii_next_excitation = np.min((len(i_excitation), 1))
  197. ii_next_refocusing = np.min((len(i_refocusing), 1))
  198. k_traj = np.zeros_like(gw)
  199. k = np.zeros((3, 1))
  200. for i in range(len(i_periods) - 1):
  201. i_period_end = i_periods[i + 1] - 1
  202. k_period = np.concatenate(
  203. (k, gw[:, i_periods[i] - 1 : i_period_end] * self.grad_raster_time),
  204. axis=1,
  205. )
  206. k_period = np.cumsum(k_period, axis=1)
  207. k_traj[:, i_periods[i] - 1 : i_period_end] = k_period[:, 1:]
  208. k = k_period[:, -1]
  209. if (
  210. ii_next_excitation > 0
  211. and i_excitation[ii_next_excitation - 1] == i_period_end
  212. ):
  213. k[:] = 0
  214. k_traj[:, i_period_end - 1] = np.nan
  215. ii_next_excitation = min(len(i_excitation), ii_next_excitation + 1)
  216. if (
  217. ii_next_refocusing > 0
  218. and i_refocusing[ii_next_refocusing - 1] == i_period_end
  219. ):
  220. k = -k
  221. ii_next_refocusing = min(len(i_refocusing), ii_next_refocusing + 1)
  222. k = k.reshape((-1, 1)) # To be compatible with np.concatenate
  223. k_traj_adc = []
  224. for _k_traj_row in k_traj:
  225. result = np.interp(
  226. xp=np.array(range(1, k_traj.shape[1] + 1)) * self.grad_raster_time,
  227. fp=_k_traj_row,
  228. x=k_time,
  229. )
  230. k_traj_adc.append(result)
  231. k_traj_adc = np.stack(k_traj_adc)
  232. t_adc = k_time
  233. return k_traj_adc, k_traj, t_excitation, t_refocusing, t_adc
  234. def calculate_kspacePP(
  235. self,
  236. trajectory_delay: Union[int, float, np.ndarray] = 0,
  237. gradient_offset: int = 0,
  238. ) -> Tuple[np.array, np.array, np.array, np.array, np.array]:
  239. """
  240. Calculates the k-space trajectory of the entire pulse sequence.
  241. Parameters
  242. ----------
  243. trajectory_delay : int, default=0
  244. Compensation factor in seconds (s) to align ADC and gradients in the reconstruction.
  245. Returns
  246. -------
  247. k_traj_adc : numpy.array
  248. K-space trajectory sampled at `t_adc` timepoints.
  249. k_traj : numpy.array
  250. K-space trajectory of the entire pulse sequence.
  251. t_excitation : numpy.array
  252. Excitation timepoints.
  253. t_refocusing : numpy.array
  254. Refocusing timepoints.
  255. t_adc : numpy.array
  256. Sampling timepoints.
  257. """
  258. if np.any(np.abs(trajectory_delay) > 100e-6):
  259. raise Warning(
  260. f"Trajectory delay of {trajectory_delay * 1e6} us is suspiciously high"
  261. )
  262. total_duration = np.sum(self.block_durations)
  263. gw_data, tfp_excitation, tfp_refocusing, t_adc, _ = self.waveforms_and_times()
  264. ng = len(gw_data)
  265. # Gradient delay handling
  266. if isinstance(trajectory_delay, (int, float)):
  267. gradient_delays = [trajectory_delay] * ng
  268. else:
  269. assert (
  270. len(trajectory_delay) == ng
  271. ) # Need to have same number of gradient channels
  272. gradient_delays = [trajectory_delay] * ng
  273. # Gradient offset handling
  274. if isinstance(gradient_offset, (int, float)):
  275. gradient_offset = [gradient_offset] * ng
  276. else:
  277. assert (
  278. len(gradient_offset) == ng
  279. ) # Need to have same number of gradient channels
  280. # Convert data to piecewise polynomials
  281. gw_pp = []
  282. gw_pp_MATLAB = []
  283. for j in range(ng):
  284. wave_cnt = gw_data[j].shape[1]
  285. if wave_cnt == 0:
  286. if np.abs(gradient_offset[j]) <= eps:
  287. continue
  288. else:
  289. gw = np.array(([0, total_duration], [0, 0]))
  290. else:
  291. gw = gw_data[j]
  292. # Now gw contains the waveform from the current axis
  293. if np.abs(gradient_delays[j]) > eps:
  294. gw[1] = gw[1] - gradient_delays[j] # Anisotropic gradient delay support
  295. if not np.all(np.isfinite(gw)):
  296. raise Warning("Not all elements of the generated waveform are finite.")
  297. teps = 1e-12
  298. if gw[0, 0] > 0 and gw[1, -1] < total_duration:
  299. # teps terms to avoid integration errors over extended periods of time
  300. _temp1 = np.array(([-teps, gw[0, 0] - teps], [0, 0]))
  301. _temp2 = np.array(([gw[0, -1] + teps, total_duration + teps], [0, 0]))
  302. gw = np.hstack((_temp1, gw, _temp2))
  303. elif gw[0, 0] > 0:
  304. _temp = np.array(([-teps, gw[0, 0] - teps], [0, 0]))
  305. gw = np.hstack((_temp, gw))
  306. elif gw[0, -1] < total_duration:
  307. _temp = np.array(([gw[0, -1] + teps, total_duration + teps], [0, 0]))
  308. gw = np.hstack((gw, _temp))
  309. if np.abs(gradient_offset[j]) > eps:
  310. gw[1:] += gradient_offset[j]
  311. gw[1][gw[1] == -0.0] = 0.0
  312. # Specify window to be same as domain prevent numpy from remapping to [-1, 1]
  313. polyfit = [
  314. np.polynomial.Polynomial.fit(
  315. gw[0, i : (i + 2)],
  316. gw[1, i : (i + 2)],
  317. deg=1,
  318. window=gw[0, i : (i + 2)],
  319. )
  320. for i in range(len(gw[0]) - 1)
  321. ]
  322. polyfit = np.stack(polyfit)
  323. gw_pp.append(polyfit)
  324. ###
  325. """
  326. Fix coefs for consistency with MATLAB:
  327. 1. MATLAB returns coefficients in descending order whereas Numpy returns coefficients in ascending order.
  328. 2. MATLAB returns local coefficients that will NOT match Numpy's outputs. Refer to the equation under the
  329. "output arguments" section of `mkpp` MATLAB docs to convert and match coefficient outputs.
  330. 3. Finally, MATLAB seems to not store any -1 < x < 1 coefficients, so we zero them out.
  331. """
  332. polyfit_MATLAB = []
  333. for i in range(len(polyfit)):
  334. polyfit_MATLAB_i = copy(polyfit[i])
  335. lower = polyfit_MATLAB_i.domain[0]
  336. co = polyfit_MATLAB_i.coef
  337. co = co[::-1] # Reverse
  338. a = co[0]
  339. b = co[1] + (a * lower)
  340. if -1 < a < 1: # to investigate
  341. a = 0
  342. if -1 < b < 1:
  343. b = 0
  344. # co = [b, a] # Un-reverse for Numpy
  345. co = [a, b]
  346. polyfit_MATLAB_i.coef = co
  347. polyfit_MATLAB.append(polyfit_MATLAB_i)
  348. gw_pp_MATLAB.append(polyfit_MATLAB)
  349. ###
  350. # Calculate slice positions. For now we entirely rely on the excitation -- ignoring complicated interleaved
  351. # refocused sequences
  352. if len(tfp_excitation) > 0:
  353. # Position in x, y, z
  354. slice_pos = np.zeros((len(gw_data), tfp_excitation.shape[1]))
  355. for j in range(len(gw_data)):
  356. if gw_pp[j] is None:
  357. slice_pos[j] = np.empty_like((1, slice_pos.shape[1]))
  358. else:
  359. slice_pos[j] = np.divide(
  360. tfp_excitation[1], self.ppval_numpy(gw_pp[j], tfp_excitation[0])
  361. )
  362. slice_pos[~np.isfinite(slice_pos)] = 0 # Reset undefined to 0
  363. t_slice_pos = tfp_excitation[0]
  364. else:
  365. slice_pos = []
  366. t_slice_pos = []
  367. # FNINT
  368. def fnint(arr_poly):
  369. pieces = len(arr_poly)
  370. breaks = np.stack([pp.domain for pp in arr_poly])
  371. breaks = np.append(breaks[:, 0], breaks[-1, 1])
  372. coefs = np.stack([pp.coef for pp in arr_poly])
  373. order = len(arr_poly[1].coef)
  374. dim = 1
  375. coefs = coefs / np.tile(range(order, 0, -1), (dim * pieces, 1))
  376. xs = np.diff(breaks[:-1])
  377. index = np.arange(pieces - 1)
  378. vv = xs * coefs[index, 0]
  379. for i in range(1, order):
  380. vv = xs * (vv + coefs[index, i])
  381. last = np.cumsum(np.insert(vv, 0, 0)).reshape((-1, 1))
  382. coefs = np.hstack((coefs[:, :order], last))
  383. arr_poly_integ = []
  384. for i in range(pieces):
  385. arr_poly_integ.append(
  386. np.polynomial.Polynomial(
  387. coefs[i],
  388. domain=[breaks[i], breaks[i + 1]],
  389. window=[breaks[i], breaks[i + 1]],
  390. )
  391. )
  392. return arr_poly_integ, coefs, breaks
  393. # =========
  394. # Integrate waveforms as PPs to produce gradient moments
  395. gm_pp = []
  396. tc = []
  397. for i in range(ng):
  398. if gw_pp[i] is None:
  399. continue
  400. res_fnint, res_coefs, res_breaks = fnint(gw_pp_MATLAB[i])
  401. gm_pp.append(res_fnint)
  402. tc.append(res_breaks)
  403. # "Sample" ramps for display purposes otherwise piecewise-linear display (plot) fails
  404. ii = np.nonzero(np.abs(res_coefs[:, 0]) > eps)[0]
  405. if len(ii) > 0:
  406. tca = []
  407. for j in range(len(ii)):
  408. res = (
  409. np.arange(
  410. np.floor(float(res_breaks[ii[j]] / self.grad_raster_time)),
  411. np.ceil(
  412. (res_breaks[ii[j] + 1] / self.grad_raster_time) + 1
  413. ),
  414. )
  415. * self.grad_raster_time
  416. )
  417. tca.extend(res)
  418. tc.append(tca)
  419. tc = np.array(list(chain(*tc)))
  420. if len(tfp_excitation) == 0:
  421. t_excitation = np.array([])
  422. else:
  423. t_excitation = tfp_excitation[0]
  424. if len(tfp_refocusing) == 0:
  425. t_refocusing = np.array([])
  426. else:
  427. t_refocusing = tfp_refocusing[0]
  428. t_acc = 1e-10 # Temporal accuracy
  429. t_acc_inv = 1 / t_acc
  430. # tc = self.__flatten_jagged_arr(tc)
  431. t_ktraj = t_acc * np.unique(
  432. np.round(
  433. t_acc_inv
  434. * np.array(
  435. [
  436. *tc,
  437. 0,
  438. *np.array(t_excitation) * 2 - self.rf_raster_time,
  439. *np.array(t_excitation) - self.rf_raster_time,
  440. *t_excitation,
  441. *np.array(t_refocusing) - self.rf_raster_time,
  442. *t_refocusing,
  443. *t_adc,
  444. total_duration,
  445. ]
  446. )
  447. )
  448. )
  449. i_excitation = np.isin(t_ktraj, t_acc * np.round(t_acc_inv * t_excitation))
  450. i_excitation = np.nonzero(i_excitation)[0] # Convert boolean array into indices
  451. i_refocusing = np.isin(t_ktraj, t_acc * np.round(t_acc_inv * t_refocusing))
  452. i_refocusing = np.nonzero(i_refocusing)[0] # Convert boolean array into indices
  453. i_adc = np.isin(t_ktraj, t_acc * np.round(t_acc_inv * t_adc))
  454. i_adc = np.nonzero(i_adc)[0] # Convert boolean array into indices
  455. i_periods = np.unique([1, *i_excitation, *i_refocusing, len(t_ktraj) - 1])
  456. if len(i_excitation) > 0:
  457. ii_next_excitation = 1
  458. else:
  459. ii_next_excitation = 0
  460. if len(i_refocusing) > 0:
  461. ii_next_refocusing = 1
  462. else:
  463. ii_next_refocusing = 0
  464. k_traj = np.zeros((3, len(t_ktraj)))
  465. for i in range(ng):
  466. if gw_pp_MATLAB[i] is None:
  467. continue
  468. it = np.logical_and(
  469. t_ktraj >= t_acc * np.round(t_acc_inv * res_breaks[0]),
  470. t_ktraj <= t_acc * np.round(t_acc_inv * res_breaks[-1]),
  471. )
  472. k_traj[i, it] = self.ppval_MATLAB(gm_pp[i], t_ktraj[it])
  473. if t_ktraj[it[-1]] < t_ktraj[-1]:
  474. k_traj[i, it[-1] + 1 :] = k_traj[i, it[-1]]
  475. # Convert gradient moments to kspace
  476. dk = -k_traj[:, 0]
  477. for i in range(len(i_periods) - 1):
  478. i_period = i_periods[i]
  479. i_period_end = i_periods[i + 1]
  480. if ii_next_excitation > 0 and i_excitation[ii_next_excitation] == i_period:
  481. if np.abs(t_ktraj[i_period] - t_excitation[ii_next_excitation]) > t_acc:
  482. raise Warning(
  483. f"np.abs(t_ktraj[i_period]-t_excitation[ii_next_excitation]) < {t_acc} failed for ii_next_excitation={ii_next_excitation} error={t_ktraj(i_period) - t_excitation(ii_next_excitation)}"
  484. )
  485. dk = -k_traj[:, i_period]
  486. if i_period > 1:
  487. # Use nans to mark the excitation points since they interrupt the plots
  488. k_traj[:, i_period - 1] = np.NaN
  489. # -1 on len(i_excitation) for 0-based indexing
  490. ii_next_excitation = np.minimum(
  491. len(i_excitation) - 1, ii_next_excitation + 1
  492. )
  493. elif (
  494. ii_next_refocusing > 0 and i_refocusing[ii_next_refocusing] == i_period
  495. ):
  496. dk = -k_traj[:, i_period]
  497. dk = -2 * k_traj[:, i_period] - dk
  498. # -1 on len(i_excitation) for 0-based indexing
  499. ii_next_refocusing = np.minimum(
  500. len(i_refocusing) - 1, ii_next_refocusing + 1
  501. )
  502. k_traj[:, i_period:i_period_end] = (
  503. k_traj[:, i_period:i_period_end] + dk[:, None]
  504. )
  505. k_traj[:, i_period_end] = k_traj[:, i_period_end] + dk
  506. k_traj_adc = k_traj[:, i_adc]
  507. return k_traj_adc, k_traj, t_excitation, t_refocusing, t_adc
  508. def check_timing(self) -> Tuple[bool, List[str]]:
  509. """
  510. Checks timing of all blocks and objects in the sequence optionally returns the detailed error log. This
  511. function also modifies the sequence object by adding the field "TotalDuration" to sequence definitions.
  512. Returns
  513. -------
  514. is_ok : bool
  515. Boolean flag indicating timing errors.
  516. error_report : str
  517. Error report in case of timing errors.
  518. """
  519. error_report = []
  520. is_ok = True
  521. num_blocks = len(self.block_events)
  522. total_duration = 0
  523. for block_counter in range(num_blocks):
  524. block = self.get_block(block_counter + 1)
  525. events = [e for e in vars(block).values() if e is not None]
  526. res, rep, duration = ext_check_timing(self.system, *events)
  527. is_ok = is_ok and res
  528. # Check the stored total block duration
  529. if np.abs(duration - self.block_durations[block_counter]) > eps:
  530. rep += "Inconsistency between the stored block duration and the duration of the block content"
  531. is_ok = False
  532. duration = self.block_durations[block_counter]
  533. # Check RF dead times
  534. if block.rf is not None:
  535. if block.rf.delay - block.rf.dead_time < -eps:
  536. rep += (
  537. f"Delay of {block.rf.delay * 1e6} us is smaller than the RF dead time "
  538. f"{block.rf.dead_time * 1e6} us"
  539. )
  540. is_ok = False
  541. if (
  542. block.rf.delay + block.rf.t[-1] + block.rf.ringdown_time - duration
  543. > eps
  544. ):
  545. rep += (
  546. f"Time between the end of the RF pulse at {block.rf.delay + block.rf.t[-1]} and the end "
  547. f"of the block at {duration * 1e6} us is shorter than rf_ringdown_time"
  548. )
  549. is_ok = False
  550. # Check ADC dead times
  551. if block.adc is not None:
  552. if block.adc.delay - self.system.adc_dead_time < -eps:
  553. rep += "adc.delay < system.adc_dead_time"
  554. is_ok = False
  555. if (
  556. block.adc.delay
  557. + block.adc.num_samples * block.adc.dwell
  558. + self.system.adc_dead_time
  559. - duration
  560. > eps
  561. ):
  562. rep += "adc: system.adc_dead_time (post-ADC) violation"
  563. is_ok = False
  564. # Update report
  565. if len(rep) != 0:
  566. error_report.append(f"Event: {block_counter} - {rep}\n")
  567. total_duration += duration
  568. # Check if all the gradients in the last block are ramped down properly
  569. if len(events) != 0 and all([isinstance(e, SimpleNamespace) for e in events]):
  570. for e in range(len(events)):
  571. if not isinstance(events[e], list) and events[e].type == "grad":
  572. if events[e].last != 0:
  573. error_report.append(
  574. f"Event {num_blocks - 1} gradients do not ramp to 0 at the end of the sequence"
  575. )
  576. self.set_definition("TotalDuration", total_duration)
  577. return is_ok, error_report
  578. def duration(self) -> Tuple[int, int, np.ndarray]:
  579. """
  580. Returns the total duration of this sequence, and the total count of blocks and events.
  581. Returns
  582. -------
  583. duration : int
  584. Duration of this sequence in seconds (s).
  585. num_blocks : int
  586. Number of blocks in this sequence.
  587. event_count : np.ndarray
  588. Number of events in this sequence.
  589. """
  590. num_blocks = len(self.block_events)
  591. event_count = np.zeros(len(self.block_events[1]))
  592. duration = 0
  593. for block_counter in range(num_blocks):
  594. event_count += self.block_events[block_counter + 1] > 0
  595. duration += self.block_durations[block_counter]
  596. return duration, num_blocks, event_count
  597. def __flatten_jagged_arr(self, arr: np.array) -> np.array:
  598. # Sanity check: we don't need to do anything if we have a flat array
  599. def __flat_check(arr: np.array) -> bool:
  600. return all([not isinstance(x, Iterable) for x in arr])
  601. if __flat_check(arr):
  602. return arr
  603. # Flatten the array simply -- 1 level deep
  604. arr_flattened = list(itertools.chain(*arr))
  605. # Not flat yet?
  606. if __flat_check(arr_flattened):
  607. return arr_flattened
  608. # Flatten the array -- 2 levels deep
  609. any_ragged = [isinstance(x, Iterable) for x in arr_flattened]
  610. if np.any(any_ragged):
  611. idx_ragged = np.array(np.where(any_ragged)[0])
  612. for i in range(len(idx_ragged)):
  613. ii = idx_ragged[i]
  614. # If we are not at the end of the list, we need to update the indices of the remaining elements
  615. # Because once we expand and insert this list element, the indices of the remaining elements
  616. # will be shifted by len(this element)
  617. if i != len(idx_ragged) - 1:
  618. idx_ragged[i + 1 :] += len(arr_flattened[ii])
  619. arr_flattened = np.insert(arr_flattened, ii, arr_flattened[ii])
  620. return arr_flattened
  621. def flip_grad_axis(self, axis: str) -> None:
  622. """
  623. Invert all gradients along the corresponding axis/channel. The function acts on all gradient objects already
  624. added to the sequence object.
  625. Parameters
  626. ----------
  627. axis : str
  628. Gradients to invert or scale. Must be one of 'x', 'y' or 'z'.
  629. """
  630. self.mod_grad_axis(axis, modifier=-1)
  631. def get_block(self, block_index: int) -> SimpleNamespace:
  632. """
  633. Return a block of the sequence specified by the index. The block is created from the sequence data with all
  634. events and shapes decompressed.
  635. See also:
  636. - `pypulseq.Sequence.sequence.Sequence.set_block()`.
  637. - `pypulseq.Sequence.sequence.Sequence.add_block()`.
  638. Parameters
  639. ----------
  640. block_index : int
  641. Index of block to be retrieved from `Sequence`.
  642. Returns
  643. -------
  644. SimpleNamespace
  645. Event identified by `block_index`.
  646. """
  647. return block.get_block(self, block_index)
  648. def get_definition(self, key: str) -> str:
  649. """
  650. Return value of the definition specified by the key. These definitions can be added manually or read from the
  651. header of a sequence file defined in the sequence header. An empty array is returned if the key is not defined.
  652. See also `pypulseq.Sequence.sequence.Sequence.set_definition()`.
  653. Parameters
  654. ----------
  655. key : str
  656. Key of definition to retrieve.
  657. Returns
  658. -------
  659. str
  660. Definition identified by `key` if found, else returns ''.
  661. """
  662. if key in self.definitions:
  663. return self.definitions[key]
  664. else:
  665. return ""
  666. def get_extension_type_ID(self, extension_string: str) -> int:
  667. """
  668. Get numeric extension ID for `extension_string`. Will automatically create a new ID if unknown.
  669. Parameters
  670. ----------
  671. extension_string : str
  672. Given string extension ID.
  673. Returns
  674. -------
  675. extension_id : int
  676. Numeric ID for given string extension ID.
  677. """
  678. if extension_string not in self.extension_string_idx:
  679. if len(self.extension_numeric_idx) == 0:
  680. extension_id = 1
  681. else:
  682. extension_id = 1 + max(self.extension_numeric_idx)
  683. self.extension_numeric_idx.append(extension_id)
  684. self.extension_string_idx.append(extension_string)
  685. assert len(self.extension_numeric_idx) == len(self.extension_string_idx)
  686. else:
  687. num = self.extension_string_idx.index(extension_string)
  688. extension_id = self.extension_numeric_idx[num]
  689. return extension_id
  690. def get_extension_type_string(self, extension_id: int) -> str:
  691. """
  692. Get string extension ID for `extension_id`.
  693. Parameters
  694. ----------
  695. extension_id : int
  696. Given numeric extension ID.
  697. Returns
  698. -------
  699. extension_str : str
  700. String ID for the given numeric extension ID.
  701. Raises
  702. ------
  703. ValueError
  704. If given numeric extension ID is unknown.
  705. """
  706. if extension_id in self.extension_numeric_idx:
  707. num = self.extension_numeric_idx.index(extension_id)
  708. else:
  709. raise ValueError(
  710. f"Extension for the given ID - {extension_id} - is unknown."
  711. )
  712. extension_str = self.extension_string_idx[num]
  713. return extension_str
  714. def mod_grad_axis(self, axis: str, modifier: int) -> None:
  715. """
  716. Invert or scale all gradients along the corresponding axis/channel. The function acts on all gradient objects
  717. already added to the sequence object.
  718. Parameters
  719. ----------
  720. axis : str
  721. Gradients to invert or scale. Must be one of 'x', 'y' or 'z'.
  722. modifier : int
  723. Scaling value.
  724. Raises
  725. ------
  726. ValueError
  727. If invalid `axis` is passed. Must be one of 'x', 'y','z'.
  728. RuntimeError
  729. If same gradient event is used on multiple axes.
  730. """
  731. if axis not in ["x", "y", "z"]:
  732. raise ValueError(
  733. f"Invalid axis. Must be one of 'x', 'y','z'. Passed: {axis}"
  734. )
  735. channel_num = ["x", "y", "z"].index(axis)
  736. other_channels = [0, 1, 2]
  737. other_channels.remove(channel_num)
  738. # Go through all event table entries and list gradient objects in the library
  739. all_grad_events = np.array(list(self.block_events.values()))
  740. all_grad_events = all_grad_events[:, 2:5]
  741. selected_events = np.unique(all_grad_events[:, channel_num])
  742. selected_events = selected_events[selected_events != 0]
  743. other_events = np.unique(all_grad_events[:, other_channels])
  744. if len(np.intersect1d(selected_events, other_events)) > 0:
  745. raise RuntimeError(
  746. "mod_grad_axis does not yet support the same gradient event used on multiple axes."
  747. )
  748. for i in range(len(selected_events)):
  749. self.grad_library.data[selected_events[i]][0] *= modifier
  750. if (
  751. self.grad_library.type[selected_events[i]] == "g"
  752. and self.grad_library.lengths[selected_events[i]] == 5
  753. ):
  754. # Need to update first and last fields
  755. self.grad_library.data[selected_events[i]][3] *= modifier
  756. self.grad_library.data[selected_events[i]][4] *= modifier
  757. def plot(
  758. self,
  759. label: str = str(),
  760. show_blocks: bool = False,
  761. save: bool = False,
  762. time_range=(0, np.inf),
  763. time_disp: str = "s",
  764. grad_disp: str = "kHz/m",
  765. plot_now: bool = True
  766. ) -> None:
  767. """
  768. Plot `Sequence`.
  769. Parameters
  770. ----------
  771. label : str, defualt=str()
  772. Plot label values for ADC events: in this example for LIN and REP labels; other valid labes are accepted as
  773. a comma-separated list.
  774. save : bool, default=False
  775. Boolean flag indicating if plots should be saved. The two figures will be saved as JPG with numerical
  776. suffixes to the filename 'seq_plot'.
  777. show_blocks : bool, default=False
  778. Boolean flag to indicate if grid and tick labels at the block boundaries are to be plotted.
  779. time_range : iterable, default=(0, np.inf)
  780. Time range (x-axis limits) for plotting the sequence. Default is 0 to infinity (entire sequence).
  781. time_disp : str, default='s'
  782. Time display type, must be one of `s`, `ms` or `us`.
  783. grad_disp : str, default='s'
  784. Gradient display unit, must be one of `kHz/m` or `mT/m`.
  785. plot_now : bool, default=True
  786. If true, function immediately shows the plots, blocking the rest of the code until plots are exited.
  787. If false, plots are shown when plt.show() is called. Useful if plots are to be modified.
  788. plot_type : str, default='Gradient'
  789. Gradients display type, must be one of either 'Gradient' or 'Kspace'.
  790. """
  791. mpl.rcParams["lines.linewidth"] = 0.75 # Set default Matplotlib linewidth
  792. valid_time_units = ["s", "ms", "us"]
  793. valid_grad_units = ["kHz/m", "mT/m"]
  794. valid_labels = get_supported_labels()
  795. if (
  796. not all([isinstance(x, (int, float)) for x in time_range])
  797. or len(time_range) != 2
  798. ):
  799. raise ValueError("Invalid time range")
  800. if time_disp not in valid_time_units:
  801. raise ValueError("Unsupported time unit")
  802. if grad_disp not in valid_grad_units:
  803. raise ValueError(
  804. "Unsupported gradient unit. Supported gradient units are: "
  805. + str(valid_grad_units)
  806. )
  807. fig1, fig2 = plt.figure(1), plt.figure(2)
  808. sp11 = fig1.add_subplot(311)
  809. sp12 = fig1.add_subplot(312, sharex=sp11)
  810. sp13 = fig1.add_subplot(313, sharex=sp11)
  811. fig2_subplots = [
  812. fig2.add_subplot(311, sharex=sp11),
  813. fig2.add_subplot(312, sharex=sp11),
  814. fig2.add_subplot(313, sharex=sp11),
  815. ]
  816. t_factor_list = [1, 1e3, 1e6]
  817. t_factor = t_factor_list[valid_time_units.index(time_disp)]
  818. g_factor_list = [1e-3, 1e3 / self.system.gamma]
  819. g_factor = g_factor_list[valid_grad_units.index(grad_disp)]
  820. t0 = 0
  821. label_defined = False
  822. label_idx_to_plot = []
  823. label_legend_to_plot = []
  824. label_store = dict()
  825. for i in range(len(valid_labels)):
  826. label_store[valid_labels[i]] = 0
  827. if valid_labels[i] in label.upper():
  828. label_idx_to_plot.append(i)
  829. label_legend_to_plot.append(valid_labels[i])
  830. if len(label_idx_to_plot) != 0:
  831. p = parula.main(len(label_idx_to_plot) + 1)
  832. label_colors_to_plot = p(np.arange(len(label_idx_to_plot)))
  833. label_colors_to_plot = np.roll(label_colors_to_plot, -1, axis=0).tolist()
  834. # Block timings
  835. block_edges = np.cumsum([0, *self.block_durations])
  836. block_edges_in_range = block_edges[
  837. (block_edges >= time_range[0]) * (block_edges <= time_range[1])
  838. ]
  839. if show_blocks:
  840. for sp in [sp11, sp12, sp13, *fig2_subplots]:
  841. sp.set_xticks(t_factor * block_edges_in_range)
  842. sp.set_xticklabels(rotation=90)
  843. for block_counter in range(len(self.block_events)):
  844. block = self.get_block(block_counter + 1)
  845. is_valid = time_range[0] <= t0 <= time_range[1]
  846. if is_valid:
  847. if getattr(block, "label", None) is not None:
  848. for i in range(len(block.label)):
  849. if block.label[i].type == "labelinc":
  850. label_store[block.label[i].label] += block.label[i].value
  851. else:
  852. label_store[block.label[i].label] = block.label[i].value
  853. label_defined = True
  854. if getattr(block, "adc", None) is not None: # ADC
  855. adc = block.adc
  856. # From Pulseq: According to the information from Klaus Scheffler and indirectly from Siemens this
  857. # is the present convention - the samples are shifted by 0.5 dwell
  858. t = adc.delay + (np.arange(int(adc.num_samples)) + 0.5) * adc.dwell
  859. sp11.plot(t_factor * (t0 + t), np.zeros(len(t)), "rx")
  860. sp13.plot(
  861. t_factor * (t0 + t),
  862. np.angle(
  863. np.exp(1j * adc.phase_offset)
  864. * np.exp(1j * 2 * np.pi * t * adc.freq_offset)
  865. ),
  866. "b.",
  867. markersize=0.25,
  868. )
  869. if label_defined and len(label_idx_to_plot) != 0:
  870. cycler = mpl.cycler(color=label_colors_to_plot)
  871. sp11.set_prop_cycle(cycler)
  872. label_colors_to_plot = np.roll(
  873. label_colors_to_plot, -1, axis=0
  874. ).tolist()
  875. arr_label_store = list(label_store.values())
  876. lbl_vals = np.take(arr_label_store, label_idx_to_plot)
  877. t = t0 + adc.delay + (adc.num_samples - 1) / 2 * adc.dwell
  878. _t = [t_factor * t] * len(lbl_vals)
  879. # Plot each label individually to retrieve each corresponding Line2D object
  880. p = itertools.chain.from_iterable(
  881. [
  882. sp11.plot(__t, _lbl_vals, ".")
  883. for __t, _lbl_vals in zip(_t, lbl_vals)
  884. ]
  885. )
  886. if len(label_legend_to_plot) != 0:
  887. sp11.legend(p, label_legend_to_plot, loc="upper left")
  888. label_legend_to_plot = []
  889. if getattr(block, "rf", None) is not None: # RF
  890. rf = block.rf
  891. tc, ic = calc_rf_center(rf)
  892. time = rf.t
  893. signal = rf.signal
  894. if np.abs(signal[0]) != 0:
  895. signal = np.insert(signal, obj=0, values=0)
  896. time = np.insert(time, obj=0, values=time[0])
  897. ic += 1
  898. if np.abs(signal[-1]) != 0:
  899. signal = np.append(signal, 0)
  900. time = np.append(time, time[-1])
  901. sp12.plot(t_factor * (t0 + time + rf.delay), np.abs(signal))
  902. sp13.plot(
  903. t_factor * (t0 + time + rf.delay),
  904. np.angle(
  905. signal
  906. * np.exp(1j * rf.phase_offset)
  907. * np.exp(1j * 2 * math.pi * time * rf.freq_offset)
  908. ),
  909. t_factor * (t0 + tc + rf.delay),
  910. np.angle(
  911. signal[ic]
  912. * np.exp(1j * rf.phase_offset)
  913. * np.exp(1j * 2 * math.pi * time[ic] * rf.freq_offset)
  914. ),
  915. "xb",
  916. )
  917. grad_channels = ["gx", "gy", "gz"]
  918. for x in range(len(grad_channels)): # Gradients
  919. if getattr(block, grad_channels[x], None) is not None:
  920. grad = getattr(block, grad_channels[x])
  921. if grad.type == "grad":
  922. # We extend the shape by adding the first and the last points in an effort of making the
  923. # display a bit less confusing...
  924. time = grad.delay + [0, *grad.tt, grad.shape_dur]
  925. waveform = g_factor * np.array(
  926. (grad.first, *grad.waveform, grad.last)
  927. )
  928. else:
  929. time = np.cumsum(
  930. [
  931. 0,
  932. grad.delay,
  933. grad.rise_time,
  934. grad.flat_time,
  935. grad.fall_time,
  936. ]
  937. )
  938. waveform = (
  939. g_factor * grad.amplitude * np.array([0, 0, 1, 1, 0])
  940. )
  941. fig2_subplots[x].plot(t_factor * (t0 + time), waveform)
  942. t0 += self.block_durations[block_counter]
  943. grad_plot_labels = ["x", "y", "z"]
  944. sp11.set_ylabel("ADC")
  945. sp12.set_ylabel("RF mag (Hz)")
  946. sp13.set_ylabel("RF/ADC phase (rad)")
  947. sp13.set_xlabel(f"t ({time_disp})")
  948. for x in range(3):
  949. _label = grad_plot_labels[x]
  950. fig2_subplots[x].set_ylabel(f"G{_label} ({grad_disp})")
  951. fig2_subplots[-1].set_xlabel(f"t ({time_disp})")
  952. # Setting display limits
  953. disp_range = t_factor * np.array([time_range[0], min(t0, time_range[1])])
  954. [x.set_xlim(disp_range) for x in [sp11, sp12, sp13, *fig2_subplots]]
  955. # Grid on
  956. for sp in [sp11, sp12, sp13, *fig2_subplots]:
  957. sp.grid()
  958. fig1.tight_layout()
  959. fig2.tight_layout()
  960. if save:
  961. fig1.savefig("seq_plot1.jpg")
  962. fig2.savefig("seq_plot2.jpg")
  963. if plot_now:
  964. plt.show()
  965. def ppval_numpy(self, arr_np_poly: np.ndarray, xq: np.ndarray):
  966. """
  967. Perform piece-wise polynomial evaluation on Numpy's polyfit objects.
  968. Parameters
  969. ==========
  970. arr_np_poly : Iterable[np.poly1d]
  971. xq : np.array
  972. Returns
  973. =======
  974. """
  975. result = []
  976. result2 = []
  977. breaks = np.array([p.domain for p in arr_np_poly])
  978. breaks = np.array([*breaks[::2].flatten(), breaks[-1, -1]])
  979. last_match = 0
  980. for x in xq:
  981. # Evaluate polynomial only when we find x's corresponding pp window
  982. for i in range(last_match, len(arr_np_poly)):
  983. corresponding_pp = x > arr_np_poly[i].domain
  984. if corresponding_pp.tolist() == [True, False]:
  985. value = np.polynomial.polynomial.polyval(x, arr_np_poly[i].coef)
  986. result.append(value)
  987. last_match = i
  988. break
  989. corresponding_pp2 = np.nonzero((x < breaks))[0][0]
  990. value2 = np.polynomial.polynomial.polyval(
  991. x, arr_np_poly[corresponding_pp2].coef
  992. )
  993. result2.append(value2)
  994. return result
  995. def ppval_MATLAB(self, arr_MATLAB_poly: np.ndarray, xq: np.ndarray) -> np.array:
  996. """
  997. Perform piece-wise polynomial evaluation on MATLAB's polyfit objects.
  998. Returns
  999. =======
  1000. """
  1001. def __ppval(c, x, b, i, n):
  1002. return c[n] * (x - b) ** i + __ppval(c, x, b, i + 1, n - 1) if n >= 0 else 0
  1003. coefs = np.array([p.coef for p in arr_MATLAB_poly])
  1004. breaks = np.array([p.domain for p in arr_MATLAB_poly])
  1005. breaks = np.array([*breaks[::2].flatten(), breaks[-1, -1]])
  1006. result = []
  1007. for x in xq:
  1008. corresponding_pp = np.nonzero((x < breaks))[0][0]
  1009. c = coefs[corresponding_pp - 1]
  1010. b = breaks[corresponding_pp - 1]
  1011. # c = coefs[2]
  1012. # res = c[0] * (x - breaks[2]) ** 2 + c[1] * (x - breaks[2]) + c[2]
  1013. res = __ppval(c, x, b, i=0, n=len(c) - 1)
  1014. result.append(res)
  1015. return np.array(result)
  1016. def read(self, file_path: str, detect_rf_use: bool = False) -> None:
  1017. """
  1018. Read `.seq` file from `file_path`.
  1019. Parameters
  1020. ----------
  1021. detect_rf_use
  1022. file_path : str
  1023. Path to `.seq` file to be read.
  1024. """
  1025. read(self, path=file_path, detect_rf_use=detect_rf_use)
  1026. def register_adc_event(self, event: EventLibrary) -> int:
  1027. return block.register_adc_event(self, event)
  1028. def register_grad_event(
  1029. self, event: SimpleNamespace
  1030. ) -> Union[int, Tuple[int, int]]:
  1031. return block.register_grad_event(self, event)
  1032. def register_label_event(self, event: SimpleNamespace) -> int:
  1033. return block.register_label_event(self, event)
  1034. def register_rf_event(self, event: SimpleNamespace) -> Tuple[int, List[int]]:
  1035. return block.register_rf_event(self, event)
  1036. def rf_from_lib_data(self, lib_data: list, use: str = str()) -> SimpleNamespace:
  1037. """
  1038. Construct RF object from `lib_data`.
  1039. Parameters
  1040. ----------
  1041. lib_data : list
  1042. RF envelope.
  1043. use : str, default=str()
  1044. RF event use.
  1045. Returns
  1046. -------
  1047. rf : SimpleNamespace
  1048. RF object constructed from `lib_data`.
  1049. """
  1050. rf = SimpleNamespace()
  1051. rf.type = "rf"
  1052. amplitude, mag_shape, phase_shape = lib_data[0], lib_data[1], lib_data[2]
  1053. shape_data = self.shape_library.data[mag_shape]
  1054. compressed = SimpleNamespace()
  1055. compressed.num_samples = shape_data[0]
  1056. compressed.data = shape_data[1:]
  1057. mag = decompress_shape(compressed)
  1058. shape_data = self.shape_library.data[phase_shape]
  1059. compressed.num_samples = shape_data[0]
  1060. compressed.data = shape_data[1:]
  1061. phase = decompress_shape(compressed)
  1062. rf.signal = amplitude * mag * np.exp(1j * 2 * np.pi * phase)
  1063. time_shape = lib_data[3]
  1064. if time_shape > 0:
  1065. shape_data = self.shape_library.data[time_shape]
  1066. compressed.num_samples = shape_data[0]
  1067. compressed.data = shape_data[1:]
  1068. rf.t = decompress_shape(compressed) * self.rf_raster_time
  1069. rf.shape_dur = (
  1070. np.ceil((rf.t[-1] - eps) / self.rf_raster_time) * self.rf_raster_time
  1071. )
  1072. else: # Generate default time raster on the fly
  1073. rf.t = (np.arange(1, len(rf.signal) + 1) - 0.5) * self.rf_raster_time
  1074. rf.shape_dur = len(rf.signal) * self.rf_raster_time
  1075. rf.delay = lib_data[4]
  1076. rf.freq_offset = lib_data[5]
  1077. rf.phase_offset = lib_data[6]
  1078. rf.dead_time = self.system.rf_dead_time
  1079. rf.ringdown_time = self.system.rf_ringdown_time
  1080. if use != "":
  1081. use_cases = {
  1082. "e": "excitation",
  1083. "r": "refocusing",
  1084. "i": "inversion",
  1085. "s": "saturation",
  1086. "p": "preparation",
  1087. }
  1088. rf.use = use_cases[use] if use in use_cases else "undefined"
  1089. return rf
  1090. def set_block(self, block_index: int, *args: SimpleNamespace) -> None:
  1091. """
  1092. Replace block at index with new block provided as block structure, add sequence block, or create a new block
  1093. from events and store at position specified by index. The block or events are provided in uncompressed form and
  1094. will be stored in the compressed, non-redundant internal libraries.
  1095. See also:
  1096. - `pypulseq.Sequence.sequence.Sequence.get_block()`
  1097. - `pypulseq.Sequence.sequence.Sequence.add_block()`
  1098. Parameters
  1099. ----------
  1100. block_index : int
  1101. Index at which block is replaced.
  1102. args : SimpleNamespace
  1103. Block or events to be replaced/added or created at `block_index`.
  1104. """
  1105. block.set_block(self, block_index, *args)
  1106. def set_definition(
  1107. self, key: str, value: Union[float, int, list, np.ndarray, str, tuple]
  1108. ) -> None:
  1109. """
  1110. Modify a custom definition of the sequence. Set the user definition 'key' to value 'value'. If the definition
  1111. does not exist it will be created.
  1112. See also `pypulseq.Sequence.sequence.Sequence.get_definition()`.
  1113. Parameters
  1114. ----------
  1115. key : str
  1116. Definition key.
  1117. value : int, list, np.ndarray, str or tuple
  1118. Definition value.
  1119. """
  1120. if key == "FOV":
  1121. if np.max(value) > 1:
  1122. text = "Definition FOV uses values exceeding 1 m. "
  1123. text += "New Pulseq interpreters expect values in units of meters."
  1124. warn(text)
  1125. self.definitions[key] = value
  1126. def set_extension_string_ID(self, extension_str: str, extension_id: int) -> None:
  1127. """
  1128. Set numeric ID for the given string extension ID.
  1129. Parameters
  1130. ----------
  1131. extension_str : str
  1132. Given string extension ID.
  1133. extension_id : int
  1134. Given numeric extension ID.
  1135. Raises
  1136. ------
  1137. ValueError
  1138. If given numeric or string extension ID is not unique.
  1139. """
  1140. if (
  1141. extension_str in self.extension_string_idx
  1142. or extension_id in self.extension_numeric_idx
  1143. ):
  1144. raise ValueError("Numeric or string ID is not unique")
  1145. self.extension_numeric_idx.append(extension_id)
  1146. self.extension_string_idx.append(extension_str)
  1147. assert len(self.extension_numeric_idx) == len(self.extension_string_idx)
  1148. def test_report(self) -> str:
  1149. """
  1150. Analyze the sequence and return a text report.
  1151. """
  1152. return ext_test_report(self)
  1153. def waveforms_and_times(
  1154. self, append_RF: bool = False
  1155. ) -> Tuple[np.ndarray, np.ndarray, np.ndarray, np.ndarray, np.ndarray]:
  1156. """
  1157. Decompress the entire gradient waveform. Returns gradient waveforms as a tuple of `np.ndarray` of
  1158. `gradient_axes` (typically 3) dimensions. Each `np.ndarray` contains timepoints and the corresponding
  1159. gradient amplitude values. Additional return values are time points of excitations, refocusings and ADC
  1160. sampling points.
  1161. Parameters
  1162. ----------
  1163. append_RF : bool, default=False
  1164. Boolean flag to indicate if RF wave shapes are to be appended after the gradients.
  1165. Returns
  1166. -------
  1167. wave_data : np.ndarray
  1168. tfp_excitation : np.ndarray
  1169. Contains time moments, frequency and phase offsets of the excitation RF pulses (similar for `
  1170. tfp_refocusing`).
  1171. tfp_refocusing : np.ndarray
  1172. t_adc: np.ndarray
  1173. Contains times of all ADC sample points.
  1174. fp_adc : np.ndarray
  1175. Contains frequency and phase offsets of each ADC object (not samples).
  1176. """
  1177. grad_channels = ["gx", "gy", "gz"]
  1178. num_blocks = len(self.block_events)
  1179. # Collect shape pieces
  1180. if append_RF:
  1181. shape_channels = len(grad_channels) + 1 # Last 'channel' is RF
  1182. else:
  1183. shape_channels = len(grad_channels)
  1184. shape_pieces = np.empty((shape_channels, num_blocks), dtype="object")
  1185. # Also collect RF and ADC timing data
  1186. # t_excitation, t_refocusing, t_adc
  1187. tfp_excitation = []
  1188. tfp_refocusing = []
  1189. t_adc = []
  1190. fp_adc = []
  1191. curr_dur = 0
  1192. out_len = np.zeros(shape_channels) # Last 'channel' is RF
  1193. for block_counter in range(num_blocks):
  1194. block = self.get_block(block_counter + 1)
  1195. for j in range(len(grad_channels)):
  1196. grad = getattr(block, grad_channels[j])
  1197. if grad is not None: # Gradients
  1198. if grad.type == "grad":
  1199. # Check if we have an extended trapezoid or an arbitrary gradient on a regular raster
  1200. tt_rast = grad.tt / self.grad_raster_time + 0.5
  1201. if np.all(
  1202. np.abs(tt_rast - np.arange(1, len(tt_rast) + 1)) < eps
  1203. ): # Arbitrary gradient
  1204. """
  1205. Arbitrary gradient: restore & recompress shape - if we had a trapezoid converted to shape we
  1206. have to find the "corners" and we can eliminate internal samples on the straight segments
  1207. but first we have to restore samples on the edges of the gradient raster intervals for that
  1208. we need the first sample.
  1209. """
  1210. max_abs = np.max(np.abs(grad.waveform))
  1211. odd_step1 = np.array([grad.first, *(2 * grad.waveform)])
  1212. odd_step2 = odd_step1 * (
  1213. np.mod(np.arange(1, len(odd_step1) + 1), 2) * 2 - 1
  1214. )
  1215. waveform_odd_rest = np.cumsum(odd_step2) * (
  1216. np.mod(np.arange(1, len(odd_step2)), 2) * 2 - 1
  1217. )
  1218. else: # Extended trapezoid
  1219. out_len[j] += len(grad.tt)
  1220. shape_pieces[j, block_counter] = np.array(
  1221. [
  1222. curr_dur + grad.delay + grad.tt,
  1223. grad.waveform,
  1224. ]
  1225. )
  1226. else:
  1227. if np.abs(grad.flat_time) > eps:
  1228. out_len[j] += 4
  1229. _temp = np.vstack(
  1230. (
  1231. curr_dur
  1232. + grad.delay
  1233. + np.cumsum(
  1234. [
  1235. 0,
  1236. grad.rise_time,
  1237. grad.flat_time,
  1238. grad.fall_time,
  1239. ]
  1240. ),
  1241. grad.amplitude * np.array([0, 1, 1, 0]),
  1242. )
  1243. )
  1244. shape_pieces[j, block_counter] = _temp
  1245. else:
  1246. out_len[j] += 3
  1247. _temp = np.vstack(
  1248. (
  1249. curr_dur
  1250. + grad.delay
  1251. + np.cumsum([0, grad.rise_time, grad.fall_time]),
  1252. grad.amplitude * np.array([0, 1, 0]),
  1253. )
  1254. )
  1255. shape_pieces[j, block_counter] = _temp
  1256. if block.rf is not None: # RF
  1257. rf = block.rf
  1258. t = rf.delay + calc_rf_center(rf)[0]
  1259. if not hasattr(rf, "use") or block.rf.use in [
  1260. "excitation",
  1261. "undefined",
  1262. ]:
  1263. tfp_excitation.append(
  1264. [curr_dur + t, block.rf.freq_offset, block.rf.phase_offset]
  1265. )
  1266. elif block.rf.use == "refocusing":
  1267. tfp_refocusing.append(
  1268. [curr_dur + t, block.rf.freq_offset, block.rf.phase_offset]
  1269. )
  1270. if append_RF:
  1271. rf_piece = np.array(
  1272. [
  1273. curr_dur + rf.delay + rf.t,
  1274. rf.signal
  1275. * np.exp(
  1276. 1j
  1277. * (rf.phase_offset + 2 * np.pi * rf.freq_offset * rf.t)
  1278. ),
  1279. ]
  1280. )
  1281. out_len[-1] += len(rf.t)
  1282. if np.abs(rf.signal[0]) > 0:
  1283. pre = np.array([[curr_dur + rf.delay + rf.t[0] - eps], [0]])
  1284. rf_piece = np.hstack((pre, rf_piece))
  1285. out_len[-1] += pre.shape[1]
  1286. if np.abs(rf.signal[-1]) > 0:
  1287. post = np.array([[curr_dur + rf.delay + rf.t[-1] + eps], [0]])
  1288. rf_piece = np.hstack((rf_piece, post))
  1289. out_len[-1] += post.shape[1]
  1290. shape_pieces[-1, block_counter] = rf_piece
  1291. if block.adc is not None: # ADC
  1292. t_adc.extend(
  1293. np.arange(block.adc.num_samples)
  1294. + 0.5 * block.adc.dwell
  1295. + block.adc.delay
  1296. + curr_dur
  1297. )
  1298. fp_adc.append([block.adc.freq_offset, block.adc.phase_offset])
  1299. curr_dur += self.block_durations[block_counter]
  1300. # Collect wave data
  1301. wave_data = np.empty(shape_channels, dtype="object")
  1302. for j in range(shape_channels):
  1303. wave_data[j] = np.zeros((2, int(out_len[j])))
  1304. # TODO: This is weird, and needs to be fixed. Time points are also complex this way, and spends 4 times more memory than necessary.
  1305. if append_RF:
  1306. wave_data[j] = np.zeros((2, int(out_len[j])), dtype=np.complex128)
  1307. wave_cnt = np.zeros(shape_channels, dtype=int)
  1308. for block_counter in range(num_blocks):
  1309. for j in range(shape_channels):
  1310. if shape_pieces[j, block_counter] is not None:
  1311. wave_data_local = shape_pieces[j, block_counter]
  1312. length = wave_data_local.shape[1]
  1313. if (
  1314. wave_cnt[j] == 0
  1315. or wave_data[j][0, wave_cnt[j] - 1] != wave_data_local[0, 0]
  1316. ):
  1317. wave_data[j][
  1318. :, wave_cnt[j] + np.arange(length)
  1319. ] = wave_data_local
  1320. wave_cnt[j] += length
  1321. else:
  1322. wave_data[j][
  1323. :, wave_cnt[j] + np.arange(length - 1)
  1324. ] = wave_data_local[:, 1:]
  1325. wave_cnt[j] += length - 1
  1326. if wave_cnt[j] != len(np.unique(wave_data[j][0, : wave_cnt[j]])):
  1327. raise Warning(
  1328. "Not all elements of the generated time vector are unique."
  1329. )
  1330. # Trim output data
  1331. for j in range(shape_channels):
  1332. if wave_cnt[j] < wave_data[j].shape[1]:
  1333. wave_data[j] = wave_data[j][:, : wave_cnt[j]]
  1334. tfp_excitation = np.array(tfp_excitation).transpose()
  1335. tfp_refocusing = np.array(tfp_refocusing)
  1336. t_adc = np.array(t_adc)
  1337. fp_adc = np.array(fp_adc)
  1338. return wave_data, tfp_excitation, tfp_refocusing, t_adc, fp_adc
  1339. def waveforms_export(self, time_range=(0, np.inf)) -> dict:
  1340. """
  1341. Plot `Sequence`.
  1342. Parameters
  1343. ----------
  1344. time_range : iterable, default=(0, np.inf)
  1345. Time range (x-axis limits) for all waveforms. Default is 0 to infinity (entire sequence).
  1346. Returns
  1347. -------
  1348. all_waveforms: dict
  1349. Dictionary containing the following sequence waveforms and time array(s):
  1350. - `t_adc` - ADC timing array [seconds]
  1351. - `t_rf` - RF timing array [seconds]
  1352. - `t_rf_centers`: `rf_t_centers`,
  1353. - `t_gx`: x gradient timing array,
  1354. - `t_gy`: y gradient timing array,
  1355. - `t_gz`: z gradient timing array,
  1356. - `adc` - ADC complex signal (amplitude=1, phase=adc phase) [a.u.]
  1357. - `rf` - RF complex signal
  1358. - `rf_centers`: RF centers array,
  1359. - `gx` - x gradient
  1360. - `gy` - y gradient
  1361. - `gz` - z gradient
  1362. - `grad_unit`: [kHz/m],
  1363. - `rf_unit`: [Hz],
  1364. - `time_unit`: [seconds],
  1365. """
  1366. # Check time range validity
  1367. if (
  1368. not all([isinstance(x, (int, float)) for x in time_range])
  1369. or len(time_range) != 2
  1370. ):
  1371. raise ValueError("Invalid time range")
  1372. t0 = 0
  1373. adc_t_all = np.array([])
  1374. adc_signal_all = np.array([], dtype=complex)
  1375. rf_t_all = np.array([])
  1376. rf_signal_all = np.array([], dtype=complex)
  1377. rf_t_centers = np.array([])
  1378. rf_signal_centers = np.array([], dtype=complex)
  1379. gx_t_all = np.array([])
  1380. gy_t_all = np.array([])
  1381. gz_t_all = np.array([])
  1382. gx_all = np.array([])
  1383. gy_all = np.array([])
  1384. gz_all = np.array([])
  1385. for block_counter in range(len(self.block_events)): # For each block
  1386. block = self.get_block(block_counter + 1) # Retrieve it
  1387. is_valid = (
  1388. time_range[0] <= t0 <= time_range[1]
  1389. ) # Check if "current time" is within requested range.
  1390. if is_valid:
  1391. # Case 1: ADC
  1392. if hasattr(block, "adc") and block.adc!=None:
  1393. adc = block.adc # Get adc info
  1394. # From Pulseq: According to the information from Klaus Scheffler and indirectly from Siemens this
  1395. # is the present convention - the samples are shifted by 0.5 dwell
  1396. t = adc.delay + (np.arange(int(adc.num_samples)) + 0.5) * adc.dwell
  1397. adc_t = t0 + t
  1398. adc_signal = np.exp(1j * adc.phase_offset) * np.exp(
  1399. 1j * 2 * np.pi * t * adc.freq_offset
  1400. )
  1401. adc_t_all = np.append(adc_t_all, adc_t)
  1402. adc_signal_all = np.append(adc_signal_all, adc_signal)
  1403. if hasattr(block, "rf") and block.rf!=None:
  1404. rf = block.rf
  1405. tc, ic = calc_rf_center(rf)
  1406. t = rf.t + rf.delay
  1407. tc = tc + rf.delay
  1408. # Debug - visualize
  1409. # sp12.plot(t_factor * (t0 + t), np.abs(rf.signal))
  1410. # plt.show()
  1411. # sp13.plot(t_factor * (t0 + t), np.angle(rf.signal * np.exp(1j * rf.phase_offset)
  1412. # * np.exp(1j * 2 * math.pi * rf.t * rf.freq_offset)),
  1413. # t_factor * (t0 + tc), np.angle(rf.signal[ic] * np.exp(1j * rf.phase_offset)
  1414. # * np.exp(1j * 2 * math.pi * rf.t[ic] * rf.freq_offset)),
  1415. # 'xb')
  1416. # plt.show()
  1417. rf_t = t0 + t
  1418. rf = (
  1419. rf.signal
  1420. * np.exp(1j * rf.phase_offset)
  1421. * np.exp(1j * 2 * math.pi * rf.t * rf.freq_offset)
  1422. )
  1423. rf_t_all = np.append(rf_t_all, rf_t)
  1424. rf_signal_all = np.append(rf_signal_all, rf)
  1425. rf_t_centers = np.append(rf_t_centers, rf_t[ic])
  1426. rf_signal_centers = np.append(rf_signal_centers, rf[ic])
  1427. grad_channels = ['gx', 'gy', 'gz']
  1428. for x in range(
  1429. len(grad_channels)
  1430. ): # Check each gradient channel: x, y, and z
  1431. if hasattr(
  1432. block, grad_channels[x]
  1433. ):
  1434. # If this channel is on in current block
  1435. grad = getattr(block, grad_channels[x])
  1436. if grad == None:
  1437. break
  1438. elif grad.type == "grad": # Arbitrary gradient option
  1439. # In place unpacking of grad.t with the starred expression
  1440. g_t = (
  1441. t0
  1442. + grad.delay
  1443. + [
  1444. 0,
  1445. *(grad.t + (grad.t[1] - grad.t[0]) / 2),
  1446. grad.t[-1] + grad.t[1] - grad.t[0],
  1447. ]
  1448. )
  1449. g = 1e-3 * np.array((grad.first, *grad.waveform, grad.last))
  1450. else: # Trapezoid gradient option
  1451. g_t = t0 + np.cumsum(
  1452. [
  1453. 0,
  1454. grad.delay,
  1455. grad.rise_time,
  1456. grad.flat_time,
  1457. grad.fall_time,
  1458. ]
  1459. )
  1460. g = 1e-3 * grad.amplitude * np.array([0, 0, 1, 1, 0])
  1461. if grad.channel == "x":
  1462. gx_t_all = np.append(gx_t_all, g_t)
  1463. gx_all = np.append(gx_all, g)
  1464. elif grad.channel == "y":
  1465. gy_t_all = np.append(gy_t_all, g_t)
  1466. gy_all = np.append(gy_all, g)
  1467. elif grad.channel == "z":
  1468. gz_t_all = np.append(gz_t_all, g_t)
  1469. gz_all = np.append(gz_all, g)
  1470. t0 += self.block_durations[
  1471. block_counter
  1472. ] # "Current time" gets updated to end of block just examined
  1473. all_waveforms = {
  1474. "t_adc": adc_t_all,
  1475. "t_rf": rf_t_all,
  1476. "t_rf_centers": rf_t_centers,
  1477. "t_gx": gx_t_all,
  1478. "t_gy": gy_t_all,
  1479. "t_gz": gz_t_all,
  1480. "adc": adc_signal_all,
  1481. "rf": rf_signal_all,
  1482. "rf_centers": rf_signal_centers,
  1483. "gx": gx_all,
  1484. "gy": gy_all,
  1485. "gz": gz_all,
  1486. "grad_unit": "[kHz/m]",
  1487. "rf_unit": "[Hz]",
  1488. "time_unit": "[seconds]",
  1489. }
  1490. return all_waveforms
  1491. def write(self, name: str, create_signature: bool = True) -> None:
  1492. """
  1493. Write the sequence data to the given filename using the open file format for MR sequences.
  1494. See also `pypulseq.Sequence.read_seq.read()`.
  1495. Parameters
  1496. ----------
  1497. name : str
  1498. Filename of `.seq` file to be written to disk.
  1499. create_signature : bool, default=True
  1500. Boolean flag to indicate if the file has to be signed.
  1501. """
  1502. write_seq(self, name, create_signature)