solvers

edrixs.solvers.ed_1v1c_fort(comm, shell_name, *, shell_level=None, v_soc=None, c_soc=0, v_noccu=1, slater=None, ext_B=None, on_which='spin', v_cfmat=None, v_othermat=None, do_ed=True, ed_solver=2, neval=1, nvector=1, ncv=3, idump=False, maxiter=500, eigval_tol=1e-08, min_ndim=1000)[source]

Perform ED for the case of one valence shell plus one Core-shell with Fortran ED solver.

The hopping and Coulomb terms of both the initial and intermediate Hamiltonians will be constructed and written to files (hopping_i.in, hopping_n.in, coulomb_i.in and coulomb_n.in). Fock basis for the initial Hamiltonian will be written to file (fock_i.in).

ED will be only performed on the initial Hamiltonian to find a few lowest eigenstates do_ed=True. Only input files will be written if do_ed=False. Due to large Hilbert space, the ed_fsolver written in Fortran will be called. mpi4py and a MPI environment (mpich or openmpi) are required to launch ed_fsolver.

If do_ed=True, it will output the eigenvalues in file (eigvals.dat) and eigenvectors in files (eigvec.n), where n means the n-th eigenvectors. The eigvec.n files will be used later as the inputs for XAS and RIXS calculations.

Parameters:
comm: MPI_comm

The MPI communicator from mpi4py.

shell_name: tuple of two strings

Names of valence and core shells. The 1st (2nd) string in the tuple is for the valence (core) shell.

  • The 1st strings can only be ‘s’, ‘p’, ‘t2g’, ‘d’, ‘f’

  • The 2nd string can be ‘s’, ‘p’, ‘p12’, ‘p32’, ‘d’, ‘d32’, ‘d52’, ‘f’, ‘f52’, ‘f72’.

For example: shell_name=(‘d’, ‘p32’) may indicate a \(L_3\) edge transition from core \(2p_{3/2}\) shell to valence \(3d\) shell for Ni.

shell_level: tuple of two float numbers

Energy level of valence (1st element) and core (2nd element) shells.

They will be set to zero if not provided.

v_soc: tuple of two float numbers

Spin-orbit coupling strength of the valence shell, v1_soc[0] for the initial Hamiltonian, and v1_soc[1] for the intermediate Hamiltonian.

They will be set to zero if not provided.

c_soc: float number

Spin-orbit coupling strength of core electrons.

v_noccu: int number

Total number of electrons in valence shells.

slater: tuple of two lists

Slater integrals for initial (1st list) and intermediate (2nd list) Hamiltonians. The order of the elements in each list should be like this:

[FX_vv, FX_vc, GX_vc, FX_cc],

where X are integers with ascending order, it can be X=0, 2, 4, 6 or X=1, 3, 5. One can ignore all the continuous zeros at the end of the list.

For example, if the full list is: [F0_dd, F2_dd, F4_dd, 0, F2_dp, 0, 0, 0, 0], one can just provide [F0_dd, F2_dd, F4_dd, 0, F2_dp]

All the Slater integrals will be set to zero if slater==None.

ext_B: tuple of three float numbers

Vector of external magnetic field with respect to global \(xyz\)-axis applied on the valence shell.

It will be set to zeros if not provided.

on_which: string

Apply Zeeman exchange field on which sector. Options are ‘spin’, ‘orbital’ or ‘both’.

v_cfmat: 2d complex array

Crystal field splitting Hamiltonian of the valence shell. The dimension and the orbital order should be consistent with the type of the valence shell.

They will be zeros if not provided.

v_othermat: 2d complex array

Other possible Hamiltonian of the valence shell. The dimension and the orbital order should be consistent with the type of the valence shell.

They will be zeros if not provided.

do_ed: logical

If do_end=True, diagonalize the Hamitlonian to find a few lowest eigenstates, return the eigenvalues and density matirx, and write the eigenvectors in files eigvec.n, otherwise, just write out the input files, do not perform the ED.

ed_solver: int

Type of ED solver, options can be 0, 1, 2

  • 0: use Lapack to fully diagonalize Hamiltonian to get all the eigenvalues.

  • 1: use standard Lanczos algorithm to find only a few lowest eigenvalues, no re-orthogonalization has been applied, so it is not very accurate.

  • 2: use parallel version of Arpack library to find a few lowest eigenvalues, it is accurate and is the recommeded choice in real calculations of XAS and RIXS.

neval: int

Number of eigenvalues to be found. For ed_solver=2, the value should not be too small, neval > 10 is usually a safe value.

nvector: int

Number of eigenvectors to be found and written into files.

ncv: int

Used for ed_solver=2, it should be at least ncv > neval + 2. Usually, set it a little bit larger than neval, for example, set ncv=200 when neval=100.

idump: logical

Whether to dump the eigenvectors to files “eigvec.n”, where n means the n-th vectors.

maxiter: int

Maximum number of iterations in finding all the eigenvalues, used for ed_solver=1, 2.

eigval_tol: float

The convergence criteria of eigenvalues, used for ed_solver=1, 2.

min_ndim: int

The minimum dimension of the Hamiltonian when the ed_solver=1, 2 can be used, otherwise, ed_solver=1 will be used.

Returns:
eval_i: 1d float array, shape=(neval, )

The eigenvalues of initial Hamiltonian.

denmat: 2d complex array, shape=(nvector, v_norb, v_norb))

The density matrix in the eigenstates.

edrixs.solvers.ed_1v1c_py(shell_name, *, shell_level=None, v_soc=None, c_soc=0, v_noccu=1, slater=None, ext_B=None, on_which='spin', v_cfmat=None, v_othermat=None, loc_axis=None, verbose=0)[source]

Perform ED for the case of two atomic shells, one valence plus one Core shell with pure Python solver. For example, for Ni-\(L_3\) edge RIXS, they are 3d valence and 2p core shells.

It will use scipy.linalag.eigh to exactly diagonalize both the initial and intermediate Hamiltonians to get all the eigenvalues and eigenvectors, and the transition operators will be built in the many-body eigenvector basis.

This solver is only suitable for small size of Hamiltonian, typically the dimension of both initial and intermediate Hamiltonian are smaller than 10,000.

Parameters:
shell_name: tuple of two strings

Names of valence and core shells. The 1st (2nd) string in the tuple is for the valence (core) shell.

  • The 1st string can only be ‘s’, ‘p’, ‘t2g’, ‘d’, ‘f’,

  • The 2nd string can be ‘s’, ‘p’, ‘p12’, ‘p32’, ‘d’, ‘d32’, ‘d52’, ‘f’, ‘f52’, ‘f72’.

For example: shell_name=(‘d’, ‘p32’) indicates a \(L_3\) edge transition from core \(p_{3/2}\) shell to valence \(d\) shell.

shell_level: tuple of two float numbers

Energy level of valence (1st element) and core (2nd element) shells.

They will be set to zero if not provided.

v_soc: tuple of two float numbers

Spin-orbit coupling strength of valence electrons, for the initial (1st element) and intermediate (2nd element) Hamiltonians.

They will be set to zero if not provided.

c_soc: a float number

Spin-orbit coupling strength of core electrons.

v_noccu: int number

Number of electrons in valence shell.

slater: tuple of two lists

Slater integrals for initial (1st list) and intermediate (2nd list) Hamiltonians. The order of the elements in each list should be like this:

[FX_vv, FX_vc, GX_vc, FX_cc],

where X are integers with ascending order, it can be X=0, 2, 4, 6 or X=1, 3, 5. One can ignore all the continuous zeros at the end of the list.

For example, if the full list is: [F0_dd, F2_dd, F4_dd, 0, F2_dp, 0, 0, 0, 0], one can just provide [F0_dd, F2_dd, F4_dd, 0, F2_dp]

All the Slater integrals will be set to zero if slater=None.

ext_B: tuple of three float numbers

Vector of external magnetic field with respect to global \(xyz\)-axis.

They will be set to zero if not provided.

on_which: string

Apply Zeeman exchange field on which sector. Options are ‘spin’, ‘orbital’ or ‘both’.

v_cfmat: 2d complex array

Crystal field splitting Hamiltonian of valence electrons. The dimension and the orbital order should be consistent with the type of valence shell.

They will be zeros if not provided.

v_othermat: 2d complex array

Other possible Hamiltonian of valence electrons. The dimension and the orbital order should be consistent with the type of valence shell.

They will be zeros if not provided.

loc_axis: 3*3 float array

The local axis with respect to which local orbitals are defined.

  • x: local_axis[:,0],

  • y: local_axis[:,1],

  • z: local_axis[:,2].

It will be an identity matrix if not provided.

verbose: int

Level of writting data to files. Hopping matrices, Coulomb tensors, eigvenvalues will be written if verbose > 0.

Returns:
eval_i: 1d float array

The eigenvalues of initial Hamiltonian.

eval_n: 1d float array

The eigenvalues of intermediate Hamiltonian.

trans_op: 3d complex array

The matrices of transition operators in the eigenvector basis. Their components are defined with respect to the global \(xyz\)-axis.

edrixs.solvers.ed_2v1c_fort(comm, shell_name, *, shell_level=None, v1_soc=None, v2_soc=None, c_soc=0, v_tot_noccu=1, slater=None, v1_ext_B=None, v2_ext_B=None, v1_on_which='spin', v2_on_which='spin', v1_cfmat=None, v2_cfmat=None, v1_othermat=None, v2_othermat=None, hopping_v1v2=None, do_ed=True, ed_solver=2, neval=1, nvector=1, ncv=3, idump=False, maxiter=500, eigval_tol=1e-08, min_ndim=1000)[source]

Perform ED for the case of two valence shell plus one core-shell with Fortran solver. For example, for Ni \(K\)-edge RIXS, \(1s\rightarrow 4p\) transition, the valence shells involved in RIXS are \(3d\) and \(4p\).

The hopping and Coulomb terms of both the initial and intermediate Hamiltonians will be constructed and written to files (hopping_i.in, hopping_n.in, coulomb_i.in and coulomb_n.in). Fock basis for the initial Hamiltonian will be written to file (fock_i.in).

ED will be only performed on the initial Hamiltonian to find a few lowest eigenstates do_ed=True. Only input files will be written if do_ed=False. Due to large Hilbert space, the ed_fsolver written in Fortran will be called. mpi4py and a MPI environment (mpich or openmpi) are required to launch ed_fsolver.

If do_ed=True, it will output the eigenvalues in file (eigvals.dat) and eigenvectors in files (eigvec.n), where n means the n-th eigenvectors. The eigvec.n files will be used later as the inputs for XAS and RIXS calculations.

Parameters:
comm: MPI_comm

The MPI communicator from mpi4py.

shell_name: tuple of three strings

Names of valence and core shells. The 1st (2nd) string in the tuple is for the 1st (2nd) valence shell, and the 3rd one is for the core shell.

  • The 1st and 2nd strings can only be ‘s’, ‘p’, ‘t2g’, ‘d’, ‘f’

  • The 3nd string can be ‘s’, ‘p’, ‘p12’, ‘p32’, ‘d’, ‘d32’, ‘d52’, ‘f’, ‘f52’, ‘f72’.

For example: shell_name=(‘d’, ‘p’, ‘s’) may indicate a \(K\) edge transition from core \(1s\) shell to valence \(3d\) and \(4p\) shell for Ni.

shell_level: tuple of three float numbers

Energy level of valence (1st and 2nd elements) and core (3nd element) shells.

They will be set to zero if not provided.

v1_soc: tuple of two float numbers

Spin-orbit coupling strength of the 1st valence shell, v1_soc[0] for the initial Hamiltonian, and v1_soc[1] for the intermediate Hamiltonian.

They will be set to zero if not provided.

v2_soc: tuple of two float numbers

Spin-orbit coupling strength of the 2nd valence shell, v2_soc[0] for the initial Hamiltonian, and v2_soc[1] for the intermediate Hamiltonian.

They will be set to zero if not provided.

c_soc: float number

Spin-orbit coupling strength of core electrons.

v_tot_noccu: int number

Total number of electrons in valence shells.

slater: tuple of two lists

Slater integrals for initial (1st list) and intermediate (2nd list) Hamiltonians. The order of the elements in each list should be like this:

[FX_v1v1, FX_v1v2, GX_v1v2, FX_v2v2, FX_v1c, GX_v1c, FX_v2c, GX_v2c],

where X are integers with ascending order, it can be X=0, 2, 4, 6 or X=1, 3, 5. One can ignore all the continuous zeros at the end of the list.

For example, if the full list is: [F0_dd, F2_dd, F4_dd, 0, F2_dp, 0, 0, 0, 0], one can just provide [F0_dd, F2_dd, F4_dd, 0, F2_dp]

All the Slater integrals will be set to zero if slater==None.

v1_ext_B: tuple of three float numbers

Vector of external magnetic field with respect to global \(xyz\)-axis applied on the 1st valence shell.

It will be set to zeros if not provided.

v2_ext_B: tuple of three float numbers

Vector of external magnetic field with respect to global \(xyz\)-axis applied on the 2nd valence shell.

It will be set to zeros if not provided.

v1_on_which: string

Apply Zeeman exchange field on which sector. Options are ‘spin’, ‘orbital’ or ‘both’. For the 1st valence shell.

v2_on_which: string

Apply Zeeman exchange field on which sector. Options are ‘spin’, ‘orbital’ or ‘both’. For the 2nd valence shell.

v1_cfmat: 2d complex array

Crystal field splitting Hamiltonian of the 1st valence shell. The dimension and the orbital order should be consistent with the type of the 1st valence shell.

They will be zeros if not provided.

v2_cfmat: 2d complex array

Crystal field splitting Hamiltonian of the 2nd valence shell. The dimension and the orbital order should be consistent with the type of the 2nd valence shell.

They will be zeros if not provided.

v1_othermat: 2d complex array

Other possible Hamiltonian of the 1st valence shell. The dimension and the orbital order should be consistent with the type of the 1st valence shell.

They will be zeros if not provided.

v2_othermat: 2d complex array

Other possible Hamiltonian of the 2nd valence shell. The dimension and the orbital order should be consistent with the type of the 2nd valence shell.

They will be zeros if not provided.

hopping_v1v2: 2d complex array

Hopping between the two valence shells. The 1st-index (2nd-index) is the 1st (2nd) valence shell.

They will be zeros if not provided.

do_ed: logical

If do_end=True, diagonalize the Hamitlonian to find a few lowest eigenstates, return the eigenvalues and density matirx, and write the eigenvectors in files eigvec.n, otherwise, just write out the input files, do not perform the ED.

ed_solver: int

Type of ED solver, options can be 0, 1, 2

  • 0: use Lapack to fully diagonalize Hamiltonian to get all the eigenvalues.

  • 1: use standard Lanczos algorithm to find only a few lowest eigenvalues, no re-orthogonalization has been applied, so it is not very accurate.

  • 2: use parallel version of Arpack library to find a few lowest eigenvalues, it is accurate and is the recommeded choice in real calculations of XAS and RIXS.

neval: int

Number of eigenvalues to be found. For ed_solver=2, the value should not be too small, neval > 10 is usually a safe value.

nvector: int

Number of eigenvectors to be found and written into files.

ncv: int

Used for ed_solver=2, it should be at least ncv > neval + 2. Usually, set it a little bit larger than neval, for example, set ncv=200 when neval=100.

idump: logical

Whether to dump the eigenvectors to files “eigvec.n”, where n means the n-th vectors.

maxiter: int

Maximum number of iterations in finding all the eigenvalues, used for ed_solver=1, 2.

eigval_tol: float

The convergence criteria of eigenvalues, used for ed_solver=1, 2.

min_ndim: int

The minimum dimension of the Hamiltonian when the ed_solver=1, 2 can be used, otherwise, ed_solver=1 will be used.

Returns:
eval_i: 1d float array, shape=(neval, )

The eigenvalues of initial Hamiltonian.

denmat: 2d complex array, shape=(nvector, v1v2_norb, v1v2_norb))

The density matrix in the eigenstates.

edrixs.solvers.ed_siam_fort(comm, shell_name, nbath, *, siam_type=0, v_noccu=1, static_core_pot=0, c_level=0, c_soc=0, trans_c2n=None, imp_mat=None, imp_mat_n=None, bath_level=None, bath_level_n=None, hyb=None, hyb_n=None, hopping=None, hopping_n=None, slater=None, ext_B=None, on_which='spin', do_ed=0, ed_solver=2, neval=1, nvector=1, ncv=3, idump=False, maxiter=1000, eigval_tol=1e-08, min_ndim=1000)[source]

Find the ground state of the initial Hamiltonian of a Single Impuirty Anderson Model (SIAM), and also prepare input files, hopping_i.in, hopping_n.in, coulomb_i.in, coulomb_n.in for following XAS and RIXS calculations.

Parameters:
comm: MPI_Comm

MPI Communicator

shell_name: tuple of two strings

Names of valence and core shells. The 1st (2nd) string in the tuple is for the valence (core) shell.

  • The 1st string can only be ‘s’, ‘p’, ‘t2g’, ‘d’, ‘f’,

  • The 2nd string can be ‘s’, ‘p’, ‘p12’, ‘p32’, ‘d’, ‘d32’, ‘d52’, ‘f’, ‘f52’, ‘f72’.

For example: shell_name=(‘d’, ‘p32’) indicates a \(L_3\) edge transition from core \(p_{3/2}\) shell to valence \(d\) shell.

nbath: int

Number of bath sites.

siam_type: int

Type of SIAM Hamiltonian,

  • 0: diagonal hybridization function, parameterized by imp_mat, bath_level and hyb

  • 1: general hybridization function, parameterized by matrix hopping

if siam_type=0, only imp_mat, bath_level and hyb are required, if siam_type=1, only hopping is required.

v_noccu: int

Number of total occupancy of impurity and baths orbitals, required when do_ed=1, 2

static_core_pot: float

Static core hole potential.

c_level: float

Energy level of core shell.

c_soc: float

Spin-orbit coupling strength of core electrons.

trans_c2n: 2d complex array

The transformation matrix from the spherical harmonics basis to the basis on which the imp_mat and hybridization function (bath_level, hyb, hopping) are defined.

imp_mat: 2d complex array

Impurity matrix for the impurity site, including CF or SOC, for siam_type=0 and the initial configurations.

imp_mat_n: 2d complex array

Impurity matrix for the impurity site, including CF or SOC, for siam_type=0 and the intermediate configurations. If imp_mat_n=None, imp_mat will be used.

bath_level: 2d complex array

Energy level of bath sites, 1st (2nd) dimension is for different bath sites (orbitals), for siam_type=0 and the initial configurations.

bath_level_n: 2d complex array

Energy level of bath sites, 1st (2nd) dimension is for different bath sites (orbitals), for siam_type=0 and the intermediate configurations. If bath_level_n=None, bath_level will be used.

hyb: 2d complex array

Hybridization strength of bath sites, 1st (2nd) dimension is for different bath sites (orbitals), for siam_type=0 and the initial configurations.

hyb_n: 2d complex array

Hybridization strength of bath sites, 1st (2nd) dimension is for different bath sites (orbitals), for siam_type=0 and the intermediate configurations. If hyb_n=None, hyb will be used.

hopping: 2d complex array

General hopping matrix when siam_type=1, including imp_mat and hybridization functions, for siam_type=1 and the initial configurations.

hopping_n: 2d complex array

General hopping matrix when siam_type=1, including imp_mat and hybridization functions, for siam_type=1 and the intermediate configurations. If hopping_n=None, hopping will be used.

slater: tuple of two lists

Slater integrals for initial (1st list) and intermediate (2nd list) Hamiltonians. The order of the elements in each list should be like this:

[FX_vv, FX_vc, GX_vc, FX_cc],

where X are integers with ascending order, it can be X=0, 2, 4, 6 or X=1, 3, 5. One can ignore all the continuous zeros at the end of the list.

For example, if the full list is: [F0_dd, F2_dd, F4_dd, 0, F2_dp, 0, 0, 0, 0], one can just provide [F0_dd, F2_dd, F4_dd, 0, F2_dp]

All the Slater integrals will be set to zero if slater=None.

ext_B: tuple of three float numbers

Vector of external magnetic field with respect to global \(xyz\)-axis.

They will be set to zero if not provided.

on_which: string

Apply Zeeman exchange field on which sector. Options are ‘spin’, ‘orbital’ or ‘both’.

do_ed: int
  • 0: First, search the ground state in different subspaces of total occupancy \(N\) with ed_solver=1, and then do a more accurate ED in the subspace \(N\) where the ground state lies to find a few lowest eigenstates, return the eigenvalues and density matirx, and write the eigenvectors in files eigvec.n

  • 1: Only do ED for given occupancy number v_noccu, return eigenvalues and density matrix, write eigenvectors to files eigvec.n

  • 2: Do not do ED, only write parameters into files: hopping_i.in, hopping_n.in, coulomb_i.in, coulomb_n.in for later XAS or RIXS calculations.

ed_solver: int

Type of ED solver, options can be 0, 1, 2

  • 0: use Lapack to fully diagonalize Hamiltonian to get all the eigenvalues.

  • 1: use standard Lanczos algorithm to find only a few lowest eigenvalues, no re-orthogonalization has been applied, so it is not very accurate.

  • 2: use parallel version of Arpack library to find a few lowest eigenvalues, it is accurate and is the recommeded choice in real calculations of XAS and RIXS.

neval: int

Number of eigenvalues to be found. For ed_solver=2, the value should not be too small, neval > 10 is usually a safe value.

nvector: int

Number of eigenvectors to be found and written into files.

ncv: int

Used for ed_solver=2, it should be at least ncv > neval + 2. Usually, set it a little bit larger than neval, for example, set ncv=200 when neval=100.

idump: logical

Whether to dump the eigenvectors to files “eigvec.n”, where n means the n-th vectors.

maxiter: int

Maximum number of iterations in finding all the eigenvalues, used for ed_solver=1, 2.

eigval_tol: float

The convergence criteria of eigenvalues, used for ed_solver=1, 2.

min_ndim: int

The minimum dimension of the Hamiltonian when the ed_solver=1, 2 can be used, otherwise, ed_solver=1 will be used.

Returns:
eval_i: 1d float array

Eigenvalues of initial Hamiltonian.

denmat: 2d complex array

Density matrix.

noccu_gs: int

Occupancy of the ground state.

edrixs.solvers.rixs_1v1c_fort(comm, shell_name, ominc, eloss, *, gamma_c=0.1, gamma_f=0.1, v_noccu=1, thin=1.0, thout=1.0, phi=0, pol_type=None, num_gs=1, nkryl=200, linsys_max=500, linsys_tol=1e-08, temperature=1.0, loc_axis=None, scatter_axis=None)[source]

Calculate RIXS for the case with one valence shell plus one core shell with Fortran solver.

Parameters:
comm: MPI_comm

MPI communicator.

shell_name: tuple of two strings

Names of valence and core shells. The 1st (2nd) string in the tuple is for the valence (core) shell.

  • The 1st string can only be ‘s’, ‘p’, ‘t2g’, ‘d’, ‘f’,

  • The 2nd string can be ‘s’, ‘p’, ‘p12’, ‘p32’, ‘d’, ‘d32’, ‘d52’, ‘f’, ‘f52’, ‘f72’.

For example: shell_name=(‘d’, ‘p32’) may indicate a \(L_3\) edge transition from core \(2p_{3/2}\) shell to valence \(3d\) shell for Ni.

ominc: 1d float array

Incident energy of photon.

eloss: 1d float array

Energy loss.

gamma_c: a float number or a 1d float array with same shape as ominc.

The core-hole life-time broadening factor. It can be a constant value or incident energy dependent.

gamma_f: a float number or a 1d float array with same shape as eloss.

The final states life-time broadening factor. It can be a constant value or energy loss dependent.

v_noccu: int

Total occupancy of valence shells.

thin: float number

The incident angle of photon (in radian).

thout: float number

The scattered angle of photon (in radian).

phi: float number

Azimuthal angle (in radian), defined with respect to the \(x\)-axis of scattering axis: scatter_axis[:,0].

pol_type: list of 4-elements-tuples

Type of polarizations. It has the following form:

(str1, alpha, str2, beta)

where, str1 (str2) can be ‘linear’, ‘left’, ‘right’, and alpha (beta) is the angle (in radians) between the linear polarization vector and the scattering plane.

It will set pol_type=[(‘linear’, 0, ‘linear’, 0)] if not provided.

num_gs: int

Number of initial states used in RIXS calculations.

nkryl: int

Maximum number of poles obtained.

linsys_max: int

Maximum iterations of solving linear equations.

linsys_tol: float

Convergence for solving linear equations.

temperature: float number

Temperature (in K) for boltzmann distribution.

loc_axis: 3*3 float array

The local axis with respect to which local orbitals are defined.

  • x: local_axis[:,0],

  • y: local_axis[:,1],

  • z: local_axis[:,2].

It will be an identity matrix if not provided.

scatter_axis: 3*3 float array

The local axis defining the scattering geometry. The scattering plane is defined in the local \(zx\)-plane.

  • local \(x\)-axis: scatter_axis[:,0]

  • local \(y\)-axis: scatter_axis[:,1]

  • local \(z\)-axis: scatter_axis[:,2]

It will be set to an identity matrix if not provided.

Returns:
rixs: 3d float array, shape=(len(ominc), len(eloss), len(pol_type))

The calculated RIXS spectra. The 1st dimension is for the incident energy, the 2nd dimension is for the energy loss and the 3rd dimension is for different polarizations.

poles: 2d list of dict, shape=(len(ominc), len(pol_type))

The calculated RIXS poles. The 1st dimension is for incident energy, and the 2nd dimension is for different polarizations.

edrixs.solvers.rixs_1v1c_py(eval_i, eval_n, trans_op, ominc, eloss, *, gamma_c=0.1, gamma_f=0.01, thin=1.0, thout=1.0, phi=0.0, pol_type=None, gs_list=None, temperature=1.0, scatter_axis=None)[source]

Calculate RIXS for the case of one valence shell plus one core shell with Python solver.

This solver is only suitable for small size of Hamiltonian, typically the dimension of both initial and intermediate Hamiltonian are smaller than 10,000.

Parameters:
eval_i: 1d float array

The eigenvalues of the initial Hamiltonian.

eval_n: 1d float array

The eigenvalues of the intermediate Hamiltonian.

trans_op: 3d complex array

The transition operators in the eigenstates basis.

ominc: 1d float array

Incident energy of photon.

eloss: 1d float array

Energy loss.

gamma_c: a float number or a 1d float array with same shape as ominc.

The core-hole life-time broadening factor. It can be a constant value or incident energy dependent.

gamma_f: a float number or a 1d float array with same shape as eloss.

The final states life-time broadening factor. It can be a constant value or energy loss dependent.

thin: float number

The incident angle of photon (in radian).

thout: float number

The scattered angle of photon (in radian).

phi: float number

Azimuthal angle (in radian), defined with respect to the \(x\)-axis of scattering axis: scatter_axis[:,0].

pol_type: list of 4-elements-tuples

Type of polarizations. It has the following form:

(str1, alpha, str2, beta)

where, str1 (str2) can be ‘linear’, ‘left’, ‘right’, and alpha (beta) is the angle (in radians) between the linear polarization vector and the scattering plane.

It will set pol_type=[(‘linear’, 0, ‘linear’, 0)] if not provided.

gs_list: 1d list of ints

The indices of initial states which will be used in RIXS calculations.

It will set gs_list=[0] if not provided.

temperature: float number

Temperature (in K) for boltzmann distribution.

scatter_axis: 3*3 float array

The local axis defining the scattering plane. The scattering plane is defined in the local \(zx\)-plane.

  • local \(x\)-axis: scatter_axis[:,0]

  • local \(y\)-axis: scatter_axis[:,1]

  • local \(z\)-axis: scatter_axis[:,2]

It will be an identity matrix if not provided.

Returns:
rixs: 3d float array

The calculated RIXS spectra. The 1st dimension is for the incident energy, the 2nd dimension is for the energy loss and the 3rd dimension is for different polarizations.

edrixs.solvers.rixs_2v1c_fort(comm, shell_name, ominc, eloss, *, gamma_c=0.1, gamma_f=0.1, v_tot_noccu=1, trans_to_which=1, thin=1.0, thout=1.0, phi=0, pol_type=None, num_gs=1, nkryl=200, linsys_max=500, linsys_tol=1e-08, temperature=1.0, loc_axis=None, scatter_axis=None)[source]

Calculate RIXS for the case with 2 valence shells plus 1 core shell.

Parameters:
comm: MPI_comm

MPI communicator.

shell_name: tuple of three strings

Names of valence and core shells. The 1st (2nd) string in the tuple is for the 1st (2nd) valence shell, and the 3rd one is for the core shell.

  • The 1st and 2nd strings can only be ‘s’, ‘p’, ‘t2g’, ‘d’, ‘f’,

  • The 3nd string can be ‘s’, ‘p’, ‘p12’, ‘p32’, ‘d’, ‘d32’, ‘d52’, ‘f’, ‘f52’, ‘f72’.

For example: shell_name=(‘d’, ‘p’, ‘s’) may indicate a \(K\) edge transition from core \(1s\) shell to valence \(3d\) and \(4p\) shell for Ni.

ominc: 1d float array

Incident energy of photon.

eloss: 1d float array

Energy loss.

gamma_c: a float number or a 1d float array with same shape as ominc.

The core-hole life-time broadening factor. It can be a constant value or incident energy dependent.

gamma_f: a float number or a 1d float array with same shape as eloss.

The final states life-time broadening factor. It can be a constant value or energy loss dependent.

v_tot_noccu: int

Total occupancy of valence shells.

trans_to_which: int

Photon transition to which valence shell.

  • 1: to 1st valence shell,

  • 2: to 2nd valence shell.

thin: float number

The incident angle of photon (in radian).

thout: float number

The scattered angle of photon (in radian).

phi: float number

Azimuthal angle (in radian), defined with respect to the \(x\)-axis of scattering axis: scatter_axis[:,0].

pol_type: list of 4-elements-tuples

Type of polarizations. It has the following form:

(str1, alpha, str2, beta)

where, str1 (str2) can be ‘linear’, ‘left’, ‘right’, and alpha (beta) is the angle (in radians) between the linear polarization vector and the scattering plane.

It will set pol_type=[(‘linear’, 0, ‘linear’, 0)] if not provided.

num_gs: int

Number of initial states used in RIXS calculations.

nkryl: int

Maximum number of poles obtained.

linsys_max: int

Maximum iterations of solving linear equations.

linsys_tol: float

Convergence for solving linear equations.

temperature: float number

Temperature (in K) for boltzmann distribution.

loc_axis: 3*3 float array

The local axis with respect to which local orbitals are defined.

  • x: local_axis[:,0],

  • y: local_axis[:,1],

  • z: local_axis[:,2].

It will be an identity matrix if not provided.

scatter_axis: 3*3 float array

The local axis defining the scattering geometry. The scattering plane is defined in the local \(zx\)-plane.

  • local \(x\)-axis: scatter_axis[:,0]

  • local \(y\)-axis: scatter_axis[:,1]

  • local \(z\)-axis: scatter_axis[:,2]

It will be set to an identity matrix if not provided.

Returns:
rixs: 3d float array, shape=(len(ominc), len(eloss), len(pol_type))

The calculated RIXS spectra. The 1st dimension is for the incident energy, the 2nd dimension is for the energy loss and the 3rd dimension is for different polarizations.

poles: 2d list of dict, shape=(len(ominc), len(pol_type))

The calculated RIXS poles. The 1st dimension is for incident energy, and the 2nd dimension is for different polarizations.

edrixs.solvers.rixs_siam_fort(comm, shell_name, nbath, ominc, eloss, *, gamma_c=0.1, gamma_f=0.1, v_noccu=1, thin=1.0, thout=1.0, phi=0, pol_type=None, num_gs=1, nkryl=200, linsys_max=1000, linsys_tol=1e-10, temperature=1.0, loc_axis=None, scatter_axis=None)[source]

Calculate RIXS for single impurity Anderson model with Fortran solver.

Parameters:
comm: MPI_comm

MPI communicator.

shell_name: tuple of two strings

Names of valence and core shells. The 1st (2nd) string in the tuple is for the valence (core) shell.

  • The 1st string can only be ‘s’, ‘p’, ‘t2g’, ‘d’, ‘f’,

  • The 2nd string can be ‘s’, ‘p’, ‘p12’, ‘p32’, ‘d’, ‘d32’, ‘d52’, ‘f’, ‘f52’, ‘f72’.

For example: shell_name=(‘d’, ‘p32’) may indicate a \(L_3\) edge transition from core \(2p_{3/2}\) shell to valence \(3d\) shell for Ni.

nbath: int

Number of bath sites.

ominc: 1d float array

Incident energy of photon.

eloss: 1d float array

Energy loss.

gamma_c: a float number or a 1d float array with same shape as ominc.

The core-hole life-time broadening factor. It can be a constant value or incident energy dependent.

gamma_f: a float number or a 1d float array with same shape as eloss.

The final states life-time broadening factor. It can be a constant value or energy loss dependent.

v_noccu: int

Total occupancy of valence shells.

thin: float number

The incident angle of photon (in radian).

thout: float number

The scattered angle of photon (in radian).

phi: float number

Azimuthal angle (in radian), defined with respect to the \(x\)-axis of scattering axis: scatter_axis[:,0].

pol_type: list of 4-elements-tuples

Type of polarizations. It has the following form:

(str1, alpha, str2, beta)

where, str1 (str2) can be ‘linear’, ‘left’, ‘right’, and alpha (beta) is the angle (in radians) between the linear polarization vector and the scattering plane.

It will set pol_type=[(‘linear’, 0, ‘linear’, 0)] if not provided.

num_gs: int

Number of initial states used in RIXS calculations.

nkryl: int

Maximum number of poles obtained.

linsys_max: int

Maximum iterations of solving linear equations.

linsys_tol: float

Convergence for solving linear equations.

temperature: float number

Temperature (in K) for boltzmann distribution.

loc_axis: 3*3 float array

The local axis with respect to which local orbitals are defined.

  • x: local_axis[:,0],

  • y: local_axis[:,1],

  • z: local_axis[:,2].

It will be an identity matrix if not provided.

scatter_axis: 3*3 float array

The local axis defining the scattering geometry. The scattering plane is defined in the local \(zx\)-plane.

  • local \(x\)-axis: scatter_axis[:,0]

  • local \(y\)-axis: scatter_axis[:,1]

  • local \(z\)-axis: scatter_axis[:,2]

It will be set to an identity matrix if not provided.

Returns:
rixs: 3d float array, shape=(len(ominc), len(eloss), len(pol_type))

The calculated RIXS spectra. The 1st dimension is for the incident energy, the 2nd dimension is for the energy loss and the 3rd dimension is for different polarizations.

poles: 2d list of dict, shape=(len(ominc), len(pol_type))

The calculated RIXS poles. The 1st dimension is for incident energy, and the 2nd dimension is for different polarizations.

edrixs.solvers.xas_1v1c_fort(comm, shell_name, ominc, *, gamma_c=0.1, v_noccu=1, thin=1.0, phi=0, pol_type=None, num_gs=1, nkryl=200, temperature=1.0, loc_axis=None, scatter_axis=None)[source]

Calculate XAS for the case with one valence shells plus one core shell with Fortran solver.

Parameters:
comm: MPI_comm

MPI communicator.

shell_name: tuple of two strings

Names of valence and core shells. The 1st (2nd) string in the tuple is for the valence (core) shell.

  • The 1st string can only be ‘s’, ‘p’, ‘t2g’, ‘d’, ‘f’,

  • The 2nd string can be ‘s’, ‘p’, ‘p12’, ‘p32’, ‘d’, ‘d32’, ‘d52’, ‘f’, ‘f52’, ‘f72’.

For example: shell_name=(‘d’, ‘p32’) may indicate a \(L_3\) edge transition from core \(2p_{3/2}\) shell to valence \(3d\) shell for Ni.

ominc: 1d float array

Incident energy of photon.

gamma_c: a float number or a 1d float array with the same shape as ominc.

The core-hole life-time broadening factor. It can be a constant value or incident energy dependent.

v_noccu: int

Total occupancy of valence shells.

thin: float number

The incident angle of photon (in radian).

phi: float number

Azimuthal angle (in radian), defined with respect to the \(x\)-axis of the local scattering axis: scatter_axis[:,0].

pol_type: list of tuples

Type of polarization, options can be:

  • (‘linear’, alpha), linear polarization, where alpha is the angle between the polarization vector and the scattering plane in radians.

  • (‘left’, 0), left circular polarization.

  • (‘right’, 0), right circular polarization.

  • (‘isotropic’, 0). isotropic polarization.

It will set pol_type=[(‘isotropic’, 0)] if not provided.

num_gs: int

Number of initial states used in XAS calculations.

nkryl: int

Maximum number of poles obtained.

temperature: float number

Temperature (in K) for boltzmann distribution.

loc_axis: 3*3 float array

The local axis with respect to which local orbitals are defined.

  • x: local_axis[:,0],

  • y: local_axis[:,1],

  • z: local_axis[:,2].

It will be an identity matrix if not provided.

scatter_axis: 3*3 float array

The local axis defining the scattering geometry. The scattering plane is defined in the local \(zx\)-plane.

  • local \(x\)-axis: scatter_axis[:,0]

  • local \(y\)-axis: scatter_axis[:,1]

  • local \(z\)-axis: scatter_axis[:,2]

It will be set to an identity matrix if not provided.

Returns:
xas: 2d array, shape=(len(ominc), len(pol_type))

The calculated XAS spectra. The first dimension is for ominc, and the second dimension if for different polarizations.

poles: list of dict, shape=(len(pol_type), )

The calculated XAS poles for different polarizations.

edrixs.solvers.xas_1v1c_py(eval_i, eval_n, trans_op, ominc, *, gamma_c=0.1, thin=1.0, phi=0, pol_type=None, gs_list=None, temperature=1.0, scatter_axis=None)[source]

Calculate XAS for the case of one valence shell plus one core shell with Python solver.

This solver is only suitable for small size of Hamiltonian, typically the dimension of both initial and intermediate Hamiltonian are smaller than 10,000.

Parameters:
eval_i: 1d float array

The eigenvalues of the initial Hamiltonian.

eval_n: 1d float array

The eigenvalues of the intermediate Hamiltonian.

trans_op: 3d complex array

The transition operators in the eigenstates basis.

ominc: 1d float array

Incident energy of photon.

gamma_c: a float number or a 1d float array with the same shape as ominc.

The core-hole life-time broadening factor. It can be a constant value or incident energy dependent.

thin: float number

The incident angle of photon (in radian).

phi: float number

Azimuthal angle (in radian), defined with respect to the \(x\)-axis of the scattering axis: scatter_axis[:,0].

pol_type: list of tuples

Type of polarization, options can be:

  • (‘linear’, alpha), linear polarization, where alpha is the angle between the polarization vector and the scattering plane in radians.

  • (‘left’, 0), left circular polarization.

  • (‘right’, 0), right circular polarization.

  • (‘isotropic’, 0). isotropic polarization.

It will set pol_type=[(‘isotropic’, 0)] if not provided.

gs_list: 1d list of ints

The indices of initial states which will be used in XAS calculations.

It will set gs_list=[0] if not provided.

temperature: float number

Temperature (in K) for boltzmann distribution.

scatter_axis: 3*3 float array

The local axis defining the scattering plane. The scattering plane is defined in the local \(zx\)-plane.

local \(x\)-axis: scatter_axis[:,0]

local \(y\)-axis: scatter_axis[:,1]

local \(z\)-axis: scatter_axis[:,2]

It will be set to an identity matrix if not provided.

Returns:
xas: 2d float array

The calculated XAS spectra. The 1st dimension is for the incident energy, and the 2nd dimension is for different polarizations.

edrixs.solvers.xas_2v1c_fort(comm, shell_name, ominc, *, gamma_c=0.1, v_tot_noccu=1, trans_to_which=1, thin=1.0, phi=0, pol_type=None, num_gs=1, nkryl=200, temperature=1.0, loc_axis=None, scatter_axis=None)[source]

Calculate XAS for the case with two valence shells plus one core shell with Fortran solver.

Parameters:
comm: MPI_comm

MPI communicator.

shell_name: tuple of three strings

Names of valence and core shells. The 1st (2nd) string in the tuple is for the 1st (2nd) valence shell, and the 3rd one is for the core shell.

  • The 1st and 2nd strings can only be ‘s’, ‘p’, ‘t2g’, ‘d’, ‘f’,

  • The 3nd string can be ‘s’, ‘p’, ‘p12’, ‘p32’, ‘d’, ‘d32’, ‘d52’, ‘f’, ‘f52’, ‘f72’.

For example: shell_name=(‘d’, ‘p’, ‘s’) may indicate a \(K\) edge transition from core \(1s\) shell to valence \(3d\) and \(4p\) shell for Ni.

ominc: 1d float array

Incident energy of photon.

gamma_c: a float number or a 1d float array with the same shape as ominc.

The core-hole life-time broadening factor. It can be a constant value or incident energy dependent.

v_tot_noccu: int

Total occupancy of valence shells.

trans_to_which: int

Photon transition to which valence shell.

  • 1: to 1st valence shell,

  • 2: to 2nd valence shell.

thin: float number

The incident angle of photon (in radian).

phi: float number

Azimuthal angle (in radian), defined with respect to the \(x\)-axis of the local scattering axis: scatter_axis[:,0].

pol_type: list of tuples

Type of polarization, options can be:

  • (‘linear’, alpha), linear polarization, where alpha is the angle between the polarization vector and the scattering plane in radians.

  • (‘left’, 0), left circular polarization.

  • (‘right’, 0), right circular polarization.

  • (‘isotropic’, 0). isotropic polarization.

It will set pol_type=[(‘isotropic’, 0)] if not provided.

num_gs: int

Number of initial states used in XAS calculations.

nkryl: int

Maximum number of poles obtained.

temperature: float number

Temperature (in K) for boltzmann distribution.

loc_axis: 3*3 float array

The local axis with respect to which local orbitals are defined.

  • x: local_axis[:,0],

  • y: local_axis[:,1],

  • z: local_axis[:,2].

It will be an identity matrix if not provided.

scatter_axis: 3*3 float array

The local axis defining the scattering geometry. The scattering plane is defined in the local \(zx\)-plane.

  • local \(x\)-axis: scatter_axis[:,0]

  • local \(y\)-axis: scatter_axis[:,1]

  • local \(z\)-axis: scatter_axis[:,2]

It will be set to an identity matrix if not provided.

Returns:
xas: 2d array, shape=(len(ominc), len(pol_type))

The calculated XAS spectra. The first dimension is for ominc, and the second dimension if for different polarizations.

poles: list of dict, shape=(len(pol_type), )

The calculated XAS poles for different polarizations.

edrixs.solvers.xas_siam_fort(comm, shell_name, nbath, ominc, *, gamma_c=0.1, v_noccu=1, thin=1.0, phi=0, pol_type=None, num_gs=1, nkryl=200, temperature=1.0, loc_axis=None, scatter_axis=None)[source]

Calculate XAS for single impurity Anderson model (SIAM) with Fortran solver.

Parameters:
comm: MPI_comm

MPI communicator.

shell_name: tuple of two strings

Names of valence and core shells. The 1st (2nd) string in the tuple is for the valence (core) shell.

  • The 1st string can only be ‘s’, ‘p’, ‘t2g’, ‘d’, ‘f’,

  • The 2nd string can be ‘s’, ‘p’, ‘p12’, ‘p32’, ‘d’, ‘d32’, ‘d52’, ‘f’, ‘f52’, ‘f72’.

For example: shell_name=(‘d’, ‘p32’) may indicate a \(L_3\) edge transition from core \(2p_{3/2}\) shell to valence \(3d\) shell for Ni.

nbath: int

Number of bath sites.

ominc: 1d float array

Incident energy of photon.

gamma_c: a float number or a 1d float array with the same shape as ominc.

The core-hole life-time broadening factor. It can be a constant value or incident energy dependent.

v_noccu: int

Total occupancy of valence shells.

thin: float number

The incident angle of photon (in radian).

phi: float number

Azimuthal angle (in radian), defined with respect to the \(x\)-axis of the local scattering axis: scatter_axis[:,0].

pol_type: list of tuples

Type of polarization, options can be:

  • (‘linear’, alpha), linear polarization, where alpha is the angle between the polarization vector and the scattering plane in radians.

  • (‘left’, 0), left circular polarization.

  • (‘right’, 0), right circular polarization.

  • (‘isotropic’, 0). isotropic polarization.

It will set pol_type=[(‘isotropic’, 0)] if not provided.

num_gs: int

Number of initial states used in XAS calculations.

nkryl: int

Maximum number of poles obtained.

temperature: float number

Temperature (in K) for boltzmann distribution.

loc_axis: 3*3 float array

The local axis with respect to which local orbitals are defined.

  • x: local_axis[:,0],

  • y: local_axis[:,1],

  • z: local_axis[:,2].

It will be an identity matrix if not provided.

scatter_axis: 3*3 float array

The local axis defining the scattering geometry. The scattering plane is defined in the local \(zx\)-plane.

  • local \(x\)-axis: scatter_axis[:,0]

  • local \(y\)-axis: scatter_axis[:,1]

  • local \(z\)-axis: scatter_axis[:,2]

It will be set to an identity matrix if not provided.

Returns:
xas: 2d array, shape=(len(ominc), len(pol_type))

The calculated XAS spectra. The first dimension is for ominc, and the second dimension if for different polarizations.

poles: list of dict, shape=(len(pol_type), )

The calculated XAS poles for different polarizations.