(git:6a2e663)
qs_operators_ao Module Reference

Functions/Subroutines

subroutine, public build_lin_mom_matrix (qs_env, matrix)
 Calculation of the linear momentum matrix <mu|∂|nu> over Cartesian Gaussian functions. More...
 
subroutine, public build_ang_mom_matrix (qs_env, matrix, rc)
 Calculation of the angular momentum matrix over Cartesian Gaussian functions. More...
 
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. More...
 
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. More...
 
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. More...
 

Detailed Description

History
created 07.2005
Author
MI (07.2005)

Function/Subroutine Documentation

◆ build_lin_mom_matrix()

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.

Parameters
qs_env...
matrix...
Date
27.02.2009
Author
VW
Version
1.0

Definition at line 73 of file qs_operators_ao.F.

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

◆ build_ang_mom_matrix()

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.

Parameters
qs_env...
matrix...
rc...
Date
27.02.2009
Author
VW
Version
1.0

Definition at line 443 of file qs_operators_ao.F.

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

◆ p_xyz_ao()

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.

Parameters
opmatrix representation of the p operator calculated in terms of the contracted basis functions
qs_envenvironment for the lists and the basis sets
minimum_imagetake into account only the first neighbors in the lists
History
06.2005 created [MI]
Author
MI

Definition at line 817 of file qs_operators_ao.F.

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

◆ rrc_xyz_ao()

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.

Parameters
opmatrix representation of the p operator calculated in terms of the contracted basis functions
qs_envenvironment for the lists and the basis sets
rcreference vector position
ordermaximum order of the momentum, for the dipole order = 1, order = -2 for quad only
minimum_imagetake into account only the first neighbors in the lists
soft...
History
03.2006 created [MI] 06.2019 added quarupole only option (A.Bussy)
Author
MI

Definition at line 1083 of file qs_operators_ao.F.

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

◆ rrc_xyz_der_ao()

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]

Parameters
opmatrix representation of the p operator calculated in terms of the contracted basis functions
op_der...
qs_envenvironment for the lists and the basis sets
rcreference vector position
ordermaximum order of the momentum, for the dipole order = 1
minimum_imagetake into account only the first neighbors in the lists
soft...
History
03.2006 created [MI]
Author
MI
Note
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

Definition at line 1349 of file qs_operators_ao.F.

Here is the call graph for this function: