XAFS Analysis in IFEFFIT consists of modeling XAFS
(k) in terms of a
sum of paths with the XAFS contribution from each path being
modeled by a calculation from FEFF. A Path is an abstract formalism for
breaking the EXAFS into manageable pieces based on the scattering path that
a photo-electron takes. Many EXAFS analysis methodologies use concepts
such as 'shell' or 'coordination sphere'. To some extent the distinction
is only conceptual, and therefore unimportant. In other ways, the path
formalism is clearly a superior way to robustly analyze complex EXAFS data
that includes the effects of multiple scattering. If you're used
to thinking of your EXAFS data as having contributions based on shells or
coordination spheres, paths aren't that big of a change.
IFEFFIT relies on FEFF for basic information about all of its paths. This means that you need to run FEFF and sort through its results (neither of which are trivial tasks) prior to defining paths with IFEFFIT. The outputs of FEFF that IFEFFIT needs is the set of feffnnnn.dat files, where each file represents the result for a single path. Some versions of FEFF can write all of the path information into a single file called feff.bin as an alternative to a large set of feffnnnn.dat files. IFEFFIT can read some versions of feff.bin, but this discussion will focus on using the more universal feffnnnn.dat files.
Once FEFF has run and a set of feffnnnn.dat files exists, using them within IFEFFIT is fairly easy. The path() command is used to define paths within IFEFFIT. The information needed to describe an XAFS path is a bit too complicated to conveniently hold as a set of IFEFFIT scalars and arrays - it could be done, but the bookkeeping would be a nightmare. Instead, paths are defined and stored internally, so that you can refer to a path by an integer index. Though a break from the idea in section 2 that all data is stored in program variables that you can access easily and directly, there are commands to convert the path data to program variables.
The path() command takes the name of a feffnnnn.dat file and a unique
IFEFFIT path index as its primary arguments. The IFEFFIT
path index is an integer (from 1 to 10000). It does not need to be the
same index as FEFF used, but that is often a convenient way to do it.
In addition, the path() command takes several optional arguments,
known as Path Parameters that represent the parameters used to
modify the EXAFS function for that path. The path parameters include such
things as
R, the change in path distance (well, half-path-length
for multiple scattering paths),
, the mean-square displacement of
the bond, and an E0 shift. The most important Path Parameters are
| Path Parameter keyword | Description |
| file | Name of feffnnnn.dat file |
| index | path index: a unique integer identification |
| label | path label: a string describing the path |
| s02 | S02, amplitude reduction factor |
| e0 | E0, energy shift |
| delr | |
| sigma2 | |
| third | C3, third cumulant |
a full list is given in the Reference Guide. A simple path definition would look like this:
Ifeffit> path(index=1,file = feff0001.dat,
s02 = 0.9, sigma2 = 0.001)
This would define path 1, and assign values for S02 and
Ifeffit> s02_001 = 0.9
Ifeffit> e0_001 = -1.
Ifeffit> sig2_001 = 0.001
Ifeffit> path(index=1,file = feff0001.dat,
s02 = s02_001,
e0 = e0_001,
sigma2 = sig2_001 )
where the parameter s02 for path 1 now depends on the definition of
the scalar s02_001. This is a fairly simple example, and the
definition for the s02 parameter for path 1 could be more
complex - say, min(1.2,s02_001) to set an upper bound on this
parameter. In addition, the s02_001 scalar here could have a more
complex definition or even be a fitting variable defined with guess
s02_001 = 0.9.
Additional paths are defined in the same way:
Ifeffit> path(index=2,file = feff0002.dat,
s02 = s02_001, e0 = e0_001,
sigma2 = sig2_001 )
The only required keyword for a path is the integer index (if the
keyword comes first in the list of arguments, the index can be
dropped, so a completely equivalent definition is
Ifeffit> path(2, file = feff0002.dat) Ifeffit> path(2, s02 = s02_001 ) Ifeffit> path(2, e0 = 2.0 ) Ifeffit> path(2, sigma2 = sig2_001 )subsequent path() commands will set or overwrite the definitions for other path parameters. This makes it fairly simple to change the definition of a path parameter half way through an analysis, say between fits, where something like
Ifeffit> path(2, s02 = 0.85)will change one path parameter and leave the rest as previously defined.
Once defined, combining paths to give a
(k) function is very easy with
the command ff2chi(). ff2chi() takes its name from the FEFF
module that combine FEFF path files into a single
(k) function.
Because IFEFFIT allows you to alter the paths through the Path
Parameters, and use paths calculated from different runs of FEFF, this
version of ff2chi() is significantly more flexible and powerful than
FEFF's own version.
ff2chi() takes a list of path indices, and creates arrays for
(k). Assuming that the paths have already been defined, a command
like this:
Ifeffit> ff2chi(1,2,3, group = ff)will add paths 1, 2, and 3 together and create arrays ff.k and ff.chi containing the sum of these paths. Each of the paths will be altered according to its own Path Parameters in the sum. The list syntax '1, 2, 3' could also be replaced by '1-3'. Lists of the form 1-3,5,7-10,23,11 are also allowed. An important note is that (currently), paths listed twice are used twice. This is an easy mistake to make, so beware4.
Besides the list of paths and group name for the output arrays,
ff2chi() has a few other arguments, most of which you won't really
need. I'll just mention a few convenient ones here. The arguments
kmin and kmax can set the k-range of the output arrays, and
the argument do_real will cause the ``real part'' of
(k) to be
written to the array $group.chi_real - the confusion over real
and imaginary parts of the complex
(k) shows up again! The
$group.chi array that ff2chi() always generates is the
imaginary part of the complex
(k), and should correspond to the
experimentally derived
(k).
Once you've combined a few paths with ff2chi() or, as we'll see shortly, done a fit with feffit(), you may want to find out some information about the path parameters for a set of paths. There are a few ways to get this information. The show() command, first discussed in section 2.5, has a ``@path'' modifier to give some information about a defined path.
Ifeffit> show @path=1
PATH 1
feff = feffcu01.dat
id = Cu metal first neighbor
reff = 2.547800, degen = 12.000000
s02 = 0.937373, e0 = 0.510790
dr = 0.000525, ss2 = 0.003496
3rd = 0.000000, 4th = 0.000000
ei = 0.000000, dphase = 0.000000
This shows all the scalar values of the path parameters (reff is
the half-path length, degen is the path degeneracy, and so on).
The argument to the @path modifier to show() can be any
valid path list as described in section 6.2
(for example: 1,2,4-9,12). To see all the defined paths
show @paths will work.
The show @path family of commands only prints the values of the path parameters. Another way to get path information is to convert them to program variables with the get_path() command. This command only takes one path at a time (not a path list), but allows you to get the path parameter values into program variables for later manipulation. The syntax for get_path() is
Ifeffit> get_path(path = 1, prefix = path1)or, equivalently and more simply
Ifeffit> get_path(1, path1)This will create scalars path1_s02, path1_e0, path1_ei, path1_delr, path1_sigma2, path1_third, path1_fourth, path1_degen, path1_reff and give the values of the path parameters for path 1. In addition, the strings $path1_file and $path1_id will contain the values of the path FEFF file name and path identification string, respectively. The prefix argument to get_path() can be any valid variable name. Of course, array variables for a given path can be formed with the ff2chi() command, as described in the previous section. Just to be clear, though
Ifeffit> ff2chi(1, group = path1)will create arrays of path1.k and path1.chi. Other arrays for the real part, magnitude, and phase of
At some point in the analysis of XAFS data, you'll want to refine a set of path parameters so that a sum of paths best matches the data. The feffit() command is the basic tool for fitting EXAFS data to a set of paths. In some sense, feffit() is just a fancy version of ff2chi(). Another view might be that the feffit() command is the main point of the IFEFFIT library, so that ff2chi() is a very simple version of feffit(). In any event, the two functions are closely related and purposefully implemented with as much overlap as possible so that once you feel comfortable with the paths() and ff2chi() commands, feffit() should be a fairly small step.
feffit() has four basic requirements:
First, feffit() needs some
(k) data to fit. As you can probably
guess, this data needs to be a named array for
(k) just as you would
use for the Fourier transform routines (ie, starting at k = 0, and
progressing in steps of
0.05 Å-1). You specify the
(k)
array with feffit(chi= data.chi,...). If the data you have is
not on the expected grid, you'll either need to interpolate it on to the
expected grid or specify the array of k values, with
Ifeffit> feffit(chi= data.chi,k=data.k, ...)Second, feffit() needs a list of paths for the sum-over-paths that makes up the model XAFS. This sum is essentially the same as for ff2chi(), comprising of a simple list of path indices. Path lists of the '1, 2, 3', '1-3', and '1-3,5,7-10,23,51' are accepted, so that the feffit() command would now look like (assuming the data is on the proper k grid)
Ifeffit> feffit(chi= data.chi,1-3, ...)Third, feffit() needs Fourier transform and fit range parameters defined. This reflects the fact that XAFS is a band-limited phenomenon, with a finite k- and R-range. As mentioned in section 5.4.1, the Fourier transform parameters can be read either from the appropriately named program variables (kmin, kmax, kweight, dk1, dk2, and so forth) or from the command argument list (using these same identifiers as keywords). In addition to the k-space Fourier transform parameters, feffit() needs the fit R-range, specified by rmin and rmax (this reflects the fact that feffit() normally does the fit in R-space - fitting in k-space is possible, as well). Like the Fourier transform parameters, the R-range can be set from program variables (rmin and rmax or in the feffit() command statement itself. We could say either
Ifeffit> kmin = 3.0, kmax = 15.0, kweight = 2, dk1 = 2, dk2= 2 Ifeffit> rmin = 1.5, rmax = 3.5 Ifeffit> feffit(chi= data.chi, 1-3 )or
Ifeffit> feffit(chi= data.chi, 1-3, rmin=1.5, rmax=3.5,
kmin=3,kmax=15,kweight=2,dk=2)
The fourth requirement for a fit may seem obvious: feffit() needs some variables to vary. Fitting Variables are like regular scalars, except that they are defined with the guess() command instead of set(), as in guess e0 = 1. . Fitting variables also carry around information about their estimated uncertainty, and correlation with other fitting variables. We'll discuss that in the next section.
Putting it all together, a simple fit of a single FEFF path to
(k)
data would look like this:
Ifeffit> read_data(file = cu_chi.dat, group= data, type = chi)
Ifeffit> set s02_001 = 0.9
Ifeffit> guess e0_001 = 1.
Ifeffit> guess dr_001 = 0.01
Ifeffit> guess sig2_001 = 0.001
Ifeffit> path(index = 1,
file = feff0001.dat,
s02 = s02_001,
e0 = e0_001,
delr = dr_001,
sigma2 = sig2_001 )
Ifeffit> kmin = 3.0, kmax = 16.0, kweight = 2, dk1 = 2
Ifeffit> rmin = 1.5, rmax = 3.2
Ifeffit> feffit(data.chi, 1, group=fit)
The feffit() command will update the values of all the fitting
variables, and also create arrays for the best-fit Ifeffit> fftf(data.chi) Ifeffit> fftf(fit.chi)for you. The best-fit values and estimated uncertainties for the variables can be seen with show @variables, as will be further discussed in the next section.
In addition to what has been discussed so far, feffit() has several optional parameters for a wide range of things such as specifying uncertainty in the input data, whether to refine background-spline parameters with the structural refinement, how to fit in k-space (with or without Fourier filtering), and an IFEFFIT macro to run for each fit iteration. In the interest of brevity, I'll resist discussing any of these fascinating topics in this brief tutorial.
A fit is not very useful without some idea of the reliability of the fit and some estimate of the uncertainty in the fitted variables. The feffit() command automatically calculates the estimated uncertainties and goodness-of-fit parameters for you. The particulars of how these values are determined is outside the scope of this tutorial - the attention here will be on how to view and manipulate the values for these statistics that feffit() determines.
The feffit() command writes several fitting statistics to help you
determine the goodness of fit. First, IFEFFIT will write the number of
fit iterations to &fit_iteration. Though not necessarily useful
as a fitting statistic, this can be valuable to see if a fit is taking far
too many iterations, probably indicating a problem with the set up or
model. The number of variables is stored in n_varys. The number
of independent points in the band-limited data is estimated as
2
k
R/
, for fit k-range
k and R-range
R and stored in n_idp. The estimated uncertainty in the
data itself is stored in the variables epsilon_r for the
uncertainty in the data
(R) and epsilon_k for the uncertainty
in the data
(k).
The main goodness-of-fit statistics are r_factor, chi_square, and chi_reduced, which are three different ways of scaling the sum-of-squares of the final misfit. The full definition of these is given elsewhere, but briefly, r_factor is the misfit scaled to the data itself, so it a relative misfit, chi_square is scaled to the estimated uncertainty in the data (epsilon_r for data fit in R-space). chi_reduced is chi_square divided by the number of ``free parameters in the fit'', given by the difference of n_idp and n_varys. More information on these terms and concepts is elsewhere in the IFEFFIT documentation.
For each fitting variable, IFEFFIT will create (or overwrite) a variable named delta_VAR, and write the estimated uncertainty for variable VAR to this variable. For example, the fit shown above would put the uncertainties in the fitting variables in delta_s02_001, delta_e0_001, delta_sig2_001, and delta_dr_001.
The estimated uncertainties are, of course, very valuable for assessing fit results. In addition, it often instructive to know the correlations between pairs of variables resulting from the fit. To extract the correlations from IFEFFIT, you use the correl() command, which simply takes the names of two fitting variables, and the name of an output scalar variable to store the resulting correlation to. Thus,
Ifeffit> correl(s02_001, sig2_001, out= cor1) Ifeffit> print cor1 0.887130839At this point, it is an error to give correl() the name of any variable that is not a fitting variable.