Fitting EPR spectra

This user guide explains how to extract magnetic parameters from experimental EPR spectra by fitting an EasySpin simulation result to the experimental spectra using least-squares fitting techniques. EasySpin's function esfit contains all the necessary machinery and is easy to use.

This tutorial contains the following topics: There are the following advanced topics:

To get the most out of the EasySpin fitting features, you should work through this entire tutorial and study the associated examples.

The fitting function

EasySpin's fitting function is esfit and can be called with up to seven parameters:

esfit('pepper',spc,Sys0,Vary,Exp)
esfit('pepper',spc,Sys0,Vary,Exp,SimOpt)
esfit('pepper',spc,Sys0,Vary,Exp,SimOpt,FitOpt)
Here is what the various parameters mean: The following two parameters are optional The next few sections explain how to set up these input parameters so that esfit can be called.
Loading the experimental spectrum

Before you can fit parameters to your experimental spectrum, you have to import the experimental data into MATLAB. There are several ways to do this, depending in which format your data are stored in.

If your spectrum was acquired by a Bruker spectrometer, it is most likely either in ESP format (having file extensions .par and .spc) or in BES3T format (with file extensions .DSC and .DTA). Both formats, as well as formats from other spectrometer brands, are handled by EasySpin's function eprload. Here is an example: Suppose your data are stored in the files dpph.DSC and dpph.DTA. To load them into MATLAB, type

[B,spc,Params] = eprload('dpph.DTA');

This returns the magnetic-field values in B, the spectral data in spc, and all the parameters in the structure Params.

Often, experimental spectral data come in text files, containing two columns of data, where the first column is the magnetic field and the second column contains the spectral intensity. To load such a file into MATLAB, use the function textread:

[B,spc] = textread('dpph.txt','%f %f');

textread is a very flexible function that can accomodate many different text file formats. E.g. it can skip header lines if told to do so, it can accomodate comment lines, etc. See the MATLAB documentation for more details.

Start set of parameters

The fitting function needs to know with which parameter values to start looking for an optimal fit. These starting parameters are given as a spin system in the third input parameter to esfit. For example, if the spectrum is that of a simple S=1/2 with rhombic g and an isotropic linewidth, one would start with

Sys0.g = [1.99, 2.07, 2.11];
Sys0.lw = 1;   % mT

Some algorithms use this set of parameters as the starting set, while others use it to define the center of the search range.

Parameter search range

Next, esfit has to be told which parameters of the spin system given in Sys0 should be fitted (and which not), and by how much they can be varied during the least-squares fitting process.

This information is given in the fourth input parameter to esfit. If only the linewidth should be fitted, and it should be varied by at most +/-0.3 mT, use

Vary.lw = 0.3;  % mT

With Sys.lw = 1, this would mean that the search range extends from 0.7 to 1.3. esfit makes every effort to stay within the specified search range.

If the second principal value of the g tensor and the linewidth should be fitted simultaneously, use

Vary.g = [0, 0.02, 0];
Vary.lw = 0.3;

In essence, all the fields in Vary must have the same names as those in Sys0, and any non-zero value indicates that that parameter should be fitted by varying it by at most plus or minus the given amount. Setting the value for any parameter to zero excludes it from the fitting.

It is advisable not to vary more than about 4 parameters at the same time, as the efficiency of essentially all fitting algorithms decreases tremendously with the number of fitting parameters.

An important note for fitting slow-motion EPR spectra with chili: Do not use tcorr or Diff for fitting the correlation time or the diffusion tensor, but rather the logarithmic forms logtcorr and logDiff.

Experimental parameters

Of course, the fitting function must also be provided with the experimental parameters for the spectral simulations. These are given as the fifth parameter, a structure Exp, which is directly forwarded to the simulation functions pepper, garlic, chili, etc. For details, see the documentation of these functions.

A minimal example would be to give the microwave frequency and magnetic field range:

Exp.mwFreq = 9.5;  % GHz
Exp.Range = [200 400];  % mT

The microwave frequency, the field range, and all other experimental parameters must correspond to the ones used in acquiring the spectral data.

Performing the least-squares fitting

Now we are all set to call the fitting function:

esfit('pepper',spc,Sys0,Vary,Exp);

The least-squares fitting is started, and a figure window pops up that contains an interactive graphical user interface (GUI), shown below. The GUI allows you to control which parameters should be varied during a fit, lets you pick another fitting algorithm, target, scaling method, and starting point parameter set. In addition, multiple fitting results from separate runs can be stored, compared and used as starting points for other fits. The fitting results can be exported to the workspace, so that they are available in the command window for later use.

This GUI is what you probably want, but if you prefer to run esfit on the command line only, then call it requesting outputs:

bestSys = esfit('pepper',spc,Sys0,Vary,Exp);

This returns the result of the fitting, the optimized spin system or list of spin systems, in bestSys.

If you want to pass settings to the simulation function, collect them in an additional stucture SimOpt and pass them as sixth parameter:

esfit('pepper',spc,Sys0,Vary,Exp,SimOpt);

It is possible to provide a structure with settings for the fitting function. If we call this structure FitOpt, the syntax is

esfit('pepper',spc,Sys0,Vary,Exp,SimOpt,FitOpt);

The possible settings in this last structure are the topic of the rest of this tutorial.

Fitting methods

Beyond a good starting parameter set or search range, the performance of the fitting depends crucially on two things: the choice of the optimization algorithm, and the choice of the target function. Let's have a look at each of them in turn.

Optimization algorithms

EasySpin provides several optimization algorithms that are in widespread use: (1) the Nelder/Mead downhill simplex method, (2) the Levenberg/Marquardt algorithm, (3) Monte Carlo random search, (4) a genetic algorithm, (5) a systematic grid search, as well as others.

The first two are local search algorithms, which start from a given starting set of parameter values and try to work their way down a nearby valley of the parameter space to find the minimum. Both methods are quite fast, although there are some differences in general performance between them: The downhill simplex is somewhat slower than Levenberg/Marquardt, but it is more robust in the sense that it does not get stuck in a local minimum as easily as Levenberg/Marquardt.

The latter three are global search methods: they do not have a single starting parameter set, but use many, distributed over the entire parameter search space. The Monte Carlo method simply performs a series of random trial simulations and picks the best one. It is very inefficient. The systematic grid search is better: It covers the parameter space with a grid and then does simulations for each knot of the grid, in random order. Thus, no point is simulated twice, and the method is more efficient than the Monte Carlo search. However, if the minimum is between two grid points, it will never be found.

The third global method is a genetic algorithm: It makes simulations for several, let's say N, parameter sets (called a population), computes the fitting error (called the fitness) for all of them and then proceeds to generate N new parameter sets from the old ones using mechansims like mutation, cross-over and reproduction. This way, a new generation of parameter sets is (pro-)created, just like in biological evolution. The benefit of this algorithm is that if a good parameter is encountered, it is likely to propagate down the generations and across the population.

To select one of the algorithms, specify it in the Method field of the fitting options

FitOpt.Method = 'simplex';     % for Nelder/Mead downhill simplex
FitOpt.Method = 'levmar';      % for Levenberg/Marquardt
FitOpt.Method = 'montecarlo';  % for Monte Carlo
FitOpt.Method = 'genetic';     % for the genetic algorithm
FitOpt.Method = 'grid';        % for grid search

and then supply this option structure as the seventh parameter to the fitting function, for example

esfit('pepper',spc,Sys0,Vary,Exp,[],FitOpt);

If you don't specify the algorithm, EasySpin uses the downhill simplex by default.

Each algorithm has some internal parameters that can be used to fine-tune its performance. In general, it is not necessary to fiddle with those parameters. For more details, see the documentation of esfit.

Target function

The second important setting is the choice of the target function. esfit computes the error of the simulated spectrum using the root-mean-square-deviation (rmsd, i.e. the square root of the mean of the square of the deviations), where the deviation is the difference between the experimental and the simulated spectrum.

Fitting speed can often be signifantly increased, however, if one used not the spectra directly, but their integrals or similar transforms. EasySpin supports several settings here: 'fcn' (use data as is), 'int' (integral), 'dint' (double integral), 'diff' (first derivative), 'fft' (Fourier transform). The target function setting string is simply appended to the 'Method' field, after a space:

FitOpt.Method = 'simplex int';   % simplex algorithm with integrals for rmsd
FitOpt.Method = 'genetic fcn';   % genetic algorithm with spectra as is
FitOpt.Method = 'levmar dint';   % Levenberg/Marquardt with double integral

Usually, 'fcn' is an excellent choice, but in the case of many lines 'int' can be better - provided the baseline in the experimental spectrum is good. The other settings ('dint', 'diff', and 'fft') have advantages in some situations.

Normally, EasySpin just scales the simulated spectrum in a way that its maximum absolute value coincides with the maximum absolute value of the experimental spectrum. However, if the experimental spectrum is very noisy, this is not ideal, and the best scaling of simulated spectrum is better determined by least-squares fitting. The scaling method can be set in the field FitOpt.Scaling. The most common settings are

FitOpt.Scaling = 'maxabs';    % scaling so that the maximum absolute values coincide
FitOpt.Scaling = 'lsq';       % determine scaling via least-squares fitting

Example

A simple mock example of this would be

% simulate a noisy spectrum as input data for fitting
Sys.g = 2; Sys.lw = 2;
Exp.mwFreq = 9.5; Exp.Range = [335 345];
[x,y] = pepper(Sys,Exp);
y = 3*addnoise(y,2);

% define parameters and settings for fitting the noisy data
Sys.g = 2.001; Sys.lw = 1.8;
Vary.g = 0.002; Vary.lw = 0.4;
FitOpt.Scaling = 'maxabs';

% start fitting
esfit('pepper',y,Sys,Vary,Exp,[],FitOpt);

Change FitOpt.Scaling to see how the two scaling methods differ. Baseline correction can be directly incorporated in the scaling method - see the documentation of esfit for more details.

For some situations, e.g. for magnetization data, it is necessary not to scale at all. In this case, use FitOpt.Scaling = 'none'.

Hybrid methods

With EasySpin it is easily possible to perform so-called hybrid least-squares fittings, where one optimization algorithm is used to locate a potential minimum, and another one is used to refine the parameters at that minimum. The first of these two steps often employs algorithms that are able to locate minima globally and do not get stuck easily in a local minimum. The disadvantage of these methods is that they are often slow. The second step closes in on the minimum by using a much faster local method. There are two ways to perform such a two-stage fitting: using the UI, or writing your own code.

If you use the UI, you have complete control over when a fitting algorithm terminates and which one you want to start next. Any sequence of fitting steps where you use the result of a previous fit as starting values for the next constitutes a 'hybrid method'. But of course, the UI lets you do much more complex operations.

Alternatively, you can write code that does two- or multistage fitting. Let's look at an example with a two-stage fitting using genetic algorithm followed by Levenberg/Marquardt. This can be set up by calling esfit twice with different settings in FitOpt.

% first stage: genetic algorithm
FitOpt.Method = 'genetic fcn';
Sys1 = esfit('pepper',y,Sys0,Vary,Exp,[],FitOpt);

% second stage: Levenberg/Marquardt
FitOpt.Method = 'levmar int';
Sys2 = esfit('pepper',y,Sys1,Vary,Exp,[],FitOpt)

Of course, you'll probably have to change some of the termination criteria for the two algorithms so that the genetic search narrows down a minimum only coarsely, and the local search can then much more efficiently finalize the fit by refining the parameters.

Termination criteria

The simplest way to terminate a fitting is by pressing the Stop button on the UI. This will interrupt the optimization process whether the result has converged or not.

Without pressing the Stop button, esfit stops the fitting when it thinks it has reached a minimum in the error function, when it has taken more than a given amount of time, or if the set number of simulations are reached. Let's have a look at these possibilities in turn.

esfit considers a local least-squares fit using the simplex or the Levenberg/Marquardt algorithm to be converged if the change in parameters from one simulation iteration to the next falls below a certain threshold and/or if the change in the error from iteration to iteration falls below another threshold. Both thresholds have pre-set values, but can be adjusted by supplying appropriate fields in the FitOpt structure:

FitOpt.TolFun = 1e-3;          % termination tolerance for error change
FitOpt.TolStep = 1e-3;         % termination tolerance for parameter step, Levenberg/Marquardt
FitOpt.TolEdgeLength = 1e-3;   % termination tolerance for parameter step, simplex

The global methods terminate also if the maximum number of simulations are reached: the Monte Carlo search does only a pre-set number of simulations (FitOpt.N), the grid search stops if all the grid points are simulated (see the option FitOpt.GridSize), and the genetic algorithm stops at the latest after a certain number of generations (see FitOpt.PopulationSize and FitOpt.Generations).

In a field FitOpt.maxTime, the fitting function can be told to terminate after a given amount of time, even if the fitting did not converge in terms of TolStep and TolFun. This can be useful when running several fittings overnight from a script.

FitOpt.maxTime = 60*8;     % maximum time, in minutes
Multiple components

So far we have looked at the fitting of a spectrum with a single spectral component. EasySpin can perform least-squares fitting of spectra with multiple components. For each component, a system structure and a variation structure must be given. For example, this is how esfit is called for a two-component fit:

esfit('pepper',spc,{Sys1,Sys2},{Vary1,Vary2},Exp,SimOpt,FitOpt);
Each spin system (Sys1 and Sys2) contains the magnetic parameters, and the corresponding variation structure (Vary1 for Sys1, and Vary2 for Sys2) contains the parameters to be varied.

In addition, the weight of each spin system has to be given in the field weight:

Sys1.weight = 0.8;
Sys2.weight = 0.4;
These weights are used as prefactors when constructing the total spectrum as a weighted sum of the component spectra. The weights need not add up to 1. They can even be included in the fit, for example:
Vary1.weight = 0.3;
User-defined simulation functions

Say you need to constrain the spin Hamiltonian values, e.g. two hyperfine tensors have to be equal. Or you want to fit a distance from which a hyperfine coupling is derived, instead of the hyperfine coupling directly. Or you want to fit the microwave phase, which is not defined in the spin system structure, but in the experimental parameters.

For these and similar cases, you can define your own custom simulation function that implements the constraints or calculations, and then use esfit with your custom function. The requirements on your function are the following:

Example: Constraint between systems

Here is an example of a function mymy.m that simulates a spectrum constraining two hyperfine tensors to be identical:

function y = mymy(Sys,Exp,Opt)
fullSys = Sys;
fullSys.A = [Sys.A; Sys.A];
fullSys.Nucs = '1H,1H';
[x,y] = pepper(fullSys,Exp,Opt)

This function takes the A from the input Sys structure, constructs a new spin system structure fullSys with two identical A's, and then performs the actual simulation.

With this function, esfit should be called in the following way (assuming Exp has been defined):

Sys.A = [2 3 5];
Sys.lwpp = 0.1;
Vary.A = [1 1 3];
esfit('mymy',data,Sys,Vary,Exp);

Example: Varying non-spin system parameters

Here is a second example. If you want to fit the microwave phase, use this custom function:

function y = mymy(Sys,Exp,Opt);
Exp.mwPhase = Sys.mwPhase;
[x,y] = pepper(Sys,Exp,Opt);

and call it as follows:

Sys.mwPhase = 0;
Vary.mwPhase = 20*pi/180;    % 20 degrees
esfit('mymy',data,Sys,Vary,Exp);

As you can see, you have to put the parameter you want to fit into the first input structure (called Sys), even if it is not a spin system parameter. This won't confuse EasySpin's simulation functions, since they disregard any field they do not use. In other words, as long as the new field name (in this case Sys.mwPhase), does not clash with any existing field names that EasySpin simulation functions use in the spin system (e.g. Sys.Nucs, Sys.A, etc.), a user-defined simulation function should work and allows for a great deal of customization. To ensure that there is no field name clashing, check the documentation for the spin system and the relevant simulation function.

Example: globally fitting multiple spectra

During the fitting process, esfit compares a vector containing the input data with a vector of simulated data and determines the error. Normally, each of these two vectors contain data from only a single spectrum (or a spectrum containing multiple components). However, by defining a custom function that combines spectra along the same dimension of the vector, i.e. concatenates the spectra, one can fit multiple spectra simultaneously. For example, to fit spectra that were acquired at X- and Q-band microwave frequencies with the same sample:

function y = globalfit(Sys,Exp,Opt);

% X-band
Exp.mwFreq = 9.5;
y1 = pepper(Sys,Exp,Opt);

% Q-band
Exp.mwFreq = 34;
y2 = pepper(Sys,Exp,Opt);

scale = Sys.scale;  % custom field for relative scaling

% Concatenate the spectra
y = [y1(:); scale*y2(:)];
Note that each spectrum might have a different number of data points, in which case two more custom fields would be required to specify the number of data points in each (e.g. Exp.nPoints = Exp.nX and Exp.nPoints = Exp.nQ for the X- and Q-band spectra, respectively).

esfit calls the custom simulation function once per iteration. If Sys is a cell array with multiple spin systems, such as Sys = {Sys1,Sys2}, then this will be passed as a cell array to the custom function.

esfit is able to handle simulation functions that do not provide the simulated data in their first output argument. If you want to use a non-conforming function, just tell esfit how it should be called and which output argument should be used in the fitting. This information is provided in FitOpt.OutArg, which should contain two numbers. The first number is the number of output arguments to expect, and the second number is the output argument to use in the fitting. Say, you have a function myfun that needs to be called with 4 output arguments ([a,b,c,d] = myfun(...), and you want to use the second one for fitting, then use

FitOpt.OutArg = [4 2];  % 4 outputs, use 2nd for fitting