(git:6a2e663)
ai_moments Module Reference

Calculation of the moment integrals over Cartesian Gaussian-type functions. More...

Functions/Subroutines

subroutine, public diff_momop_velocity (la_max, npgfa, zeta, rpgfa, la_min, lb_max, npgfb, zetb, rpgfb, lb_min, order, rac, rbc, difmab, lambda, iatom, jatom)
 This returns the derivative of the moment integrals [a|\mu|b], with respect to the primitive on the right difmab(:, :, beta, alpha) = < a | r_beta | ∂_alpha b > * (iatom - jatom) More...
 
subroutine, public diff_momop2 (la_max, npgfa, zeta, rpgfa, la_min, lb_max, npgfb, zetb, rpgfb, lb_min, order, rac, rbc, difmab, mab_ext, deltaR, iatom, jatom)
 This returns the derivative of the moment integrals [a|\mu|b], with respect to the position of the primitive on the left and right, i.e. [da/dR_ai|\mu|b] + [a|\mu|d/dR_bi] [da/dR_ai|\mu|b] = 2*zeta*[a+1i|\mu|b] - Ni(a)[a-1i|\mu|b] [a|\mu|d/dR_bi] = 2*zetb*[a|\mu|b+1i] - Ni(b)[a|\mu|b-1i] order indicates the max order of the moment operator to be calculated 1: dipole 2: quadrupole ... More...
 
subroutine, public contract_cossin (cos_block, sin_block, iatom, ncoa, nsgfa, sgfa, sphi_a, ldsa, jatom, ncob, nsgfb, sgfb, sphi_b, ldsb, cosab, sinab, ldab, work, ldwork)
 ... More...
 
subroutine, public cossin (la_max_set, npgfa, zeta, rpgfa, la_min_set, lb_max, npgfb, zetb, rpgfb, lb_min, rac, rbc, kvec, cosab, sinab, dcosab, dsinab)
 ... More...
 
subroutine, public moment (la_max, npgfa, zeta, rpgfa, la_min, lb_max, npgfb, zetb, rpgfb, lc_max, rac, rbc, mab)
 ... More...
 
subroutine, public diffop (la_max, npgfa, zeta, rpgfa, la_min, lb_max, npgfb, zetb, rpgfb, lb_min, rab, difab)
 This returns the derivative of the overlap integrals [a|b], with respect to the position of the primitive on the left, i.e. [da/dR_ai|b] = 2*zeta*[a+1i|b] - Ni(a)[a-1i|b]. More...
 
subroutine, public diff_momop (la_max, npgfa, zeta, rpgfa, la_min, lb_max, npgfb, zetb, rpgfb, lb_min, order, rac, rbc, difmab, mab_ext)
 This returns the derivative of the moment integrals [a|\mu|b], with respect to the position of the primitive on the left, i.e. [da/dR_ai|\mu|b] = 2*zeta*[a+1i|\mu|b] - Ni(a)[a-1i|\mu|b] order indicates the max order of the moment operator to be calculated 1: dipole 2: quadrupole ... More...
 
subroutine, public dipole_force (la_max, npgfa, zeta, rpgfa, la_min, lb_max, npgfb, zetb, rpgfb, lb_min, order, rac, rbc, pab, forcea, forceb)
 This returns the derivative of the dipole integrals [a|x|b], with respect to the position of the primitive on the left and right, i.e. [da/dR_ai|\mu|b] = 2*zeta*[a+1i|\mu|b] - Ni(a)[a-1i|\mu|b]. More...
 

Detailed Description

Calculation of the moment integrals over Cartesian Gaussian-type functions.

Literature
S. Obara and A. Saika, J. Chem. Phys. 84, 3963 (1986)
History
none
Author
J. Hutter (16.02.2005)

Function/Subroutine Documentation

◆ diff_momop_velocity()

subroutine, public ai_moments::diff_momop_velocity ( integer, intent(in)  la_max,
integer, intent(in)  npgfa,
real(kind=dp), dimension(:), intent(in)  zeta,
real(kind=dp), dimension(:), intent(in)  rpgfa,
integer, intent(in)  la_min,
integer, intent(in)  lb_max,
integer, intent(in)  npgfb,
real(kind=dp), dimension(:), intent(in)  zetb,
real(kind=dp), dimension(:), intent(in)  rpgfb,
integer, intent(in)  lb_min,
integer, intent(in)  order,
real(kind=dp), dimension(3), intent(in)  rac,
real(kind=dp), dimension(3), intent(in)  rbc,
real(kind=dp), dimension(:, :, :, :), intent(out)  difmab,
integer, intent(in)  lambda,
integer, intent(in), optional  iatom,
integer, intent(in), optional  jatom 
)

This returns the derivative of the moment integrals [a|\mu|b], with respect to the primitive on the right difmab(:, :, beta, alpha) = < a | r_beta | ∂_alpha b > * (iatom - jatom)

Parameters
la_max...
npgfa...
zeta...
rpgfa...
la_min...
lb_max...
npgfb...
zetb...
rpgfb...
lb_min...
order...
rac...
rbc...
difmab...
lambdaThe atom on which we take the derivative
iatom...
jatom...
Author
Edward Ditler

Definition at line 76 of file ai_moments.F.

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

◆ diff_momop2()

subroutine, public ai_moments::diff_momop2 ( integer, intent(in)  la_max,
integer, intent(in)  npgfa,
real(kind=dp), dimension(:), intent(in)  zeta,
real(kind=dp), dimension(:), intent(in)  rpgfa,
integer, intent(in)  la_min,
integer, intent(in)  lb_max,
integer, intent(in)  npgfb,
real(kind=dp), dimension(:), intent(in)  zetb,
real(kind=dp), dimension(:), intent(in)  rpgfb,
integer, intent(in)  lb_min,
integer, intent(in)  order,
real(kind=dp), dimension(3), intent(in)  rac,
real(kind=dp), dimension(3), intent(in)  rbc,
real(kind=dp), dimension(:, :, :, :), intent(out)  difmab,
real(kind=dp), dimension(:, :, :), optional, pointer  mab_ext,
real(kind=dp), dimension(:, :), intent(in), optional, pointer  deltaR,
integer, intent(in), optional  iatom,
integer, intent(in), optional  jatom 
)

This returns the derivative of the moment integrals [a|\mu|b], with respect to the position of the primitive on the left and right, i.e. [da/dR_ai|\mu|b] + [a|\mu|d/dR_bi] [da/dR_ai|\mu|b] = 2*zeta*[a+1i|\mu|b] - Ni(a)[a-1i|\mu|b] [a|\mu|d/dR_bi] = 2*zetb*[a|\mu|b+1i] - Ni(b)[a|\mu|b-1i] order indicates the max order of the moment operator to be calculated 1: dipole 2: quadrupole ...

Parameters
la_max...
npgfa...
zeta...
rpgfa...
la_min...
lb_max...
npgfb...
zetb...
rpgfb...
lb_min...
order...
rac...
rbc...
difmab...
mab_ext...
deltaRneeded for weighted derivative
iatom...
jatom... SL August 2015, ED 2021

Definition at line 161 of file ai_moments.F.

Here is the call graph for this function:

◆ contract_cossin()

subroutine, public ai_moments::contract_cossin ( real(dp), dimension(:, :), pointer  cos_block,
real(dp), dimension(:, :), pointer  sin_block,
integer, intent(in)  iatom,
integer, intent(in)  ncoa,
integer, intent(in)  nsgfa,
integer, intent(in)  sgfa,
real(dp), dimension(:, :), intent(in)  sphi_a,
integer, intent(in)  ldsa,
integer, intent(in)  jatom,
integer, intent(in)  ncob,
integer, intent(in)  nsgfb,
integer, intent(in)  sgfb,
real(dp), dimension(:, :), intent(in)  sphi_b,
integer, intent(in)  ldsb,
real(dp), dimension(:, :), intent(in)  cosab,
real(dp), dimension(:, :), intent(in)  sinab,
integer, intent(in)  ldab,
real(dp), dimension(:, :)  work,
integer, intent(in)  ldwork 
)

...

Parameters
cos_block...
sin_block...
iatom...
ncoa...
nsgfa...
sgfa...
sphi_a...
ldsa...
jatom...
ncob...
nsgfb...
sgfb...
sphi_b...
ldsb...
cosab...
sinab...
ldab...
work...
ldwork...

Definition at line 257 of file ai_moments.F.

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

◆ cossin()

subroutine, public ai_moments::cossin ( integer, intent(in)  la_max_set,
integer, intent(in)  npgfa,
real(kind=dp), dimension(:), intent(in)  zeta,
real(kind=dp), dimension(:), intent(in)  rpgfa,
integer, intent(in)  la_min_set,
integer, intent(in)  lb_max,
integer, intent(in)  npgfb,
real(kind=dp), dimension(:), intent(in)  zetb,
real(kind=dp), dimension(:), intent(in)  rpgfb,
integer, intent(in)  lb_min,
real(kind=dp), dimension(3), intent(in)  rac,
real(kind=dp), dimension(3), intent(in)  rbc,
real(kind=dp), dimension(3), intent(in)  kvec,
real(kind=dp), dimension(:, :), intent(inout)  cosab,
real(kind=dp), dimension(:, :), intent(inout)  sinab,
real(kind=dp), dimension(:, :, :), intent(inout), optional  dcosab,
real(kind=dp), dimension(:, :, :), intent(inout), optional  dsinab 
)

...

Parameters
la_max_set...
npgfa...
zeta...
rpgfa...
la_min_set...
lb_max...
npgfb...
zetb...
rpgfb...
lb_min...
rac...
rbc...
kvec...
cosab...
sinab...
dcosab...
dsinab...

Definition at line 336 of file ai_moments.F.

Here is the caller graph for this function:

◆ moment()

subroutine, public ai_moments::moment ( integer, intent(in)  la_max,
integer, intent(in)  npgfa,
real(kind=dp), dimension(:), intent(in)  zeta,
real(kind=dp), dimension(:), intent(in)  rpgfa,
integer, intent(in)  la_min,
integer, intent(in)  lb_max,
integer, intent(in)  npgfb,
real(kind=dp), dimension(:), intent(in)  zetb,
real(kind=dp), dimension(:), intent(in)  rpgfb,
integer, intent(in)  lc_max,
real(kind=dp), dimension(3), intent(in)  rac,
real(kind=dp), dimension(3), intent(in)  rbc,
real(kind=dp), dimension(:, :, :), intent(inout)  mab 
)

...

Parameters
la_max...
npgfa...
zeta...
rpgfa...
la_min...
lb_max...
npgfb...
zetb...
rpgfb...
lc_max...
rac...
rbc...
mab...

Definition at line 983 of file ai_moments.F.

Here is the caller graph for this function:

◆ diffop()

subroutine, public ai_moments::diffop ( integer, intent(in)  la_max,
integer, intent(in)  npgfa,
real(kind=dp), dimension(:), intent(in)  zeta,
real(kind=dp), dimension(:), intent(in)  rpgfa,
integer, intent(in)  la_min,
integer, intent(in)  lb_max,
integer, intent(in)  npgfb,
real(kind=dp), dimension(:), intent(in)  zetb,
real(kind=dp), dimension(:), intent(in)  rpgfb,
integer, intent(in)  lb_min,
real(kind=dp), dimension(3), intent(in)  rab,
real(kind=dp), dimension(:, :, :), intent(out)  difab 
)

This returns the derivative of the overlap integrals [a|b], with respect to the position of the primitive on the left, i.e. [da/dR_ai|b] = 2*zeta*[a+1i|b] - Ni(a)[a-1i|b].

Parameters
la_max...
npgfa...
zeta...
rpgfa...
la_min...
lb_max...
npgfb...
zetb...
rpgfb...
lb_min...
rab...
difab...

Definition at line 1515 of file ai_moments.F.

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

◆ diff_momop()

subroutine, public ai_moments::diff_momop ( integer, intent(in)  la_max,
integer, intent(in)  npgfa,
real(kind=dp), dimension(:), intent(in)  zeta,
real(kind=dp), dimension(:), intent(in)  rpgfa,
integer, intent(in)  la_min,
integer, intent(in)  lb_max,
integer, intent(in)  npgfb,
real(kind=dp), dimension(:), intent(in)  zetb,
real(kind=dp), dimension(:), intent(in)  rpgfb,
integer, intent(in)  lb_min,
integer, intent(in)  order,
real(kind=dp), dimension(3), intent(in)  rac,
real(kind=dp), dimension(3), intent(in)  rbc,
real(kind=dp), dimension(:, :, :, :), intent(out)  difmab,
real(kind=dp), dimension(:, :, :), optional, pointer  mab_ext 
)

This returns the derivative of the moment integrals [a|\mu|b], with respect to the position of the primitive on the left, i.e. [da/dR_ai|\mu|b] = 2*zeta*[a+1i|\mu|b] - Ni(a)[a-1i|\mu|b] order indicates the max order of the moment operator to be calculated 1: dipole 2: quadrupole ...

Parameters
la_max...
npgfa...
zeta...
rpgfa...
la_min...
lb_max...
npgfb...
zetb...
rpgfb...
lb_min...
order...
rac...
rbc...
difmab...
mab_ext...
Note

Definition at line 1580 of file ai_moments.F.

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

◆ dipole_force()

subroutine, public ai_moments::dipole_force ( integer, intent(in)  la_max,
integer, intent(in)  npgfa,
real(kind=dp), dimension(:), intent(in)  zeta,
real(kind=dp), dimension(:), intent(in)  rpgfa,
integer, intent(in)  la_min,
integer, intent(in)  lb_max,
integer, intent(in)  npgfb,
real(kind=dp), dimension(:), intent(in)  zetb,
real(kind=dp), dimension(:), intent(in)  rpgfb,
integer, intent(in)  lb_min,
integer, intent(in)  order,
real(kind=dp), dimension(3), intent(in)  rac,
real(kind=dp), dimension(3), intent(in)  rbc,
real(kind=dp), dimension(:, :), intent(in)  pab,
real(kind=dp), dimension(:, :), intent(inout)  forcea,
real(kind=dp), dimension(:, :), intent(inout)  forceb 
)

This returns the derivative of the dipole integrals [a|x|b], with respect to the position of the primitive on the left and right, i.e. [da/dR_ai|\mu|b] = 2*zeta*[a+1i|\mu|b] - Ni(a)[a-1i|\mu|b].

Parameters
la_max...
npgfa...
zeta...
rpgfa...
la_min...
lb_max...
npgfb...
zetb...
rpgfb...
lb_min...
order...
rac...
rbc...
pab...
forcea...
forceb...
Note

Definition at line 1664 of file ai_moments.F.

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