/
eels.py
638 lines (551 loc) · 25.4 KB
/
eels.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
import numpy as np
from numbers import Number
from scipy import constants
# hyperspy dependency
import hyperspy.api as hs
from hyperspy.external.progressbar import progressbar
from hyperspy.defaults_parser import preferences
# own dependency
from signals import SignalMixin
class ModifiedEELS(hs.signals.EELSSpectrum, SignalMixin):
def __init__(self, *args, **kwargs):
"""
Modified hypespy EELS signal. Input can be array or EELSSpectrum type.
In the last case, additional *args and **kwargs are discarded.
"""
if len(args) > 0 and isinstance(args[0], hs.signals.EELSSpectrum):
# Pretend it is a hs signal, copy axes and metadata
sdict = args[0]._to_dictionary()
hs.signals.EELSSpectrum.__init__(self, **sdict)
else:
hs.signals.EELSSpectrum.__init__(self, *args, **kwargs)
def remove_negative_intensity(self, inplace=False):
'''
By definition, electron energy-loss spectral intensity is positive but,
sometimes, our spectral treatments overlook this important fact.
Parameters
----------
inplace : bool
Performs the operation in-place.
Returns
-------
spc : EELSSpectrum
Negative intensity is set to 0.
'''
if inplace:
spc = self
else:
spc = self.deepcopy()
spc.data[spc.data < 0.] = 0
return spc
def remove_intensity(self, left_value=None, right_value=None, inplace=False):
'''
Remove all intensity in a range.
Parameters
----------
threshold : {int, float}
Integer index or energy-loss units.
inplace : bool
Performs the operation in-place.
Returns
-------
spc : EELSSpectrum
Spectral intensity in this range is set to 0.
'''
self._check_signal_dimension_equals_one()
try:
signal_range_from_roi = hs.hyperspy.misc.utils.signal_range_from_roi
left_value, right_value = signal_range_from_roi(left_value)
except TypeError:
# It was not a ROI, we carry on
pass
eax = self.axes_manager.signal_axes[0]
i1, i2 = eax._get_index(left_value), eax._get_index(right_value)
if inplace:
spc = self
else:
spc = self.deepcopy()
spc.data[(slice(None),)*eax.index_in_array+(slice(i1, i2),Ellipsis)]=0.
return spc
def get_effective_collection_angle(self):
"""
Calculates the effective collection angle for the whole energy axis.
The beam energy, convergence and collection angles need to be set in the
metadata.
Returns
-------
bstar : Signal1D
Contains the effective collection energy for each energy-loss in
mrad. It is a signal with the same length as the signal dimension.
References
----------
The calculation is done following Egerton (3rd edition) page 276.
"""
try:
e0 = self.metadata.Acquisition_instrument.TEM.beam_energy
except BaseException:
raise AttributeError("Please define the beam energy."
"You can do this e.g. by using the "
"set_microscope_parameters method")
try:
alpha = self.metadata.Acquisition_instrument.TEM.convergence_angle
except BaseException:
raise AttributeError("Please define the convergence semi-angle. "
"You can do this e.g. by using the "
"set_microscope_parameters method")
try:
beta = self.metadata.Acquisition_instrument.TEM.Detector.\
EELS.collection_angle
except BaseException:
raise AttributeError("Please define the collection semi-angle. "
"You can do this e.g. by using the "
"set_microscope_parameters method")
energy = self.axes_manager.signal_axes[0].axis
tgt = e0 * (1.+e0/1022.) / (1.+e0/511.)
thetae = (energy+1e-6) / tgt
# A2,B2,T2 ARE SQUARES OF ANGLES IN RADIANS**2
a2 = alpha * alpha * 1e-6 + 1e-10
b2 = beta * beta * 1e-6
t2 = thetae * thetae * 1e-6
eta1 = np.sqrt((a2+b2+t2)**2 - 4.*a2*b2) - a2 - b2 - t2
eta2 = 2.*b2 * np.log(0.5/t2*(np.sqrt((a2+t2-b2)**2 + 4.*b2*t2) \
+ a2 + t2 - b2))
eta3 = 2.*a2 * np.log(0.5/t2*(np.sqrt((b2+t2-a2)**2 + 4.*a2*t2) \
+ b2 + t2 - a2))
eta = (eta1+eta2+eta3) / a2 / np.log(4./t2)
f1 = (eta1+eta2+eta3) / 2. / a2 / np.log(1.+b2/t2)
f2 = f1
if (alpha/beta > 1.):
f2 = f1 * a2/b2
bstar = thetae * np.sqrt(np.exp(f2*np.log(1.+b2/t2))-1.)
bstar = self._get_signal_signal(bstar)
bstar.set_signal_type('Signal1D')
return bstar
def model_zero_loss_peak_mirror(self, hanning=True):
'''
Model the zero-loss peak as a mirror using the negative energy-losses.
This model only makes sense if the ZLP is symmetric around zero.
Parameters
----------
hanning : bool
Optionally apply a Hanning taper to the border, so that they intensity
decays smoothily to 0.
Returns
-------
zlp : EELSSignal
Mirrored ZLP model.
'''
zlp = self.deepcopy()
eax = zlp.axes_manager.signal_axes[0]
if eax.high_index*0.5 > eax.value2index(0.):
eaxshift = - eax.scale*0.5
elif eax.high_index*0.5 < eax.value2index(0.):
eaxshift = eax.scale*0.5
else:
eaxshift = 0
eax.offset = eax.offset + eaxshift
plz = zlp.deepcopy().isig[0.::-1]
xae = plz.axes_manager.signal_axes[0]
xae.scale *= -1
xae.offset *= -1
plz.get_dimensions_from_data()
Ea = xae.low_value
Eb = xae.high_value
Ec = eax.high_value
i1 = eax.value2index(Ea)
if Eb<Ec:
i2 = eax.value2index(Eb)
i3 = -1
elif Eb>Ec:
i2 = eax.high_index
i3 = i2-xae.high_index-1
else:
i2 = eax.high_index
i3 = -1
zlp.data[(slice(None),)*eax.index_in_array+(slice(i1, i2), Ellipsis)]= \
plz.data[(slice(None),)*xae.index_in_array+(slice(0, i3), Ellipsis)]
zlp.data[(slice(None),)*eax.index_in_array + \
(slice(i2, None), Ellipsis)] = 0.
eax = zlp.axes_manager.signal_axes[0]
eax.offset = eax.offset-eaxshift
if hanning:
zlp_h = zlp.isig[:Eb-eaxshift]
zlp_h.hanning_taper('both')
return zlp
def model_zero_loss_peak_tail(self, signal_range, show_progressbar=None,
*args, **kwargs):
'''
Model the zero-loss peak tail using a (power-law) model fit. The fit
window is set using the signal_range tuple in axis units (not indices).
Spectral intensity at energies above the signal_range is substituted by
the model tail. The fit is performed using `remove_background`,
*args and **kwargs are passed to this method.
Parameters
----------
signal_range : tuple
Initial and final position of the fit window. Given in axis units. The
components can be single or multidimensional. If multidimensional, an
array or a signal with the same dimensions as the navigation dimension
should be used.
Returns
-------
zlp : EELSSignal
Modeled tail ZLP model.
Examples
--------
>>> s = hs.load('some_eels.h5')
>>> zlp = s.model_zero_loss_peak_tail((0.5, 1.),
... fast=False,
... show_progressbar=False)
'''
self._check_signal_dimension_equals_one()
if not isinstance(signal_range, tuple):
raise AttributeError('signal_range not recognized:'
'must be a tuple!')
if len(signal_range) != 2:
raise AttributeError('signal_range not recognized,'
'must be len = 2')
axis = self.axes_manager.signal_axes[0]
if isinstance(signal_range[0], Number) and (
isinstance(signal_range[1], Number)):
zlp = self - self.remove_background(signal_range,
show_progressbar=show_progressbar,
*args, **kwargs)
I2 = axis.value2index(signal_range[1])
ids = (slice(None),)*axis.index_in_array + (slice(None, I2), Ellipsis)
zlp.data[ids] = self.data[ids]
return zlp
if isinstance(signal_range[0], (np.ndarray, hs.signals.BaseSignal)) or (
isinstance(signal_range[1], (np.ndarray, hs.signals.BaseSignal))):
signal_range_ini = self._check_adapt_map_input(signal_range[0])
signal_range_fin = self._check_adapt_map_input(signal_range[1])
for name in ['signal_range_ini', 'signal_range_fin']:
parameter = eval(name)
if isinstance(parameter, ValueError):
parameter.args = (parameter.args[0]+name,)
raise parameter
zlp = self.deepcopy()
for si in progressbar(self, disable=not show_progressbar):
indices = self.axes_manager.indices
E1 = signal_range_ini.inav[indices].data[0]
E2 = signal_range_fin.inav[indices].data[0]
ri = si.remove_background((E1, E2), show_progressbar=False,
*args, **kwargs)
I2 = axis.value2index(E2)
ids = (*indices, slice(I2, None), Ellipsis)
zlp.data[ids] = si.data[I2:] - ri.data[I2:]
return zlp
def power_law_extrapolation_until(self, window_size=20, total_size=1024,
hanning=True, *args, **kwargs):
'''
Usefullity mathod. Extends the spectrum using `power_law_extrapolation`
until the resulting spectrum has the given total number of pixels.
Parameters
----------
window_size : {int, float}
Window size. Alternatively, this parameter can be a float specifying
the size in energy-loss units.
total_size : {int, float}
Total number of pixels. Alternatively, this parameter can be a float
specifying the size in energy-loss units.
Returns
-------
spc : EELSSpectrum
The extended spectrum.
'''
eax = self.axes_manager.signal_axes[0]
if type(total_size) is float:
offset = self.axes_manager[-1].offset - self.axes_manager[-1].scale
total_size = int(round((total_size-offset)/eax.scale))
if type(window_size) is float:
window_size = int(round(window_size/eax.scale))
extrapolation_size = total_size-eax.size
if extrapolation_size < 0:
raise AttributeError("total_size is less than spectral axis size")
spc = self.power_law_extrapolation(window_size,
extrapolation_size,
*args,
**kwargs)
if hanning:
spc.hanning_taper('both')
return spc
def relativistic_kramers_kronig(self,
zlp=None,
n=None,
t=None,
delta=0.9,
fsmooth=None,
iterations=20,
chi2_target=1e-4,
average=False,
full_output=True,
show_progressbar=None,
*args,
**kwargs):
r"""Calculate the complex dielectric function from a single scattering
distribution (SSD) using the Kramers-Kronig relations and a relativistic
correction for thin slab geometry.
The input SSD should be and EELSSpectrum instance, containing only
inelastic scattering information (elastic and plural scattering
deconvolved). The dielectric information is obtained by normalization of
the inelastic scattering using the elastic scattering intensity and
either refractive index or thickness information.
A full complex dielectric function (CDF) is obtained by Kramers-Kronig
transform, solved using FFT as in `kramers_kronig_analysis`. This inital
guess for the CDF is improved in an iterative loop, devised to
approximately subtract the relativistic contribution supposing an
unoxidized planar surface.
The loop runs until a chi-square target has been achieved or for a
maximum number of iterations. This behavior can be modified using the
parameters below. This method does not account for instrumental and
finite-size effects.
Note: either refractive index or thickness (`n` or `t`) are required.
If both are None or if both are provided an exception is raised. Many
input types are accepted for zlp, n and t parameters, which are parsed
using `self._check_adapt_map_input`, see the documentation therein for
more information.
Parameters
----------
zlp: {None, number, ndarray, Signal}
ZLP intensity. It is optional (can be None) if t is given,
full_output is False and no iterations are run. In any other case,
the ZLP is required either to perform the normalization step,
to calculate the thickness and/or to calculate the relativistic
correction.
n: {None, number, ndarray, Signal}
The medium refractive index. Used for normalization of the
SSD to obtain the energy loss function. If given the thickness
is estimated and returned. It is only required when `t` is None.
t: {None, number, ndarray, Signal}
The sample thickness in nm. Used for normalization of the
SSD to obtain the energy loss function. It is only required when
`n` is None.
delta : {None, float}
Optionally apply a fractional limit to the relativistic correction
in order to improve stability. Can be None, if no limit is desired.
A value of around 0.9 ensures the correction is never larger than
the original EELS signal, producing a negative spectral region.
fsmooth : {None, float}
Optionally apply a gaussian filter to the relativistic correction
in order to eliminate high-frequency noise. The cut-off is set in
the energy-loss scale, e.g. fsmooth = 1.5 (eV).
iterations: {None, int}
Number of the iterations for the internal loop to remove the
relativistic contribution. If None, the loop runs until a chi-square
target has been achieved (see below).
chi2_target : float
The average chi-square test score is measured in each iteration, and
the reconstruction loop terminates when the target score is reached.
See `_chi2_score` for more information.
average : bool
If True, use the average of the obtained dielectric functions over
the navigation dimensions to calculate the relativistic correction.
False by default, should only be used when analyzing spectra from a
homogenous sample, as only one dielectric function is retrieved.
This switch has no effect if only one spectrum is being analyzed.
full_output : bool
If True, return a dictionary that contains the estimated
thickness if `t` is None and the estimated relativistic correction
if `iterations` > 1.
Returns
-------
eps: DielectricFunction instance
The complex dielectric function results,
.. math::
\epsilon = \epsilon_1 + i*\epsilon_2,
contained in an DielectricFunction instance.
output: Dictionary (optional)
A dictionary of optional outputs with the following keys:
``thickness``
The estimated thickness in nm calculated by normalization of
the corrected spectrum (only when `t` is None).
``relativistic correction``
The estimated relativistic correction at the final iteration.
Raises
------
ValueError
If both `n` and `t` are undefined (None).
AttributeError
If the beam_energy or the collection semi-angle are not defined in
metadata.
See also
--------
get_relativistic_spectrum, _check_adapt_map_input
"""
# prepare data arrays
if iterations == 1:
# In this case s.data is not modified so there is no need to make
# a deep copy.
s = self.isig[0.:]
else:
s = self.isig[0.:].deepcopy()
sorig = self.isig[0.:]
# Avoid singularity at 0
if s.axes_manager.signal_axes[0].axis[0] == 0:
s = s.isig[1:]
sorig = self.isig[1:]
axis = s.axes_manager.signal_axes[0]
eaxis = axis.axis.copy()
# Constants and units, electron mass, beam energy and collection angle
me = constants.value(
'electron mass energy equivalent in MeV') * 1e3 # keV
try:
e0 = s.metadata.Acquisition_instrument.TEM.beam_energy
except BaseException:
raise AttributeError("Please define the beam energy."
"You can do this e.g. by using the "
"set_microscope_parameters method")
try:
beta = s.metadata.Acquisition_instrument.TEM.Detector.\
EELS.collection_angle
except BaseException:
raise AttributeError("Please define the collection semi-angle. "
"You can do this e.g. by using the "
"set_microscope_parameters method")
# Mapped parameters, zlp, n and t
if isinstance(zlp, hs.signals.Signal1D):
if (zlp.axes_manager.signal_dimension == 1) and (
zlp.axes_manager.navigation_shape ==
self.axes_manager.navigation_shape):
zlp = zlp.integrate1D(axis.index_in_axes_manager)
elif zlp is None and (full_output or iterations > 1):
raise AttributeError("Please define the zlp parameter when "
"full output or iterations > 1 are selected.")
zlp = self._check_adapt_map_input(zlp)
n = self._check_adapt_map_input(n)
t = self._check_adapt_map_input(t)
for name in ['zlp', 'n', 't']:
parameter = eval(name)
if isinstance(parameter, ValueError):
parameter.args = (parameter.args[0]+name,)
raise parameter
# select refractive or thickness loop
if n is None and t is None:
raise ValueError('Thickness and refractive index undefined.'
'Please provide one of them.')
elif n is not None and t is not None:
raise ValueError('Thickness and refractive index both defined.'
'Please provide only one of them.')
elif n is not None:
refractive_loop = True
if (zlp is not None) and (full_output is True or iterations > 1):
t = self._get_navigation_signal().T
elif t is not None:
refractive_loop = False
if zlp is None:
raise ValueError('Zero-loss intensity is needed for thickness '
'normalization. Provide also parameter zlp')
# Slicer to get the signal data from 0 to axis.size
slicer = s.axes_manager._get_data_slice(
[(axis.index_in_array, slice(None, axis.size)), ])
# Kinetic definitions
ke = e0 * (1 + e0 / 2. / me) / (1 + e0 / me) ** 2
tgt = e0 * (2 * me + e0) / (me + e0)
rk0 = 2590 * (1 + e0 / me) * np.sqrt(2 * ke / me)
# prepare the output dielectric function
eps = s._deepcopy_with_new_data(np.zeros_like(s.data, np.complex128))
eps.set_signal_type("DielectricFunction")
eps.metadata.General.title = (self.metadata.General.title +
'KKA dielectric function')
if eps.tmp_parameters.has_item('filename'):
eps.tmp_parameters.filename = (
self.tmp_parameters.filename +
'_CDF_after_Kramers_Kronig_transform')
from dielectric import ModifiedCDF
eps = ModifiedCDF(eps)
eps_corr = eps.deepcopy()
# progressbar support
if show_progressbar is None:
show_progressbar = preferences.General.show_progressbar
pbar = progressbar(total = iterations,
desc = '1.00e+30',
disable=not show_progressbar)
# initialize iteration control
io = 0
chi2 = chi2_target*1e3
while (io < iterations) and (chi2 > chi2_target):
# Calculation of the ELF by normalization of the SSD
Im = s.data / (np.log(1 + (beta * tgt / eaxis) ** 2)) / axis.scale
if refractive_loop:
# normalize using the refractive index.
K = (Im / eaxis).sum(axis=axis.index_in_array) * axis.scale
K = (K / (np.pi / 2) / (1 - 1. / n.data ** 2))
# Calculate the thickness only if possible and required
if full_output or iterations > 1:
te = (332.5 * K * ke / zlp.data)
t.data = te.squeeze()
else:
# normalize using the thickness
K = t.data * zlp.data / (332.5 * ke)
Im = Im / K[..., None] if len(self) != 1 else Im / K
# Kramers-Kronig transform
esize = 2 * axis.size
q = -2 * np.fft.fft(Im, esize, axis.index_in_array).imag / esize
q[slicer] *= -1
q = np.fft.fft(q, axis=axis.index_in_array)
Re = q[slicer].real + 1
epsabs = (Re ** 2 + Im ** 2)
eps.data = Re / epsabs + 1j* Im / epsabs
del Im, Re, q, epsabs
if average and (eps.axes_manager.navigation_dimension > 0):
eps_corr.data[:] = eps.data.mean(
eps.axes_manager.navigation_indices_in_array,
keepdims=True)
else:
eps_corr.data = eps.data.copy()
if full_output or iterations > 1:
# Relativistic correction
# Calculates relativistic correction from the Kroeger equation
# The difference with the relativistic DCS is subtracted
scorr = eps_corr.get_relativistic_spectrum(zlp=zlp, t=t,
output='diff',
*args, **kwargs)
# Limit the fractional correction
if delta is not None:
fcorr = np.clip(scorr.data / sorig.data, -delta, delta)
scorr.data = fcorr * sorig.data
# smooth
if fsmooth is not None:
scorr.gaussian_filter(fsmooth)
# Apply correction
s.data = sorig.data - scorr.data
s.data[s.data < 0.] = 0.
if io > 0:
#chi2 = ((scorr.data-smemory)**2/smemory**2).sum()
chi2 = smemory._chi2_score(scorr)
chi2str = '{:0.2e}'.format(chi2)
pbar.set_description(chi2str)
smemory = scorr.deepcopy()
io += 1
pbar.update(1)
pbar.close()
if full_output:
output = {}
tstr = self.metadata.General.title
if refractive_loop:
t.metadata.General.title = tstr + ', r-KKA thickness'
output['thickness'] = t
scorr.metadata.General.title = tstr + ', r-KKA correction'
output['relativistic correction'] = scorr
return eps, output
else:
return eps
def _chi2_score(self, sig, p=0.5):
"""
Calculate the mean chi-2 score removing outliers in a pedestrian way.
Parameters
----------
sig : signal
Must be of equal dimensions to self or this won't work.
p : float
Remove outliers coefficient: percentage of the data left out above and
below.
Returns
-------
chi2 : float
Average of the chi-square test score.
"""
chi2 = (sig.data - self.data)**2 / self.data**2
chi2 = np.nan_to_num(chi2)
a_min, a_max = np.percentile(chi2, p), np.percentile(chi2, 100-p)
return np.clip(chi2, a_min, a_max).mean()