API
NMRBase
NMRTools.NMRBase.AbstractNMRData
— TypeAbstractNMRData <: DimensionalData.AbstractDimArray
Abstract supertype for objects that wrap an array of NMR data, and metadata about its contents.
AbstractNMRData
s inherit from AbstractDimArray
from DimensionalData.jl. They can be indexed as regular Julia arrays or with DimensionalData.jl Dimension
s.
NMRTools.NMRBase.Coherence
— TypeNMRTools.NMRBase.ComplexLineshape
— TypeComplexLineshape
Return a complex-valued lineshape when used in calculations
NMRTools.NMRBase.CosWindow
— TypeCosWindow(tmax)
Apodization by a pure cosine function. Acquisition time is tmax
.
See also Cos²Window
, SineWindow
.
NMRTools.NMRBase.Cos²Window
— TypeCos²Window(tmax)
Apodization by a pure cosine squared function. Acquisition time is tmax
.
See also CosWindow
, SineWindow
.
NMRTools.NMRBase.ExponentialWindow
— TypeExponentialWindow(lb, tmax)
Exponential window function, with a line broadening of lb
Hz. Acquisition time is tmax
.
NMRTools.NMRBase.FrequencyDimension
— TypeFrequencyDimension <: NMRDimension
Abstract supertype for frequency dimensions used in NMRData objects. Concrete types F1Dim
, F2Dim
, F3Dim
and F4Dim
are generated for use in creating objects.
See also NonFrequencyDimension
.
NMRTools.NMRBase.GaussWindow
— TypeGaussWindow(expHz, gaussHz, center, tmax)
Abstract representation of Lorentz-to-Gauss window functions, applying an inverse exponential of expHz
Hz, and a gaussian broadening of gaussHz
Hz, with maximum at center
(between 0 and 1). Acquisition time is tmax
.
Specialises to LorentzToGaussWindow
when center
is zero, otherwise GeneralGaussWindow
.
NMRTools.NMRBase.GradientDimension
— TypeGradientDimension <: NonFrequencyDimension <: NMRDimension
Abstract supertype for gradient-encoded dimensions used in NMRData objects. Concrete types G1Dim
, G2Dim
, G3Dim
and G4Dim
are generated for use in creating objects.
NMRTools.NMRBase.HasNonFrequencyDimension
— Type@traitdef HasNonFrequencyDimension{D}
A trait indicating whether the data object has a non-frequency domain dimension.
Example
@traitfn f(x::X) where {X; HasNonFrequencyDimension{X}} = "This spectrum has a non-frequency domain dimension!"
@traitfn f(x::X) where {X; Not{HasNonFrequencyDimension{X}}} = "This is a pure frequency-domain spectrum!"
NMRTools.NMRBase.LineshapeComplexity
— TypeLineshapeComplexity
Abstract type to specify calculation of RealLineshape
or ComplexLineshape
in function calls.
NMRTools.NMRBase.MQ
— TypeMQ(coherences, label=="")
Representation of a multiple-quantum coherence. Coherences are specified as a tuple of tuples, of the form (nucleus, coherenceorder)
Examples
julia> MQ(((H1,1), (C13,-1)), "ZQ")
MQ(((H1, 1), (C13, -1)), "ZQ")
julia> MQ(((H1,3), (C13,1)), "QQ")
MQ(((H1, 3), (C13, 1)), "QQ")
NMRTools.NMRBase.NMRData
— TypeNMRData <: AbstractNMRData
NMRData(A::AbstractArray{T,N}, dims; kw...)
NMRData(A::AbstractNMRData; kw...)
A generic AbstractNMRData
for NMR array data. It holds memory-backed arrays.
Keywords
dims
:Tuple
ofNMRDimension
s for the array.name
:Symbol
name for the array, which will also retreive named layers ifNMRData
is used on a multi-layered file like a NetCDF.missingval
: value representing missing data, normally detected from the file. Set manually when you know the value is not specified or is incorrect. This will not change any values in the NMRData, it simply assigns which value is treated as missing.
Internal Keywords
In some cases it is possible to set these keywords as well.
data
: can replace the data in anAbstractNMRData
refdims
:Tuple of
positionDimension
s the array was sliced from, defaulting to()
.
NMRTools.NMRBase.NMRDimension
— TypeNMRDimension
Abstract supertype for all axes used in NMRData objects.
See also FrequencyDimension
and NonFrequencyDimension
.
NMRTools.NMRBase.NMRToolsError
— TypeNMRToolsError(message)
An error arising in NMRTools.
NMRTools.NMRBase.NonFrequencyDimension
— TypeNonFrequencyDimension <: NMRDimension
Abstract supertype for non-frequency dimensions used in NMRData objects. Sub-types include TimeDimension
, GradientDimension
, and UnknownDimension
.
See also FrequencyDimension
.
NMRTools.NMRBase.Nucleus
— TypeNucleus
Enumeration of common nuclei associated with biomolecular NMR. Nuclei are named e.g. H1
, C13
.
Defined nuclei: H1
, H2
, C12
, C13
, N14
, N15
, F19
, P31
.
See also spin
, gyromagneticratio
, Coherence
.
NMRTools.NMRBase.NullWindow
— TypeNullWindow(tmax)
No apodization applied. Acquisition time is tmax
.
NMRTools.NMRBase.RealLineshape
— TypeRealLineshape
Return a real-valued lineshape when used in calculations
NMRTools.NMRBase.SQ
— TypeSQ(nucleus::Nucleus, label=="")
Representation of a single quantum coherence on a given nucleus.
NMRTools.NMRBase.SineWindow
— TypeSineWindow(offset, endpoint, power, tmax)
Abstract window function representing multiplication by sine/cosine functions. Acquisition time is tmax
.
\[\left[\sin\left( \pi\cdot\mathrm{offset} + \frac{\left(\mathrm{end} - \mathrm{offset}\right)\pi t}{\mathrm{tmax}} \right)\right]^\mathrm{power}\]
Specialises to CosWindow
, Cos²Window
or GeneralSineWindow
.
Arguments
offset
: initial value is $\sin(\mathrm{offset}\cdot\pi)$ (0 to 1)endpoint
: initial value is $\sin(\mathrm{endpoint}\cdot\pi)$ (0 to 1)pow
: sine exponent
NMRTools.NMRBase.TimeDimension
— TypeTimeDimension <: NonFrequencyDimension <: NMRDimension
Abstract supertype for time dimensions used in NMRData objects. Concrete types T1Dim
, T2Dim
, T3Dim
and T4Dim
are generated for time-domains representing frequency evolution, and TrelaxDim
and TkinDim
are generated for representing relaxation and real-time kinetics.
NMRTools.NMRBase.UnknownDimension
— TypeUnknownDimension <: NonFrequencyDimension <: NMRDimension
Abstract supertype for unknown, non-frequency dimensions used in NMRData objects. Concrete types X1Dim
, X2Dim
, X3Dim
and X4Dim
are generated for use in creating objects.
NMRTools.NMRBase.UnknownWindow
— TypeUnknownWindow(tmax)
Unknown apodization applied. Acquisition time is tmax
.
NMRTools.NMRBase.WindowFunction
— TypeWindowFunction
Abstract type to represent apodization functions.
Window functions are represented by subtypes of the abstract type WindowFunction
, each of which contain appropriate parameters to specify the particular function applied. In addition, the acquisition time tmax
is also stored (calculated at the point the window function is applied, i.e. after linear prediction but before zero filling).
Base.stack
— Methodstack(expts::Vector{NMRData})
Combine a collection of equally-sized NMRData into one larger array, by arranging them along a new dimension, of type UnknownDimension
.
Throws a DimensionMismatch
if data are not of compatible shapes.
DimensionalData.Dimensions.LookupArrays.metadata
— Functionmetadata(nmrdata, key)
metadata(nmrdata, dim, key)
metadata(nmrdimension, key)
Return the metadata for specified key, or nothing
if not found. Keys are passed as symbols.
Examples (spectrum metadata)
:ns
: number of scans:ds
: number of dummy scans:rg
: receiver gain:ndim
: number of dimensions:title
: spectrum title (contents of title pdata file):filename
: spectrum filename:pulseprogram
: title of pulse program used for acquisition:experimentfolder
: path to experiment:noise
: RMS noise level
Examples (dimension metadata)
:pseudodim
: flag indicating non-frequency domain data:npoints
: final number of (real) data points in dimension (after extraction):td
: number of complex points acquired:tdzf
: number of complex points when FT executed, including LP and ZF:bf
: base frequency, in MHz:sf
: carrier frequency, in MHz:offsethz
: carrier offset from bf, in Hz:offsetppm
: carrier offset from bf, in ppm:swhz
: spectrum width, in Hz:swppm
: spectrum width, in ppm:region
: extracted region, expressed as a range in points, otherwise missing:window
:WindowFunction
indicating applied apodization:referenceoffset
: referencing (in ppm) applied to the dimension
See also estimatenoise!
.
DimensionalData.Dimensions.LookupArrays.units
— Functionunits(nmrdata)
units(nmrdata, dim)
units(nmrdimension)
Return the physical units associated with an NMRData
structure or an NMRDimension
.
DimensionalData.Dimensions.label
— Functionlabel(nmrdata)
label(nmrdata, dim)
label(nmrdimension)
Return a short label associated with an NMRData
structure or an NMRDimension
. By default, for a spectrum this is obtained from the first line of the title file. For a frequency dimension, this is normally something of the form 1H chemical shift (ppm)
.
See also label!
.
NMRTools.NMRBase._lineshape
— Function_lineshape(ω, R2, ωaxis, window, complexity)
Internal function to calculate a resonance lineshape with frequency ω
and relaxation rate R2
, calculated at frequencies ωaxis
and with apodization according to the specified window function.
NMRTools.NMRBase.acqus
— Methodacqus(nmrdata)
acqus(nmrdata, key)
acqus(nmrdata, key, index)
Return data from a Bruker acqus file, or nothing
if it does not exist. Keys can be passed as symbols or strings. If no key is specified, a dictionary is returned representing the entire acqus file.
If present, the contents of auxilliary files such as vclist
and vdlist
can be accessed using this function.
Examples
julia> acqus(expt, :pulprog)
"zgesgp"
julia> acqus(expt, "TE")
276.9988
julia> acqus(expt, :p, 1)
9.2
julia> acqus(expt, "D", 1)
0.1
julia> acqus(expt, :vclist)
11-element Vector{Int64}:
[...]
See also metadata
.
NMRTools.NMRBase.add_offset
— Methodadd_offset!(data::NMRData, dim_ref, offset)
Add an offset to a frequency dimension in an NMRData object. The dimension can be specified as a numerical index or an object like F1Dim
. The metadata is copied using replacedimension
, and an entry is added or updated in the dimension metadata to record the offset change.
NMRTools.NMRBase.apod
— Functionapod(spec::NMRData, dimension, zerofill=true)
Return the time-domain apodization function for the specified axis, as a vector of values.
NMRTools.NMRBase.coherenceorder
— Functioncoherenceorder(coherence)
Calculate the total coherence order.
Examples
julia> coherenceorder(SQ(H1))
1
julia> coherenceorder(MQ(((H1,1),(C13,1))))
2
julia> coherenceorder(MQ(((H1,1),(C13,-1))))
0
julia> coherenceorder(MQ(((H1,3),(C13,1))))
4
julia> coherenceorder(MQ(((H1,0),)))
0
NMRTools.NMRBase.data
— Methoddata(nmrdata, dim)
Return the numerical data associated with the specified dimension.
NMRTools.NMRBase.data
— Methoddata(nmrdata)
Return the numerical data associated with the specified NMRData.
NMRTools.NMRBase.data
— Methoddata(nmrdimension)
Return the numerical data associated with an NMR dimension.
NMRTools.NMRBase.decimate
— Functiondecimate(data, n dims=1)
Decimate NMR data into n-point averages along the specified dimension. Note that data are averaged and not summed. Noise metadata is not updated.
NMRTools.NMRBase.estimatenoise!
— Methodestimatenoise!(nmrdata)
Estimate the rms noise level in the data and update :noise
metadata.
If called on an Array
of data, each item will be updated.
Algorithm
Data are sorted into numerical order, and the highest and lowest 12.5% of data are discarded (so that 75% of the data remain). These values are then fitted to a truncated gaussian distribution via maximum likelihood analysis.
The log-likelihood function is:
\[\log L(\mu, \sigma) = \sum_i{\log P(y_i, \mu, \sigma)}\]
where the likelihood of an individual data point is:
\[\log P(y,\mu,\sigma) = \log\frac{ \phi\left(\frac{x-\mu}{\sigma}\right) }{ \sigma \cdot \left[\Phi\left(\frac{b-\mu}{\sigma}\right) - \Phi\left(\frac{a-\mu}{\sigma}\right)\right]}\]
and $\phi(x)$ and $\Phi(x)$ are the standard normal pdf and cdf functions.
NMRTools.NMRBase.getω
— Methodgetω(axis, δ)
Return the offset (in rad/s) for a chemical shift (or list of shifts) on a frequency axis.
NMRTools.NMRBase.getω
— Methodgetω(axis)
Return the offsets (in rad/s) for points along a frequency axis.
NMRTools.NMRBase.gyromagneticratio
— Functiongyromagneticratio(n::Nucleus)
gyromagneticratio(c::Coherence)
Return the gyromagnetic ratio in Hz/T of a nucleus, or calculate the effective gyromagnetic ratio of a coherence. This is equal to the product of the individual gyromagnetic ratios with their coherence orders.
Returns nothing
if not defined.
Examples
julia> gyromagneticratio(H1)
2.6752218744e8
julia> gyromagneticratio(SQ(H1))
2.6752218744e8
julia> gyromagneticratio(MQ(((H1,1),(C13,1))))
3.3480498744e8
julia> gyromagneticratio(MQ(((H1,0),)))
0.0
NMRTools.NMRBase.hasnonfrequencydimension
— Methodhasnonfrequencydimension(spectrum)
Return true if the spectrum contains a non-frequency domain dimension.
Example
julia> y2=loadnmr("exampledata/2D_HN/test.ft2");
julia> hasnonfrequencydimension(y2)
false
julia> y3=loadnmr("exampledata/pseudo3D_HN_R2/ft/test%03d.ft2");
julia> hasnonfrequencydimension(y3)
true
NMRTools.NMRBase.label!
— Functionlabel!(nmrdata, labeltext)
label!(nmrdata, dim, labeltext)
label!(nmrdimension, labeltext)
Set the label associated with an NMRData
structure or an NMRDimension
.
See also label
.
NMRTools.NMRBase.lineshape
— Methodlineshape(axis, δ, R2, complexity=RealLineshape())
Return a simulated real- or complex-valued spectrum for a resonance with chemical shift δ
and relaxation rate R2
, using the parameters and window function associated with the specified axis.
NMRTools.NMRBase.missingval
— Functionmissingval(x)
Returns the value representing missing data in the dataset.
NMRTools.NMRBase.pcls
— Methodpcls(A, y)
Compute the phase-constrained least squares solution:
\[y = A x e^{i\phi}\]
following the algorithm of Bydder (2010) Lin Alg & Apps.
Returns the tuple (x, ϕ)
, containing the component amplitudes and global phase.
If passed a matrix Y
, the function will return a matrix of component amplitudes and a list of global phases corresponding to each row.
Arguments
A
: (m,n) complex matrix with component spectray
: (m,) complex vector containing the observed spectrum
NMRTools.NMRBase.replacedimension
— Methodreplacedimension(nmrdata, olddimnumber, newdim)
Return a new NMRData, in which the numbered axis is replaced by a new Dimension
.
NMRTools.NMRBase.scale
— Methodscale(d::AbstractNMRData)
Return a scaling factor for the data combining the number of scans, receiver gain, and, if specified, the sample concentration.
\[\mathrm{scale} = \mathrm{ns} \cdot \mathrm{rg} \cdot \mathrm{conc}\]
NMRTools.NMRBase.setgradientlist
— Methodsetgradientlist(A::NMRData, [dimnumber], relativegradientlist, Gmax=nothing)
Return a new NMRData with a gradient axis containing the passed values. If no maximum strength is specified, a default gradient strength of 0.55 T m⁻¹ will be set, but a warning raised for the user.
If a dimension number is specified, that dimension will be replaced. If not, the function will search for a unique non-frequency dimension, and replace that. If there are multiple non-frequency dimensions, the dimension number must be specified explicitly, and the function will throw an error.
NMRTools.NMRBase.setkinetictimes
— Methodsetkinetictimes(A::NMRData, [dimnumber], tvals, units="")
Return a new NMRData with a kinetic time axis containing the passed values (and optionally, units). If a dimension number is specified, that dimension will be replaced with a TkinDim
. If not, the function will search for a unique non-frequency dimension, and replace that. If there are multiple non-frequency dimensions, the dimension number must be specified explicitly, and the function will throw an error.
NMRTools.NMRBase.setrelaxtimes
— Methodsetrelaxtimes(A::NMRData, [dimnumber], tvals, units="")
Return a new NMRData with a relaxation time axis containing the passed values (and optionally, units). If a dimension number is specified, that dimension will be replaced with a TrelaxDim
. If not, the function will search for a unique non-frequency dimension, and replace that. If there are multiple non-frequency dimensions, the dimension number must be specified explicitly, and the function will throw an error.
NMRTools.NMRBase.spin
— Methodspin(n::Nucleus)
Return the spin quantum number of nucleus n
, or nothing
if not defined.
Examples
julia> spin(H1)
1//2
See also Coherence
.
NMRIO
NMRTools.NMRIO.FQList
— TypeFQList(values, unit::Symbol, relative::Bool)
Represents a frequency list. unit
can be :Hz
or :ppm
, and relative
indicates whether the frequency is given relative to SFO (true) or BF (false).
Raw values can be extracted using the data
function, or (better) as absolute chemical shifts (in ppm) or relative offsets (in Hz) using getppm
and getoffset
functions.
NMRTools.NMRIO.getformat
— Methodgetformat(filename)
Take an input filename and return either :ucsf, :nmrpipe, :pdata (bruker processed), or :unknown after checking whether the filename matches any known format.
NMRTools.NMRIO.getoffset
— Methodgetoffset(f::FQList, ax::FrequencyDimension)
Return frequency list values as offsets relative to the spectrometer frequency, in Hz.
See also: getppm
NMRTools.NMRIO.getppm
— Methodgetppm(f::FQList, ax::FrequencyDimension)
Return frequency list values in ppm (in absolute terms, i.e. relative to 0 ppm).
See also: getoffset
NMRTools.NMRIO.loadnmr
— Methodloadnmr(filename, experimentfolder=nothing)
Main function for loading NMR data. experimentfolder
contains the path to an experiment directory, for identification of metadata, if the filename is not directly within an experiment.
Returns an NMRData
structure, or throws an NMRToolsError
is there is a problem.
Examples
nmrPipe import:
loadnmr("exampledata/1D_1H/1/test.ft1");
loadnmr("exampledata/1D_19F/1/test.ft1");
loadnmr("exampledata/2D_HN/1/test.ft2");
loadnmr("exampledata/pseudo2D_XSTE/1/test.ft1");
loadnmr("exampledata/pseudo3D_HN_R2/1/ft/test%03d.ft2");
Bruker pdata import:
loadnmr("exampledata/1D_19F/1");
loadnmr("exampledata/1D_19F/1/");
loadnmr("exampledata/1D_19F/1/pdata/1");
loadnmr("exampledata/1D_19F/1/pdata/1/");
ucsf (Sparky) import:
loadnmr("exampledata/2D_HN/1/hmqc.ucsf");
loadnmr("exampledata/1D_19F/1/");
loadnmr("exampledata/1D_19F/1/pdata/1");
loadnmr("exampledata/1D_19F/1/pdata/1/");
NMRTools.NMRIO.loadnmrpipe1d
— Methodloadnmrpipe1d(filename, md, mdax)
Return an NMRData array containing spectrum and associated metadata.
NMRTools.NMRIO.loadnmrpipe2d
— Methodloadnmrpipe2d(filename, md, mdax)
Return NMRData containing spectrum and associated metadata.
NMRTools.NMRIO.loadnmrpipe3d
— Methodloadnmrpipe3d(filename, md, mdax)
Return NMRData containing spectrum and associated metadata.
NMRTools.NMRIO.loadpdata
— Functionloadpdata(filename, allcomponents=false)
Filename will be a reference to a pdata folder.
NMRTools.NMRIO.parsefqlist
— Methodparsefqlist(filename)
Return contents of the specified fqlist.
fqlists can have several different formats:
first line | reference | unit
| sfo | Hz
sfo hz | sfo | Hz sfo ppm | sfo | ppm bf hz | bf | Hz bf ppm | bf | ppm p | sfo | ppm P | bf | ppm
NMRTools.NMRIO.parsenmrpipeheader
— Methodparsenmrpipeheader(header)
Pass a 512 x 4 byte array containing an nmrPipe header file, and returns dictionaries of metadata. The nmrPipe header format is defined in fdatap.h
.
Return values
md
: dictionary of spectrum metadatamdax
: array of dictionaries containing axis metadata
Examples
md, mdax = parsenmrpipeheader(header)
NMRTools.NMRIO.parseucsfaxisheader
— Methodparse ucsf axis header into a dictionary
NMRTools.NMRIO.parseucsfheader
— Methodparse ucsf file header into a dictionary
NMRTools.NMRIO.parsevalist
— Methodreturn valist contents in dB
NMRTools.NMRIO.parsevdlist
— Methodreturn vdlist contents in seconds
NMRTools.NMRIO.parsevplist
— Methodreturn vplist contents in seconds
NMRTools.NMRIO.sumexpts
— Methodsumexpts(outputexpt, inputexpts...; weights=[])
Sum a collection of Bruker NMR experiments, with optional weighting factors. Data will be truncated to fit the shortest data file.
Arguments
outputexpt
: Output experiment name (string or integer)inputexpts...
: Input experiment names (strings or integers)weights
: Optional vector of weights for each input experiment
Details
- Experiment names should be either strings with filenames or integers converted to strings
- If weights are supplied, the length must match the number of input experiments
- First input experiment is used as a template for the output experiment
- For 1D experiments, 'fid' files are added; for nD experiments, 'ser' files are added
- Processed data files in pdata/1 are removed, other pdata subdirectories are deleted
- Updates the title file with information about the summed experiments
- Updates the number of scans in the acqus file to be the sum of individual experiments
- Prompts before overwriting any existing experiment
NMRTools.NMRIO.topspinversion
— Methodparsetopspinversion(acqusfilename)
Return the TopSpin version as a VersionNumber (e.g. v"4.2.0").
This is obtained from the end of the first line of the acqus file, e.g.
##TITLE= Parameter file, TopSpin 4.2.0
##TITLE= Parameter file, TOPSPIN Version 2.1
##TITLE= Parameter file, TOPSPIN Version 3.2