make_sigpy_pulse.py 8.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268
  1. import math
  2. from types import SimpleNamespace
  3. from typing import Tuple, Union
  4. import numpy as np
  5. import sigpy.mri.rf as rf
  6. import sigpy.plot as pl
  7. from seqgen.pypulseq.make_trapezoid import make_trapezoid
  8. from seqgen.pypulseq.opts import Opts
  9. from seqgen.pypulseq.sigpy_pulse_opts import SigpyPulseOpts
  10. def sigpy_n_seq(
  11. flip_angle: float,
  12. delay: float = 0,
  13. duration: float = 4e-3,
  14. freq_offset: float = 0,
  15. center_pos: float = 0.5,
  16. max_grad: float = 0,
  17. max_slew: float = 0,
  18. phase_offset: float = 0,
  19. return_gz: bool = True,
  20. slice_thickness: float = 0,
  21. system: Opts = Opts(),
  22. time_bw_product: float = 4,
  23. pulse_cfg: SigpyPulseOpts = SigpyPulseOpts(),
  24. use: str = str(),
  25. ) -> Union[SimpleNamespace, Tuple[SimpleNamespace, SimpleNamespace, SimpleNamespace]]:
  26. """
  27. Creates a radio-frequency sinc pulse event using the sigpy rf pulse library and optionally accompanying slice select, slice select rephasing
  28. trapezoidal gradient events.
  29. Parameters
  30. ----------
  31. flip_angle : float
  32. Flip angle in radians.
  33. apodization : float, optional, default=0
  34. Apodization.
  35. center_pos : float, optional, default=0.5
  36. Position of peak.5 (midway).
  37. delay : float, optional, default=0
  38. Delay in seconds (s).
  39. duration : float, optional, default=4e-3
  40. Duration in seconds (s).
  41. freq_offset : float, optional, default=0
  42. Frequency offset in Hertz (Hz).
  43. max_grad : float, optional, default=0
  44. Maximum gradient strength of accompanying slice select trapezoidal event.
  45. max_slew : float, optional, default=0
  46. Maximum slew rate of accompanying slice select trapezoidal event.
  47. phase_offset : float, optional, default=0
  48. Phase offset in Hertz (Hz).
  49. return_gz:bool, default=False
  50. Boolean flag to indicate if slice-selective gradient has to be returned.
  51. slice_thickness : float, optional, default=0
  52. Slice thickness of accompanying slice select trapezoidal event. The slice thickness determines the area of the
  53. slice select event.
  54. system : Opts, optional
  55. System limits. Default is a system limits object initialised to default values.
  56. time_bw_product : float, optional, default=4
  57. Time-bandwidth product.
  58. use : str, optional, default=str()
  59. Use of radio-frequency sinc pulse. Must be one of 'excitation', 'refocusing' or 'inversion'.
  60. Returns
  61. -------
  62. rf : SimpleNamespace
  63. Radio-frequency sinc pulse event.
  64. gz : SimpleNamespace, optional
  65. Accompanying slice select trapezoidal gradient event. Returned only if `slice_thickness` is provided.
  66. gzr : SimpleNamespace, optional
  67. Accompanying slice select rephasing trapezoidal gradient event. Returned only if `slice_thickness` is provided.
  68. Raises
  69. ------
  70. ValueError
  71. If invalid `use` parameter was passed. Must be one of 'excitation', 'refocusing' or 'inversion'.
  72. If `return_gz=True` and `slice_thickness` was not provided.
  73. """
  74. valid_use_pulses = ["excitation", "refocusing", "inversion"]
  75. if use != "" and use not in valid_use_pulses:
  76. raise ValueError(
  77. f"Invalid use parameter. Must be one of 'excitation', 'refocusing' or 'inversion'. Passed: {use}"
  78. )
  79. if pulse_cfg.pulse_type == "slr":
  80. [signal, t, pulse] = make_slr(
  81. flip_angle=flip_angle,
  82. time_bw_product=time_bw_product,
  83. duration=duration,
  84. system=system,
  85. pulse_cfg=pulse_cfg,
  86. disp=True,
  87. )
  88. if pulse_cfg.pulse_type == "sms":
  89. [signal, t, pulse] = make_sms(
  90. flip_angle=flip_angle,
  91. time_bw_product=time_bw_product,
  92. duration=duration,
  93. system=system,
  94. pulse_cfg=pulse_cfg,
  95. disp=True,
  96. )
  97. rfp = SimpleNamespace()
  98. rfp.type = "rf"
  99. rfp.signal = signal
  100. rfp.t = t
  101. rfp.freq_offset = freq_offset
  102. rfp.phase_offset = phase_offset
  103. rfp.dead_time = system.rf_dead_time
  104. rfp.ringdown_time = system.rf_ringdown_time
  105. rfp.delay = delay
  106. if use != "":
  107. rfp.use = use
  108. if rfp.dead_time > rfp.delay:
  109. rfp.delay = rfp.dead_time
  110. if return_gz:
  111. if slice_thickness == 0:
  112. raise ValueError("Slice thickness must be provided")
  113. if max_grad > 0:
  114. system.max_grad = max_grad
  115. if max_slew > 0:
  116. system.max_slew = max_slew
  117. BW = time_bw_product / duration
  118. amplitude = BW / slice_thickness
  119. area = amplitude * duration
  120. gz = make_trapezoid(
  121. channel="z", system=system, flat_time=duration, flat_area=area
  122. )
  123. gzr = make_trapezoid(
  124. channel="z",
  125. system=system,
  126. area=-area * (1 - center_pos) - 0.5 * (gz.area - area),
  127. )
  128. if rfp.delay > gz.rise_time:
  129. gz.delay = (
  130. math.ceil((rfp.delay - gz.rise_time) / system.grad_raster_time)
  131. * system.grad_raster_time
  132. )
  133. if rfp.delay < (gz.rise_time + gz.delay):
  134. rfp.delay = gz.rise_time + gz.delay
  135. if rfp.ringdown_time > 0:
  136. t_fill = np.arange(1, round(rfp.ringdown_time / 1e-6) + 1) * 1e-6
  137. rfp.t = np.concatenate((rfp.t, rfp.t[-1] + t_fill))
  138. rfp.signal = np.concatenate((rfp.signal, np.zeros(len(t_fill))))
  139. # Following 2 lines of code are workarounds for numpy returning 3.14... for np.angle(-0.00...)
  140. negative_zero_indices = np.where(rfp.signal == -0.0)
  141. rfp.signal[negative_zero_indices] = 0
  142. if return_gz:
  143. return rfp, gz, gzr, pulse
  144. else:
  145. return rfp
  146. def make_slr(
  147. flip_angle: float,
  148. time_bw_product: float = 4,
  149. duration: float = 0,
  150. system: Opts = Opts(),
  151. pulse_cfg: SigpyPulseOpts = SigpyPulseOpts(),
  152. disp: bool = False,
  153. ):
  154. N = int(round(duration / 1e-6))
  155. t = np.arange(1, N + 1) * system.rf_raster_time
  156. # Insert sigpy
  157. ptype = pulse_cfg.ptype
  158. ftype = pulse_cfg.ftype
  159. d1 = pulse_cfg.d1
  160. d2 = pulse_cfg.d2
  161. cancel_alpha_phs = pulse_cfg.cancel_alpha_phs
  162. pulse = rf.slr.dzrf(
  163. n=N,
  164. tb=time_bw_product,
  165. ptype=ptype,
  166. ftype=ftype,
  167. d1=d1,
  168. d2=d2,
  169. cancel_alpha_phs=cancel_alpha_phs,
  170. )
  171. flip = np.sum(pulse) * system.rf_raster_time * 2 * np.pi
  172. signal = pulse * flip_angle / flip
  173. if disp:
  174. pl.LinePlot(pulse)
  175. pl.LinePlot(signal)
  176. # Simulate it
  177. [a, b] = rf.sim.abrm(
  178. pulse,
  179. np.arange(
  180. -20 * time_bw_product, 20 * time_bw_product, 40 * time_bw_product / 2000
  181. ),
  182. True,
  183. )
  184. Mxy = 2 * np.multiply(np.conj(a), b)
  185. pl.LinePlot(Mxy)
  186. return signal, t, pulse
  187. def make_sms(
  188. flip_angle: float,
  189. time_bw_product: float = 4,
  190. duration: float = 0,
  191. system: Opts = Opts(),
  192. pulse_cfg: SigpyPulseOpts = SigpyPulseOpts(),
  193. disp: bool = False,
  194. ):
  195. N = int(round(duration / 1e-6))
  196. t = np.arange(1, N + 1) * system.rf_raster_time
  197. # Insert sigpy
  198. ptype = pulse_cfg.ptype
  199. ftype = pulse_cfg.ftype
  200. d1 = pulse_cfg.d1
  201. d2 = pulse_cfg.d2
  202. cancel_alpha_phs = pulse_cfg.cancel_alpha_phs
  203. n_bands = pulse_cfg.n_bands
  204. band_sep = pulse_cfg.band_sep
  205. phs_0_pt = pulse_cfg.phs_0_pt
  206. pulse_in = rf.slr.dzrf(
  207. n=N,
  208. tb=time_bw_product,
  209. ptype=ptype,
  210. ftype=ftype,
  211. d1=d1,
  212. d2=d2,
  213. cancel_alpha_phs=cancel_alpha_phs,
  214. )
  215. pulse = rf.multiband.mb_rf(
  216. pulse_in, n_bands=n_bands, band_sep=band_sep, phs_0_pt=phs_0_pt
  217. )
  218. flip = np.sum(pulse) * system.rf_raster_time * 2 * np.pi
  219. signal = pulse * flip_angle / flip
  220. if disp:
  221. pl.LinePlot(pulse_in)
  222. pl.LinePlot(pulse)
  223. pl.LinePlot(signal)
  224. # Simulate it
  225. [a, b] = rf.sim.abrm(
  226. pulse,
  227. np.arange(
  228. -20 * time_bw_product, 20 * time_bw_product, 40 * time_bw_product / 2000
  229. ),
  230. True,
  231. )
  232. Mxy = 2 * np.multiply(np.conj(a), b)
  233. pl.LinePlot(Mxy)
  234. return signal, t, pulse