nmrsim package

Submodules

nmrsim.discrete module

Non-quantum mechanical solutions for specific second-order patterns.

These are adapted from the routines from WINDNMR [1] by Hans Reich, U. Wisconsin, and equations from Pople, Schneider and Bernstein [2]. Note that many of the names for arguments, etc. are non-Pythonic but chosen to match the WINDNMR interface and source code (for now).

The partials module provides the following functions:

  • AB: simulates an AB quartet.
  • AB2: simulates an AB2 system.
  • ABX: simulates an ABX system.
  • ABX3: simulates an ABX3 system.
  • AAXX: simulates an AA’XX’ system.
  • AABB: simulates an AA’BB’ system.

References

[1]WINDNMR-Pro home page: https://www.chem.wisc.edu/areas/reich/plt/windnmr.htm
[2]Pople, J.A.; Schneider, W.G.; Bernstein, H.J. High-Resolution Nuclear Magnetic Resonance. New York: McGraw-Hill, 1959.
nmrsim.discrete.AABB(Vab, Jaa, Jbb, Jab, Jab_prime, Vcentr, normalize=True, **kwargs)[source]

A wrapper for a second-order AA’BB’ calculation, but using the same arguments as WINDNMR.

Parameters:
  • Vab (float) – the difference in frequency (Hz) between Ha and Hb in the absence of coupling. A positive number indicates vb > va.
  • Jbb, Jab, Jab_prime (Jaa,) – Jaa is the JAA’ coupling constant; Jxx the JXX’; Jax the JAX; and JAX_prime the JAX’.
  • Vcentr (float) – the frequency for the center of the signal.
  • normalize (bool) – whether the signal intensity should be normalized.
Returns:

a list of (frequency, intensity) tuples.

Return type:

[(float, float)..]

nmrsim.discrete.AAXX(Jaa, Jxx, Jax, Jax_prime, Vcentr, normalize=True)[source]

Simulates one half (‘A’ part) of an AA’XX’ spin system.

All frequencies are in Hz.

Parameters:
  • Jax, Jax, Jax_prime (Jaa,) – Jaa is the JAA’ coupling constant; Jxx the JXX’; Jax the JAX; and JAX_prime the JAX’.
  • Vcentr (float) – the frequency for the center of the signal.
  • normalize (bool) – whether the signal intensity should be normalized (to 2).
Returns:

a list of (frequency, intensity) tuples.

Return type:

[(float, float)..]

nmrsim.discrete.AB(Jab, Vab, Vcentr, normalize=True)[source]

Calculates the signal frequencies and intensities for two strongly coupled protons (Ha and Hb).

Parameters:
  • Jab (float) – The coupling constant (Hz) between Ha and Hb
  • Vab (float) – The chemical shift difference (Hz) between Ha and Hb in the absence of coupling.
  • Vcentr (float) – The frequency (Hz) for the center of the AB quartet.
  • normalize (bool) – Whether the signal intensity should be normalized.
Returns:

A list of four (frequency, intensity) tuples.

Return type:

[(float, float)..]

nmrsim.discrete.AB2(Jab, Vab, Vcentr, normalize=True)[source]

Calculates signal frequencies and intensities for an AB2 spin system.

Parameters:
  • Jab (float) – the Ha-Hb coupling constant (Hz).
  • Vab (float) – the difference in the frequencies (Hz). A positive value means vb > va; negative means va > vb.
  • Vcentr (float) – the frequency (Hz) for the center of the AB2 signal.
  • normalize (bool) – whether the signal intensity should be normalized.
Returns:

a list of (frequency, intensity) tuples.

Return type:

[(float, float)..]

nmrsim.discrete.ABX(Jab, Jax, Jbx, Vab, Vcentr, vx, normalize=True)[source]

Non-QM approximation for an ABX spin system. The approximation assumes that Hx is very far away in chemical shift from Ha/Hb.

Parameters:
  • Jab (float) – The Ha-Hb coupling constant (Hz).
  • Jax (float) – The Ha-Hx coupling constant (Hz).
  • Jbx (float) – The Hb-Hx coupling constant (Hz).
  • Vab (float) – The difference in the frequencies (in the absence of coupling) of Ha and Hb (Hz).
  • Vcentr (float) – The frequency (Hz) for the center of the AB signal.
  • vx (float) – The frequency (Hz) for Hx in the absence of coupling.
  • normalize (bool (optional)) – whether the signal intensity should be normalized. If false, the total signal intensity happens to be ~12.
Returns:

a list of (frequency, intensity) tuples.

Return type:

[(float, float)..]

nmrsim.discrete.ABX3(Jab, Jax, Jbx, Vab, Vcentr)[source]

Simulation of the AB part of an ABX3 spin system.

Parameters:
  • Jab (float) – the Ha-Hb coupling constant (Hz).
  • Jax (float) – the Ha-Hb coupling constant (Hz).
  • Jbx (float) – the Ha-Hb coupling constant (Hz).
  • Vab (float) – the difference in the frequencies (Hz) of Ha and Hb in the absence of coupling. Positive when vb > va.
  • Vcentr (float) – the frequency (Hz) for the center of the AB signal.
Returns:

a list of (frequency, intensity) tuples.

Return type:

[(float, float)..]

nmrsim.dnmr module

The dnmr module provides functions for calculating DNMR line shapes, and classes to describe DNMR systems.

The dnmr module provides the following classes:

  • DnmrTwoSinglets: a sumulation of the lineshape for two uncoupled nuclei undergoing exchange.
  • DnmrAB: a simulation of the lineshape for two coupled nuclei undergoing exchange (i.e. an AB (or AX) pattern at the slow exchange limit).

The dnmr module provides the following functions:

  • dnmr_two_singlets: for simulating the lineshape for two uncoupled nuclei undergoing exchange [3].
  • dnmr_AB : for simulating the lineshape for two coupled nuclei undergoing exchange (i.e. an AB (or AX) pattern at the slow exchange limit) [4].

References

[3]Sandström, J. Dynamic NMR Spectroscopy; Academic Press: New York, 1982.
[4]
  1. Brown, K.C.; Tyson, R.L.; Weil, J.A. J. Chem. Educ. 1998, 75, 1632.
  2. an important math correction to the previous reference:
TODO: add reference to correction
class nmrsim.dnmr.DnmrAB(va=165.0, vb=135.0, J=12.0, k=12.0, w=0.5, limits=None, points=800)[source]

Bases: object

Simulate the DNMR lineshape for two coupled nuclei undergoing exchange (AB or AX pattern at the slow-exchange limit).

Parameters:
  • vb (va,) – frequencies of a and b nuclei (at the slow exchange limit, in the absence of coupling)
  • J (int or float) – the coupling constant between the two nuclei.
  • k (int or float) – rate constant for state A–> state B
  • w (int or float) – peak widths at half height (at the slow-exchange limit).
  • limits ((int or float, int or float), optional) – The minimum and maximum frequencies (in any order) for the simulation.
  • points (int) – The length of the returned arrays (i.e. the number of points plotted).

See also

DnmrAB
A class representation for this simulation.

References

See the documentation for the nmrsim.dnmr module.

J

The coupling constant (Hz) between the two nuclei.

Returns:
Return type:int or float
k

The rate constant (Hz) for state A–> state B (must be >0).

Returns:
Return type:int or float
limits

Give minimum and maximum frequencies for the simulated lineshape.

Returns:
Return type:(int or float, int or float)
lineshape()[source]

Return the x, y lineshape data for the simulation.

Returns:x, y – Arrays for the x (frequency) and y (intensity) lineshape data points.
Return type:numpy.array, numpy.array
points

Give the length of the returned arrays (i.e. the number of points plotted).

Returns:
Return type:int
va

The frequency of nucleus “a” (Hz) at the slow-exchange limit, in the absence of coupling.

Returns:
Return type:int or float
vb

The frequency of nucleus “b” (Hz) at the slow-exchange limit, in the absence of coupling.

Returns:
Return type:int or float
w

The peak width (Hz) at half height (at the slow-exchange limit).

Returns:
Return type:int or float
class nmrsim.dnmr.DnmrTwoSinglets(va=1, vb=0, k=0.01, wa=0.5, wb=0.5, pa=0.5, limits=None, points=800)[source]

Bases: object

A DNMR simulation for two uncoupled nuclei undergoing exchange.

Parameters:
  • vb (va,) – The frequencies (Hz) of nuclei ‘a’ and ‘b’ at the slow exchange limit.
  • k (int or float) – The rate constant (Hz) for state a–> state b
  • wb (wa,) – The peak widths at half height for the ‘a’ and ‘b’ singlets at the slow-exchange limit.
  • pa (float (0 <= pa <= 1)) – The fraction of the population in state a
  • limits ((int or float, int or float), optional) – The minimum and maximum frequencies (in any order) for the simulation.
  • points (int) – The length of the returned arrays (i.e. the number of points plotted).

See also

DnmrTwoSinglets
A class representation for this simulation
k

The rate constant (Hz) for state A–> state B (must be >0).

Returns:
Return type:int or float
limits

The minimum and maximum frequencies for the simulated lineshape.

Returns:
Return type:(int or float, int or float)
lineshape()[source]

Calculate and return the lineshape for the DNMR spectrum.

Returns:x, y – Arrays for the x (frequency) and y (intensity) lineshape data points.
Return type:numpy.array, numpy.array
pa

The fraction of the population in state a. Must be >=0 and <=1.

Returns:
Return type:float
points

The length of the returned arrays (i.e. the number of points plotted).

Returns:
Return type:int
va

The frequency of nucleus “a” (Hz) at the slow-exchange limit.

Returns:
Return type:int or float
vb

The frequency of nucleus “b” (Hz) at the slow-exchange limit.

Returns:
Return type:int or float
wa

The peak width at half height (Hz) for the ‘a’ singlet at the slow-exchange limit.

Returns:
Return type:int or float
wb

The peak width at half height (Hz) for the ‘b’ singlet at the slow-exchange limit.

Returns:
Return type:int or float
nmrsim.dnmr.dnmr_AB(va, vb, J, k, w, limits=None, points=800)[source]

Simulate the DNMR lineshape for two coupled nuclei undergoing exchange (AB or AX pattern at the slow-exchange limit).

Parameters:
  • vb (va,) – frequencies of a and b nuclei (at the slow exchange limit, in the absence of coupling)
  • J (float) – the coupling constant between the two nuclei.
  • k (float) – rate constant for state A–> state B
  • w (float) – peak widths at half height (at the slow-exchange limit).
  • limits ((int or float, int or float), optional) – The minimum and maximum frequencies (in any order) for the simulation.
  • points (int) – The length of the returned arrays (i.e. the number of points plotted).
Returns:

x, y – Arrays for the x (frequency) and y (intensity) lineshape data points.

Return type:

numpy.array, numpy.array

See also

DnmrAB()
A class representation for this simulation.

References

See the documentation for the nmrsim.dnmr module.

nmrsim.dnmr.dnmr_two_singlets(va, vb, ka, wa, wb, pa, limits=None, points=800)[source]

Create a the lineshape for a DNMR spectrum of two uncoupled spin-half nuclei.

Parameters:
  • vb (va,) – The frequencies (Hz) of nuclei ‘a’ and ‘b’ at the slow exchange limit.
  • ka (int or float) – The rate constant (Hz) for state a–> state b
  • wb (wa,) – The peak widths at half height for the ‘a’ and ‘b’ singlets at the slow-exchange limit.
  • pa (float (0 <= pa <= 1)) – The fraction of the population in state a
  • limits ((int or float, int or float), optional) – The minimum and maximum frequencies (in any order) for the simulation.
  • points (int) – The length of the returned arrays (i.e. the number of points plotted).
Returns:

x, y – Arrays for the x (frequency) and y (intensity) lineshape data points.

Return type:

numpy.array, numpy.array

See also

DnmrTwoSinglets()
A class representation for this simulation.

References

See the documentation for the nmrsim.dnmr module.

nmrsim.firstorder module

“Functions for calculating first-order spectra.

The nmrsim.firstorder module provides the following functions:

  • multiplet: performs first-order splitting of a signal into multiple signals.
  • first_order_spin_system: provides a peaklist for several nuclei, using the
    same v/J parameters that are used for second-order spin systems. See nmrsim.qm for details on these parameters.
nmrsim.firstorder.first_order_spin_system(v, J)[source]

Create a first-order peaklist of several multiplets from the same v/J arguments used for qm calculations.

This allows a user to model several multiplets at once, rather than creating each multiplet individually. It also provides a “toggle” where the user, or a higher-level function/class (such as nmrsim.SpinSystem) can decide whether a spin system is modeled as first order or second order.

Parameters:
  • v (array-like [float...]) – an array of frequencies
  • J (2D array-like (square)) – a matrix of J coupling constants
Returns:

a combined peaklist of signals for all the multiplets in the spin system.

Return type:

[(float, float)..]

nmrsim.firstorder.multiplet(signal, couplings)[source]

Splits a set of signals into first-order multiplets.

Parameters:
  • signal ((float, float)) – a (frequency (Hz), intensity) tuple;
  • couplings ([(float, int)..]) – A list of (J, # of nuclei) tuples. The order of the tuples in couplings does not matter. e.g. to split a signal into a dt, J = 8, 5 Hz, use: couplings = [(8, 2), (5, 3)]
Returns:

a sorted peaklist for the multiplet that results from splitting the signal by each J.

Return type:

[(float, float)..]

nmrsim.math module

A collection of functions for processing simulated NMR spectra.

Terms used: signal: a pair (e.g. tuple) of frequency, intensity values peaklist: a list (or 1D-array-like) of signals.

Provides the following functions:

  • add_peaks: combines a list of signals into one signal of average frequency and summed intensity.
  • reduce_peaks: processes a peaklist so that signals within a frequency tolerance are added together.
  • normalize_peaklist: scales a peaklist so that intensities add to a specific value.
  • lorentz: given a frequency, a signal and a linewidth, calculates an intensity. Used to calculate Lorentzian lineshapes for signals.
  • get_intensity: given a lineshape and a frequency, find the intensity at the datapoint closest to that frequency.
nmrsim.math.add_lorentzians(linspace, peaklist, w)[source]

Given a numpy linspace, a peaklist of (frequency, intensity) tuples, and a linewidth, returns an array of y coordinates for the total line shape.

Parameters:
  • linspace (array-like) – Normally a numpy.linspace of x coordinates corresponding to frequency in Hz.
  • peaklist ([(float, float)..]) – A list of (frequency, intensity) tuples.
  • w (float) – Peak width at half maximum intensity.
Returns:

an array of y coordinates corresponding to intensity.

Return type:

[float…]

nmrsim.math.add_peaks(plist)[source]

Reduces a list of (frequency, intensity) tuples to an (average frequency, total intensity) tuple.

Parameters:plist ([(float, float)..]) – a list of (frequency, intensity) tuples
Returns:a tuple of (average frequency, total intensity)
Return type:(float, float)
nmrsim.math.get_intensity(lineshape, x)[source]

A crude method to find the intensity of data point closest to frequency x. Better: interpolate between two data points if match isn’t exact (TODO?)

Parameters:
  • lineshape (tuple of (x, y) arrays for frequency, intensity data) –
  • x (frequency lookup) –
Returns:

float

Return type:

the intensity at that frequency

nmrsim.math.get_maxima(lineshape)[source]

Crude function that returns maxima in the lineshape.

Parameters:lineshape (tuple of frequency, intensity arrays) –
Returns:
Return type:a list of (frequency, intensity) tuples for individual maxima.
nmrsim.math.lorentz(v, v0, I, w)[source]

A lorentz function that takes linewidth at half intensity (w) as a parameter.

When v = v0, and w = 0.5 (Hz), the function returns intensity I.

Parameters:
  • v (float) – The frequency (x coordinate) in Hz at which to evaluate intensity (y coordinate).
  • v0 (float) – The center of the distribution.
  • I (float) – the relative intensity of the signal
  • w (float) – the peak width at half maximum intensity
Returns:

the intensity (y coordinate) for the Lorentzian distribution evaluated at frequency v.

Return type:

float

nmrsim.math.normalize_peaklist(peaklist, n=1)[source]

Normalize the intensities in a peaklist so that total intensity equals value n (nominally the number of nuclei giving rise to the signal).

Parameters:
  • peaklist ([(float, float)..]) – a list of (frequency, intensity) tuples.
  • n (int or float (optional)) – total intensity to normalize to (default = 1).
nmrsim.math.reduce_peaks(plist_, tolerance=0)[source]

Takes a list of (x, y) tuples and adds together tuples whose values are within a certain tolerance limit.

Parameters:
  • plist ([(float, float)..]) – A list of (x, y) tuples
  • tolerance (float) – tuples that differ in x by <= tolerance are combined using add_peaks
Returns:

a list of (x, y) tuples where all x values differ by > tolerance

Return type:

[(float, float)..]

nmrsim.plt module

The plt module provides convenience functions for creating matplotlib plots, plus applying Lorentzian distributions about signals.

The plt module provides the following functions:

  • add_lorentzians: Creates lineshape data from a provided linspace (array of x coordinates) and peaklist).
  • mplplot: Creates a lineshape plot from a peaklist and returns the x, y plot data.
  • mplplot_stick: Creates a “stick” (matplotlib “stem” plot) plot from a peaklist and returns the x, y plot data.
  • mplplot_lineshape: Creates a lineshape plot from provided x, y lineshape data and returns the x, y plot data.
nmrsim.plt.mplplot(peaklist, w=1, y_min=-0.01, y_max=1, points=800, limits=None, hidden=False)[source]

A matplotlib plot of the simulated lineshape for a peaklist.

Parameters:
  • peaklist ([(float, float)..]) – A list of (frequency, intensity) tuples.
  • w (float) – Peak width at half height
  • y_min (float or int) – Minimum intensity for the plot.
  • y_max (float or int) – Maximum intensity for the plot.
  • points (int) – Number of data points.
  • limits ((float, float)) – Frequency limits for the plot.
  • hidden (bool) – Whether showing the plot should be omitted (e.g. to not block CI tests)
Returns:

x, y – Arrays for frequency (x) and intensity (y) for the simulated lineshape.

Return type:

numpy.array

nmrsim.plt.mplplot_lineshape(x, y, y_min=None, y_max=None, limits=None, hidden=False)[source]

A matplotlib plot that accepts arrays of x and y coordinates.

Parameters:
  • x (array-like) – The list of x coordinates.
  • y (array-like) – The list of y coordinates.
  • y_min (float or int) – Minimum intensity for the plot. Default is -10% max y.
  • y_max (float or int) – Maximum intensity for the plot. Default is 110% max y.
  • limits ((float, float)) – Frequency limits for the plot.
  • hidden (bool) – Whether showing the plot should be omitted (e.g. to not block CI tests)
Returns:

x, y

Return type:

The original x, y arguments.

nmrsim.plt.mplplot_stick(peaklist, y_min=-0.01, y_max=1, limits=None, hidden=False)[source]

A matplotlib plot of a spectrum in “stick” (stem) style.

Parameters:
  • peaklist ([(float, float)..]) – A list of (frequency, intensity) tuples.
  • y_min (float or int) – Minimum intensity for the plot.
  • y_max (float or int) – Maximum intensity for the plot.
  • limits ((float, float)) – Frequency limits for the plot.
  • hidden (bool) – Whether showing the plot should be omitted (e.g. to not block CI tests)
Returns:

The arrays of x and y coordinates used for the plot.

Return type:

numpy.array, numpy.array

nmrsim.qm module

qm contains functions for the quantum-mechanical (second-order) calculation of NMR spectra.

The qm module provides the following attributes:

  • CACHE : bool (default True)
    Whether saving to disk of partial solutions is allowed.
  • SPARSE : bool (default True)
    Whether the sparse library can be used.

The qm module provides the following functions:

  • qm_spinsystem: The high-level function for computing a second-order simulation from frequency and J-coupling data.
  • hamiltonian_dense: Calculate a spin Hamiltonian using dense arrays (slower).
  • hamiltonian_sparse: Calculate a spin Hamiltonian using cached sparse arrays (faster).
  • solve_hamiltonian: Calculate a peaklist from a spin Hamiltonian.
  • secondorder_dense: Calculate a peaklist for a second-order spin system, using dense arrays (slower).
  • secondorder_sparse: Calculate a peaklist for a second-order spin system, using cached sparse arrays (faster).

Notes

Because numpy.matrix is marked as deprecated, starting with Version 0.2.0 the qm code was refactored to a) accommodate this deprecation and b) speed up the calculations. The fastest calculations rely on:

1. the pydata/sparse library. SciPy’s sparse depends on numpy.matrix, and they currently recommend that pydata/sparse be used for now.

2. Caching partial solutions for spin operators and transition matrices as .npz files.

If the pydata/sparse package is no longer available, and/or if distributing the library with .npz files via PyPI is problematic, then a backup is required. The qm module for now provides two sets of functions for calculating second-order spectra: one using pydata/sparse and caching, and the other using neither.

nmrsim.qm.hamiltonian_dense(v, J)[source]

Calculate the spin Hamiltonian as a dense array.

Parameters:
  • v (array-like) – list of frequencies in Hz (in the absence of splitting) for each nucleus.
  • J (2D array-like) – matrix of coupling constants. J[m, n] is the coupling constant between v[m] and v[n].
Returns:

H – a sparse spin Hamiltonian.

Return type:

numpy.ndarray

nmrsim.qm.hamiltonian_sparse(v, J)[source]

Calculate the spin Hamiltonian as a sparse array.

Parameters:
  • v (array-like) – list of frequencies in Hz (in the absence of splitting) for each nucleus.
  • J (2D array-like) – matrix of coupling constants. J[m, n] is the coupling constant between v[m] and v[n].
Returns:

H – a sparse spin Hamiltonian.

Return type:

sparse.COO

nmrsim.qm.qm_spinsystem(*args, cache=True, sparse=True, **kwargs)[source]

Calculates second-order spectral data (frequency and intensity of signals) for n spin-half nuclei.

Currently, n is capped at 11 spins.

Parameters:
  • freqs ([float...]) – a list of n nuclei frequencies in Hz.
  • couplings (array-like) – An n, n array of couplings in Hz. The order of nuclei in the list corresponds to the column and row order in the matrix, e.g. couplings[0][1] and [1]0] are the J coupling between the nuclei of freqs[0] and freqs[1].
  • normalize (bool (optional keyword argument; default = True)) – True if the intensities should be normalized so that total intensity equals the total number of nuclei.
Returns:

peaklist – of [frequency, intensity] pairs.

Return type:

[[float, float]..] numpy 2D array

Other Parameters:
 
  • cache (bool (default = nmrsim.qm.CACHE)) – Whether caching of partial solutions (for acceleration) is allowed. Currently CACHE = True, but this provides a hook to modify nmrsim for platforms such as Raspberry Pi where storage space is a concern.
  • sparse (bool (default = nmrsim.qm.SPARSE)) – Whether the pydata sparse library for sparse matrices is available. Currently SPARSE = True, but this provides a hook to modify nmrsim should the sparse library become unavailable (see notes).
  • cutoff (float) – The intensity cutoff for reporting signals (default is 0.001).

Notes

With numpy.matrix marked for deprecation, the scipy sparse array functionality is on shaky ground, and the current recommendation is to use the pydata sparse library. In case a problem arises in the numpy/scipy/ sparse ecosystem, SPARSE provides a hook to use a non-sparse-dependent alternative.

nmrsim.qm.secondorder_dense(freqs, couplings, normalize=True, **kwargs)[source]

Calculates second-order spectral data (freqency and intensity of signals) for n spin-half nuclei.

Parameters:
  • freqs ([float...]) – a list of n nuclei frequencies in Hz
  • couplings (array-like) – an n, n array of couplings in Hz. The order of nuclei in the list corresponds to the column and row order in the matrix, e.g. couplings[0][1] and [1]0] are the J coupling between the nuclei of freqs[0] and freqs[1].
  • normalize (bool) – True if the intensities should be normalized so that total intensity equals the total number of nuclei.
Returns:

peaklist – numpy 2D array of [frequency, intensity] pairs.

Return type:

[[float, float]..]

Other Parameters:
 

cutoff (float) – The intensity cutoff for reporting signals (default is 0.001).

nmrsim.qm.secondorder_sparse(freqs, couplings, normalize=True, **kwargs)[source]

Calculates second-order spectral data (frequency and intensity of signals) for n spin-half nuclei.

Parameters:
  • freqs ([float...]) – a list of n nuclei frequencies in Hz
  • couplings (array-like) – an n, n array of couplings in Hz. The order of nuclei in the list corresponds to the column and row order in the matrix, e.g. couplings[0][1] and [1]0] are the J coupling between the nuclei of freqs[0] and freqs[1].
  • normalize (bool) – True if the intensities should be normalized so that total intensity equals the total number of nuclei.
Returns:

peaklist – of [frequency, intensity] pairs.

Return type:

[[float, float]..] numpy 2D array

Other Parameters:
 

cutoff (float) – The intensity cutoff for reporting signals (default is 0.001).

nmrsim.qm.solve_hamiltonian(H, nspins, **kwargs)[source]

Calculates frequencies and intensities of signals from a spin Hamiltonian and number of spins.

Parameters:
  • H (numpy.ndarray (2D)) – The spin Hamiltonian
  • nspins (int) – The number of spins in the system
Returns:

Return type:

[[float, float]..] numpy 2D array of frequency, intensity pairs.

Other Parameters:
 

cutoff (float) – The intensity cutoff for reporting signals (default is 0.001).

Module contents

nmrsim

The nmrsim package provides tools for simulating nuclear magnetic resonance (NMR) spectra.

The API is still in flux. Currently, it includes the following modules:

  • dnmr: for modeling Dynamic NMR systems
  • firstorder: for modeling first-order spectra
  • math: core math routines for handling NMR data
  • partial: uses non-quantum-mechanical solutions for common second-order NMR patterns
  • plt: convenience plotting routines for NMR results
  • qm: quantum-mechanical second-order simulation of spin systems (currently capped at 11 nuclei)

Currently, only spin-1/2 nuclei are accommodated.

The top-level nmrsim namespace provides the following classes:

  • Multiplet: a representation of a first-order multiplet (e.g. quartet; doublet of triplets).
  • SpinSystem: a representation of a set of coupled nuclei (modeled as either first-order or second-order).
  • Spectrum: a collection of components such as Multiplets or SpinSystems that contribute to a total NMR spectrum simulation.

Definitions of Terms Used

In naming classes, functions, methods, data types etc. certain phrases, taken from NMR nomenclature, have the following interpretations:

  • multiplet (e.g. the nmrsim.Multiplet class): a first-order simulation for one signal (i.e. one or more chemical shift-equivalent nuclei). Examples: doublet, triplet, doublet of triplets, but not an AB quartet (which is a second-order pattern for two nuclei).
  • spin system (e.g. the SpinSystem class): a simulation of a set of coupled nuclei.
  • spectrum (e.g. the Spectrum class): a complete collection of first- and/or second-order components for simulating a total NMR spectrum. ‘Spectrum’ can also refer in general to the simulation results for the system, e.g a peaklist or lineshape (see below).
  • peak: a pair of frequency (Hz), intensity values corresponding to a resonance in an NMR spectrum. For example, a 1H triplet centered at 100 Hz with J = 10 Hz would have the following peaks: (110, 0.25), (100, 0.5), (90, 0.25).
  • peaklist: a list of peaks (e.g. [(110, 0.25), (100, 0.5), (90, 0.25)] for the above triplet).
  • lineshape: a pair of [x coordinates…], [y coordinates] arrays for plotting the lineshape of a spectrum.

The following idioms are used for arguments: * v for a frequency or list of frequencies (similar to the Greek lowercase “nu” character). * I for a signal intensity (despite being a PEP8 naming violation) * J for coupling constant data (exact format depends on the implementation).