Package specific_analyses :: Package relax_disp :: Module data
[hide private]
[frames] | no frames]

Module data

source code

Module for handling relaxation dispersion data within the relax data store.

Ordering of data

The dispersion data model is based on the following concepts, in order of importance:

Indices

The data structures used in this module consist of many different index types which follow the data ordering above. These are abbreviated as:

Functions [hide private]
float
average_intensity(spin=None, exp_type=None, frq=None, offset=None, point=None, time=None, sim_index=None, error=False)
Return the average peak intensity for the spectrometer frequency, dispersion point, and relaxation time.
source code
List of dict()
calc_rotating_frame_params(spin=None, spin_id=None, fields=None, verbosity=0)
Calculates and rotating frame parameters, calculated from:
source code
 
check_intensity_errors(spin_id=None)
Check if intensity errors have already been calculated by the user.
source code
int
count_exp()
Count the number of experiments present.
source code
int
count_frq()
Count the number of spectrometer frequencies present.
source code
int
count_relax_times(exp_type=None, frq=None, offset=None, point=None, ei=None)
Count the number of relaxation times present.
source code
 
count_spins(spins=None)
Count the number of selected spins in the spin cluster.
source code
 
cpmg_setup(spectrum_id=None, cpmg_frq=None, ncyc_even=True)
Set the CPMG frequency associated with a given spectrum.
source code
str, float
decompose_r20_key(key=None)
Decompose the unique R20 key into the experiment type and spectrometer frequency.
source code
list of str
find_intensity_keys(exp_type=None, frq=None, offset=None, point=None, time=None, raise_error=True)
Return the key corresponding to the spectrometer frequency, dispersion point, and relaxation time.
source code
str
generate_r20_key(exp_type=None, frq=None)
Generate the unique R20 key from the experiment type and spectrometer frequency.
source code
str
get_curve_type(id=None)
Return the unique curve type.
source code
str
get_exp_type(id=None)
Return the experiment type for the given ID.
source code
bool
has_cpmg_exp_type()
Determine if the current data pipe contains CPMG experiment types.
source code
bool
has_disp_data(spins=None, spin_ids=None, exp_type=None, frq=None, offset=None, point=None)
Determine if dispersion data exists for the given data combination.
source code
bool
has_exponential_exp_type()
Determine if the current data pipe contains exponential curves.
source code
bool
has_fixed_time_exp_type()
Determine if the current data pipe contains fixed time data.
source code
bool
has_proton_mmq_cpmg()
Determine if the current data pipe contains either proton SQ or MQ (MMQ) CPMG data.
source code
bool
has_proton_mq_cpmg()
Determine if the current data pipe contains proton MQ CPMG data.
source code
bool
has_proton_sq_cpmg()
Determine if the current data pipe contains proton SQ CPMG data.
source code
bool
has_r1rho_exp_type()
Determine if the current data pipe contains R1rho experiment types.
source code
 
insignificance(level=0.0)
Deselect all spins with insignificant dispersion profiles.
source code
boolean, rank-4 list of numpy rank-1 float arrays, rank-4 list of numpy rank-1 float arrays, rank-3 list of numpy rank-1 float arrays, rank-4 list of numpy rank-1 float arrays, rank-2 list of numpy rank-1 float arrays, rank-4 list of numpy rank-1 float arrays, rank-4 list of numpy rank-1 float arrays, rank-4 list of numpy rank-1 float arrays
interpolate_disp(spin=None, spin_id=None, si=None, num_points=None, extend_hz=None, relax_times=None)
Interpolate function for 2D Grace plotting function for the dispersion curves.
source code
boolean, rank-4 list of numpy rank-1 float arrays, rank-4 list of numpy rank-1 float arrays, rank-3 list of numpy rank-1 float arrays, rank-4 list of numpy rank-1 float arrays, rank-2 list of numpy rank-1 float arrays, rank-4 list of numpy rank-1 float arrays, rank-4 list of numpy rank-1 float arrays, rank-4 list of numpy rank-1 float arrays
interpolate_offset(spin=None, spin_id=None, si=None, num_points=None, extend_ppm=None, relax_times=None)
Interpolate function for 2D Grace plotting function for the dispersion curves, interpolating through spin-lock offset in rad/s.
source code
bool
is_cpmg_exp_type(id=None)
Determine if the given spectrum ID corresponds to a CPMG experiment type.
source code
bool
is_r1_optimised(model=None)
Should R1 values be optimised?
source code
bool
is_r1rho_exp_type(id=None)
Determine if the given spectrum ID corresponds to a R1rho experiment type.
source code
list of str
loop_cluster(skip_desel=True)
Loop over the spin groupings for one model applied to multiple spins.
source code
str, (int)
loop_exp(return_indices=False)
Generator method for looping over all experiment types.
source code
str, float, (int, int)
loop_exp_frq(return_indices=False)
Generator method for looping over the exp and frq data.
source code
str, float, float, (int, int, int)
loop_exp_frq_offset(return_indices=False)
Generator method for looping over the exp, frq, and offset data.
source code
str, float, float, float, (int, int, int, int)
loop_exp_frq_offset_point(return_indices=False)
Generator method for looping over the exp, frq, offset, and point data.
source code
str, float, float, float, (int, int, int, int)
loop_exp_frq_offset_point_time(return_indices=False)
Generator method for looping over the exp, frq, offset, and point data.
source code
str, float, float, (int, int, int)
loop_exp_frq_point(return_indices=False)
Generator method for looping over the exp, frq, and point data.
source code
str, float, float, float, (int, int, int, int)
loop_exp_frq_point_time(return_indices=False)
Generator method for looping over the exp, frq, point, and time data.
source code
float, (int)
loop_frq(return_indices=False)
Generator method for looping over all spectrometer frequencies.
source code
float, float, (int, int)
loop_frq_offset(exp_type=None, return_indices=False)
Generator method for looping over the spectrometer frequencies and dispersion points.
source code
float, float, (int, int)
loop_frq_point(exp_type=None, return_indices=False)
Generator method for looping over the spectrometer frequencies and dispersion points.
source code
str
loop_frq_offset_point_key(exp_type=None)
Generator method for looping over the spectrometer frequencies, spin-lock offsets and dispersion points (returning the key).
source code
float, float, float
loop_frq_point_time(exp_type=None, return_indices=False)
Generator method for looping over the spectrometer frequencies, dispersion points, and relaxation times.
source code
float, (int)
loop_offset(exp_type=None, frq=None, return_indices=False)
Generator method for looping over the spin-lock offset values.
source code
float, float, (int, int)
loop_offset_point(exp_type=None, frq=None, skip_ref=True, return_indices=False)
Generator method for looping over the offsets and dispersion points.
source code
float, (int)
loop_point(exp_type=None, frq=None, offset=None, time=None, skip_ref=True, return_indices=False)
Generator method for looping over the dispersion points.
source code
str
loop_spectrum_ids(exp_type=None, frq=None, offset=None, point=None, time=None)
Generator method for selectively looping over the spectrum IDs.
source code
float
loop_time(exp_type=None, frq=None, offset=None, point=None, return_indices=False)
Generator method for looping over the relaxation times.
source code
int
num_exp_types()
Count the number of experiment types present.
source code
 
pack_back_calc_r2eff(spin=None, spin_id=None, si=None, back_calc=None, proton_mmq_flag=False)
Store the back calculated R2eff data for the given spin.
source code
 
plot_disp_curves(dir=None, y_axis='r2_eff', x_axis='disp', num_points=1000, extend_hz=500.0, extend_ppm=500.0, interpolate='disp', force=False) source code
 
plot_disp_curves_to_file(file_name_ini=None, dir=None, y_axis=None, x_axis=None, interpolate=None, num_points=None, extend_hz=None, extend_ppm=None, force=None, proton_mmq_flag=None) source code
 
plot_exp_curves(file=None, dir=None, force=None, norm=None)
Custom 2D Grace plotting function for the exponential curves.
source code
 
r20_from_min_r2eff(force=True, verbosity=1)
Set the R20 values to the minimum R2eff values.
source code
 
r2eff_read(id=None, file=None, dir=None, disp_frq=None, offset=None, spin_id_col=None, mol_name_col=None, res_num_col=None, res_name_col=None, spin_num_col=None, spin_name_col=None, data_col=None, error_col=None, sep=None)
Read R2eff/R1rho values directly from a file whereby each row corresponds to a different spin.
source code
 
r2eff_read_spin(id=None, spin_id=None, file=None, dir=None, disp_point_col=None, offset_col=None, data_col=None, error_col=None, sep=None)
Read R2eff/R1rho values from file whereby each row is a different dispersion point.
source code
 
randomise_R1(spin=None, ri_id=None, N=None)
Randomise the R1 data for the given spin for use in the Monte Carlo simulations.
source code
 
relax_time(time=0.0, spectrum_id=None)
Set the relaxation time period associated with a given spectrum.
source code
rank-2 list of numpy rank-1 float64 arrays
return_cpmg_frqs(ref_flag=True)
Return the list of nu_CPMG frequencies.
source code
numpy rank-1 float64 array
return_cpmg_frqs_single(exp_type=None, frq=None, offset=None, time=None, ref_flag=True)
Return the list of nu_CPMG frequencies.
source code
 
return_grace_data_vs_disp(y_axis=None, x_axis=None, interpolate=None, exp_type=None, ei=None, current_spin=None, spin_id=None, si=None, back_calc=None, cpmg_frqs_new=None, spin_lock_nu1_new=None, chemical_shifts=None, offsets_inter=None, tilt_angles_inter=None, Delta_omega_inter=None, w_eff_inter=None, interpolated_flag=None, graph_index=None, data=None, set_labels=None, set_colours=None, x_axis_type_zero=None, symbols=None, symbol_sizes=None, linetype=None, linestyle=None, axis_labels=None) source code
 
return_grace_data_vs_offset(y_axis=None, x_axis=None, interpolate=None, exp_type=None, ei=None, current_spin=None, spin_id=None, si=None, back_calc=None, cpmg_frqs_new=None, spin_lock_nu1_new=None, chemical_shifts=None, offsets_inter=None, tilt_angles_inter=None, Delta_omega_inter=None, w_eff_inter=None, interpolated_flag=None, graph_index=None, data=None, set_labels=None, set_colours=None, x_axis_type_zero=None, symbols=None, symbol_sizes=None, linetype=None, linestyle=None, axis_labels=None) source code
 
return_grace_file_name_ini(y_axis=None, x_axis=None, interpolate=None) source code
 
return_grace_x_y_axis_labels(y_axis=None, x_axis=None, exp_type=None, interpolate=None) source code
str, int, float, int, int
return_x_y_data_labels_settings(spin=None, data_type=None, y_axis=None, exp_type=None, frq=None, offset=None, point=None, interpolated_flag=None)
Return the X and Y labels and plot settings, according to selected axis to plot for.
source code
 
return_grace_x_y_point(data_type=None, y_axis=None, x_axis=None, interpolate=None, data_key=None, spin=None, back_calc=None, offset=None, point=None, r1=None, r1_err=None, w_eff=None, theta=None, err=False) source code
int
return_index_from_disp_point(value, exp_type=None)
Convert the dispersion point data into the corresponding index.
source code
int
return_index_from_exp_type(exp_type=None)
Convert the experiment type into the corresponding index.
source code
int
return_index_from_frq(value)
Convert the dispersion point data into the corresponding index.
source code
int
return_index_from_disp_point_key(key, exp_type=None)
Convert the dispersion point key into the corresponding index.
source code
str
return_key_from_di(mi=None, di=None)
Convert the dispersion point index into the corresponding key.
source code
rank-3 list of numpy rank-1 float arrays, rank-4 list of numpy rank-1 float arrays, rank-2 list of numpy rank-1 float arrays, rank-4 list of numpy rank-1 float arrays, rank-4 list of numpy rank-1 float arrays, rank-4 list of numpy rank-1 float arrays
return_offset_data(spins=None, spin_ids=None, field_count=None, spin_lock_offset=None, fields=None)
Return numpy arrays of the chemical shifts, offsets and tilt angles.
source code
str
return_param_key_from_data(exp_type=None, frq=0.0, offset=0.0, point=0.0)
Generate the unique key from the spectrometer frequency and dispersion point.
source code
numpy rank-2 float array
return_r1_data(spins=None, spin_ids=None, field_count=None, sim_index=None)
Return the R1 data structures for off-resonance R1rho experiments.
source code
numpy rank-2 float array
return_r1_err_data(spins=None, spin_ids=None, field_count=None, sim_index=None)
Return the R1 error data structures for off-resonance R1rho experiments.
source code
lists of numpy float arrays, lists of numpy float arrays, lists of numpy float arrays, numpy rank-2 int array
return_r2eff_arrays(spins=None, spin_ids=None, fields=None, field_count=None, sim_index=None)
Return numpy arrays of the R2eff/R1rho values and errors.
source code
numpy rank-2 float64 array
return_relax_times()
Return the list of relaxation times.
source code
rank-2 list of numpy rank-1 float64 arrays
return_spin_lock_nu1(ref_flag=True)
Return the list of spin-lock field strengths.
source code
numpy rank-1 float64 array
return_spin_lock_nu1_single(exp_type=None, frq=None, offset=None, ref_flag=True)
Return the list of spin-lock field strengths.
source code
float
return_value_from_frq_index(mi=None)
Return the spectrometer frequency corresponding to the frequency index.
source code
float
return_value_from_offset_index(ei=None, mi=None, oi=None)
Return the offset corresponding to the offset index.
source code
 
set_exp_type(spectrum_id=None, exp_type=None)
Select the relaxation dispersion experiment type performed.
source code
bool
spin_has_frq_data(spin=None, frq=None)
Determine if the spin has intensity data for the given spectrometer frequency.
source code
list of SpinContainer instances
spin_ids_to_containers(spin_ids)
Take the list of spin IDs and return the corresponding spin containers.
source code
 
spin_lock_field(spectrum_id=None, field=None)
Set the spin-lock field strength (nu1) for the given spectrum.
source code
 
spin_lock_offset(spectrum_id=None, offset=None)
Set the spin-lock offset (omega_rf) for the given spectrum.
source code
 
write_disp_curves(dir=None, force=None)
Write out the dispersion curves to text files.
source code
Variables [hide private]
  R20_KEY_FORMAT = '%s - %.8f MHz'
  Y_AXIS_R2_EFF = 'r2_eff'
  Y_AXIS_R2_R1RHO = 'r2_r1rho'
  X_AXIS_DISP = 'disp'
  X_AXIS_W_EFF = 'w_eff'
  X_AXIS_THETA = 'theta'
  INTERPOLATE_DISP = 'disp'
  INTERPOLATE_OFFSET = 'offset'
  COLOUR_ORDER = [4, 15, 2, 13, 11, 1, 3, 5, 6, 7, 8, 9, 10, 12,...
  __package__ = 'specific_analyses.relax_disp'

Imports: cos, pi, sin, sqrt, array, concatenate, float64, int32, max, ones, unique, zeros, F_OK, access, expanduser, gauss, search, sys, warn, EXP_TYPE_CPMG_DQ, EXP_TYPE_CPMG_MQ, EXP_TYPE_CPMG_PROTON_MQ, EXP_TYPE_CPMG_PROTON_SQ, EXP_TYPE_CPMG_SQ, EXP_TYPE_CPMG_ZQ, EXP_TYPE_DESC_CPMG_DQ, EXP_TYPE_DESC_CPMG_MQ, EXP_TYPE_DESC_CPMG_PROTON_MQ, EXP_TYPE_DESC_CPMG_PROTON_SQ, EXP_TYPE_DESC_CPMG_SQ, EXP_TYPE_DESC_CPMG_ZQ, EXP_TYPE_DESC_R1RHO, EXP_TYPE_LIST, EXP_TYPE_LIST_CPMG, EXP_TYPE_LIST_R1RHO, EXP_TYPE_R1RHO, MODEL_B14, MODEL_B14_FULL, MODEL_DPL94, MODEL_LIST_FIT_R1, MODEL_LIST_MMQ, MODEL_LIST_NUMERIC_CPMG, MODEL_LIST_R1RHO_FULL, MODEL_LIST_R1RHO_ON_RES, MODEL_MP05, MODEL_NOREX, MODEL_NS_R1RHO_2SITE, MODEL_PARAMS, MODEL_R2EFF, MODEL_TAP03, MODEL_TP02, PARAMS_R20, RelaxError, RelaxNoSpectraError, RelaxNoSpinError, RelaxSpinTypeError, isNaN, extract_data, get_file_path, open_write_file, strip, write_data, frequency_to_ppm, frequency_to_ppm_from_rad, frequency_to_rad_per_s, rotating_frame_params, periodic_table, write_xy_data, write_xy_header, script_grace2images, read_spin_data, write_spin_data, section, RelaxWarning, RelaxNoSpinWarning, check_mol_res_spin_data, exists_mol_res_spin_data, generate_spin_id_unique, generate_spin_string, return_spin, spin_loop, check_pipe, add_result_file, desel_spin, return_attached_protons, add_spectrum_id, error_analysis, check_frequency, get_frequency, value, specific_analyses, check_exp_type, check_interpolate_offset_cpmg_model, check_missing_r1, check_mixed_curve_types, S_IRWXU, S_IRGRP, S_IROTH, chmod, sep


Function Details [hide private]

average_intensity(spin=None, exp_type=None, frq=None, offset=None, point=None, time=None, sim_index=None, error=False)

source code 

Return the average peak intensity for the spectrometer frequency, dispersion point, and relaxation time.

This is for handling replicate peak intensity data.

Parameters:
  • spin (SpinContainer instance) - The spin container to average the peak intensities for.
  • exp_type (str) - The experiment type.
  • frq (float) - The spectrometer frequency.
  • offset (float) - The spin-lock or hard pulse offset.
  • point (float) - The dispersion point data (either the spin-lock field strength in Hz or the nu_CPMG frequency in Hz).
  • time (float) - The relaxation time period.
  • sim_index (None or int) - The simulation index. This should be None for the measured intensity values.
  • error (bool) - A flag which if True will average and return the peak intensity errors.
Returns: float
The average peak intensity value.

calc_rotating_frame_params(spin=None, spin_id=None, fields=None, verbosity=0)

source code 

Calculates and rotating frame parameters, calculated from:

  • The spectrometer frequency.
  • The spin-lock or hard pulse offset.
  • The dispersion point data (the spin-lock field strength in Hz).

    The return will be for each spin,

  • Rotating frame tilt angle ( theta = arctan(w_1 / Omega) ) [rad]
  • The average resonance offset in the rotating frame ( Domega = w_{pop_ave} - w_rf ) [rad/s]
  • Effective field in rotating frame ( w_eff = sqrt( Omega^2 + w_1^2 ) ) [rad/s]

    Calculations are mentioned in the manual

Parameters:
  • spin (SpinContainer instance) - The spin system specific data container
  • spin_id (None or str) - The spin ID string.
  • fields (rank-2 list of floats) - The spin-lock field strengths to use instead of the user loaded values - to enable interpolation. The dimensions are {Ei, Mi}.
  • verbosity (int) - A flag specifying to print calculations.
Returns: List of dict()
List with dict() of theta, Domega, w_eff and list of dict() keys.

check_intensity_errors(spin_id=None)

source code 

Check if intensity errors have already been calculated by the user.

Parameters:
  • spin_id (str) - The spin identification string.

count_exp()

source code 

Count the number of experiments present.

Returns: int
The experiment count

count_frq()

source code 

Count the number of spectrometer frequencies present.

Returns: int
The spectrometer frequency count

count_relax_times(exp_type=None, frq=None, offset=None, point=None, ei=None)

source code 

Count the number of relaxation times present.

Parameters:
  • exp_type (str) - The experiment type.
  • frq (float) - The spectrometer frequency in Hz.
  • offset (None or float) - The spin-lock or hard pulse offset value in ppm.
  • point (float) - The dispersion point data (either the spin-lock field strength in Hz or the nu_CPMG frequency in Hz).
  • ei (str) - The experiment type index.
Returns: int
The relaxation time count for the given experiment.

cpmg_setup(spectrum_id=None, cpmg_frq=None, ncyc_even=True)

source code 

Set the CPMG frequency associated with a given spectrum.

Parameters:
  • spectrum_id (str) - The spectrum identification string.
  • cpmg_frq (float) - The frequency, in Hz, of the CPMG pulse train.
  • ncyc_even (bool) - A flag which if True means that the number of CPMG blocks must be even. This is pulse sequence dependant.

decompose_r20_key(key=None)

source code 

Decompose the unique R20 key into the experiment type and spectrometer frequency.

Parameters:
  • key (str) - The unique R20 key.
Returns: str, float
The experiment and the spectrometer frequency in Hz.

find_intensity_keys(exp_type=None, frq=None, offset=None, point=None, time=None, raise_error=True)

source code 

Return the key corresponding to the spectrometer frequency, dispersion point, and relaxation time.

Parameters:
  • exp_type (str) - The experiment type.
  • frq (float) - The spectrometer frequency.
  • offset (None or float) - The optional offset value for off-resonance R1rho-type data.
  • point (float) - The dispersion point data (either the spin-lock field strength in Hz or the nu_CPMG frequency in Hz).
  • time (float) - The relaxation time period.
  • raise_error (bool) - A flag which if True will cause a RelaxError to be raised if no keys could be found.
Returns: list of str
The keys corresponding to the spectrometer frequency, dispersion point, and relaxation time.

generate_r20_key(exp_type=None, frq=None)

source code 

Generate the unique R20 key from the experiment type and spectrometer frequency.

Parameters:
  • exp_type (str) - The experiment type.
  • frq (float) - The spectrometer frequency in Hz.
Returns: str
The unique R20 key.

get_curve_type(id=None)

source code 

Return the unique curve type.

Parameters:
  • id (str) - The spectrum ID. If not supplied, then all data will be assumed.
Returns: str
The curve type - either 'fixed time' or 'exponential'.

get_exp_type(id=None)

source code 

Return the experiment type for the given ID.

Parameters:
  • id (str) - The spectrum ID.
Returns: str
The experiment type corresponding to the ID.

has_cpmg_exp_type()

source code 

Determine if the current data pipe contains CPMG experiment types.

Returns: bool
True if CPMG experiment types exist, False otherwise.

has_disp_data(spins=None, spin_ids=None, exp_type=None, frq=None, offset=None, point=None)

source code 

Determine if dispersion data exists for the given data combination.

Parameters:
  • spins (list of SpinContainer instances) - The list of spin containers in the cluster.
  • spin_ids (list of str) - The list of spin IDs for the cluster.
  • exp_type (str) - The experiment type.
  • frq (float) - The spectrometer frequency.
  • offset (None or float) - For R1rho-type data, the spin-lock offset value in ppm.
  • point (float) - The dispersion point data (either the spin-lock field strength in Hz or the nu_CPMG frequency in Hz).
Returns: bool
True if dispersion data exists, False otherwise.

has_exponential_exp_type()

source code 

Determine if the current data pipe contains exponential curves.

Returns: bool
True if spectral data for exponential curves exist, False otherwise.

has_fixed_time_exp_type()

source code 

Determine if the current data pipe contains fixed time data.

Returns: bool
True if spectral data for fixed time data exists, False otherwise.

has_proton_mmq_cpmg()

source code 

Determine if the current data pipe contains either proton SQ or MQ (MMQ) CPMG data.

This is only for the MMQ models.

Returns: bool
True if either proton SQ or MQ CPMG data exists, False otherwise.

has_proton_mq_cpmg()

source code 

Determine if the current data pipe contains proton MQ CPMG data.

This is only for the MMQ models.

Returns: bool
True if proton MQ CPMG data exists, False otherwise.

has_proton_sq_cpmg()

source code 

Determine if the current data pipe contains proton SQ CPMG data.

This is only for the MMQ models.

Returns: bool
True if proton SQ CPMG data exists, False otherwise.

has_r1rho_exp_type()

source code 

Determine if the current data pipe contains R1rho experiment types.

Returns: bool
True if R1rho experiment types exist, False otherwise.

insignificance(level=0.0)

source code 

Deselect all spins with insignificant dispersion profiles.

Parameters:
  • level (float) - The R2eff/R1rho value in rad/s by which to judge insignificance. If the maximum difference between two points on all dispersion curves for a spin is less than this value, that spin will be deselected.

interpolate_disp(spin=None, spin_id=None, si=None, num_points=None, extend_hz=None, relax_times=None)

source code 

Interpolate function for 2D Grace plotting function for the dispersion curves.

Parameters:
  • spin (SpinContainer instance.) - The specific spin data container.
  • spin_id (str) - The spin ID string.
  • si (int) - The index of the given spin in the cluster.
  • num_points (int) - The number of points to generate the interpolated fitted curves with.
  • extend_hz (float) - How far to extend the interpolated fitted curves to (in Hz).
  • relax_times (rank-4 list of floats) - The experiment specific fixed time period for relaxation (in seconds). The dimensions are {Ei, Mi, Oi, Di, Ti}.
Returns: boolean, rank-4 list of numpy rank-1 float arrays, rank-4 list of numpy rank-1 float arrays, rank-3 list of numpy rank-1 float arrays, rank-4 list of numpy rank-1 float arrays, rank-2 list of numpy rank-1 float arrays, rank-4 list of numpy rank-1 float arrays, rank-4 list of numpy rank-1 float arrays, rank-4 list of numpy rank-1 float arrays
The interpolated_flag, list of back calculated R2eff/R1rho values in rad/s {Ei, Si, Mi, Oi, Di}, list of interpolated frequencies for cpmg_frqs in Hz {Ei, Si, Mi, Oi, Di}, interpolated spin-lock offsets in rad/s {Ei, Si, Mi, Oi}, list of interpolated spin-lock field strength frequencies for spin_lock_nu1_new in Hz {Ei, Si, Mi, Oi, Di}, chemical shifts in rad/s {Ei, Si, Mi}, interpolated rotating frame tilt angles theta {Ei, Si, Mi, Oi, Di}, interpolated average resonance offset in the rotating frame Omega in rad/s {Ei, Si, Mi, Oi, Di} and the interpolated effective field in rotating frame w_eff in rad/s {Ei, Si, Mi, Oi, Di}.

interpolate_offset(spin=None, spin_id=None, si=None, num_points=None, extend_ppm=None, relax_times=None)

source code 

Interpolate function for 2D Grace plotting function for the dispersion curves, interpolating through spin-lock offset in rad/s.

Parameters:
  • spin (SpinContainer instance.) - The specific spin data container.
  • spin_id (str) - The spin ID string.
  • si (int) - The index of the given spin in the cluster.
  • num_points (int) - The number of points to generate the interpolated fitted curves with.
  • extend_ppm (float) - How far to extend the interpolated fitted curves to in offset ppm.
  • relax_times (rank-4 list of floats) - The experiment specific fixed time period for relaxation (in seconds). The dimensions are {Ei, Mi, Oi, Di, Ti}.
Returns: boolean, rank-4 list of numpy rank-1 float arrays, rank-4 list of numpy rank-1 float arrays, rank-3 list of numpy rank-1 float arrays, rank-4 list of numpy rank-1 float arrays, rank-2 list of numpy rank-1 float arrays, rank-4 list of numpy rank-1 float arrays, rank-4 list of numpy rank-1 float arrays, rank-4 list of numpy rank-1 float arrays
The interpolated_flag, list of back calculated R2eff/R1rho values in rad/s {Ei, Si, Mi, Oi, Di}, list of interpolated frequencies for cpmg_frqs in Hz {Ei, Si, Mi, Oi, Di}, interpolated spin-lock offsets in rad/s {Ei, Si, Mi, Oi}, list of interpolated spin-lock field strength frequencies for spin_lock_nu1_new in Hz {Ei, Si, Mi, Oi, Di}, chemical shifts in rad/s {Ei, Si, Mi}, interpolated rotating frame tilt angles theta {Ei, Si, Mi, Oi, Di}, interpolated average resonance offset in the rotating frame Omega in rad/s {Ei, Si, Mi, Oi, Di} and the interpolated effective field in rotating frame w_eff in rad/s {Ei, Si, Mi, Oi, Di}.

is_cpmg_exp_type(id=None)

source code 

Determine if the given spectrum ID corresponds to a CPMG experiment type.

Parameters:
  • id (str) - The spectrum ID string.
Returns: bool
True if the spectrum ID corresponds to a CPMG experiment type, False otherwise.

is_r1_optimised(model=None)

source code 

Should R1 values be optimised?

Parameters:
  • model (str) - The model to test for.
Returns: bool
True if the R1 values should be optimised, False if loaded values should be used instead.

is_r1rho_exp_type(id=None)

source code 

Determine if the given spectrum ID corresponds to a R1rho experiment type.

Parameters:
  • id (str) - The spectrum ID string.
Returns: bool
True if the spectrum ID corresponds to a R1rho experiment type, False otherwise.

loop_cluster(skip_desel=True)

source code 

Loop over the spin groupings for one model applied to multiple spins.

Parameters:
  • skip_desel (bool) - A flag which if True will cause deselected spins or spin clusters to be skipped.
Returns: list of str
The list of spin IDs per block will be yielded.

loop_exp(return_indices=False)

source code 

Generator method for looping over all experiment types.

Parameters:
  • return_indices (bool) - A flag which if True will cause the experiment type index to be returned as well.
Returns: str, (int)
The experiment type, and the index if asked.

loop_exp_frq(return_indices=False)

source code 

Generator method for looping over the exp and frq data.

These are the experiment types and spectrometer frequencies.

Parameters:
  • return_indices (bool) - A flag which if True will cause the experiment type and spectrometer frequency indices to be returned as well.
Returns: str, float, (int, int)
The experiment type and spectrometer frequency in Hz, and the indices if asked.

loop_exp_frq_offset(return_indices=False)

source code 

Generator method for looping over the exp, frq, and offset data.

These are the experiment types, spectrometer frequencies and spin-lock offset data.

Parameters:
  • return_indices (bool) - A flag which if True will cause the experiment type, spectrometer frequency and spin-lock offset indices to be returned as well.
Returns: str, float, float, (int, int, int)
The experiment type, spectrometer frequency in Hz and spin-lock offset data, and the indices if asked.

loop_exp_frq_offset_point(return_indices=False)

source code 

Generator method for looping over the exp, frq, offset, and point data.

These are the experiment types, spectrometer frequencies, spin-lock offset data, and dispersion points.

Parameters:
  • return_indices (bool) - A flag which if True will cause the experiment type, spectrometer frequency, spin-lock offset and dispersion point indices to be returned as well.
Returns: str, float, float, float, (int, int, int, int)
The experiment type, spectrometer frequency in Hz, spin-lock offset data and dispersion point data (either the spin-lock field strength in Hz or the nu_CPMG frequency in Hz), and the indices if asked.

loop_exp_frq_offset_point_time(return_indices=False)

source code 

Generator method for looping over the exp, frq, offset, and point data.

These are the experiment types, spectrometer frequencies, spin-lock offset data, and dispersion points.

Parameters:
  • return_indices (bool) - A flag which if True will cause the experiment type, spectrometer frequency, spin-lock offset and dispersion point indices to be returned as well.
Returns: str, float, float, float, (int, int, int, int)
The experiment type, spectrometer frequency in Hz, spin-lock offset data and dispersion point data (either the spin-lock field strength in Hz or the nu_CPMG frequency in Hz), and the indices if asked.

loop_exp_frq_point(return_indices=False)

source code 

Generator method for looping over the exp, frq, and point data.

These are the experiment types, spectrometer frequencies and dispersion points.

Parameters:
  • return_indices (bool) - A flag which if True will cause the experiment type, spectrometer frequency and dispersion point indices to be returned as well.
Returns: str, float, float, (int, int, int)
The experiment type, spectrometer frequency in Hz and dispersion point data (either the spin-lock field strength in Hz or the nu_CPMG frequency in Hz), and the indices if asked.

loop_exp_frq_point_time(return_indices=False)

source code 

Generator method for looping over the exp, frq, point, and time data.

These are the experiment types, spectrometer frequencies, dispersion points, and relaxation times.

Parameters:
  • return_indices (bool) - A flag which if True will cause the experiment type, spectrometer frequency, dispersion point, and relaxation time indices to be returned as well.
Returns: str, float, float, float, (int, int, int, int)
The experiment type, spectrometer frequency in Hz, dispersion point data (either the spin-lock field strength in Hz or the nu_CPMG frequency in Hz), the relaxation time, and the indices if asked.

loop_frq(return_indices=False)

source code 

Generator method for looping over all spectrometer frequencies.

Parameters:
  • return_indices (bool) - A flag which if True will cause the spectrometer frequency index to be returned as well.
Returns: float, (int)
The spectrometer frequency in Hz, and the index if asked.

loop_frq_offset(exp_type=None, return_indices=False)

source code 

Generator method for looping over the spectrometer frequencies and dispersion points.

Parameters:
  • exp_type (str) - The experiment type.
  • return_indices (bool) - A flag which if True will cause the spectrometer frequency and dispersion point indices to be returned as well.
Returns: float, float, (int, int)
The spectrometer frequency in Hz and dispersion point data (either the spin-lock field strength in Hz or the nu_CPMG frequency in Hz).

loop_frq_point(exp_type=None, return_indices=False)

source code 

Generator method for looping over the spectrometer frequencies and dispersion points.

Parameters:
  • exp_type (str) - The experiment type.
  • return_indices (bool) - A flag which if True will cause the spectrometer frequency and dispersion point indices to be returned as well.
Returns: float, float, (int, int)
The spectrometer frequency in Hz and dispersion point data (either the spin-lock field strength in Hz or the nu_CPMG frequency in Hz).

loop_frq_offset_point_key(exp_type=None)

source code 

Generator method for looping over the spectrometer frequencies, spin-lock offsets and dispersion points (returning the key).

Parameters:
  • exp_type (str) - The experiment type.
Returns: str
The key corresponding to the spectrometer frequency, offset and dispersion point.

loop_frq_point_time(exp_type=None, return_indices=False)

source code 

Generator method for looping over the spectrometer frequencies, dispersion points, and relaxation times.

Parameters:
  • exp_type (str) - The experiment type.
  • return_indices (bool) - A flag which if True will cause the spectrometer frequency, dispersion point, and relaxation time indices to be returned as well.
Returns: float, float, float
The spectrometer frequency in Hz, dispersion point data (either the spin-lock field strength in Hz or the nu_CPMG frequency in Hz), and the relaxation time.

loop_offset(exp_type=None, frq=None, return_indices=False)

source code 

Generator method for looping over the spin-lock offset values.

Parameters:
  • exp_type (str) - The experiment type.
  • frq (float) - The spectrometer frequency.
  • return_indices (bool) - A flag which if True will cause the offset index to be returned as well.
Returns: float, (int)
The spin-lock offset value and the index if asked.

loop_offset_point(exp_type=None, frq=None, skip_ref=True, return_indices=False)

source code 

Generator method for looping over the offsets and dispersion points.

Parameters:
  • exp_type (str) - The experiment type.
  • frq (float) - The spectrometer frequency.
  • skip_ref (bool) - A flag which if True will cause the reference point to be skipped.
  • return_indices (bool) - A flag which if True will cause the offset and dispersion point indices to be returned as well.
Returns: float, float, (int, int)
The offsets in ppm and the dispersion point data (either the spin-lock field strength in Hz or the nu_CPMG frequency in Hz), and the index if asked.

loop_point(exp_type=None, frq=None, offset=None, time=None, skip_ref=True, return_indices=False)

source code 

Generator method for looping over the dispersion points.

Parameters:
  • exp_type (str) - The experiment type.
  • frq (float) - The spectrometer frequency.
  • offset (None or float) - The spin-lock or hard pulse offset value in ppm.
  • time (float) - The relaxation time period.
  • skip_ref (bool) - A flag which if True will cause the reference point to be skipped.
  • return_indices (bool) - A flag which if True will cause the experiment type index to be returned as well.
Returns: float, (int)
Dispersion point data for the given indices (either the spin-lock field strength in Hz or the nu_CPMG frequency in Hz), and the index if asked.

loop_spectrum_ids(exp_type=None, frq=None, offset=None, point=None, time=None)

source code 

Generator method for selectively looping over the spectrum IDs.

Parameters:
  • exp_type (str) - The experiment type.
  • frq (float) - The spectrometer frequency.
  • offset (None or float) - For R1rho-type data, the spin-lock offset value in ppm.
  • point (float) - The dispersion point data (either the spin-lock field strength in Hz or the nu_CPMG frequency in Hz).
  • time (float) - The relaxation time period.
Returns: str
The spectrum ID.

loop_time(exp_type=None, frq=None, offset=None, point=None, return_indices=False)

source code 

Generator method for looping over the relaxation times.

Parameters:
  • exp_type (str) - The experiment type.
  • frq (float) - The spectrometer frequency in Hz.
  • offset (None or float) - The spin-lock or hard pulse offset value in ppm.
  • point (float) - The dispersion point data (either the spin-lock field strength in Hz or the nu_CPMG frequency in Hz).
  • return_indices (bool) - A flag which if True will cause the relaxation time index to be returned as well.
Returns: float
The relaxation time.

num_exp_types()

source code 

Count the number of experiment types present.

Returns: int
The number of experiment types.

pack_back_calc_r2eff(spin=None, spin_id=None, si=None, back_calc=None, proton_mmq_flag=False)

source code 

Store the back calculated R2eff data for the given spin.

Parameters:
  • spin (SpinContainer instance) - The spin data container to store the data in.
  • spin_id (str) - The spin ID string.
  • si (int) - The index of the given spin in the cluster.
  • back_calc (list of lists of lists of lists of float) - The back calculated data. The first index corresponds to the experiment type, the second is the spin of the cluster, the third is the magnetic field strength, and the fourth is the dispersion point.
  • proton_mmq_flag (bool) - The flag specifying if proton SQ or MQ CPMG data exists for the spin.

plot_exp_curves(file=None, dir=None, force=None, norm=None)

source code 

Custom 2D Grace plotting function for the exponential curves.

Parameters:
  • file (str) - The name of the Grace file to create.
  • dir (str) - The optional directory to place the file into.
  • force (bool) - Boolean argument which if True causes the file to be overwritten if it already exists.
  • norm (bool) - The normalisation flag which if set to True will cause all graphs to be normalised to a starting value of 1.

r20_from_min_r2eff(force=True, verbosity=1)

source code 

Set the R20 values to the minimum R2eff values.

For a 2 field cpmg experiment with model CR72, that would drop number of uniform grid search point from gridNr^5 to gridNr^3. For standard 21 grid Nr, it would make the grid search 441 times faster.

Parameters:
  • force (bool) - A flag forcing the overwriting of current values.
  • verbosity (int) - A flag specifying to print the setting of values.

r2eff_read(id=None, file=None, dir=None, disp_frq=None, offset=None, spin_id_col=None, mol_name_col=None, res_num_col=None, res_name_col=None, spin_num_col=None, spin_name_col=None, data_col=None, error_col=None, sep=None)

source code 

Read R2eff/R1rho values directly from a file whereby each row corresponds to a different spin. If the spin-container already contain r2eff values with the 'frequency of the CPMG pulse' or 'spin-lock field strength', the frequency will be changed by a infinitesimal small value of + 0.001 Hz. This allow for duplicates or more of the same frequency.

Parameters:
  • id (str) - The experiment ID string to associate the data with.
  • file (str) - The name of the file to open.
  • dir (str or None) - The directory containing the file (defaults to the current directory if None).
  • disp_frq (float) - For CPMG-type data, the frequency of the CPMG pulse train. For R1rho-type data, the spin-lock field strength (nu1). The units must be Hertz.
  • offset (None or float) - For R1rho-type data, the spin-lock offset value in ppm.
  • spin_id_col (int or None) - The column containing the spin ID strings. If supplied, the mol_name_col, res_name_col, res_num_col, spin_name_col, and spin_num_col arguments must be none.
  • mol_name_col (int or None) - The column containing the molecule name information. If supplied, spin_id_col must be None.
  • res_name_col (int or None) - The column containing the residue name information. If supplied, spin_id_col must be None.
  • res_num_col (int or None) - The column containing the residue number information. If supplied, spin_id_col must be None.
  • spin_name_col (int or None) - The column containing the spin name information. If supplied, spin_id_col must be None.
  • spin_num_col (int or None) - The column containing the spin number information. If supplied, spin_id_col must be None.
  • data_col (int or None) - The column containing the R2eff/R1rho data in Hz.
  • error_col (int or None) - The column containing the R2eff/R1rho errors.
  • sep (str or None) - The column separator which, if None, defaults to whitespace.

r2eff_read_spin(id=None, spin_id=None, file=None, dir=None, disp_point_col=None, offset_col=None, data_col=None, error_col=None, sep=None)

source code 

Read R2eff/R1rho values from file whereby each row is a different dispersion point.

Parameters:
  • id (str) - The experiment ID string to associate the data with. This will be modified to include the dispersion point data as "%s_%s" % (id, disp_point).
  • spin_id (str) - The spin ID string.
  • file (str) - The name of the file to open.
  • dir (str or None) - The directory containing the file (defaults to the current directory if None).
  • disp_point_col (int) - The column containing the dispersion point information. For CPMG-type data, this is the frequency of the CPMG pulse train. For R1rho-type data, this is the spin-lock field strength (nu1). The units must be Hertz.
  • offset_col (None or int) - This is for R1rho data - the dispersion point column can be substituted for the offset values in Hertz.
  • data_col (int) - The column containing the R2eff/R1rho data in Hz.
  • error_col (int) - The column containing the R2eff/R1rho errors.
  • sep (str or None) - The column separator which, if None, defaults to whitespace.

randomise_R1(spin=None, ri_id=None, N=None)

source code 

Randomise the R1 data for the given spin for use in the Monte Carlo simulations.

Parameters:
  • spin (SpinContainer instance) - The spin container to randomise the data for.
  • ri_id (str) - The relaxation data ID string.
  • N (int) - The number of randomisations to perform.

relax_time(time=0.0, spectrum_id=None)

source code 

Set the relaxation time period associated with a given spectrum.

Parameters:
  • time (float) - The time, in seconds, of the relaxation period.
  • spectrum_id (str) - The spectrum identification string.

return_cpmg_frqs(ref_flag=True)

source code 

Return the list of nu_CPMG frequencies.

Parameters:
  • ref_flag (bool) - A flag which if False will cause the reference spectrum frequency of None to be removed from the list.
Returns: rank-2 list of numpy rank-1 float64 arrays
The list of nu_CPMG frequencies in Hz. It has the dimensions {Ei, Mi, Oi}.

return_cpmg_frqs_single(exp_type=None, frq=None, offset=None, time=None, ref_flag=True)

source code 

Return the list of nu_CPMG frequencies.

Parameters:
  • exp_type (str) - The experiment type.
  • frq (float) - The spectrometer frequency in Hz.
  • offset (None or float) - The hard pulse offset, if desired.
  • time (float) - The relaxation time period.
  • ref_flag (bool) - A flag which if False will cause the reference spectrum frequency of None to be removed from the list.
Returns: numpy rank-1 float64 array
The list of nu_CPMG frequencies in Hz.

return_x_y_data_labels_settings(spin=None, data_type=None, y_axis=None, exp_type=None, frq=None, offset=None, point=None, interpolated_flag=None)

source code 

Return the X and Y labels and plot settings, according to selected axis to plot for.

Parameters:
  • spin (SpinContainer instance) - The specific spin data container.
  • data_type (str) - String flag to tell which data type to return for. Option can be either "data", "back_calculated", "interpolated" or "residual".
  • exp_type (str) - The experiment type.
  • y_axis (str) - String flag to tell which data on Y axis to plot for. Option can be either "%s" which plot 'r2eff' for CPMG experiments or 'r1rho' for R1rho experiments or option can be "%s", which for R1rho experiments plot R2 = R1rho / sin^2(theta) - R_1 / tan^2(theta) = (R1rho - R_1 * cos^2(theta) ) / sin^2(theta).
  • frq (float) - The spectrometer frequency in Hz.
  • offset (None or float) - The spin-lock offset.
  • point (float) - The Spin-lock field strength (Hz).
  • interpolated_flag (bool) - Flag telling if the graph should be interpolated.
Returns: str, int, float, int, int
The data label, the data symbol, the data symbol size, the data line type, the data line style.

return_index_from_disp_point(value, exp_type=None)

source code 

Convert the dispersion point data into the corresponding index.

Parameters:
  • value (float) - The dispersion point data (either the spin-lock field strength in Hz or the nu_CPMG frequency in Hz).
  • exp_type (str) - The experiment type.
Returns: int
The corresponding index.

return_index_from_exp_type(exp_type=None)

source code 

Convert the experiment type into the corresponding index.

Parameters:
  • exp_type (str) - The experiment type.
Returns: int
The corresponding index.

return_index_from_frq(value)

source code 

Convert the dispersion point data into the corresponding index.

Parameters:
  • value (float) - The spectrometer frequency in Hz.
Returns: int
The corresponding index.

return_index_from_disp_point_key(key, exp_type=None)

source code 

Convert the dispersion point key into the corresponding index.

Parameters:
  • exp_type (str) - The experiment type.
  • key (str) - The dispersion point or R2eff/R1rho key.
Returns: int
The corresponding index.

return_key_from_di(mi=None, di=None)

source code 

Convert the dispersion point index into the corresponding key.

Parameters:
  • mi (int) - The spectrometer frequency index.
  • di (int) - The dispersion point or R2eff/R1rho index.
Returns: str
The corresponding key.

return_offset_data(spins=None, spin_ids=None, field_count=None, spin_lock_offset=None, fields=None)

source code 

Return numpy arrays of the chemical shifts, offsets and tilt angles.

Indices

The data structures consist of many different index types. These are:

  • Ei: The index for each experiment type.
  • Si: The index for each spin of the spin cluster.
  • Mi: The index for each magnetic field strength.
  • Oi: The index for each spin-lock offset. In the case of CPMG-type data, this index is always zero.
  • Di: The index for each dispersion point (either the spin-lock field strength or the nu_CPMG frequency).
Parameters:
  • spins (list of SpinContainer instances) - The list of spin containers in the cluster.
  • spin_ids (list of str) - The list of spin IDs for the cluster.
  • field_count (int) - The number of spectrometer field strengths. This may not be equal to the length of the fields list as the user may not have set the field strength.
  • spin_lock_offset (list of lists of numpy rank-1 float arrays) - The spin-lock offsets to use instead of the user loaded values - to enable interpolation.
  • fields (rank-2 list of floats) - The spin-lock field strengths to use instead of the user loaded values - to enable interpolation. The dimensions are {Ei, Mi}.
Returns: rank-3 list of numpy rank-1 float arrays, rank-4 list of numpy rank-1 float arrays, rank-2 list of numpy rank-1 float arrays, rank-4 list of numpy rank-1 float arrays, rank-4 list of numpy rank-1 float arrays, rank-4 list of numpy rank-1 float arrays
interpolated spin-lock offsets in rad/s {Ei, Si, Mi, Oi}, interpolated spin-lock field strength frequencies in Hz {Ei, Si, Mi, Oi, Di}, the chemical shifts in rad/s {Ei, Si, Mi}, the interpolated rotating frame tilt angles {Ei, Si, Mi, Oi, Di}, interpolated average resonance offset in the rotating frame in rad/s {Ei, Si, Mi, Oi, Di} and the interpolated effective field in rotating frame in rad/s {Ei, Si, Mi, Oi, Di}.

return_param_key_from_data(exp_type=None, frq=0.0, offset=0.0, point=0.0)

source code 

Generate the unique key from the spectrometer frequency and dispersion point.

Parameters:
  • exp_type (str) - The experiment type.
  • frq (float) - The spectrometer frequency in Hz.
  • offset (None or float) - The optional offset value for off-resonance R1rho-type data.
  • point (float) - The dispersion point data (either the spin-lock field strength in Hz or the nu_CPMG frequency in Hz).
Returns: str
The unique key.

return_r1_data(spins=None, spin_ids=None, field_count=None, sim_index=None)

source code 

Return the R1 data structures for off-resonance R1rho experiments.

Parameters:
  • spins (list of SpinContainer instances) - The list of spin containers in the cluster.
  • spin_ids (list of str) - The list of spin IDs for the cluster.
  • field_count (int) - The number of spectrometer field strengths. This may not be equal to the length of the fields list as the user may not have set the field strength.
  • sim_index (None or int) - The index of the simulation to return the R1 data of. This should be None if the normal data is required.
Returns: numpy rank-2 float array
The R1 relaxation data.

return_r1_err_data(spins=None, spin_ids=None, field_count=None, sim_index=None)

source code 

Return the R1 error data structures for off-resonance R1rho experiments.

Parameters:
  • spins (list of SpinContainer instances) - The list of spin containers in the cluster.
  • spin_ids (list of str) - The list of spin IDs for the cluster.
  • field_count (int) - The number of spectrometer field strengths. This may not be equal to the length of the fields list as the user may not have set the field strength.
  • sim_index (None or int) - The index of the simulation to return the R1 data of. This should be None if the normal data is required.
Returns: numpy rank-2 float array
The R1 relaxation error data.

return_r2eff_arrays(spins=None, spin_ids=None, fields=None, field_count=None, sim_index=None)

source code 

Return numpy arrays of the R2eff/R1rho values and errors.

Parameters:
  • spins (list of SpinContainer instances) - The list of spin containers in the cluster.
  • spin_ids (list of str) - The list of spin IDs for the cluster. In the case of multi-quantum systems, these will be different to the spins argument and instead refer to the second spin of the pair.
  • fields (list of float) - The list of spectrometer field strengths.
  • field_count (int) - The number of spectrometer field strengths. This may not be equal to the length of the fields list as the user may not have set the field strength.
  • sim_index (None or int) - The index of the simulation to return the data of. This should be None if the normal data is required.
Returns: lists of numpy float arrays, lists of numpy float arrays, lists of numpy float arrays, numpy rank-2 int array
The numpy array structures of the R2eff/R1rho values, errors, missing data, and corresponding Larmor frequencies. For each structure, the first dimension corresponds to the experiment types, the second the spins of a spin block, the third to the spectrometer field strength, and the fourth is the dispersion points. For the Larmor frequency structure, the fourth dimension is omitted. For R1rho-type data, an offset dimension is inserted between the spectrometer field strength and the dispersion points.

return_relax_times()

source code 

Return the list of relaxation times.

Returns: numpy rank-2 float64 array
The list of relaxation times in s.

return_spin_lock_nu1(ref_flag=True)

source code 

Return the list of spin-lock field strengths.

Parameters:
  • ref_flag (bool) - A flag which if False will cause the reference spectrum frequency of None to be removed from the list.
Returns: rank-2 list of numpy rank-1 float64 arrays
The list of spin-lock field strengths in Hz. It has the dimensions {Ei, Mi, Oi}.

return_spin_lock_nu1_single(exp_type=None, frq=None, offset=None, ref_flag=True)

source code 

Return the list of spin-lock field strengths.

Parameters:
  • exp_type (str) - The experiment type.
  • frq (float) - The spectrometer frequency in Hz.
  • offset (None or float) - The spin-lock offset.
  • ref_flag (bool) - A flag which if False will cause the reference spectrum frequency of None to be removed from the list.
Returns: numpy rank-1 float64 array
The list of spin-lock field strengths in Hz.

return_value_from_frq_index(mi=None)

source code 

Return the spectrometer frequency corresponding to the frequency index.

Parameters:
  • mi (int) - The spectrometer frequency index.
Returns: float
The spectrometer frequency in Hertz or None if no information is present.

return_value_from_offset_index(ei=None, mi=None, oi=None)

source code 

Return the offset corresponding to the offset index.

Parameters:
  • ei (int) - The experiment type index.
  • mi (int) - The spectrometer frequency index.
  • oi (int) - The offset index.
Returns: float
The offset in Hertz or None if no information is present.

set_exp_type(spectrum_id=None, exp_type=None)

source code 

Select the relaxation dispersion experiment type performed.

Parameters:
  • spectrum_id (str) - The spectrum ID string.
  • exp (str) - The relaxation dispersion experiment type.

spin_has_frq_data(spin=None, frq=None)

source code 

Determine if the spin has intensity data for the given spectrometer frequency.

Parameters:
  • spin (SpinContainer instance) - The specific spin data container.
  • frq (float) - The spectrometer frequency.
Returns: bool
True if data for that spectrometer frequency is present, False otherwise.

spin_ids_to_containers(spin_ids)

source code 

Take the list of spin IDs and return the corresponding spin containers.

This is useful for handling the data from the model_loop() method.

Parameters:
  • spin_ids (list of str) - The list of spin ID strings.
Returns: list of SpinContainer instances
The list of spin containers.

spin_lock_field(spectrum_id=None, field=None)

source code 

Set the spin-lock field strength (nu1) for the given spectrum.

Parameters:
  • spectrum_id (str) - The spectrum ID string.
  • field (int or float) - The spin-lock field strength (nu1) in Hz.

spin_lock_offset(spectrum_id=None, offset=None)

source code 

Set the spin-lock offset (omega_rf) for the given spectrum.

Parameters:
  • spectrum_id (str) - The spectrum ID string.
  • offset (int or float) - The spin-lock offset (omega_rf) in ppm.

write_disp_curves(dir=None, force=None)

source code 

Write out the dispersion curves to text files.

One file will be created per spin system.

Parameters:
  • dir (str) - The optional directory to place the file into.
  • force (bool) - If True, the files will be overwritten if they already exists.

Variables Details [hide private]

COLOUR_ORDER

Value:
[4,
 15,
 2,
 13,
 11,
 1,
 3,
 5,
...