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

class schrodinger.application.jaguar.spin_simulator_v2.Parameters(offset: float, sweep: float, spectrum_type: str, axis_units: str = 'ppm', invert_axis: bool = True, max_subgraph_size: int = 14, timestep_npoints: int = 8192, spectrum_npoints: int = 16536)

Bases: object

offset: float
sweep: float
spectrum_type: str
axis_units: str = 'ppm'
invert_axis: bool = True
max_subgraph_size: int = 14
timestep_npoints: int = 8192
spectrum_npoints: int = 16536
__init__(offset: float, sweep: float, spectrum_type: str, axis_units: str = 'ppm', invert_axis: bool = True, max_subgraph_size: int = 14, timestep_npoints: int = 8192, spectrum_npoints: int = 16536) None
schrodinger.application.jaguar.spin_simulator_v2.get_spectral_params(spectrum_spin: str, isotopes: List[str], shifts: numpy.ndarray, 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.

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, timestep_npoints: int = 8192, spectrum_npoints: int = 16536, sym_groups: Optional[List[str]] = None, sym_spins: Optional[List[List[int]]] = None, max_subgraph_size: Optional[int] = 14, debug: bool = False, debug_filename: str = 'debug')

Bases: object

Core Driver class for the simulation.

SPEC_BUFFER = 0.9
LIOUV_ZERO = 1e-07
PROP_ZERO = 1e-08
PROP_NORM = 1e-08
PROP_CHOP = 1e-10
KRYLOV_SWITCHOVER = 50000
ZTE_TOL = 1e-30
ZTE_NSTEPS = 10
ZTE_MAXDEN = 0.5
INTER_CUTOFF = 1e-05
PATH_TRACE = 1e-07
PATH_DROP = 1e-30
IRREP_DROP = 1e-30
SMALL_MATRIX = 200
DENSE_MATRIX = 0.25
__init__(chemical_shifts: List[float], couplings: List[Tuple[int, int, float]], magnet: float, isotopes: List[str], spectrum_type: str, timestep_npoints: int = 8192, spectrum_npoints: int = 16536, sym_groups: Optional[List[str]] = None, sym_spins: Optional[List[List[int]]] = None, max_subgraph_size: Optional[int] = 14, debug: bool = False, debug_filename: str = 'debug')

Parameter and option preprocessing.

Symmetry information is currently unused as the build-and-prune strategy currently implmented requires too much memory to be practical. Instead a constructive scheme should be developed, if it becomes apparent that symmetry optimizations are required.

simulate_spectrum(plot=False) Tuple[numpy.ndarray, numpy.ndarray]

Driver function to compute the NMR spectrum of this spin system.

Parameters

plot – If True then a plot of the spectrum will be produced. Primarily for debugging purposes.

Returns

Spectrum X (delta in ppm) and Y (signal) arrays.

pulse_acquire() numpy.ndarray

The standard 90-acqure pulse sequence.

Returns

fid, numpy array containing the time propagated signal.

operator(operators: Union[List[str], str], spins: Union[List[int], str]) scipy.sparse._arrays.csr_array

Generates commutation superoperators from their user-friendly descriptions.

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.

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.

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.

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.

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

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.

Parameters

mult – multiplicity

Returns

Left and right product tensors

irr_sph_ten(mult: int, k: Optional[int] = None) List[scipy.sparse._arrays.coo_array]

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_superop() 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.

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.

WE ARE NOT USING THIS IN OUR SIMPLE NMR. Kept as placeholder for future development and parity with spinach.

Returns

Redfield superoperator array, presently all zeros.

k_superop() scipy.sparse._arrays.csr_array

Chemical kinetics superoperator.

WE ARE NOT USING THIS IN OUR SIMPLE NMR Kept as placeholder for future development and parity with spinach.

Returns

Chemical kinetics superoperator array, presently all zeros.

state(states: List[str], spins: Union[List[int], str]) scipy.sparse._arrays.lil_array

Generates state vectors from their user-friendly descriptions.

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.

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).

Returns

State vectors.

statevec(opspec: List[int]) scipy.sparse._arrays.lil_array

User-specified state vector. Converts a Spinach operator specification into a Liouville space state vector.

If the requested state is not in the basis, a ValueError is raised.

Parameters

opspec – Operator specification constructed by human2opspec.

Returns

The resulting state vector

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: scipy.sparse._base.spmatrix, rho: scipy.sparse._base.spmatrix, time_step: float)

A propagation step function using Krylov propagation for large matrices and scipy’s expm for small matrices.

Parameters
  • 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.

Returns

The stepped state vector.

evolution(L: scipy.sparse._arrays.csr_array, coil: scipy.sparse._arrays.lil_array, rho: scipy.sparse._base.spmatrix, 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: scipy.sparse._arrays.csr_array, timestep: float) scipy.sparse._arrays.csr_array

Calculates exponential propagators and propagator derivatives.

Parameters
  • L – the Liouvillian matrix to propagate.

  • timestep – time step value (dt)

Return P

The propagator matrix.

reduce(L: scipy.sparse._arrays.csr_array, rho: scipy.sparse._arrays.lil_array) List[scipy.sparse._arrays.csr_array]

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

Parameters
  • L – Liouvillian matrix.

  • rho – Initial state vector.

Returns

List of projector matrices.

zte(L: scipy.sparse._arrays.csr_array, rho: scipy.sparse._arrays.lil_array) scipy.sparse._arrays.csr_array

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.

Parameters
  • L – The Liouvillian to be used for time propagation

  • rho – The initial state to be used for time propagation

Returns

projector matrix into the reduced space, to be used as follows: L = P.T @ L @ P, rho = P.T @ rho

path_trace(L: scipy.sparse._base.spmatrix) List[scipy.sparse._arrays.csr_array]

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.

Parameters

L – User-supplied Liouvillian

Returns

List of projectors into the independently evolving subspaces.

normest(A: scipy.sparse._base.spmatrix) float

Estimate the norm of sparse matrix A.

Currently wraps scipy.sparse.linalg.norm as a placeholder for an optimized norm estimation.

Parameters

A – Sparse matrix to estimate the norm of.

Returns

Estimated norm of matrix A.

clean_up(A: scipy.sparse._base.spmatrix, nonzero_tol: float) scipy.sparse._base.spmatrix

Sparse matrix clean-up utility.

Parameters
  • 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, molname: str = 'TEST', runtime: Optional[float] = None, data_only=False) Tuple[numpy.ndarray, numpy.ndarray]

Simple plotting function for debugging purposes. Builds the axes from spin information.

Parameters
  • spectrum – spectrum data returned by running the simulation.

  • molname – Label to annotate the spectrum with.

  • runtime – The time taken to run the simulation, if present this will be annotated on the plot.

  • data_only – flag that only the x, y data should be returned, no plot should be made

Returns

The x and y data that make up the plot.

plot_subgraphs(subgraphs: Optional[numpy.ndarray] = None, fname: str = '')

Debugging function plotting the subgraph topology of the spin system.

Parameters
  • subgraphs – Numpy array of boolean vectors flagging which spins are involved.

  • fname – filename that the figure should be written to. If empty the figure will show interactively instead.