epsproc.MFPAD module
ePSproc MFPAD functions
- 29/10/20 v2 Updated with Wigner delay routine.
Tidied up To Do list items. Tidied up docstring.
- 05/08/19 v1 Initial python version.
Based on original Matlab code ePS_MFPAD.m
Structure
Calculate MFPAD on a grid from input ePS matrix elements. Use fast functions, pre-calculate if possible. Data in Xarray, use selection functions and multiplications based on relevant quantum numbers, other axes summed over.
- Choices for functions…
Moble’s spherical functions (quaternion based)
Provides fast wignerD, 3j and Ylm functions, uses Numba.
Install with:
>>> conda install -c conda-forge spherical_functions
To Do
More efficient computation? Use Xarray group by?
Formalism
MF-PADs numerical expansion
This is implemented numerically in epsproc.MFPAD.mfpad()
. For direct computation of $beta_{LM}$ parameters (rather than full observable expansion) see epsproc.MFBLM()
(original looped version) or epsproc.geomFunc.mfblmGeom()
(geometric tensor version).
From ePSproc: Post-processing suite for ePolyScat electron-molecule scattering calculations.
In this formalism:
\(I_{l,m,\mu}^{p_{i}\mu_{i},p_{f}\mu_{f}}(E)\) is the radial part of the dipole matrix element, determined from the initial and final state electronic wavefunctions \(\Psi_{i}^{p_{i},\mu_{i}}\) and \(\Psi_{f}^{p_{f},\mu_{f}}\), photoelectron wavefunction \(\varphi_{klm}^{(-)}\) and dipole operator \(\hat{d_{\mu}}\). Here the wavefunctions are indexed by irreducible representation (i.e. symmetry) by the labels \(p_{i}\) and \(p_{f}\), with components \(\mu_{i}\) and \(\mu_{f}\) respectively; \(l,m\) are angular momentum components, \(\mu\) is the projection of the polarization into the MF (from a value \(\mu_{0}\) in the LF). Each energy and irreducible representation corresponds to a calculation in ePolyScat.
\(T_{\mu_{0}}^{p_{i}\mu_{i},p_{f}\mu_{f}}(\theta_{\hat{k}},\phi_{\hat{k}},\theta_{\hat{n}},\phi_{\hat{n}})\) is the full matrix element (expanded in polar coordinates) in the MF, where \(\hat{k}\) denotes the direction of the photoelectron \(\mathbf{k}\)-vector, and \(\hat{n}\) the direction of the polarization vector \(\mathbf{n}\) of the ionizing light. Note that the summation over components \(\{l,m,\mu\}\) is coherent, and hence phase sensitive.
\(Y_{lm}^{*}(\theta_{\hat{k}},\phi_{\hat{k}})\) is a spherical harmonic.
\(D_{\mu,\mu_{0}}^{1}(R_{\hat{n}})\) is a Wigner rotation matrix element, with a set of Euler angles \(R_{\hat{n}}=(\phi_{\hat{n}},\theta_{\hat{n}},\chi_{\hat{n}})\), which rotates/projects the polarization into the MF.
\(I_{\mu_{0}}(\theta_{\hat{k}},\phi_{\hat{k}},\theta_{\hat{n}},\phi_{\hat{n}})\) is the final (observable) MFPAD, for a polarization \(\mu_{0}\) and summed over all symmetry components of the initial and final states, \(\mu_{i}\) and \(\mu_{f}\). Note that this sum can be expressed as an incoherent summation, since these components are (by definition) orthogonal.
\(g_{p_{i}}\) is the degeneracy of the state \(p_{i}\).
For more details see: Toffoli, D., Lucchese, R. R., Lebech, M., Houver, J. C., & Dowek, D. (2007). Molecular frame and recoil frame photoelectron angular distributions from dissociative photoionization of NO2. The Journal of Chemical Physics, 126(5), 054307. https://doi.org/10.1063/1.2432124
(MF) Wigner Delays
This is implemented numerically in epsproc.MFPAD.mfWignerDelay()
.
The partial wave, or channel resolved, Wigner delay can be given as the energy-derivative of the continuum wavefunction phase, $arg(psi_{lm}^{*})$
And the total, or group, delay from the (coherent) sum over these channels.
Here the full wavefunction $sum_{l,m}psi_{lm}^{*}$ is essentially equivalent to the phase of the (conjugate) $T_{mu_{0}}^{p_{i}mu_{i},p_{f}mu_{f}*}(theta_{hat{k}},phi_{hat{k}},theta_{hat{n}},phi_{hat{n}})$ functions given above (and by epsproc.MFPAD.mfpad()
).
For more details see, for example:
Hockett, Paul, Eugene Frumker, David M Villeneuve, and Paul B Corkum. “Time Delay in Molecular Photoionization.” Journal of Physics B: Atomic, Molecular and Optical Physics 49, no. 9 (May 2016): 095602. https://doi.org/10.1088/0953-4075/49/9/095602.
Carvalho, C.A.A. de, and H.M. Nussenzveig. “Time Delay.” Physics Reports 364, no. 2 (June 2002): 83–174. https://doi.org/10.1016/S0370-1573(01)00092-8.
Wigner, Eugene. “Lower Limit for the Energy Derivative of the Scattering Phase Shift.” Physical Review 98, no. 1 (April 1955): 145–47. https://doi.org/10.1103/PhysRev.98.145.
- epsproc.MFPAD.mfWignerDelay(dataIn, pType='phaseUW')[source]
Compute MF Wigner Delay as phase derivative of continuum wavefunction.
- Parameters
dataIn (Xarray) – Contains set(s) of wavefunctions to use, as output by
epsproc.MFPAD.mfpad()
.pType (str, optional, default = 'phaseUW') – Used to set data conversion options, as implemented in
epsproc.plotTypeSelector()
- ‘phase’ use np.angle() - ‘phaseUW’ unwapped with np.unwrap(np.angle())
- epsproc.MFPAD.mfpad(dataIn, thres=0.01, inds={'Type': 'L', 'it': 1}, res=50, R=None, p=0)[source]
- Parameters
dataIn (Xarray) – Contains set(s) of matrix elements to use, as output by epsproc.readMatEle().
thres (float, optional, default 1e-2) – Threshold value for matrix elements to use in calculation.
ind (dictionary, optional.) – Used for sub-selection of matrix elements from Xarrays. Default set for length gauage, single it component only, inds = {‘Type’:’L’,’it’:’1’}.
res (int, optional, default 50) – Resolution for output (theta,phi) grids.
R (list or Xarray of Euler angles or quaternions, optional.) – Define LF > MF polarization geometry/rotations. For default case (R = None), calls
epsproc.setPolGeoms()
, where 3 geometries are calculated, corresponding to z-pol, x-pol and y-pol cases. Defined by Euler angles (p,t,c) = [0 0 0] for z-pol, [0 pi/2 0] for x-pol, [pi/2 pi/2 0] for y-pol.p (int, optional.) – Defines LF polarization state, p = -1…1, default p = 0 (linearly pol light along z-axis). TODO: add summation over p for multiple pol states in LF.
- Returns
Ta – Xarray (theta, phi, E, Sym) of MFPADs, summed over (l,m)
Tlm – Xarray (theta, phi, E, Sym, lm) of MFPAD components, expanded over all (l,m)