(git:6a2e663)
qs_moments Module Reference

Calculates the moment integrals <a|r^m|b> and <a|r x d/dr|b> More...

Functions/Subroutines

subroutine, public dipole_velocity_deriv (qs_env, difdip, order, lambda, rc)
 This returns the derivative of the moment integrals [a|\mu|b], with respect to the basis function on the right difdip(beta, alpha) = < mu | r_beta | ∂_alpha nu > * (mu - nu) More...
 
subroutine, public build_dsdv_moments (qs_env, moments, nmoments, ref_point, ref_points, basis_type)
 Builds the moments for the derivative of the overlap with respect to nuclear velocities. More...
 
subroutine, public build_local_moment_matrix (qs_env, moments, nmoments, ref_point, ref_points, basis_type)
 ... More...
 
subroutine, public build_local_moments_der_matrix (qs_env, moments_der, nmoments_der, nmoments, ref_point, moments)
 Calculate right-hand sided derivatives of multipole moments, e. g. < a | xy d/dz | b > Optionally stores the multipole moments themselves for free. Note that the multipole moments are symmetric while their derivatives are anti-symmetric Only first derivatives are performed, e. g. x d/dy. More...
 
subroutine, public build_local_magmom_matrix (qs_env, magmom, nmoments, ref_point, ref_points, basis_type)
 ... More...
 
subroutine, public build_berry_moment_matrix (qs_env, cosmat, sinmat, kvec, sab_orb_external, basis_type)
 ... More...
 
subroutine, public build_berry_kpoint_matrix (qs_env, cosmat, sinmat, kvec, basis_type)
 ... More...
 
subroutine, public qs_moment_berry_phase (qs_env, magnetic, nmoments, reference, ref_point, unit_number)
 ... More...
 
subroutine, public qs_moment_locop (qs_env, magnetic, nmoments, reference, ref_point, unit_number, vel_reprs, com_nl)
 ... More...
 
subroutine, public dipole_deriv_ao (qs_env, difdip, deltaR, order, rcc)
 ... More...
 

Detailed Description

Calculates the moment integrals <a|r^m|b> and <a|r x d/dr|b>

History
added angular moments (JGH 11.2012)
Author
JGH (20.07.2006)

Function/Subroutine Documentation

◆ dipole_velocity_deriv()

subroutine, public qs_moments::dipole_velocity_deriv ( type(qs_environment_type), pointer  qs_env,
type(dbcsr_p_type), dimension(:, :), intent(inout)  difdip,
integer, intent(in)  order,
integer, intent(in)  lambda,
real(kind=dp), dimension(3)  rc 
)

This returns the derivative of the moment integrals [a|\mu|b], with respect to the basis function on the right difdip(beta, alpha) = < mu | r_beta | ∂_alpha nu > * (mu - nu)

Parameters
qs_env...
difdip...
orderThe order of the derivative (1 for dipole moment)
lambdaThe atom on which we take the derivative
rc...
Author
Edward Ditler

Definition at line 125 of file qs_moments.F.

Here is the call graph for this function:

◆ build_dsdv_moments()

subroutine, public qs_moments::build_dsdv_moments ( type(qs_environment_type), pointer  qs_env,
type(dbcsr_p_type), dimension(:)  moments,
integer, intent(in)  nmoments,
real(kind=dp), dimension(:), intent(in), optional  ref_point,
real(kind=dp), dimension(:, :), intent(in), optional  ref_points,
character(len=*), optional  basis_type 
)

Builds the moments for the derivative of the overlap with respect to nuclear velocities.

Parameters
qs_env...
moments...
nmoments...
ref_point...
ref_points...
basis_type...
Author
Edward Ditler

Definition at line 328 of file qs_moments.F.

Here is the call graph for this function:

◆ build_local_moment_matrix()

subroutine, public qs_moments::build_local_moment_matrix ( type(qs_environment_type), pointer  qs_env,
type(dbcsr_p_type), dimension(:), pointer  moments,
integer, intent(in)  nmoments,
real(kind=dp), dimension(:), intent(in), optional  ref_point,
real(kind=dp), dimension(:, :), intent(in), optional  ref_points,
character(len=*), optional  basis_type 
)

...

Parameters
qs_env...
moments...
nmoments...
ref_point...
ref_points...
basis_type...

Definition at line 557 of file qs_moments.F.

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

◆ build_local_moments_der_matrix()

subroutine, public qs_moments::build_local_moments_der_matrix ( type(qs_environment_type), pointer  qs_env,
type(dbcsr_p_type), dimension(:, :), intent(inout), pointer  moments_der,
integer, intent(in)  nmoments_der,
integer, intent(in)  nmoments,
real(kind=dp), dimension(:), intent(in), optional  ref_point,
type(dbcsr_p_type), dimension(:), intent(inout), optional, pointer  moments 
)

Calculate right-hand sided derivatives of multipole moments, e. g. < a | xy d/dz | b > Optionally stores the multipole moments themselves for free. Note that the multipole moments are symmetric while their derivatives are anti-symmetric Only first derivatives are performed, e. g. x d/dy.

Parameters
qs_env...
moments_derwill contain the derivatives of the multipole moments
nmoments_derorder of the moments with derivatives
nmomentsorder of the multipole moments (no derivatives, same output as build_local_moment_matrix, needs moments as arguments to store results)
ref_point...
momentscontains the multipole moments, optionally for free, up to order nmoments
Note
Adapted from rRc_xyz_der_ao in qs_operators_ao

Definition at line 788 of file qs_moments.F.

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

◆ build_local_magmom_matrix()

subroutine, public qs_moments::build_local_magmom_matrix ( type(qs_environment_type), pointer  qs_env,
type(dbcsr_p_type), dimension(:), pointer  magmom,
integer, intent(in)  nmoments,
real(kind=dp), dimension(:), intent(in), optional  ref_point,
real(kind=dp), dimension(:, :), intent(in), optional  ref_points,
character(len=*), optional  basis_type 
)

...

Parameters
qs_env...
magmom...
nmoments...
ref_point...
ref_points...
basis_type...

Definition at line 1122 of file qs_moments.F.

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

◆ build_berry_moment_matrix()

subroutine, public qs_moments::build_berry_moment_matrix ( type(qs_environment_type), pointer  qs_env,
type(dbcsr_type), pointer  cosmat,
type(dbcsr_type), pointer  sinmat,
real(kind=dp), dimension(3), intent(in)  kvec,
type(neighbor_list_set_p_type), dimension(:), optional, pointer  sab_orb_external,
character(len=*), optional  basis_type 
)

...

Parameters
qs_env...
cosmat...
sinmat...
kvec...
sab_orb_external...
basis_type...

Definition at line 1335 of file qs_moments.F.

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

◆ build_berry_kpoint_matrix()

subroutine, public qs_moments::build_berry_kpoint_matrix ( type(qs_environment_type), pointer  qs_env,
type(dbcsr_p_type), dimension(:, :), pointer  cosmat,
type(dbcsr_p_type), dimension(:, :), pointer  sinmat,
real(kind=dp), dimension(3), intent(in)  kvec,
character(len=*), optional  basis_type 
)

...

Parameters
qs_env...
cosmat...
sinmat...
kvec...
basis_type...

Definition at line 1502 of file qs_moments.F.

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

◆ qs_moment_berry_phase()

subroutine, public qs_moments::qs_moment_berry_phase ( type(qs_environment_type), pointer  qs_env,
logical, intent(in)  magnetic,
integer, intent(in)  nmoments,
integer, intent(in)  reference,
real(dp), dimension(:), pointer  ref_point,
integer, intent(in)  unit_number 
)

...

Parameters
qs_env...
magnetic...
nmoments...
reference...
ref_point...
unit_number...

Definition at line 1711 of file qs_moments.F.

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

◆ qs_moment_locop()

subroutine, public qs_moments::qs_moment_locop ( type(qs_environment_type), pointer  qs_env,
logical, intent(in)  magnetic,
integer, intent(in)  nmoments,
integer, intent(in)  reference,
real(dp), dimension(:), intent(in), pointer  ref_point,
integer, intent(in)  unit_number,
logical, intent(in), optional  vel_reprs,
logical, intent(in), optional  com_nl 
)

...

Parameters
qs_env...
magnetic...
nmoments...
reference...
ref_point...
unit_number...
vel_reprs...
com_nl...

Definition at line 2232 of file qs_moments.F.

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

◆ dipole_deriv_ao()

subroutine, public qs_moments::dipole_deriv_ao ( type(qs_environment_type), pointer  qs_env,
type(dbcsr_p_type), dimension(:, :), intent(inout), pointer  difdip,
real(kind=dp), dimension(:, :), intent(in), pointer  deltaR,
integer, intent(in)  order,
real(kind=dp), dimension(3), optional  rcc 
)

...

Parameters
qs_env...
difdip...
deltaR...
order...
rcc...
Note
calculate matrix elements <a|r_beta |db/dR_alpha > + <da/dR_alpha | r_beta | b > be aware: < a | r_beta| db/dR_alpha > = - < da/dR_alpha | r_beta | b > only valid if alpha .neq.beta if alpha=beta: < a | r_beta| db/dR_alpha > = - < da/dR_alpha | r_beta | b > - < a | b > modified from qs_efield_mo_derivatives SL July 2015

Definition at line 3093 of file qs_moments.F.

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