salt

Calculation of powder and single-crystal ENDOR spectra.

Syntax
salt(Sys,Exp)
salt(Sys,Exp,Opt)
spec = salt(...)
[rf,spec] = salt(...)
[rf,spec,info] = salt(...)
Description

This function calculates powder and single-crystal ENDOR spectra. Its calling syntax is identical to that of pepper, many of its options are equal to those of endorfrq, which is used by salt to compute line positions and amplitudes.

There are up to three output arguments

The three input arguments to the function are

Sys is a spin system structure. Fields available in Sys include all needed for the construction of a Hamiltonian and the ENDOR line width parameter lwEndor. If lwEndor is not given, it is assumed to be zero, and stick spectra are returned. The field HStrain is included in excitation window computations (see Exp.ExciteWidth).

For simulating a multi-component mixture, Sys should be a cell array of spin systems, e.g. {Sys1,Sys2} for a two-component mixture. Each of the component spin systems should have a field weight that specifies the weight of the corresponding component in the final spectrum.

Exp contains experimental parameters such as the microwave frequency, the magnetic field range and temperature. Here is a full list of its fields.

Field

Magnetic field (in mT) at which the experiment is performed.

Range

Two-element array with lower and upper limit of the rf frequency range [rfmin rfmax]. The values are in MHz. Example: Exp.Range = [1 30].

If omitted, EasySpin tries to determine the frequency range automatically.

CenterSweep

Two-element array with center and sweep width of the rf frequency range [rfcenter rfwidth]. The values are in MHz. Example: Exp.CenterSweep = [51 10].

If omitted, EasySpin tries to determine the frequency range automatically.

nPoints

Number of points on the rf frequency axis. If not given, EasySpin sets this to 1024.

mwFreq

EPR spectrometer frequency in GHz. Only needed if orientation selection should be taken into account, i.e. when Exp.ExciteWidth is given.

ExciteWidth

The excitation width of the microwave in MHz (responsible for orientation selection). The excitation profile is assumed to be Gaussian, and ExciteWidth is its FWHM (full width at half maximum). The default is Inf (infinity). To obtain the full excitation with for a given orientation, ExciteWidth is combined with HStrain from the spin system structure.

Temperature

Temperature at which the experiment is performed, in kelvin. If omitted or set to NaN, no temperature effects are computed.

SampleFrame

An Nx3 array that specifies the sample orientations for which the EPR spectrum should be computed. Each row of SampleFrame contains the three Euler rotation angles that transform the lab frame to the sample/crystal frame.

Exp.SampleFrame = [0 0 0];                   % sample/crystal frame aligned with lab frame
Exp.SampleFrame = [0 pi/2 0];                % sample/crystal frame tilted relative to lab frame
Exp.SampleFrame = [0 pi/2 pi/4];             % sample/crystal frame tilted relative to lab frame
Exp.SampleFrame = [0 0 0; 0 pi/2 pi/4];      % two samples/crystals

SampleFrame is only used for crystals and partially ordered samples (films etc). It is ignored for disordered samples (powders and frozen solutions).

CrystalSymmetry

Specifies the symmetry of the crystal. You can give either the number of the space group (between 1 and 230), the symbol of the space group, or the symbol for the point group of the space group (in either Schönflies or Hermann-Mauguin notation).

Exp.CrystalSymmetry = 'P21/c'; % space group symbol
Exp.CrystalSymmetry = 11;      % space group number (between 1 and 230)
Exp.CrystalSymmetry = 'C2h';   % point group, Schönflies notation
Exp.CrystalSymmetry = '2/m';   % point group, Hermann-Mauguin notation

When CrystalSymmetry is given, all symmetry-related sites in the crystal are included in the calculation. If CrystalSymmetry is not given, space group 1 (P1, point group C1, one site per unit cell) is assumed.

Ordering scalar (default: zero) or function handle

If a number is given in this field, it specifies the orientational distribution of the paramagnetic molecules in the sample. If not given or set to zero, the distribution is isotropic, i.e. all orientations occur with the same probability.

If a number is given, the orientational distribution is non-isotropic and computed according to the formula P(β) = exp(-U(β)) with U(β) = -λ(3 cos2β - 1)/2, where β is the angle between the sample z axis (z axis of the sample frame) and the molecular z axis ( axis of the molecular frame), and λ is the number specified in Exp.Ordering.

Typical values for λ are between about -10 and +10. For negative values, the orientational distribution function P(β) is maximal at β = 90° (preferential orientation of molecular z axis in the sample xy plane), for positive values it is maximal at β = 0° and β = 180° (preferential alignment of the molecular z axis with the sample z axis). The larger the magnitude of λ, the sharper the distributions.

To plot the orientational distribution for a given value of λ, use

lambda = 2;
beta = linspace(0,pi,1001);
U = -lambda*plegendre(2,0,cos(beta));
%U = -lambda*(3*cos(beta).^2-1)/2;  % equivalent
P = exp(-U);
plot(beta*180/pi,P);

If Exp.Ordering is a function handle, pepper calls the function to obtain the orientational distribution. The function must accept either one (beta) or two input arguments (beta and gamma), in radians. The function must accept vector arguments and return a vector P containing probabilities for each orientation, that is P(k) is the probability of finding the spin centers with orientation (relative to the sample frame) specified by beta(k) and gamma(k). Here is an example with an anonymous function:

Exp.Ordering = @(beta,gamma) gaussian(beta,0,15/180*pi);

If the function cannot be written as a singe line, define a separate function and provide it via its function handle. For example, here is a function oridist defining a von Mises-Fisher distribution

function P = oridist(beta,gamma)
gamma0 = 0;
beta0 = pi/4;
kappa = 10;
mu = ang2vec(gamma0,beta0);
x = ang2vec(gamma(:),beta(:));
P = kappa/(4*pi)/sinh(kappa)*exp(kappa*mu'*x);
P = reshape(P,size(beta));
end

Include this in the simulation via Exp.Ordering = @oridist.

When providing a function for Exp.Ordering, make sure that the symmetry used in the simulation is the same or lower than the symmetry of the distribution. Otherwise, incorrect spectra are obtained. In cases of doubt, set Opt.GridSymmetry='C1'. This always gives correct results.

If the orientational distribution provided via Exp.Ordering is very narrow, increase the number of knots in Opt.GridSize.

There are defaults for all fields except Range and Field, which have to be specified when invoking salt.

The structure Opt collects a number of parameters allowing to tune speed and accuracy of the simulation. Opt is optional, if it is omitted, pre-set values for the parameters are used. The field names and their possible values are

Method
'matrix' (default) or 'perturb1' or 'perturb2'

Specifies the algorithm used for the ENDOR simulation. By default, matrix diagonalization is used. This is an exact algorithm, but becomes slow when many nuclei are present. In such cases, large speed-ups at the cost of small losses in accuracy are possible using 'perturb1', first-order perturbation theory, or 'perturb2', second-order perturbation theory. Perturbation theory does not work for electron spins with S>1/2. It is accurate only if the hyperfine coupling is much smaller than the microwave frequency, and the nuclear quadrupole interaction is much smaller than both the hyperfine and the nuclear Zeeman interaction.

Opt.Method = 'perturb1';    % first-order perturbation theory
Opt.Method = 'matrix';      % matrix diagonalization
Nuclei

List of nuclei to include in the computation. Nuclei is a list of indices selecting those nuclei for which ENDOR peaks should be computed. 1 is the first nucleus, etc. E.g. the following specifies ENDOR of only the protons only in a copper complex.

Sys.Nucs = '63Cu,1H,1H';
Opt.Nuclei = [2 3];      % only protons
By default, all nuclei are included in the simulation.
Verbosity

Determines how much information salt prints to the screen. If Opt.Verbosity=0, salt is completely silent. 1 logs relevant information, 2 gives more details.

GridSize
[N1] or [N1 N2]

Determines the number of orientations (knots) in a powder simulation for which spectra are calculated.

Opt.GridSize = 91;       % 1° increments, no interpolation
Opt.GridSize = [46 0];   % 2° increments, no interpolation
Opt.GridSize = [31 6];   % 3° increments, 6-fold interpolation (giving 0.5° increments)
GridSymmetry
'' (default), 'Dinfh', 'D2h', 'C2h' or 'Ci'

Determines the grid symmetry used for the powder simulation. 'Dinfh' corresponds to a line from θ=0° to &theta=90° (with φ=0°), 'D2h' to one octant, 'C2h' to two octants, and 'Ci' to one hemisphere (four octants). If not given or '', then salt determines the appropriate grid symmetry automatically from the given spin system. With any other setting, pepper is forced into using the specified symmetry, even if it is incorrect for the spin system. See also hamsymm.

Sites

In crystal simulations, this gives a list of crystal sites to include in the simulation.

If Opt.Sites is empty or not given, all sites are included. If given, it must be a list of site numbers. The number of sites depends on the space group given in Exp.CrystalSymmetry. E.g. the following set limits the simulation to sites 1 and 3 of the 4 magnetically distinct sites in crystal of space group no. 47.

Exp.CrystalSymmetry = 47;  % space group Pmmm
Opt.Sites = [1 3];

In powder simulations, Opt.Sites is ignored.

separate '' (default), 'components', 'transitions', 'sites', 'orientations'

Determines whether to return the total spectrum or a list of subspectra. If set to '', the total spectrum is returned. If set to 'components', spec is a matrix with the subspectra of all components (including isotopologues). Each row in spec is one subspectrum. If 'transitions' is specified, transitions subspectra are returned (for powders only). If 'sites' is specified, site subspectra are returned (for crystals only). If 'orientations' is specified, orientation subspectra are returned (for crystals only).

The following options are only available for matrix diagonalization (Opt.Method='matrix'), but not for perturbation theory (Opt.Method='perturb1' or 'perturb2').

Transitions

Determine manually the level pairs (transitions) which are used in the spectrum calculation. If given, salt uses them and skips its automatic transition selection scheme. Level pairs are specified in Transitions(k,:) by the level numbers which start with 1 for the lowest-energy level.

Opt.Transitions = [1 2];      % transition between levels 1 and 2
Opt.Transitions = [1 2; 5 6]; % 2 transitions, 1<->2 and 5<->6
Enhancement
'off' (default) or 'on'

Switches hyperfine enhancement in the intensity computation on or off. See the same option in endorfrq.

Intensity
'on' (default) or 'off'

Switches all intensity computations on or off. Intensity computations include the quantum-mechanical transition probability, the orientation selectivity and the Boltzmann factor.

Threshold

Specifies the threshold for transition pre-selection. Only transitions with an estimated relative average intensity larger than this number are included. The relative average intensity of the strongest transition is set to 1. The default value for the threshold is 1e-4. The pre-selection is an approximate procedure, and it might miss transitions for complicated spin systems. In these cases, setting it to zero will include all transitions in the simulation. After the intensities of all included transitions are computed, the transition are screened again against this threshold in a post-selection step. If transitions are specified manually in Opt.Transitions, Opt.Threshold is ignored.

OriPreSelect
0 or 1 (default)

Specifies whether salt uses automatic orientational pre-selection to speed-up simulations. This speed-up is most noticeable for large spin systems and field/frequency settings that lead to single-crystal like spectra.

Algorithms

salt computes line positions and intensities for a set of orientations using either matrix diagonalization or perturbation theory.

The matrix diagonalization approach used in salt is identical to that of pepper, with the obvious exception of the calculation of line intensities, which is similar to that used in MAGRES (Keijzers et al, J.Chem.Soc. Faraday Trans. 1 83, 3493-3503, 1984).

First- and second-order perturbation theory is based on expressions by Iwasaki (J.Magn.Reson. 16, 417-423, 1974). It includes pseudosecular contributions. No transition moments are computed, that is, the intensities of all resonances are equal. Currently, the perturbation-theory algorithm is limited to systems with one electron spin S=1/2 (but an arbitrary number of nuclei with arbitrary spins).

For powder simulations, salt uses the same methods as pepper, orientational interpolation and interpolative projection, to construct the powder spectrum.

Examples

A full simulation of the powder ENDOR spectrum of a radical with a proton is

clear
Sys.g = 2;
Sys.Nucs = '1H';
Sys.A = [-2 1 4];
Sys.lwEndor = 0.1;
Exp.Range = [8 18];
Exp.Field = 308.46;
salt(Sys,Exp);
See also

resfields_eig, endorfrq, pepper, resfields, resfields_perturb