epsproc.geomFunc.geomCalc module¶
ePSproc geometric terms/functions
Collection of codes for geometric functions and tensors.
26/02/20 v1 Initial implementation.
-
epsproc.geomFunc.geomCalc.
CG
(QNs, dlist=['l', 'lp', 'L', 'm', 'mp', 'M'], form='xarray')[source]¶ Basic Clebsch-Gordan from 3j calculation, from table of input QNs (corresponding to CG term defn.).
This implements numerical defn. from Moble’s Spherical Functions, https://github.com/moble/spherical_functions/blob/master/spherical_functions/recursions/wigner3j.py
def clebsch_gordan(j_1, m_1, j_2, m_2, j_3, m_3)
(-1.)**(j_1-j_2+m_3) * math.sqrt(2*j_3+1) * Wigner3j(j_1, j_2, j_3, m_1, m_2, -m_3)
22/06/20 - barebones version for quick testing, should upgrade as per w3jTable (which is the back-end here in any case).
Parameters: - QNs (np.array) – List of QNs [l, lp, L, m, mp, M] to compute 3j terms for.
- form (string, optional, default = 'xarray') –
Defines return format. Options are:
- 2d, return 2D np.array, rows [l, lp, L, m, mp, M, 3j]
- xarray, return xarray
- This is nice for easy selection/indexing, but may be problematic for large Lmax if unstacked (essentailly similar to nd case).
- nd, return ND np.array, dims indexed as [l, lp, L, l+m, lp+mp, L+M], with values 3j.
- This is suitable for direct indexing, but will have a lot of zero entries and may be large.
- ndsparse, return ND sparse array, dims indexed as [l, lp, L, l+m, lp+mp, L+M], with values 3j.
Additional options are set via
remapllpL()
. This additionally sorts values by (l,lp,L) triples, which is useful in some cases.- ’dict’ : dictionary with keys (l,lp,L), coordinate tables
- ’3d’ : dictionary with keys (l,lp,L), 3D arrays indexed by [l+m, lp+mp, L+M]; this case also sets (0,0,0) term as ‘w3j0’.
- ’xdaLM’ : Xarray dataarray, with stacked dims [‘lSet’,’mSet’]
- ’xds’ : Xarray dataset, with one array per (l,lp,L)
- ’xdalist’ : List of Xarray dataarrays, one per (l,lp,L)
- dlist (list of labels, optional, default ['l','lp','L','m','mp','M']) – Used to label array for Xarray output case.
-
epsproc.geomFunc.geomCalc.
EPR
(QNs=None, p=None, ep=None, nonzeroFlag=True, form='2d', dlist=['l', 'lp', 'P', 'p', 'R-p', 'R'], phaseConvention='S')[source]¶ Define polarization tensor (LF) for 1-photon case.
Define field terms (from QM book, corrected vs. original S&U version - see
beta-general-forms\_rewrite\_290917.lyx
):\[\begin{split}\begin{equation} E_{PR}(\hat{e})=[e\otimes e^{*}]_{R}^{P}=[P]^{\frac{1}{2}}\sum_{p}(-1)^{R}\left(\begin{array}{ccc} 1 & 1 & P\\ p & R-p & -R \end{array}\right)e_{p}e_{R-p}^{*} \end{equation}\end{split}\]Parameters: - QNs (np.array, optional, default = None) – List of QNs [l, lp, P, m, mp, R] to compute 3j terms for. If not supplied all allowed terms for the one-photon case, l=lp=1, will be set.
- p (list or array, optional, default = None) – Specify polarization terms p. If set to None, all allowed values for the one-photon case will be set, p=[-1,0,1]
- ep (list or array, optional, default = None) – Relative strengths for the fields ep. If set to None, all terms will be set to unity, ep = 1
- nonzeroFlag (bool, optional, default = True) – Drop null terms before returning values if true.
- form (string, optional, default = '2d') – For options see
ep.w3jTable()
- phaseConvention (optional, str, default = 'S') –
Set phase conventions:
- ’S’ : Standard derivation.
- ’R’ : Reduced form geometric tensor derivation.
- ’E’ : ePolyScat, may have additional changes in numerics, e.g. conjugate Wigner D.
See
setPhaseConventions()
for more details.
Examples
# Generate full EPR list with defaults >>> EPRtable = EPR()
# Return as Xarray >>> EPRtable = EPR(form = ‘xarray’)
Note
Currently not handling ep correctly! Should implement as passed Xarray for correct (p,p’) assignment.
-
epsproc.geomFunc.geomCalc.
MFproj
(QNs=None, RX=None, nonzeroFlag=True, form='2d', dlist=['l', 'lp', 'P', 'mu', 'mup', 'Rp', 'R'], phaseConvention='S')[source]¶ Define MF projection term, \(\Lambda_{R',R}(R_{\hat{n}})\):
\[\begin{split}\begin{equation} \Lambda_{R',R}(R_{\hat{n}})=(-1)^{(R')}\left(\begin{array}{ccc} 1 & 1 & P\\ \mu & -\mu' & R' \end{array}\right)D_{-R',-R}^{P}(R_{\hat{n}}) \end{equation}\end{split}\]Then…
\[\begin{eqnarray} \beta_{L,-M}^{\mu_{i},\mu_{f}} & = & \sum_{P,R',R}{\color{red}E_{P-R}(\hat{e};\mu_{0})}\sum_{l,m,\mu}\sum_{l',m',\mu'}(-1)^{(\mu'-\mu_{0})}{\color{red}\Lambda_{R',R}(R_{\hat{n}};\mu,P,R,R')B_{L,-M}(l,l',m,m')}I_{l,m,\mu}^{p_{i}\mu_{i},p_{f}\mu_{f}}(E)I_{l',m',\mu'}^{p_{i}\mu_{i},p_{f}\mu_{f}*}(E) \end{eqnarray}\]Parameters: phaseConvention (optional, str, default = 'S') – Set phase conventions: - ‘S’ : Standard derivation. - ‘R’ : Reduced form geometric tensor derivation. - ‘E’ : ePolyScat, conjugate Wigner D. See setPhaseConventions()
for more details.Notes
This is very similar to \(E_{PR}\) term.
Examples
>>> lTerm, lambdaTable, lambdaD, QNs = MFproj(form = 'xarray')
-
epsproc.geomFunc.geomCalc.
betaTerm
(QNs=None, Lmin=0, Lmax=10, nonzeroFlag=True, form='2d', dlist=['l', 'lp', 'L', 'm', 'mp', 'M'], phaseConvention='S')[source]¶ Define BLM coupling tensor
Define field terms (from QM book, corrected vs. original S&U version - see
beta-general-forms\_rewrite\_290917.lyx
):\[\begin{split}\begin{equation} B_{L,M}=(-1)^{m}\left(\frac{(2l+1)(2l'+1)(2L+1)}{4\pi}\right)^{1/2}\left(\begin{array}{ccc} l & l' & L\\ 0 & 0 & 0 \end{array}\right)\left(\begin{array}{ccc} l & l' & L\\ -m & m' & M \end{array}\right) \end{equation}\end{split}\]Parameters: - QNs (np.array, optional, default = None) – List of QNs [l, lp, L, m, mp, M] to compute 3j terms for. If not supplied all allowed terms for {Lmin, Lmax} will be set. (If supplied, values for Lmin, Lmax and mFlag are not used.)
- Lmax (Lmin,) – Integer values for Lmin and Lmax respectively.
- mFlag (bool, optional, default = True) – m, mp take all values -l…+l if mFlag=True, or =0 only if mFlag=False
- nonzeroFlag (bool, optional, default = True) – Drop null terms before returning values if true.
- form (string, optional, default = '2d') – For options see
ep.w3jTable()
- dlist (list of labels, optional, default ['l','lp','L','m','mp','M']) – Used to label array for Xarray output case.
- phaseConvention (optional, str, default = 'S') – Set phase conventions:
- ‘S’ : Standard derivation.
- ‘R’ : Reduced form geometric tensor derivation.
- ‘E’ : ePolyScat, may have additional changes in numerics, e.g. conjugate Wigner D.
See
setPhaseConventions()
for more details.
Examples
>>> Lmax = 2 >>> BLMtable = betaTerm(Lmax = Lmax, form = 'xds') >>> BLMtable = betaTerm(Lmax = Lmax, form = 'xdaLM')
-
epsproc.geomFunc.geomCalc.
deltaLMKQS
(EPRX, AKQS, phaseConvention='S')[source]¶ Calculate aligned-frame “alignment” term:
\[\begin{equation} \sum_{K,Q,S}\Delta_{L,M}(K,Q,S)A_{Q,S}^{K}(t) \end{equation}\]\[\begin{split}\begin{equation} \Delta_{L,M}(K,Q,S)=(2K+1)^{1/2}(-1)^{K+Q}\left(\begin{array}{ccc} P & K & L\\ R & -Q & -M \end{array}\right)\left(\begin{array}{ccc} P & K & L\\ R' & -S & S-R' \end{array}\right) \end{equation}\end{split}\]15/06/20 IN PROGRESS
Parameters: - EPRX (Xarray) – Polarization terms in an Xarray, as set by
epsproc.geomCalc.EPR()
- AKQS (Xarray) – Alignement terms in an Xarray, as set by
epsproc.setADMs()
Returns: - AFterm (Xarray) – Full term, including multiplication and sum over (K,Q,S) (note S-Rp term is retained).
- DeltaKQS (Xarray) – Alignment term \(\Delta_{L,M}(K,Q,S)\).
- To do
- —–
- - Add optional inputs.
- - Add error checks.
- See other similar functions for schemes.
- EPRX (Xarray) – Polarization terms in an Xarray, as set by
-
epsproc.geomFunc.geomCalc.
deltalmKQSwf
(matE, AKQS, phaseConvention='S', dlist1=['l', 'E', 'K', 'm', 'mu', 'Q'], dlist2=['l', 'E', 'K', 'Lambda', 'mu0', 'S'])[source]¶ Calculate aligned-frame “alignment” term, version for AF wavefunction expansion.
\[\begin{split}\begin{equation} ^{AF}\Delta_{l,m}(K,Q,S)=(-1)^{m-\Lambda}(-1)^{\mu-\mu_{0}}(-1)^{Q-S}\left(\begin{array}{ccc} l & 1 & K\\ -m & -\mu & -Q \end{array}\right)\left(\begin{array}{ccc} l & 1 & K\\ -\Lambda & -\mu_{0} & -S \end{array}\right) \end{equation}\end{split}\]13/01/21 IN PROGRESS, adapted from existing deltaLMKQS() function.
NOTE: photon QN currently labelled as (E,mu,mu0), but may want to change to avoid confusion with EPR term.
Parameters: - matE (Xarray) – Xarray containing matrix elements, with QNs (l,m), as created by
readMatEle()
- AKQS (Xarray) – Alignement terms in an Xarray, as set by
epsproc.setADMs()
- phaseConvention (optional, str, default = 'S') – Set phase conventions with
epsproc.geomCalc.setPhaseConventions()
. To use preset phase conventions, pass existing dictionary. If matE.attrs[‘phaseCons’] is already set, this will be used instead of passed args. - dlist2 (dlist1,) – Labels for output QNs.
Returns: - AFterm (Xarray) – Full term, including multiplication and sum over (K,Q,S).
- DeltaKQS (Xarray) – Alignment term \(\Delta_{L,M}(K,Q,S)\).
- To do
- —–
- - Add optional inputs.
- - Add error checks.
- See other similar functions for schemes.
- matE (Xarray) – Xarray containing matrix elements, with QNs (l,m), as created by
-
epsproc.geomFunc.geomCalc.
remapllpL
(dataIn, QNs, form='dict', method='sel', dlist=['l', 'lp', 'L', 'm', 'mp', 'M'], verbose=0)[source]¶ Remap Wigner 3j table, with QNs (l,lp,L,m,mp,M) > tensor forms.
Tensors are stored by (l,lp,L) triples, with corresponding arrays [m,mp,M], as a dictionary, or Xarray dataarray or dataset.
Parameters: - dataIn (np.array) – Array of data values corresponding to rows (coords) in QNs.
- QNs (np.array) – List of QNs [l, lp, L, m, mp, M] to compute 3j terms for. If not supplied all allowed terms for {Lmin, Lmax} will be set. (If supplied, values for Lmin, Lmax and mFlag are not used.)
- form (str, optional, default = 'dict') – Type of output structure. - ‘dict’ : dictionary with keys (l,lp,L), coordinate tables - ‘3d’ : dictionary with keys (l,lp,L), 3D arrays indexed by [l+m, lp+mp, L+M]; this case also sets (0,0,0) term as ‘w3j0’. - ‘xdaLM’ : Xarray dataarray, with stacked dims [‘lSet’,’mSet’] - ‘xds’ : Xarray dataset, with one array per (l,lp,L) - ‘xdalist’ : List of Xarray dataarrays, one per (l,lp,L)
- method (str, optional, default = 'sel') – Method for selection from input array.
- ‘sel’ : use full selection routine,
epsproc.selQNsRow()
, should be most robust. - ‘n’ : sort based on np.unique reindexing. Should be faster, but possibly not robust…? - dlist (list, optional, default = ['l','lp','L','m','mp','M']) – Full list of dimension names to use/set.
Returns: w3j – Wigner3j(l,lp,L,m,mp,M) values corresponding to rows (coords) in QNs, type according to input form. Data structure sorted by (l,lp,L) triples.
Return type: Xarray, dictionary
-
epsproc.geomFunc.geomCalc.
setPhaseConventions
(phaseConvention='S', typeList=False)[source]¶ Set phase convention/choices for geometric functions.
20/03/20 - first attempt. Aim to centralise all phase choices here to keep things clean and easy to debug/change.
Set as dictionary for each term, to be appended to results Xarray.
Parameters: - phaseConvention (optional, str, default = 'S') – Set phase conventions: - ‘S’ : Standard derivation. - ‘R’ : Reduced form geometric tensor derivation. - ‘E’ : ePolyScat, may have additional changes in numerics, e.g. conjugate Wigner D. If a dict of phaseConventions is passed they will simply be returned - this is for transparency/consistency over multiple fns which call setPhaseConventions()… although may be an issue in some cases.
- typeList (optional, bool, default = False) – If true, return list of supported options instead of list of phase choices.
Note
If a dict of phaseConventions is passed they will simply be returned - this is for transparency/consistency over multiple fns which call setPhaseConventions()… although may be an issue in some cases.
-
epsproc.geomFunc.geomCalc.
w3jTable
(Lmin=0, Lmax=10, QNs=None, mFlag=True, nonzeroFlag=False, form='2d', dlist=['l', 'lp', 'L', 'm', 'mp', 'M'], backend='par', verbose=0)[source]¶ Calculate/tabulate all wigner 3j terms for a given problem/set of QNs.
\[\begin{split}\begin{equation} \begin{array}{ccc} l & l' & L\\ m & m' & M \end{array} \end{equation}\end{split}\]Where l, l’ take values Lmin…Lmax (default 0…10). \(\l-lp\<=L<=l+lp\) m, mp take values -l…+l if mFlag=True, or =0 only if mFlag=False
Parameters: - Lmax (Lmin,) – Integer values for Lmin and Lmax respectively.
- QNs (np.array, optional, default = None) – List of QNs [l, lp, L, m, mp, M] to compute 3j terms for.
If not supplied all allowed terms for {Lmin, Lmax} will be set.
(If supplied, values for Lmin, Lmax and mFlag are not used.)
NOTE: some return types will convert QNs to int when constructing Xarray output, unless half-int values present.
Functions using
epsproc.geomFunc.geomCalc.remapllpL()
support half-int values in Xarray. - mFlag (bool, optional, default = True) – m, mp take all values -l…+l if mFlag=True, or =0 only if mFlag=False
- nonzeroFlag (bool, optional, default = False) – Drop null terms before returning values if true.
- form (string, optional, default = '2d') –
- Defines return format. Options are:
- 2d, return 2D np.array, rows [l, lp, L, m, mp, M, 3j]
- xarray, return xarray
- This is nice for easy selection/indexing, but may be problematic for large Lmax if unstacked (essentailly similar to nd case).
- nd, return ND np.array, dims indexed as [l, lp, L, l+m, lp+mp, L+M], with values 3j.
- This is suitable for direct indexing, but will have a lot of zero entries and may be large.
- ndsparse, return ND sparse array, dims indexed as [l, lp, L, l+m, lp+mp, L+M], with values 3j.
- Additional options are set via
remapllpL()
. This additionally sorts values by (l,lp,L) triples, which is useful in some cases. - ’dict’ : dictionary with keys (l,lp,L), coordinate tables
- ’3d’ : dictionary with keys (l,lp,L), 3D arrays indexed by [l+m, lp+mp, L+M]; this case also sets (0,0,0) term as ‘w3j0’.
- ’xdaLM’ : Xarray dataarray, with stacked dims [‘lSet’,’mSet’]
- ’xds’ : Xarray dataset, with one array per (l,lp,L)
- ’xdalist’ : List of Xarray dataarrays, one per (l,lp,L)
- dlist (list of labels, optional, default ['l','lp','L','m','mp','M']) – Used to label array for Xarray output case.
- backend (str, optional, default = 'par') – See Implementation note below.
Returns: w3j – Wigner3j(l,lp,L,m,mp,M) values corresponding to rows (coords) in QNs, type according to input form.
Return type: np.array, Xarray, dictionary
- Currently set to run:
- ‘vec’:
w3jguVecCPU()
, which uses sf.Wigner3j on the back-end, with additional vectorisation over supplied QNs via Numba’s @guvectorize. - ‘par’:
w3jprange()
, which uses sf.Wigner3j on the back-end, with parallelization over QNs via Numba’s @njit with a prange loop.
- ‘vec’: