sequence.py 76 KB

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