![]() |
(git:b77b4be)
|
Functions/Subroutines | |
subroutine, public | build_lin_mom_matrix (qs_env, matrix) |
Calculation of the linear momentum matrix <mu|∂|nu> over Cartesian Gaussian functions. | |
subroutine, public | build_ang_mom_matrix (qs_env, matrix, rc) |
Calculation of the angular momentum matrix over Cartesian Gaussian functions. | |
subroutine, public | p_xyz_ao (op, qs_env, minimum_image) |
Calculation of the components of the dipole operator in the velocity form The elements of the sparse matrices are the integrals in the basis functions. | |
subroutine, public | rrc_xyz_ao (op, qs_env, rc, order, minimum_image, soft) |
Calculation of the components of the dipole operator in the length form by taking the relative position operator r-Rc, with respect a reference point Rc Probably it does not work for PBC, or maybe yes if the wfn are sufficiently localized The elements of the sparse matrices are the integrals in the basis functions. | |
subroutine, public | rrc_xyz_der_ao (op, op_der, qs_env, rc, order, minimum_image, soft) |
Calculation of the multipole operators integrals and of its derivatives of the type [\mu | op | d(\nu)/dR(\nu)]-[d(\mu)/dR(\mu)| op | \nu] by taking the relative position operator r-Rc, with respect a reference point Rc The derivative are with respect to the primitive position, The multipole operator is symmetric and if it does not depend on R(\mu) or R(\nu) therefore [\mu | op | d(\nu)/dR(\nu)] = -[d(\mu)/dR(\mu)| op | \nu] [\mu|op|d(\nu)/dR]-[d(\mu)/dR|op|\nu]=2[\mu|op|d(\nu)/dR] When it is not the case a correction term is needed. | |
subroutine, public qs_operators_ao::build_lin_mom_matrix | ( | type(qs_environment_type), pointer | qs_env, |
type(dbcsr_p_type), dimension(:), pointer | matrix | ||
) |
Calculation of the linear momentum matrix <mu|∂|nu> over Cartesian Gaussian functions.
qs_env | ... |
matrix | ... |
Definition at line 73 of file qs_operators_ao.F.
subroutine, public qs_operators_ao::build_ang_mom_matrix | ( | type(qs_environment_type), pointer | qs_env, |
type(dbcsr_p_type), dimension(:), pointer | matrix, | ||
real(dp), dimension(:), intent(in) | rc | ||
) |
Calculation of the angular momentum matrix over Cartesian Gaussian functions.
qs_env | ... |
matrix | ... |
rc | ... |
Definition at line 443 of file qs_operators_ao.F.
subroutine, public qs_operators_ao::p_xyz_ao | ( | type(dbcsr_p_type), dimension(:), pointer | op, |
type(qs_environment_type), pointer | qs_env, | ||
logical, intent(in), optional | minimum_image | ||
) |
Calculation of the components of the dipole operator in the velocity form The elements of the sparse matrices are the integrals in the basis functions.
op | matrix representation of the p operator calculated in terms of the contracted basis functions |
qs_env | environment for the lists and the basis sets |
minimum_image | take into account only the first neighbors in the lists |
Definition at line 817 of file qs_operators_ao.F.
subroutine, public qs_operators_ao::rrc_xyz_ao | ( | type(dbcsr_p_type), dimension(:), pointer | op, |
type(qs_environment_type), pointer | qs_env, | ||
real(dp), dimension(3) | rc, | ||
integer, intent(in) | order, | ||
logical, intent(in), optional | minimum_image, | ||
logical, intent(in), optional | soft | ||
) |
Calculation of the components of the dipole operator in the length form by taking the relative position operator r-Rc, with respect a reference point Rc Probably it does not work for PBC, or maybe yes if the wfn are sufficiently localized The elements of the sparse matrices are the integrals in the basis functions.
op | matrix representation of the p operator calculated in terms of the contracted basis functions |
qs_env | environment for the lists and the basis sets |
rc | reference vector position |
order | maximum order of the momentum, for the dipole order = 1, order = -2 for quad only |
minimum_image | take into account only the first neighbors in the lists |
soft | ... |
Definition at line 1083 of file qs_operators_ao.F.
subroutine, public qs_operators_ao::rrc_xyz_der_ao | ( | type(dbcsr_p_type), dimension(:), pointer | op, |
type(dbcsr_p_type), dimension(:, :), pointer | op_der, | ||
type(qs_environment_type), pointer | qs_env, | ||
real(dp), dimension(3) | rc, | ||
integer, intent(in) | order, | ||
logical, intent(in), optional | minimum_image, | ||
logical, intent(in), optional | soft | ||
) |
Calculation of the multipole operators integrals and of its derivatives of the type [\mu | op | d(\nu)/dR(\nu)]-[d(\mu)/dR(\mu)| op | \nu] by taking the relative position operator r-Rc, with respect a reference point Rc The derivative are with respect to the primitive position, The multipole operator is symmetric and if it does not depend on R(\mu) or R(\nu) therefore [\mu | op | d(\nu)/dR(\nu)] = -[d(\mu)/dR(\mu)| op | \nu] [\mu|op|d(\nu)/dR]-[d(\mu)/dR|op|\nu]=2[\mu|op|d(\nu)/dR] When it is not the case a correction term is needed.
The momentum operator [\mu|M|\nu] is symmetric, the number of components is determined by the order: 3 for order 1 (x,y,x), 9 for order 2(xx,xy,xz,yy,yz,zz) The derivative of the type [\mu | op | d(\nu)/dR_i(\nu)], where i indicates the cartesian direction, is antisymmetric only when the no component M =(r_i) or (r_i r_j) is in the same cartesian direction of the derivative, indeed d([\mu|M|\nu])/dr_i = [d(\mu)/dr_i|M|\nu] + [\mu|M|d(\nu)/dr_i] + [\mu |d(M)/dr_i|\nu] d([\mu|M|\nu])/dr_i = -[d(\mu)/dR_i(\mu)|M|\nu] -[\mu|M|d(\nu)/dR_i(\nu)] + [\mu |d(M)/dr_i|\nu] Therefore we cannot use an antisymmetric matrix
The same holds for the derivative with respect to the electronic position r taking into account that [\mu|op|d(\nu)/dR] = -[\mu|op|d(\nu)/dr]
op | matrix representation of the p operator calculated in terms of the contracted basis functions |
op_der | ... |
qs_env | environment for the lists and the basis sets |
rc | reference vector position |
order | maximum order of the momentum, for the dipole order = 1 |
minimum_image | take into account only the first neighbors in the lists |
soft | ... |
Definition at line 1349 of file qs_operators_ao.F.