API

GrapeMR.RF_pulse_analysisMethod
RF_pulse_analysis(cf::ControlField; attenuation_ref = 0.0, B1_ref = 1.0, power_ref = 500.0)

Calculate calibration analysis of shaped pulse

Arguments

  • cf::ControlField: The control field object containing the RF waveform.
  • attenuation_ref::Float64: Reference attenuation in dB (default is 0.0 dB).
  • B1_ref::Float64: Reference RF field strength in Tesla (default is 1.0 T).
  • power_ref::Float64: Reference RF power in Watts (default is 500.0 W).

Outputs:

A tuple of calculated values:

  • max_amp: Maximum RF amplitude in Hertz (Hz).
  • max_amp_tesla: Maximum RF amplitude in microtesla (µT).
  • attenuation_B1: Attenuation in decibels (dB).
  • power_max_B1: Maximum power in Watts (W).
  • power_average: Average power in Watts (W).
  • pulse_energy: Pulse energy in Joules (J).
source
GrapeMR.add_arc_to_trajectoryMethod
add_arc_to_trajectory(p, Mxy, Mz, i1::Int, i2::Int)

Add an arc between two specified quivers on an existing magnetization trajectory plot.

Arguments

  • p: Existing plot
  • Mxy: Vector of transverse magnetization values
  • Mz: Vector of longitudinal magnetization values
  • i1: Index of first vector
  • i2: Index of second vector

Returns

  • Updated plot with arc
source
GrapeMR.backward_propagation!Method
backward_propagation!(χ::AbstractMatrix, cf::ControlField, iso::Isochromat)

In-place backward propagation for the adjoint state matrix, calculating gradients for the control fields.

Arguments

  • χ::AbstractMatrix: Matrix to store the adjoint state (4xN).
  • cf::ControlField: Struct containing control field parameters.
  • iso::Isochromat: Isochromat containing the magnetization data for backward propagation.

Returns

  • χ: Updated adjoint state matrix (4xN) after backward propagation.
source
GrapeMR.backward_propagationMethod
backward_propagation(cost_grad::AbstractVector, cf::ControlField, iso::Isochromat)

Performs backward propagation for the adjoint state matrix, calculating gradients for the control fields.

Arguments

  • cost_grad::AbstractVector: Gradient of the cost function for the initial adjoint state.
  • cf::ControlField: Struct containing control field parameters.
  • iso::Isochromat: Isochromat containing the magnetization data for backward propagation.

Returns

  • χ::Matrix{Float64}: Adjoint state matrix (4xN) after backward propagation.
source
GrapeMR.bloch_matrixMethod
bloch_matrix(B1x::Float64, B1y::Float64, Bz::Float64, Γ1::Float64, Γ2::Float64)

Calculates the Bloch matrix for spin dynamics.

Arguments

  • B1x::Float64: x-component of the B1 field.
  • B1y::Float64: y-component of the B1 field.
  • Bz::Float64: z-component of the magnetic field.
  • Γ1::Float64: Longitudinal relaxation rate.
  • Γ2::Float64: Transverse relaxation rate.

Returns

  • A 4x4 Bloch matrix based on the given field components and relaxation rates.
source
GrapeMR.bohb_hyperoptMethod
bohb_hyperopt(spins::Vector{<:Spins}, 
              gp::GrapeParams, 
              Tc::LinRange, 
              max_iter::Int; 
              i::Int=5, 
              poly_start::Vector{Float64} = [5e-1, 1e-1, 1e-2], 
              poly_degree::Vector{Int} = [1, 2], 
              B1ref::Float64 = 1.0)

Performs hyperparameter optimization using Bayesian Optimization and Hyperband (BOHB).

Arguments

  • spins::Vector{<:Spins}: Vector of spin systems for the optimization process.
  • gp::GrapeParams: GRAPE algorithm parameters, including cost function and fields to optimize.
  • Tc::LinRange: Range for time control points for spline interpolation.
  • max_iter::Int: Maximum number of optimization iterations.
  • i::Int=5: Number of optimization evaluations.
  • poly_start::Vector{Float64}: Initial values for polynomial learning rate.
  • poly_degree::Vector{Int}: Degrees for the polynomial learning rate.
  • B1ref::Float64=1.0: Reference B1 field amplitude.

Returns

  • A BOHB optimization object with optimized hyperparameter configurations and costs.
source
GrapeMR.bruker_normalized_amplitudes_and_phasesMethod
bruker_normalized_amplitudes_and_phases(cf::ControlField)

Calculates amplitudes and phases normalized to 100 and 180 deg for Bruker implementation on TopSpin. All negative phase values are added a 360deg phase

Arguments

  • cf::ControlField: The control field object containing the RF waveform.

Outputs:

A tuple of normalized amplitudes and phases:

  • norm_amp: Maximum RF amplitude in Hertz (Hz).
  • norm_phase: Normalize phases in degrees.
source
GrapeMR.color_paletteMethod
color_palette(var::Union{AbstractArray, Int})

Generates a color palette using the rainbow color scheme for the specified array of indices or integer count.

Arguments

  • var::Union{AbstractArray, Int}: Array or integer specifying the number of colors required.

Returns

  • A ColorScheme object with colors matching the length or count of var.
source
GrapeMR.cost_offsetsMethod
cost_offsets(control_field::ControlField, spin_system::Spin, offsets::Vector{Float64})

Calculates the cost values for a given control field (control_field) and spin system (spin_system) over a range of B0 inhomogeneities (offsets).

Arguments

  • control_field::ControlField: The control field object that contains the RF pulse sequence.
  • spin_system::Spin: The spin system for which the cost function is to be calculated.
  • offsets::Vector{Float64}: A vector of offset frequencies (in Hz) representing B0 inhomogeneities.

Outputs

  • cost_values::Vector{Float64}: A vector containing the cost values corresponding to each offset in B0_values.
source
GrapeMR.countour_costMethod
countour_cost(cost_matrix::Matrix{Float64})

Creates a contour plot of the cost function values over the range of B0 and B1 inhomogeneities.

Arguments

  • cost_matrix::Matrix{Float64}: A 2D matrix representing the cost values for combinations of B0 and B1 inhomogeneities.

Outputs

  • c::Plot: A contour plot of the cost function map.
source
GrapeMR.create_cost_matrixMethod
create_cost_matrix(control_field::ControlField, spin_system::Spin, offsets::Vector{Float64}, b1_inhomogeneities::Vector{Float64})

Generates a cost matrix representing the cost function values across a range of B0 and B1 inhomogeneities.

Arguments

  • control_field::ControlField: The control field object that contains the RF pulse sequence
  • spin_system::Spin: The spin system for which the cost function is to be calculated.
  • offsets::Vector{Float64}: A vector of offset frequencies (in Hz) representing B0 inhomogeneities.
  • b1_inhomogeneities::Vector{Float64}: A vector representing B1 inhomogeneity percentages.

Outputs

  • cost_matrix::Matrix{Float64}: A 2D matrix where each element represents the cost value for a particular combination of B0 and B1 values. The matrix is normalized by its maximum value.
source
GrapeMR.create_folderMethod
create_folder(path::String)

Creates a folder at the specified path if it doesn't already exist.

Arguments

  • path::String: The path where the folder will be created.

Output

Returns the path to the created folder. If the folder already exists, it simply returns the path.

source
GrapeMR.create_splineMethod
create_spline(spline_time::AbstractArray, control_time_vals::AbstractArray, B1_random_vals::AbstractArray)

Creates a cubic spline interpolation for the control field based on the specified time and amplitude values.

Arguments

  • spline_time::AbstractArray: Array of time points for the spline.
  • control_time_vals::AbstractArray: Array of control time points for evaluation.
  • B1_random_vals::AbstractArray: Array of B1 amplitude values for spline generation.

Returns

  • Array of interpolated control field values at each control time.
source
GrapeMR.dynamicsMethod
dynamics(cf::ControlField, spin::Spins)

Calculates and returns the Isochromat object representing the spin dynamics under a given control field.

Arguments

  • cf::ControlField: Control field affecting the spin dynamics.
  • spin::Spin: Spin object for which dynamics are computed.

Returns

  • iso::Isochromat: Isochromat containing the calculated magnetization dynamics.
source
GrapeMR.export_brukerMethod
export_bruker(go::GrapeOutput; folder_path = pwd())

Exports the GRAPE optimization results to a Bruker-compatible .exc file format for use in TopSpin software.

Arguments

  • go::GrapeOutput: The output struct from GRAPE optimization containing control field data.
  • folder_path::String=pwd(): The directory path where the .exc file will be saved. Defaults to the current working directory.

File Structure

  • The generated file includes metadata such as date, time, amplitude, and phase ranges.
  • Amplitude and phase data points are formatted for compatibility with Bruker TopSpin.

Outputs

  • Saves a .exc file named based on the GRAPE output and current configuration in the specified folder.
source
GrapeMR.forward_propagation!Method
forward_propagation!(M::AbstractMatrix, cf::ControlField, s::Spins)

In-place version of forward propagation that updates the provided magnetization matrix.

Arguments

  • M::AbstractMatrix: Matrix to store the forward-propagated magnetization (4xN).
  • cf::ControlField: Struct containing control field parameters.
  • s::Spins: Spin struct with initial magnetization and relaxation parameters.

Returns

  • M: Updated magnetization matrix (4xN) with forward propagation results.
source
GrapeMR.forward_propagationMethod
forward_propagation(cf::ControlField, s::Spins)

Performs forward propagation for the magnetization vector under the influence of control fields.

Arguments

  • cf::ControlField: Struct containing control field parameters.
  • s::Spins: Spin struct with initial magnetization and relaxation parameters.

Returns

  • M::Matrix{Float64}: Magnetization matrix (4xN) at each time step.
source
GrapeMR.gaussian_RFMethod
gaussian_RF(N::Int, t_c::Float64, B1ref::Float64)

Generates a Gaussian-shaped RF pulse.

Arguments

  • N::Int: Number of time points.
  • t_c::Float64: Duration of the pulse in seconds.
  • B1ref::Float64: Reference amplitude for scaling the pulse.

Returns

  • A ControlField struct with Gaussian-distributed B1x and B1y components.
source
GrapeMR.get_resources_configurations_hbandMethod
get_resources_configurations_hband(η::Int, R::Float64)

Calculates the resources and configurations for each stage of the Hyperband algorithm.

Arguments

  • η::Int: Downsampling factor, typically used in Hyperband.
  • R::Float64: Maximum resource budget (e.g., number of iterations).

Returns

  • Prints the configurations n and resources r for each Hyperband outer loop.
source
GrapeMR.get_target_propertiesMethod
get_target_properties(s, colors)

Determines the color and label for a given spin target ("max", "min", or default).

Arguments

  • s::Spin: A spin object with target and label properties.
  • colors::Vector: A vector of colors from the color palette.

Returns

  • A tuple (color, label) with the color and label for the given spin target.
source
GrapeMR.gradientMethod
gradient(χ::Matrix{Float64}, M::Matrix{Float64}, H::Matrix)

Calculates the gradient of the cost function with respect to the Hamiltonian for each time step.

Arguments

  • χ::Matrix{Float64}: Adjoint state matrix.
  • M::Matrix{Float64}: Forward propagation matrix.
  • H::Matrix: Hamiltonian matrix.

Returns

  • grad::Matrix{Float64}: Gradient of the cost function, as a 1xN matrix.
source
GrapeMR.grapeMethod
grape(p::Parameters, cf::ControlField, spins::Vector{<:Spins})

Executes the GRAPE algorithm to optimize control fields for spin dynamics in an NMR/MRI system.

Arguments

  • p::Parameters: Struct containing optimization parameters, including maximum iterations and cost function.
  • cf::ControlField: Initial control field, typically as a spline function.
  • spins::Vector{<:Spins}: Vector of spins included in the optimization process.

Returns

  • grape_output::GrapeOutput: Struct containing optimization results, including optimized control fields, spin dynamics, and cost function values.
source
GrapeMR.hard_RFMethod
hard_RF(N::Int, t_c::Float64, B1ref::Float64)

Generates a hard RF pulse with constant amplitude in the x-axis and zero amplitude in the y-axis.

Arguments

  • N::Int: Number of time points.
  • t_c::Float64: Duration of the pulse in seconds.
  • B1ref::Float64: Amplitude of the RF pulse.

Returns

  • A ControlField struct with constant B1x and zero B1y components.
source
GrapeMR.hband_hyperoptMethod
hband_hyperopt(spins::Vector{<:Spins}, 
               gp::GrapeParams, 
               Tc::LinRange, 
               max_iter::Int; 
               B1ref::Float64 = 1.0)

Executes Hyperband optimization for selecting hyperparameters.

Arguments

  • spins::Vector{<:Spins}: Vector of spin systems for the optimization process.
  • gp::GrapeParams: GRAPE algorithm parameters, including cost function and fields to optimize.
  • Tc::LinRange: Range for time control points for spline interpolation.
  • max_iter::Int: Maximum number of optimization iterations.
  • B1ref::Float64=1.0: Reference B1 field amplitude.

Returns

  • A Hyperband optimization object with configurations, costs, and history.
source
GrapeMR.initialize_plotMethod
initialize_plot(title, xlabel, ylabel; zlabel=nothing)

Creates a plot object with common settings and title and axis labels.

Arguments

  • title::String: Title of the plot.
  • xlabel::String: Label for the x-axis.
  • ylabel::String: Label for the y-axis.
  • zlabel::Union{String, Nothing}: Label for the z-axis (optional, only for 3D plots).

Returns

  • p: A plot object with preconfigured labels and settings.
source
GrapeMR.plot_configMethod
plot_config()

Sets up a default plot configuration for consistency in appearance across plots.

Returns

  • p: A plot object with font sizes, frame style, and grid style preconfigured.
source
GrapeMR.plot_control_fieldsMethod
plot_control_fields(cf::ControlField; unit::String = "Hz")

Plot the control fields of an RF pulse in different units.

Arguments

  • cf::ControlField: The control field object containing the RF waveform data (B1x, B1y components).
  • unit::String: The unit for the plot. Supported units are:
    • "Hz": Plots B1x and B1y in Hertz.
    • "rad/s": Converts the control fields to rads/sec and plots the amplitude and phase.
    • "Tesla": Converts the control fields to Tesla using the gyromagnetic ratio for ¹H and plots in microtesla (µT).

Returns

  • A plot object displaying the control fields in the specified units with amplitude and phase information.
source
GrapeMR.plot_control_fields_phase_shiftMethod
plot_control_fields_phase_shift(cf::ControlField; ψ::Float64 = π)

Plot the control fields of an RF pulse after applying a phase shift.

Arguments

  • cf::ControlField: The control field object containing the RF waveform data (B1x, B1y components).
  • ψ::Float64: The phase shift to apply to the control fields (default is π radians).

Returns

  • A plot object displaying the control fields with the specified phase shift.
source
GrapeMR.plot_cost_offsetMethod
plot_cost_offset(cost_profile::Vector{Float64}, B0_values::Vector{Float64})

Creates a plot of the cost values as a function of B0 offset frequencies.

Arguments

  • cost_profile::Vector{Float64}: A vector containing the cost values to be plotted.
  • B0_values::Vector{Float64}: A vector of offset frequencies (in Hz) corresponding to the cost values.

Outputs

  • p::Plot: A plot object displaying the cost function offset profile.
source
GrapeMR.plot_cost_valuesMethod
plot_cost_values(cost::Vector{Float64}, gp::GrapeParams)

Plot the cost function convergence over iterations during the GRAPE optimization process.

Arguments

  • cost::Vector{Float64}: A vector containing the cost values at each iteration.
  • gp::GrapeParams: The parameters of the GRAPE optimization, including the cost function name.

Returns

  • A plot object displaying the convergence of the cost function over the iterations.
source
GrapeMR.plot_longitudinal_timeMethod
plot_longitudinal_time(isos::Vector{Isochromat}, t::Float64)

Plots the longitudinal magnetization (Mz) as a function of time for a set of isochromats.

Arguments

  • isos::Vector{Isochromat}: A vector of Isochromat objects containing magnetization and spin data.
  • t::Float64: Total time duration for the plot.

Returns

  • pLongTime: A plot object displaying Mz as a function of time for each isochromat.
source
GrapeMR.plot_magnetization_2DMethod
plot_magnetization_2D(isos::Vector{Isochromat})

Plots the 2D magnetization, showing transverse (Mxy) vs longitudinal (Mz) magnetization.

Arguments

  • isos::Vector{Isochromat}: A vector of Isochromat objects containing magnetization and spin data.

Returns

  • pMag2D: A plot object displaying the transverse (Mxy) vs longitudinal (Mz) magnetization.
source
GrapeMR.plot_magnetization_3DMethod
plot_magnetization_3D(isos::Vector{Isochromat})

Plots the 3D magnetization (Mx, My, Mz) for a set of isochromats.

Arguments

  • isos::Vector{Isochromat}: A vector of Isochromat objects containing magnetization and spin data.

Returns

  • pMag3D: A plot object displaying the 3D magnetization components (Mx, My, Mz) for each isochromat.
source
GrapeMR.plot_magnetization_control_fieldMethod
plot_magnetization_control_field(cf::ControlField, isos::Vector{Isochromat})

Plot the magnetization trajectory for a set of isochromats and the corresponding control field.

Arguments

  • cf::ControlField: The control field object containing the RF waveform data (B1x, B1y components).
  • isos::Vector{Isochromat}: A vector of isochromat objects representing different spin systems.

Returns

  • A plot object displaying:
    1. The magnetization trajectory in the transverse plane.
    2. The control field amplitude and phase over time.
source
GrapeMR.plot_magnetization_timeMethod
plot_magnetization_time(iso::Isochromat, t::Float64)

Plots the time evolution of magnetization components (Mx, My, Mz) for a single isochromat.

Arguments

  • iso::Isochromat: An Isochromat object containing magnetization dynamics.
  • t::Float64: Total time duration of the plot.

Returns

  • pMagTime: A plot object displaying Mx, My, and Mz as functions of time.
source
GrapeMR.plot_ss_offset_profileMethod
plot_ss_offset_profile(ss_spin::Vector{GrapeMR.SteadyState})

Plot the bSSFP offset frequency profile showing magnitude and phase.

Arguments

  • ss_spin: Vector of steady-state spins at different offset frequencies

Returns

  • Tuple of three plots: (magnitude profile, phase profile, combined plot)

Example

plots = plotssoffsetprofile(ssspins) display(plots[1]) # Show

source
GrapeMR.plot_transverse_magnetizationMethod
plot_transverse_magnetization(isos::Vector{Isochromat})

Plots the transverse magnetization (Mx and My) for a set of isochromats, with colors indicating different spin targets.

Arguments

  • isos::Vector{Isochromat}: A vector of Isochromat objects containing magnetization and spin information.

Returns

  • pTrans: A plot object displaying the transverse magnetization (Mx vs My) for each isochromat.
source
GrapeMR.plot_transverse_timeMethod
plot_transverse_time(isos::Vector{Isochromat}, t::Float64)

Plots the transverse magnetization (Mxy) as a function of time for a set of isochromats.

Arguments

  • isos::Vector{Isochromat}: A vector of Isochromat objects containing magnetization and spin data.
  • t::Float64: Total time duration for the plot.

Returns

  • pTransTime: A plot object displaying Mxy as a function of time for each isochromat.
source
GrapeMR.random_hyperoptMethod
random_hyperopt(spins::Vector{<:Spins}, 
                gp::GrapeParams, 
                Tc::LinRange, 
                max_iter::StepRange; 
                i::Int = 50, 
                poly_start::Vector{Float64} = [5e-1, 1e-1, 1e-2], 
                poly_degree::Vector{Int} = [1, 2], 
                B1ref::Float64 = 1.0)

Performs random sampling for hyperparameter optimization.

Arguments

  • spins::Vector{<:Spins}: A vector of spin systems for the optimization process.
  • gp::GrapeParams: GRAPE algorithm parameters, including cost function and fields to optimize.
  • Tc::LinRange: Range for time control points for spline interpolation.
  • max_iter::StepRange: Range for maximum optimization iterations.
  • i::Int=50: Number of random samples to evaluate.
  • poly_start::Vector{Float64}: Initial values for polynomial learning rate.
  • poly_degree::Vector{Int}: Degrees for the polynomial learning rate.
  • B1ref::Float64=1.0: Reference B1 field amplitude.

Returns

  • A hyperparameter optimization object with randomly sampled configurations and their costs.
source
GrapeMR.run_cost_analysisMethod
run_cost_analysis(grape_output::GrapeOutput, offset::Float64, b1_inhomogeneity_percent::Int)

Runs a complete cost analysis for the specified GrapeOutput object over a range of B0 and B1 inhomogeneities, and generates plots to visualize the cost function.

Arguments

  • control_field::ControlField: The control field object that contains the RF pulse sequence.
  • spin_system::Spin: The spin system for which the cost function is to be calculated.
  • offset::Float64: The range of offset frequencies (in Hz) for the B0 inhomogeneities.
  • b1_inhomogeneity_percent::Int: The percentage of B1 inhomogeneity to consider for the analysis.

Workflow

  1. Generates a range of B0 values from -offset to offset.
  2. Computes cost values for the specified B0 and B1 inhomogeneities using cost_offsets and create_cost_matrix.
  3. Plots the cost function offset profile using plot_cost_offset.
  4. Generates a heatmap of the cost matrix using heatmap_cost.
  5. Creates a contour plot of the cost matrix using countour_cost.

Example

runcostanalysis(grape_output, 50.0, 30)

source
GrapeMR.run_grape_optimizationMethod
run_grape_optimization(config_path::String)

Runs the GRAPE optimization process based on a TOML configuration file. Initializes spins, sets up parameters, executes optimization, and optionally saves and plots the results.

Arguments

  • config_path::String: Path to the TOML configuration file containing parameters for spins, optimization, and control fields.

Configuration File Structure

The TOML configuration file should include sections like: - spins: Defines spin properties. - grape_parameters: Parameters for the GRAPE optimization. - optimization_parameters: Parameters for the hyperparameter optimization. - control_field: Control field specifications. - save_files: Settings for saving outputs. - plot: Plot settings.

Returns

  • Produces and saves results depending on configuration settings, including optimized control fields, cost values, optional Bruker export, and plots.

Example

```julia rungrapeoptimization("path/to/config.toml")

source
GrapeMR.save_grape_dataMethod
save_output_data(go::GrapeMR.GrapeOutput; folder_path = pwd())

Save data related to GRAPE optimization into a folder organized by date.

Arguments

  • go::GrapeMR.GrapeOutput: The output from a GRAPE optimization process.
  • folder_path::String = pwd(): The folder path where data will be saved. Defaults to the current working directory.

Output

Returns the full path to the folder where the data was saved.

Saved Files

  • grape_output.jld2: Contains the entire GrapeOutput struct in JLD2 format.
  • dict_cost_values.csv: CSV file containing the cost values from the optimization process.
  • dict_control_field.csv: CSV file containing control field values (B1x, B1y, Bz, and RF time).
  • dict_iso_spins.csv: CSV file containing isochromat spin parameters.

If no path is provided, it saves the files inside the folder where the package was installed folder name format : yyyy-mm-dd

source
GrapeMR.save_hyperopt_dataMethod
save_output_data(bohb::GrapeMR.GrapeOutput; folder_path = pwd())

Save data related to BOHB optimization into a folder organized by date.

Arguments

  • bohb::GrapeMR.GrapeOutput: The output from a BOHB optimization process.
  • folder_path::String = pwd(): The folder path where data will be saved. Defaults to the current working directory.

Output

Returns the full path to the folder where the data was saved.

Saved Files

  • bohb.jld2: Contains BOHB result in JLD2 format.
source
GrapeMR.sinc_RFMethod
sinc_RF(N::Int, t_c::Float64, B1ref::Float64; α=π/2)

Generates a sinc RF pulse with a specified flip angle. Bandwidth hardcoded to 100 Hz. B1x and B1y have a π/2 phase difference

Arguments

  • N::Int: Number of time points.
  • t_c::Float64: Duration of the pulse in seconds.
  • B1ref::Float64: Reference amplitude for scaling the pulse.
  • α::Float64=π/2: Flip angle in radians.

Returns

  • A ControlField struct with B1x, B1y, and Bz components generated as sinc functions.
source
GrapeMR.spline_RFMethod
spline_RF(N::Int, t_c::Float64, B1ref::Float64)

Generates a cubic spline-based RF pulse.

Arguments

  • N::Int: Number of time points.
  • t_c::Float64: Duration of the pulse in seconds.
  • B1ref::Float64: Reference amplitude for scaling the pulse.

Returns

  • A ControlField struct with B1x, B1y, and Bz components generated using spline interpolation.
source
GrapeMR.steady_stateMethod
steady_state(s::GrapeMR.SteadyState)

Calculates the steady-state signal of a bSSFP sequence for a given set of parameters. The function calculates the magnetization evolution using repeated RF excitations and free precession steps until the steady state is reached.

Arguments

  • s::GrapeMR.SteadyState: A struct containing TR (Repetition Time), TE (Echo Time), T1 (Spin-lattice relaxation time), T2 (Spin-spin relaxation time), α (Flip angle in radians), Δϕ (Phase cycling increment), M_init (Initial magnetization vector), and B0inho (B0 inhomogeneity).

Outputs

  • signal::Vector{ComplexF64}: The complex steady-state signal of the spin system.
source
GrapeMR.steady_state_geometricMethod
steady_state_geometric(s::GrapeMR.SteadyState)

Calculates the transverse steady-state magnetization using a geometric solution for a bSSFP sequence. Uses a geometric approach to derive the transverse component of the steady-state magnetization.

Arguments

  • s::GrapeMR.SteadyState: Struct containing the sequence parameters and spin system properties.

Outputs

  • Mxy::Vector{Float64}: The steady-state transverse magnetization magnitude Mxy
source
GrapeMR.steady_state_geometric_MzMethod
steady_state_geometric_Mz(s::GrapeMR.SteadyState)

Calculates the longitudinal steady-state magnetization using a geometric solution for a bSSFP sequence. Uses a geometric approach to derive the longitudinal component of the steady-state magnetization.

Arguments

  • s::GrapeMR.SteadyState: Struct containing the sequence parameters and spin system properties.

Outputs

  • Mz::Vector{Float64}: The steady-state longitudinal magnetization magnitude Mz
source
GrapeMR.update!Method
update!(cf::ControlField, ∇xy::Tuple, ϵ::Float64)

Updates the control fields based on the calculated gradient and a learning rate.

Arguments

  • cf::ControlField: Control field struct to be updated.
  • ∇xy::Tuple{Matrix{Float64}, Matrix{Float64}}: Gradients for the x and y components of the field.
  • ϵ::Float64: Learning rate for gradient descent.

Returns

  • (u1x, u1y): Updated x and y control fields.
source
GrapeMR.ControlFieldType
ControlField{T, M1, Mz}

Represents the RF control field parameters for an NMR/MRI sequence.

Fields

  • B1x::M1: Matrix for the x-component of the RF field.
  • B1y::M1: Matrix for the y-component of the RF field.
  • B1_ref::T: Reference amplitude of the RF field.
  • Bz::Mz: Matrix for the z-component of the magnetic field.
  • t_control::T: Total control time for the sequence.
source
GrapeMR.GrapeOutputType
GrapeOutput{T, M1, Mz, F}

Stores the output of GRAPE optimization.

Fields

  • isochromats::Vector{Isochromat}: Vector of isochromats containing the spin dynamics.
  • control_field::ControlField{T, M1, Mz}: Optimized control field after GRAPE optimization.
  • cost_values::Vector{Float64}: Sequence of cost function values at each iteration.
  • params::Parameters{F}: Struct containing GRAPE parameters and optimization parameters.
source
GrapeMR.GrapeParamsType
GrapeParams{F}

Encapsulates key parameters for the GRAPE optimization process.

Fields

  • N::Int64: Number of time steps in the control sequence.
  • cost_function::F: Function used to compute the optimization cost.
  • fields_opt::Dict{String, Bool}: Dictionary indicating which fields to optimize (e.g., B1x, B1y).
source
GrapeMR.IsochromatType
Isochromat{S}

Represents an isochromat, combining magnetization data with a specific spin configuration.

Fields

  • magnetization::Magnetization: Magnetization data associated with a collection or single isochromat.
  • spin::S: Spin configuration for a collection or single isochromat.
source
GrapeMR.MagnetizationType
Magnetization{T, M}

Represents the magnetization dynamics of a spin system.

Fields

  • dynamics::M: Time-evolution data for magnetization, either as a 4xN matrix.
source
GrapeMR.OptimizationParamsType
OptimizationParams

Defines parameters for polynomial-based decay for ϵ in the GRAPE algorithm.

Fields

  • poly_start::Float64: Initial polynomial coefficient.
  • poly_degree::Int: Degree of the polynomial.
  • max_iter::Int: Maximum number of optimization iterations.
source
GrapeMR.ParametersType
Parameters{F}

Combines GRAPE and optimization parameters for the full optimization process.

Fields

  • grape_params::GrapeParams{F}: GRAPE-specific optimization parameters.
  • opt_params::OptimizationParams: General optimization parameters.
source
GrapeMR.SpinType
Spin <: Spins

Gets relaxation parameters and inhomogeneity effects of the spin system.

Fields

  • M_init::Vector{Float64}: Initial magnetization vector.
  • T1::Float64: Longitudinal relaxation time.
  • T2::Float64: Transverse relaxation time.
  • B0inho::Float64: B0 inhomogeneity.
  • B1inho::Float64: B1 inhomogeneity.
  • target::String: Target magnetization state.
  • label::String: Label for identifying the spin.
  • Nspins::Float64: Number of spins in this configuration.
source
GrapeMR.SpinMethod
Spin(M_ini, T1, T2, B0, B1, targets, labels)

Constructs a collection of Spin instances with specified parameters.

Arguments

  • M_ini::Vector{Float64}: Initial magnetization vector.
  • T1::Vector{Float64}: Longitudinal relaxation times.
  • T2::Vector{Float64}: Transverse relaxation times.
  • B0::Vector{Float64}: Array of B0 inhomogeneity values.
  • B1::Vector{Float64}: Array of B1 inhomogeneity values.
  • targets::Vector{String}: Target states for each spin configuration.
  • labels::Vector{String}: Labels for each spin configuration.

Returns

  • A vector of Spin instances covering all combinations of the provided parameters.
source