schrodinger.application.jaguar.spin_simulator_v2 module¶
This simulation follows the simulation algorithm used in the Spinach package, which can be found here: https://doi.org/10.1016/j.jmr.2010.11.008
- schrodinger.application.jaguar.spin_simulator_v2.print_proclamation(msg: str)¶
Utility function to print title-formatted message
- schrodinger.application.jaguar.spin_simulator_v2.where_sorted_array(a, row)¶
- schrodinger.application.jaguar.spin_simulator_v2.is_larger(a: numpy.ndarray, b: numpy.ndarray) int ¶
- schrodinger.application.jaguar.spin_simulator_v2.PauliMatrices¶
alias of
schrodinger.application.jaguar.spin_simulator_v2.PauliMats
- schrodinger.application.jaguar.spin_simulator_v2.build_pauli(mult: int) schrodinger.application.jaguar.spin_simulator_v2.PauliMats ¶
Build Pauli matrices for the given multiplicity.
We use these so often it’s worth it to cache these for the two values we might concievably allow (mult 2 and 3).
- Parameters
mult – multiplicity to create the Pauli matrices of
- Returns
namedtuple with the x, y, z, plus, and minus Pauli matrices
- schrodinger.application.jaguar.spin_simulator_v2.get_spectral_params(spectrum_spin: str, isotopes: List[str], shifts: List[float], magnet: float) Tuple[float, float] ¶
Determine sweep and offset (in Hz) for the particular system under consideration.
- Parameters
spectrum_spin – isotope that 1D NMR is probing, e.g. 1H
isotopes – list of isotopes for each spin in the system
shifts – list of chemical shifts in ppm for each spin in the system
magnet – strength of NMR magnet in Tesla
- Returns
width and offset of spectral window
- schrodinger.application.jaguar.spin_simulator_v2.apodization(fid: numpy.ndarray) numpy.ndarray ¶
Post-processing cleanup of FID to remove unwanted ringing near peaks. Conventionally applied to NMR signals before FT. Typically NMR signals are truncated for time and signal/noise constraints. This trunctation introduces artifacts which are removed by a window function.
In this case an exponential window is used.
- Parameters
fid – input FID signal.
- Returns
FID signal apodized with exponential window.
- schrodinger.application.jaguar.spin_simulator_v2.dfpt(connectivity_matrix: numpy.ndarray, max_subgraph_size: int) numpy.ndarray ¶
Graph partitioning function.
- Parameters
connectivity_matrix – Input connectivity matrix.
max_subgraph_size – Maximum size of the subgraph.
- Returns
Matrix of subgraphs.
- class schrodinger.application.jaguar.spin_simulator_v2.SpinSystem(chemical_shifts: List[float], couplings: List[Tuple[int, int, float]], magnet: float, isotopes: List[str], spectrum_type: str, sym_groups: Optional[List[str]] = None, sym_spins: Optional[List[List[int]]] = None, min_ppm: float = - 0.5, max_ppm: float = 9.0, debug=False)¶
Bases:
object
- SPEC_BUFFER = 0.9¶
- GLOBAL_IT = 0¶
- __init__(chemical_shifts: List[float], couplings: List[Tuple[int, int, float]], magnet: float, isotopes: List[str], spectrum_type: str, sym_groups: Optional[List[str]] = None, sym_spins: Optional[List[List[int]]] = None, min_ppm: float = - 0.5, max_ppm: float = 9.0, debug=False)¶
Parameter and option preprocessing.
- simulate_spectrum(plot=False) Tuple[numpy.ndarray, numpy.ndarray] ¶
- tolerances()¶
- basis()¶
Generate basis set for calculations,
- pulse_acquire(parameters: Dict) numpy.ndarray ¶
The standard 90-acqure pulse sequence.
The following parameters are currently accepted:
parameters.spins - ‘1H’, ‘13C’, etc parameters.sweep - the width of the spectral window (Hz) parameters.npoints - number time steps in the simulation parameters.zerofill - number of points to zerofill to parameters.offset - transmitter offset (Hz) parameters.axis_units :return: fid, numpy array containing
- operator(operators: List[str], spins: Union[List[int], str]) scipy.sparse._arrays.csr_array ¶
Generates commutation superoperators from their user-friendly descriptions.
- Parameters
operators – a list of strings. Each string may be ‘E’ (identity operator), ‘Lz’, ‘L+’, ‘L-’ or ‘Tl,m’ (higher spin orders as irreducible spherical tensors; l and m are both integers).
spins –
EITHER a list of integers specifying the numbers of spins on which the operators given in the ‘operators’ argument operate (this produces the corresponding multi-spin superoperator)
OR an isotope specification: ‘13C’, ‘15N’, ‘all’ (this produces a sum of single-spin operators on all the spins of the specified isotope).
- Returns
The operator sparse matrix.
Example: LzSp = spin_system.operator([‘Lz’, ‘L+’], [1, 2]) This would return the [LzS+, ] commutation superoperator with Lz on spin 1 and L+ on spin 2.
Example: Sum_Lz = spin_system.operator(‘Lz’, ‘13C’) This would return the sum of [Lz, ] commutation superoperators on all carbons in the system.
- human2opspec(operators: List[str], spins: List[int]) Tuple[List[int], float] ¶
Converts user-friendly descriptions of spin states and operators into the formal description (opspec) used by Spinach kernel.
- Parameters
operators – a list of strings. Each string can be ‘E’ (identity operator), ‘Lz’, ‘L+’, ‘L-’ (single spin orders) or ‘Tl,m’ (higher spin orders as irreducible spherical tensors; l and m are both integers).
spins – a list of integers specifying the numbers of spins on which the operators given in the ‘operators’ argument operate.
- Returns
operator specification and coefficient
Example: opspec, coefficient = spin_system.human2opspec([‘Lz’, ‘L+’], [1, 2]) This would return the required opspec and coefficient to generate an LzL+ operator or state vector with Lz on the spins 1 and L+ on spin 2.
- c_superop(opspec: List[int]) scipy.sparse._arrays.csr_array ¶
Computes commutation superoperators.
- Parameters
opspec – Operator specification
- Returns
Sparse matrix representation of operator.
- p_superop(opspec: List[int], side: str) scipy.sparse._arrays.csr_array ¶
Returns superoperators corresponding to right or left multiplication of a density matrix by a user-specified operator.
- Parameters
opspec – operator specification.
side – ‘left’ or ‘right’
- Returns
sparse matrix representation of super-operator
- ist_product_table(mult: int) Tuple[List[scipy.sparse._arrays.csr_array], List[scipy.sparse._arrays.csr_array]] ¶
Multiplication tables for irreducible spherical tensors.
NOTE: Caching on disk here is probably unneccsary. So far I have only seen multiplicity = 2, and that is shape (4, 4, 4). So we could just keep in memory after first calculation. multiplicity = 6 would be required for oxygen, but currently Pauli matrix generation is only up to 3, so no worry there yet.
- Parameters
mult – multiplicity
- Returns
Left and right product tensors
- irr_sph_ten(mult: int, k: Optional[int] = None) List[int] ¶
Returns an array of single-spin irreducible spherical tensor operators T(k,m).
- Parameters
mult – multiplicity
k – tensor rank
- Returns
List if spherical tensor operators
- offset(spins: List[str], offset_value: float) numpy.ndarray ¶
Generates the offset superoperator.
When added to the main Liouvillian, the offset superoperator shifts the master frame of a given set of spins by a given amount (in Hz).
- Parameters
spins – May be set to ‘all’, ‘1H’, ‘13C’ etc.
offset_value – Offset value in Hz.
- Returns
Offset superoperator.
- secularity()¶
Sets the interaction secularity assumptions for the coherent Liouvillian. Simulations of different experiments require different rotating frames and different secular approximations - those are stored in this function.
Modifies the SpinSystem in-place.
- h_superop() scipy.sparse._arrays.csr_array ¶
Hamiltonian superoperator and its rotational decomposition for Liouville space simulations.
- Returns
H superoperator.
- r_superop() scipy.sparse._arrays.csr_array ¶
The relaxation superoperator. Computes the Redfield superoperator by numerical evaluation of the integral in the Bloch-Redfield-Wangsness master equation in Liouville space. All options are set during the call to create function.
WE ARE NOT USING THIS IN OUR SIMPLE NMR
- k_superop() scipy.sparse._arrays.csr_array ¶
Chemical kinetics superoperator. All adjustable parameters are specified in the call to create function.
WE ARE NOT USING THIS IN OUR SIMPLE NMR
- state(states: List[str], spins: Union[List[int], str]) scipy.sparse._arrays.lil_array ¶
Generates state vectors from their user-friendly descriptions.
- Parameters
states – Each string may be ‘E’ (identity state), ‘Lz’, ‘L+’, ‘L-’ or ‘Tl,m’ (higher spin states as irreducible spherical tensors; l and m are both integers).
spins – EITHER a list of integers specifying the numbers of spins to be generated in states given in the ‘states’ argument (this produces the corresponding multi-spin state) OR an isotope specification: ‘13C’, ‘15N’, ‘all’ (this produces a sum of single-spin states on all the spins of the specified isotope)
Example:
LzSp=state(spin_system,{‘Lz’,’L+’},{1,2});
would return the LzS+ state with Lz on spin 1 and L+ on spin 2.
Example:
Sum_Lz=state(spin_system,’Lz’,’13C’);
would return the sum of Lz states on all carbons in the system.
- statevec(opspec)¶
User-specified state vector. Converts a Spinach operator specification into a Liouville space state vector.
Args: - opspec : Spinach operator specification (refer to the Spinach manual)
Returns: - rho : the resulting state vector
If the requested state is not in the basis, a ValueError is raised.
- equilibrium()¶
Thermal equilibrium state.
If the temperature is set to zero during initialization, returns the high-field, high-temperature approximation to the thermal equilibrium state. If the temperature is specified, returns the thermal equilibrium state of the Zeeman Hamiltonian, ignoring all couplings and offsets.
If the requested state is not in the basis, a ValueError is raised.
- step(L, rho, time_step)¶
A propagation step function using Krylov propagation for large matrices and scipy’s expm for small matrices.
- Arguments:
L : the Liouvillian to be used for propagation rho : the state vector (or a horizontal stack thereof) to be propagated time_step : the length of the time step to take
- evolution(L: scipy.sparse._arrays.csr_array, coil, rho, timestep: float, nsteps: int) numpy.ndarray ¶
Time evolution function. Performs all types of time propagation with automatic trajectory level state space restriction.
- Parameters
L – the Liouvillian to be used during evolution
rho – the initial state vector or a horizontal stack thereof
coil – the detection state, used when ‘observable’ is specified as the output option
timestep – the length of each time step
nsteps – number of steps in the evolution
- Returns
the resulting state, observable, or trajectory based on the specified output
- propagator(L, timestep)¶
Calculates exponential propagators and propagator derivatives.
- Args:
L : the Liouvillian timestep : time step value
- Returns:
P : the propagator matrix
- reduce(L, rho)¶
Symmetry and trajectory-level state space restriction for a user-specified Liouvillian L and initial state rho. Applies all reduction methods and returns a list of projectors into a set of independent reduced subspaces. The projectors are to be used as follows:
L = P.T @ L @ P rho = P.T @ rho
- zte(L, rho)¶
Zero track elimination function. Inspects the first few steps in the system trajectory and drops the states that did not get populated to a user-specified tolerance.
L : the Liouvillian to be used for time propagation rho : the initial state to be used for time propagation
- Returns:
- projectorprojector matrix into the reduced space,
to be used as follows: L = P.T @ L @ P, rho = P.T @ rho
- stepsize(L)¶
Estimate the Larmor time step based on the Liouvillian L.
L : Liouvillian for which the time step is estimated
- Returns:
timestep : Estimated Larmor time step
- path_trace(L)¶
Liouvillian path tracing. Treats the given Liouvillian as the adjacency matrix of a graph, computes the weakly connected subgraphs of that graph, and returns a list of projectors into the corresponding independently evolving subspaces.
L: User-supplied Liouvillian
- Returns:
projectors: List of projectors into the independently evolving subspaces
- negligible(value, mode, threshold)¶
- significant(something, type_, tolerance)¶
Determines whether a given object deserves attention given the tolerance specified. Used in the internal decision-making conducted by spin_system kernel functions.
- Parameters
something – object to check
type – type of the object (e.g., ‘tensor’, ‘vector’, etc.)
tolerance – the tolerance level for significance
- Returns
True if the object is significant, False otherwise
- normest(A, tolerance)¶
- clean_up(A: Union[scipy.sparse._arrays.csr_array, scipy.sparse._csr.csr_matrix, numpy.ndarray], nonzero_tol: float) Union[scipy.sparse._arrays.csr_array, scipy.sparse._csr.csr_matrix, numpy.ndarray] ¶
Sparse matrix clean-up utility.
Args: A: The matrix to clean up. nonzero_tol: Tolerance below which values will be set to zero.
Returns: A cleaned-up version of the input matrix.
- plot_1d(spectrum, parameters, molname: str = 'TEST', runtime: Optional[float] = None, data_only=False) Tuple[numpy.ndarray, numpy.ndarray] ¶
Build the axis and apply transformations.
Args: spin_system (dict): Dictionary containing the spin system parameters. parameters (dict): Dictionary containing various parameters needed for the axis computation.
Returns: numpy.ndarray: Computed axis.
- schrodinger.application.jaguar.spin_simulator_v2.spin(name: str)¶
Get gamma and multiplicity values for given isotope name.
- Parameters:
name (str): Name of the isotope.
- Returns:
tuple: gamma and multiplicity values.
- schrodinger.application.jaguar.spin_simulator_v2.symmetry_group(group_name)¶
Symmetry group database
- schrodinger.application.jaguar.spin_simulator_v2.remove_subset_clusters(mat: scipy.sparse._csr.csr_matrix) scipy.sparse._csr.csr_matrix ¶
Remove clusters that are subsets of other clusters.
If a cluster is purely a subset of another cluster, it does not add any value to the simulation.
- Parameters
mat – matrix with rows that may be subsets
- Returns
matrix with subset rows removed
- schrodinger.application.jaguar.spin_simulator_v2.remove_vals_below_thresh(mat: scipy.sparse._arrays.csr_array, thresh: float)¶
Floor small values to 0 in sparse matrices so that they are not stored explicitly.
- Parameters
mat – matrix possibly with some small values
thresh – threshold for small values to be floored to zeros
- schrodinger.application.jaguar.spin_simulator_v2.load_class_array(fname: str) List[numpy.ndarray] ¶