utils

edrixs.utils.CT_imp_bath(U_dd, Delta, n)[source]

Compute energies of the impurity and bath for an Anderson impurity or charge-transfer model appropriate for a \(d\)-shell transition metal compound. Core hole potential effects are omitted from this calculation because one typically performs the ground state diagonalization step leaving out these states.

Parameters:
U_dd: float

Coulomb interaction \(U_{dd}\)

Delta: float

Charge-transfer energy \(\Delta\)

ninteger

Number of electrons in the \(d\)-shell

Returns:
E_dfloat

Energy of the impurity states \(E_d\)

E_Lfloat

Energy of the bath states \(E_L\)

Notes

We emit core electrons for initial state calculation which is done by setting $n_p=0$. It is crucial to properly set the energy splitting of $E_d$ and $E_L$ as defined by $Delta$. The initial state energy is ill-defined to within an additive factor, so we arbitrarily set the energy with no ligand holes to zero. From there, obtain energies by solving the relevant linear equations.

edrixs.utils.CT_imp_bath_core_hole(U_dd, U_pd, Delta, n)[source]

Compute energies of the impurity and bath for an Anderson impurity or charge-transfer model appropriate for a \(d\)-shell transition metal compound with a core hole.

Parameters:
U_dd: float

Coulomb interaction \(U_{dd}\)

U_pd: float

Coulomb interaction \(U_{pd}\)

Delta: float

Charge-transfer energy \(\Delta\)

ninteger

Number of electrons in the \(d\)-shell \(n\)

Returns:
E_dcfloat

Energy of the impurity states \(E_{dc}\) with a core hole

E_Lcfloat

Energy of the bath states \(E_{Lc}\) with a core hole

E_pfloat

Energy of the core hole state \(E_p\)

Notes

For the XAS states, we need to introduce core states into the Hamiltonian and account for their presence when setting energies. It is crucial to define the charge transfer energy in the same way as the intial states, as this affects the way the d and ligand states interact. We do this by imposing the same condition as before just accounting for the presence of core electrons.

We now have two arbitrary choices. The overall energy is ill-defined to within an additive factor. Since the core states are not hybridized, we can also set those to whatever energy we like (eventually, the core energy will be shifted to match the experiment). Arbitrarily, we pick the same condition as before with $d^n$ no core hole as zero, and we also set the core state energy to make $d^{n+1}$ one core hole zero energy. The relevant energies can be obtained by solving the resulting linear equations.

edrixs.utils.F0F2F4F6_to_UdJH(F0, F2, F4, F6)[source]

Given \(F_0\), \(F_2\), \(F_4\) and \(F_6\), return \(U_d\) and \(J_H\), for \(f\)-orbitals.

Parameters:
F0: float

Slater integral \(F_0\).

F2: float

Slater integral \(F_2\).

F4: float

Slater integral \(F_4\).

F6: float

Slater integral \(F_6\).

Returns:
Ud: float

Coulomb interaction \(U_d\).

JH: float

Hund’s coupling \(J_H\).

edrixs.utils.F0F2F4_to_UJ(F0, F2, F4)[source]

Given \(F_0\), \(F_2\) and \(F_4\), return \(U\) and \(J\), for \(t2g\)-orbitals.

Parameters:
F0: float

Slater integral \(F_0\).

F2: float

Slater integral \(F_2\).

F4: float

Slater integral \(F_4\).

Returns:
U: float

Coulomb interaction \(U\).

J: float

Hund’s coupling \(J\).

edrixs.utils.F0F2F4_to_UdJH(F0, F2, F4)[source]

Given \(F_0\), \(F_2\) and \(F_4\), return \(U_d\) and \(J_H\), for \(d\)-orbitals.

Parameters:
F0: float

Slater integral \(F_0\).

F2: float

Slater integral \(F_2\).

F4: float

Slater integral \(F_4\).

Returns:
Ud: float

Coulomb interaction \(U_d\).

JH: float

Hund’s coupling \(J_H\).

edrixs.utils.UJ_to_UdJH(U, J)[source]

Given Kanamori \(U\) and \(J\), return \(U_d\) and \(J_H\), for \(t2g\)-orbitals.

Parameters:
U: float

Coulomb interaction \(U\).

J: float

Hund’s coupling \(J\).

Returns:
Ud: float

Coulomb interaction \(U_d\).

JH: float

Hund’s coupling \(J_{H}\).

edrixs.utils.UdJH_to_F0F2F4(Ud, JH)[source]

Given \(U_d\) and \(J_H\), return \(F_0\), \(F_2\) and \(F_4\), for \(d\)-orbitals.

Parameters:
Ud: float

Coulomb interaction \(U_d\).

JH: float

Hund’s coupling \(J_H\).

Returns:
F0: float

Slater integral \(F_0\).

F2: float

Slater integral \(F_2\).

F4: float

Slater integral \(F_4\).

edrixs.utils.UdJH_to_F0F2F4F6(Ud, JH)[source]

Given \(U_d\) and \(J_H\), return \(F_0\), \(F_2\), \(F_4\) and \(F_6\), for \(f\)-orbitals.

Parameters:
Ud: float

Coulomb interaction \(U_d\).

JH: float

Hund’s coupling \(J_H\).

Returns:
F0: float

Slater integral \(F_0\).

F2: float

Slater integral \(F_2\).

F4: float

Slater integral \(F_4\).

F6: float

Slater integral \(F_6\).

edrixs.utils.UdJH_to_UJ(Ud, JH)[source]

Given \(U_d\) and \(J_H\), return Kanamori \(U\) and \(J\), for \(t2g\)-orbitals.

Parameters:
Ud: float

Coulomb interaction \(U_d\).

JH: float

Hund’s coupling \(J_H\).

Returns:
U: float

Coulomb interaction \(U\) in Kanamori form.

J: float

Hund’s coupling \(J\) in Kanamori form.

edrixs.utils.beta_to_kelvin(beta)[source]

Convert \(\beta\) to Kelvin.

Parameters:
beta: float

Inversion temperature.

Returns:
T: float

Temperature (K).

edrixs.utils.boltz_dist(gs, T)[source]

Return Boltzmann distributition.

Parameters:
gs: 1d float array

Energy levels.

T: float

Temperature in Kelvin.

Returns:
res: 1d float array

The Boltzmann distributition.

edrixs.utils.case_to_shell_name(case)[source]

Return the shell names for different cases.

Parameters:
case: string

A string describing the shells included.

Returns:
shell_name: tuple of one or two strings

The name of shells.

Examples

>>> import edrixs
>>> edrixs.case_to_shell_name('d')
('d',)
>>> edrixs.case_to_shell_name('t2gp32')
('t2g', 'p32')
edrixs.utils.edge_to_shell_name(edge_name, with_main_qn=False)[source]

Given edge name, return shell name. If one needs to include both spin-orbit split edges of one shell, one can use string, for example, ‘L23’ means both L2 and L3 edges are considered in the calculations, and the shell name will be ‘p’.

Parameters:
edge_name: string

Standard edge name.

with_main_qn: logical

If true, the shell name will include the main quantum number.

Returns
——-
shell_name: string

Shell name.

  • with_main_qn=True, shell_name will include the main quantum number, for example, 2p

  • with_main_qn=False, shell_name will not include the main quantum number, for example, p

edrixs.utils.get_atom_data(atom, v_name, v_noccu, edge=None, trans_to_which=1, label=None)[source]

Return Slater integrals, spin-orbit coupling strength, edge energy and core hole life-time broadening for an atom with given type of valence shells, occupancy of valence shells and the x-ray resonant edge.

The Slater integrals and spin-orbit coupling are calculated by Cowan’s code (https://github.com/mretegan/atomic-parameters) based on Hartree-Fock approximation. These numbers are just initial guess and serve as a start point, you need to rescale them to reproduce your XAS or RIXS spectra.

NOTE: F0 is not calculated by Cowan’s code, they are all set to be zero, you need to set F0 by yourself.

Parameters:
atom: string

Name of atom, for example, ‘Ni’, ‘Cu’, ‘Ir’, ‘U’.

v_name: a string or a tuple of one or two strings

Names of one or two valence shells. Set it to a single string if there is only one valence shell, or set it to a tuple of two strings if there are two valence shells. For example, v_name=’3d’, v_name=(‘3d’, ‘4p’)

v_noccu: a int or a tuple of one or two ints

Occupancy of valence shells before x-ray absorption (without a core hole). Set it to a single int if there is only one valence shell, or set it to a tuple of two ints if there are two valence shells. The order should be consistent with that in v_name. For example, v_name=’3d’, v_noccu=8, v_name=(‘3d’, ‘4p’), v_noccu=(8, 0).

edge: a string

X-ray resonant edge. If edge is not None, both the information of the initial (without a core hole) and intermediate (with a core hole) Hamiltonians will be returned, otherwise, only the information of the initial Hamiltonian will be returned. It can be,

  • K, L1, L2, L3, M1, M2, M3, M4, M5, N1, N2, N3, N4, N5, N6, N7, O1, O2, O3, O4, O5, P1, P2, P3

  • L23 (both L2 and L3), M23 (both M2 and M3), M45 (both M4 and M5) N23 (both N2 and N3), N45 (both N4 and N5), N67 (both N6 and N7) O23 (both O2 and O3), O45 (both O4 and O5), P23 (both P2 and P3)

trans_to_which: int

If there are two valence shells, this variable is used to indicate which valence shell the photon transition happens.

  • 1: to the first valence shell

  • 2: to the second valence shell

label: a string or tuple of strings

User-defined symbols to label the names of Slater integrals.

  • one element, for the valence shell

  • two elements, 1st (2nd)-element for the 1st (2nd) valence shells or 1st-element for the 1st valence shell and the 2nd one for the core shell.

  • three elements, the 1st (2nd)-element is for the 1st (2nd) valence shell, the 3rd one is for the core shell

Returns:
res: dict

Atomic data. A dict likes this if edge is None:

res={
    'slater_i': [(name of slater integrals, vaule of slater integrals), (), ...],
    'v_soc': [list of SOC for each atomic shell in v_name],
}

otherwise,

res={
    'slater_i': [(name of slater integrals, value of slater integrals), (), ...],
    'slater_n': [(name of slater integrals, value of slater integrals), (), ...],
    'v_soc': [list of SOC for each atomic shell in v_name],
    'c_soc': SOC of core shell,
    'edge_ene': [list of edge energy, 2 elements if edge is any of
    "L23, M23, M45, N23, N45, N67, O23, O45, P23", othewise only 1 element],
    'gamma_c': [list of core hole life-time broadening, 2 elements if edge is any of
    "L23, M23, M45, N23, N45, N67, O23, O45, P23", othewise only 1 element],
}

Examples

>>> import edrixs
>>> res = edrixs.get_atom_data('Ni', v_name='3d', v_noccu=8)
>>> res
{'slater_i': [('F0_11', 0.0), ('F2_11', 13.886), ('F4_11', 8.67)],
 'v_soc_i': [0.112]}
>>> name, slater = [list(i) for i in zip(*res['slater_i'])]
>>> name
['F0_11', 'F2_11', 'F4_11']
>>> slater
[0.0, 13.886, 8.67]
>>> slater = [i * 0.8 for i in slater]
>>> slater[0] = edrixs.get_F0('d', slater[1], slater[2])
>>> slater
[0.5728507936507936, 11.1088, 6.936]
>>> import edrixs
>>> res=edrixs.get_atom_data('Ni', v_name='3d', v_noccu=8, edge='L3', label=('d', 'p'))
>>> res
{'slater_i': [('F0_dd', 0.0), ('F2_dd', 12.234), ('F4_dd', 7.598)],
 'v_soc_i': [0.083],
 'slater_n': [('F0_dd', 0.0), ('F2_dd', 12.234), ('F4_dd', 7.598),
  ('F0_dp', 0.0), ('F2_dp', 7.721), ('G1_dp', 5.787), ('G3_dp', 3.291),
  ('F0_pp', 0.0), ('F2_pp', 0.0)],
 'v_soc_n': [0.102],
 'c_soc': 11.507,
 'edge_ene': [852.7],
 'gamma_c': [0.275]}
>>> name_i, slater_i = [list(i) for i in zip(*res['slater_i'])]
>>> name_n, slater_n = [list(i) for i in zip(*res['slater_n'])]
>>> name_i
['F0_dd', 'F2_dd', 'F4_dd']
>>> slater_i
[0.0, 12.234, 7.598]
>>> name_n
['F0_dd', 'F2_dd', 'F4_dd', 'F0_dp', 'F2_dp', 'G1_dp', 'G3_dp', 'F0_pp', 'F2_pp']
>>> slater_n
[0.0, 12.234, 7.598, 0.0, 7.721, 5.787, 3.291, 0.0, 0.0]
>>> slater_n[0] = edrixs.get_F0('d', slater_n[1], slater_n[2])
>>> slater_n[3] = edrixs.get_F0('dp', slater_n[5], slater_n[6])
>>> slater_n
[0.6295873015873016, 12.234, 7.598, 0.5268428571428572, 7.721, 5.787, 3.291, 0.0, 0.0]
>>> import edrixs
>>> import collections
>>> res=edrixs.get_atom_data('Ni', v_name=('3d', '4p'), v_noccu=(8, 0),
...     edge='K', trans_to_which=2, label=('d', 'p', 's'))
>>> res
{'slater_i': [('F0_dd', 0.0), ('F2_dd', 12.234), ('F4_dd', 7.598),
  ('F0_dp', 0.0), ('F2_dp', 0.0), ('G1_dp', 0.0), ('G3_dp', 0.0),
  ('F0_pp', 0.0), ('F2_pp', 0.0)],
 'v_soc_i': [0.083, 0.0],
 'slater_n': [('F0_dd', 0.0), ('F2_dd', 13.851), ('F4_dd', 8.643),
  ('F0_dp', 0.0), ('F2_dp', 2.299), ('G1_dp', 0.828), ('G3_dp', 0.713),
  ('F0_pp', 0.0), ('F2_pp', 0.0),
  ('F0_ds', 0.0), ('G2_ds', 0.079),
  ('F0_ps', 0.0), ('G1_ps', 0.194),
  ('F0_ss', 0.0)],
 'v_soc_n': [0.113, 0.093],
 'c_soc': 0.0,
 'edge_ene': [8333.0],
 'gamma_c': [0.81]}
>>> slat_i = collections.OrderedDict(res['slater_i'])
>>> slat_n = collections.OrderedDict(res['slater_n'])
>>> list(slat_i.keys())
['F0_dd', 'F2_dd', 'F4_dd', 'F0_dp', 'F2_dp', 'G1_dp', 'G3_dp', 'F0_pp', 'F2_pp']
>>> list(slat_i.values())
[0.0, 12.234, 7.598, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0]
>>> list(slat_n.keys())
['F0_dd', 'F2_dd', 'F4_dd', 'F0_dp', 'F2_dp', 'G1_dp', 'G3_dp', 'F0_pp', 'F2_pp',
 'F0_ds', 'G2_ds', 'F0_ps', 'G1_ps', 'F0_ss']
>>> list(slat_n.values())
[0.0, 13.851, 8.643, 0.0, 2.299, 0.828, 0.713, 0.0, 0.0,
 0.0, 0.079, 0.0, 0.194, 0.0]
edrixs.utils.info_atomic_shell()[source]

Return a dict to describe the information of atomic shell. The key is a string of the shell’s name, the value is a 2-elments tuple, the first value is the orbital angular moment number and the second value is the dimension of the Hilbert space.

Returns:
info: dict

The dict describing info of atomic shell.

edrixs.utils.kelvin_to_beta(k)[source]

Convert temperature from Kelvin to \(\beta\).

Parameters:
k: float

Temperature in Kelvin.

Returns:
beta: float

Inversion temperature.

edrixs.utils.rescale(old_list, scale=None)[source]

Rescale a 1d list.

Parameters:
old_list: 1d list of numerical values

Input list.

scale: a list of tuples

The scaling factors.

Returns:
new_list: 1d list

The rescaled list.

edrixs.utils.slater_integrals_name(shell_name, label=None)[source]

Given shell names, return the names of required Slater integrals. The order of these names in the return list is according to the convention:

  • For 1 shell: [FX_11]

  • For 2 shells: [FX_11, FX_12, GX_12, FX_22]

  • For 3 shells: [FX_11, FX_12, GX_12, FX_22, FX_13, GX_13, FX_23, GX_23, FX_33]

where, X=0, 2, 4, 6, or X=1, 3, 5, X shoule be in ascending order.

Parameters:
shell_name: tuple of strings

Name of shells. Its length should be less and equal than 3.

label: tuple of strings

Label of shells, same shape as shell_name.

If not provided, label will be set to

  • label=(1,) for one shell, or

  • label=(1,2) for two shells, or

  • label=(1,2,3) for three shells

Returns:
res: list of strings

Names of Slater integrals.