(git:3add494)
xas_tdp_utils Module Reference

Utilities for X-ray absorption spectroscopy using TDDFPT. More...

Functions/Subroutines

subroutine, public setup_xas_tdp_prob (donor_state, qs_env, xas_tdp_env, xas_tdp_control)
 Builds the matrix that defines the XAS TDDFPT generalized eigenvalue problem to be solved for excitation energies omega. The problem has the form omega*G*C = M*C, where C contains the response orbitals coefficients. The matrix M and the metric G are stored in the given donor_state. More...
 
subroutine, public solve_xas_tdp_prob (donor_state, xas_tdp_control, xas_tdp_env, qs_env, ex_type)
 Solves the XAS TDP generalized eigenvalue problem omega*C = matrix_tdp*C using standard full diagonalization methods. The problem is Hermitian (made that way even if not TDA) More...
 
subroutine, public include_os_soc (donor_state, xas_tdp_env, xas_tdp_control, qs_env)
 Includes the SOC effects on the precomputed spin-conserving and spin-flip excitations from an open-shell calculation (UKS or ROKS). This is a perturbative treatment. More...
 
subroutine, public include_rcs_soc (donor_state, xas_tdp_env, xas_tdp_control, qs_env)
 Includes the SOC effects on the precomputed restricted closed-shell singlet and triplet excitations. This is a perturbative treatmnent. More...
 
subroutine, public rcs_amew_soc_elements (amew_soc, lr_soc, lr_overlap, domo_soc, pref_trace, pref_overall, pref_diags, symmetric)
 Computes the rcs SOC matrix elements between excited states AMEWs based on the LR orbitals. More...
 

Detailed Description

Utilities for X-ray absorption spectroscopy using TDDFPT.

Author
AB (01.2018)

Function/Subroutine Documentation

◆ setup_xas_tdp_prob()

subroutine, public xas_tdp_utils::setup_xas_tdp_prob ( type(donor_state_type), pointer  donor_state,
type(qs_environment_type), pointer  qs_env,
type(xas_tdp_env_type), pointer  xas_tdp_env,
type(xas_tdp_control_type), pointer  xas_tdp_control 
)

Builds the matrix that defines the XAS TDDFPT generalized eigenvalue problem to be solved for excitation energies omega. The problem has the form omega*G*C = M*C, where C contains the response orbitals coefficients. The matrix M and the metric G are stored in the given donor_state.

Parameters
donor_statethe donor_state for which the problem is restricted
qs_env...
xas_tdp_env...
xas_tdp_control...
Note
the matrix M is symmetric and has the form | M_d M_o | | M_o M_d |, -In the SPIN-RESTRICTED case: depending on whether we consider singlet or triplet excitation, the diagonal (M_d) and off-diagonal (M_o) parts of M differ:
  • For singlet: M_d = A + 2B + C_aa + C_ab - D M_o = 2B + C_aa + C_ab - E
  • For triplet: M_d = A + C_aa - C_ab - D M_o = C_aa - C_ab - E where other subroutines computes the matrices A, B, E, D and G, which are:
  • A: the ground-state contribution: F_pq*delta_IJ - epsilon_IJ*S_pq
  • B: the Coulomb kernel ~(pI|Jq)
  • C: the xc kernel c_aa (double derivatibe wrt to n_alpha) and C_ab (wrt n_alpha and n_beta)
  • D: the on-digonal exact exchange kernel ~(pq|IJ)
  • E: the off-diagonal exact exchange kernel ~(pJ|Iq)
  • G: the metric S_pq*delta_IJ For the xc functionals, C_aa + C_ab or C_aa - C_ab are stored in the same matrix In the above definitions, I,J label the donnor MOs and p,q the sgfs of the basis

-In the SPIN-UNRESTRICTED, spin-conserving case: the on- and off-diagonal elements of M are: M_d = A + B + C -D M_o = B + C - E where the submatrices A, B, C, D and E are:

  • A: the groun-state contribution: (F_pq*delta_IJ - epsilon_IJ*S_pq) * delta_ab
  • B: the Coulomb kernel: (pI_a|J_b q)
  • C: the xc kernel: (pI_a|fxc_ab|J_b q)
  • D: the on-diagonal exact-exchange kernel: (pq|I_a J_b) delta_ab
  • E: the off-diagonal exact-exchange kernel: (pJ_b|I_a q) delta_ab
  • G: the metric S_pq*delta_IJ*delta_ab p,q label the sgfs, I,J the donro MOs and a,b the spins

-In both above cases, the matrix M is always projected onto the unperturbed unoccupied ground state: M <= Q * M * Q^T = (1 - SP) * M * (1 - PS)

-In the SPIN-FLIP case: Only the TDA is implemented, that is, there are only on-diagonal elements: M_d = A + C - D where the submatrices A, C and D are:

  • A: the ground state-contribution: (F_pq*delta_IJ - epsilon_IJ*S_pq) * delta_ab, but here, the alph-alpha quadrant has the beta Fock matrix and the beta-beta quadrant has the alpha Fock matrix
  • C: the SF xc kernel: (pI_a|fxc|J_bq), fxc = 1/m * (vxc_a -vxc_b)
  • D: the on-diagonal exact-exchange kernel: (pq|I_a J_b) delta_ab To ensure that all excitation start from a given spin to the opposite, we then multiply by a Q projector where we swap the alpha-alpha and beta-beta spin-quadrants

All possibilities: TDA or full-TDDFT, singlet or triplet, xc or hybrid, etc are treated in the same routine to avoid recomputing stuff Under TDA, only the on-diagonal elements of M are computed In the case of non-TDA, one turns the problem Hermitian

Definition at line 167 of file xas_tdp_utils.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ solve_xas_tdp_prob()

subroutine, public xas_tdp_utils::solve_xas_tdp_prob ( type(donor_state_type), pointer  donor_state,
type(xas_tdp_control_type), pointer  xas_tdp_control,
type(xas_tdp_env_type), pointer  xas_tdp_env,
type(qs_environment_type), pointer  qs_env,
integer, intent(in)  ex_type 
)

Solves the XAS TDP generalized eigenvalue problem omega*C = matrix_tdp*C using standard full diagonalization methods. The problem is Hermitian (made that way even if not TDA)

Parameters
donor_state...
xas_tdp_control...
xas_tdp_env...
qs_env...
ex_typewhether we deal with singlets, triplets, spin-conserving open-shell or spin-flip
Note
The computed eigenvalues and eigenvectors are stored in the donor_state The eigenvectors are the LR-coefficients. In case of TDA, c^- is stored. In the general case, the sum c^+ + c^- is stored.
  • Spin-restricted: In case both singlets and triplets are considered, this routine must be called twice. This is the choice that was made because the body of the routine is exactly the same in both cases Note that for singlet we solve for u = 1/sqrt(2)*(c_alpha + c_beta) = sqrt(2)*c and that for triplets we solve for v = 1/sqrt(2)*(c_alpha - c_beta) = sqrt(2)*c
  • Spin-unrestricted: The problem is solved for the LR coefficients c_pIa as they are (not linear combination) The routine might be called twice (once for spin-conservign, one for spin-flip)

Definition at line 457 of file xas_tdp_utils.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ include_os_soc()

subroutine, public xas_tdp_utils::include_os_soc ( type(donor_state_type), pointer  donor_state,
type(xas_tdp_env_type), pointer  xas_tdp_env,
type(xas_tdp_control_type), pointer  xas_tdp_control,
type(qs_environment_type), pointer  qs_env 
)

Includes the SOC effects on the precomputed spin-conserving and spin-flip excitations from an open-shell calculation (UKS or ROKS). This is a perturbative treatment.

Parameters
donor_state...
xas_tdp_env...
xas_tdp_control...
qs_env...
Note
Using AMEWs, build an hermitian matrix with all excited states SOC coupling + the excitation energies on the diagonal. Then diagonalize it to get the new excitation energies and corresponding linear combinations of lr_coeffs. The AMEWs are normalized Only for open-shell calculations

Definition at line 1501 of file xas_tdp_utils.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ include_rcs_soc()

subroutine, public xas_tdp_utils::include_rcs_soc ( type(donor_state_type), pointer  donor_state,
type(xas_tdp_env_type), pointer  xas_tdp_env,
type(xas_tdp_control_type), pointer  xas_tdp_control,
type(qs_environment_type), pointer  qs_env 
)

Includes the SOC effects on the precomputed restricted closed-shell singlet and triplet excitations. This is a perturbative treatmnent.

Parameters
donor_state...
xas_tdp_env...
xas_tdp_control...
qs_env...
Note
Using AMEWs, build an hermitian matrix with all excited states SOC coupling + the excitation energies on the diagonal. Then diagonalize it to get the new excitation energies and corresponding linear combinations of lr_coeffs. The AMEWs are normalized Only for spin-restricted calculations The ms=-1,+1 triplets are not explicitely computed in the first place. Assume they have the same energy as the ms=0 triplets and apply the spin raising and lowering operators on the latter to get their AMEWs => this is the qusi-degenerate perturbation theory approach by Neese (QDPT)

Definition at line 1938 of file xas_tdp_utils.F.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ rcs_amew_soc_elements()

subroutine, public xas_tdp_utils::rcs_amew_soc_elements ( type(dbcsr_type)  amew_soc,
type(dbcsr_type)  lr_soc,
type(dbcsr_type)  lr_overlap,
real(dp), dimension(:, :)  domo_soc,
real(dp)  pref_trace,
real(dp)  pref_overall,
real(dp), optional  pref_diags,
logical, optional  symmetric 
)

Computes the rcs SOC matrix elements between excited states AMEWs based on the LR orbitals.

Parameters
amew_socoutput dbcsr matrix with the SOC in the AMEW basis (needs to be fully resereved)
lr_socdbcsr matrix with the SOC wrt the LR orbitals
lr_overlapdbcsr matrix with the excited states LR orbital overlap
domo_socthe SOC in the basis of the donor MOs
pref_tracesee notes
pref_overallsee notes
pref_diagssee notes
symmetricif the outcome is known to be symmetric, only elements with iex <= jex are done
Note
For an excited states pair i,j, the AMEW SOC matrix element is: soc_ij = pref_overall*(SUM(diag(lr_soc_ij)) + pref_trace*SUM(lr_overlap_ij*TRANSPOSE(domo_soc))) optionally, the value pref_diags*SUM(diag(lr_overlap_ij)) can be added (before pref_overall)

Definition at line 2881 of file xas_tdp_utils.F.

Here is the call graph for this function:
Here is the caller graph for this function: