gaussian Module

read Gaussian Basis set and wavefunction interfaces to libcint to calculate integrals



Subroutines

public subroutine print_gto_info(bas, atm, env, output)

prints info on Gaussian basis set

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: bas(:,:)
integer, intent(in) :: atm(:,:)
real(kind=long), intent(in) :: env(:)
integer, intent(in) :: output

public subroutine construct_gaussian_basis_set(a, env, atm, bas, filename, filename2, shellinset)

constructs Gaussian basis set

Arguments

Type IntentOptional Attributes Name
type(atom_set), intent(in) :: a(:)
real(kind=long), intent(inout), allocatable :: env(:)
integer, intent(inout), allocatable :: atm(:,:)
integer, intent(inout), allocatable :: bas(:,:)
character(len=*), intent(in) :: filename
character(len=*), intent(in), optional :: filename2
integer, intent(out), optional, allocatable :: shellinset(:)

public subroutine calculate_hkin_gaussian_basis_set(env, atm, bas, first_bf_of_shell, hkin)

calculates the kinetic energy matrix

Arguments

Type IntentOptional Attributes Name
real(kind=long), intent(in) :: env(:)
integer, intent(in) :: atm(:,:)
integer, intent(in) :: bas(:,:)
integer, intent(in) :: first_bf_of_shell(:)
real(kind=long), intent(out) :: hkin(:,:)

public subroutine calculate_vcoul_gaussian_basis_set(env, atm, bas, first_bf_of_shell, pos, vcoul)

calculates the external potential energy matrix

Arguments

Type IntentOptional Attributes Name
real(kind=long), intent(in) :: env(:)
integer, intent(in) :: atm(:,:)
integer, intent(in) :: bas(:,:)
integer, intent(in) :: first_bf_of_shell(:)
real(kind=long), intent(in) :: pos(3)
real(kind=long), intent(out) :: vcoul(:,:)

public subroutine testlibcint()

make some simple tests for the library it checks that the correct order of p basis functions is used

Arguments

None

public subroutine calculate_vext_gaussian_basis_set(env, atm, bas, first_bf_of_shell, vext)

computes the interaction matrix elements of the electrons with nuclei potential

Arguments

Type IntentOptional Attributes Name
real(kind=long), intent(in) :: env(:)
integer, intent(in) :: atm(:,:)
integer, intent(in) :: bas(:,:)
integer, intent(in) :: first_bf_of_shell(:)
real(kind=long), intent(out) :: vext(:,:)

public subroutine calculate_overlap_gaussian_basis_set(env, atm, bas, first_bf_of_shell, s)

computes the overlap matrix

Arguments

Type IntentOptional Attributes Name
real(kind=long), intent(in) :: env(:)
integer, intent(in) :: atm(:,:)
integer, intent(in) :: bas(:,:)
integer, intent(in) :: first_bf_of_shell(:)
real(kind=long), intent(out) :: s(:,:)

public subroutine calculate_teints_gaussian_basis_set(env, atm, bas, first_bf_of_shell, n_basis, teint, accuracy)

computes all two-electron integrals

Arguments

Type IntentOptional Attributes Name
real(kind=long), intent(in) :: env(:)
integer, intent(in) :: atm(:,:)
integer, intent(in) :: bas(:,:)
integer, intent(in) :: first_bf_of_shell(:)
integer, intent(in) :: n_basis
real(kind=long), intent(out), allocatable :: teint(:)
real(kind=long), optional :: accuracy

public subroutine calculate_1contracted_teints_direct(env, atm, bas, first_bf_of_shell, n_basis, lc, c_teint, accuracy)

returns two-electron integrals, where the first index has been contracted with LCAO orbital coefficients LC

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=long), intent(in) :: env(:)
integer, intent(in) :: atm(:,:)
integer, intent(in) :: bas(:,:)
integer, intent(in) :: first_bf_of_shell(:)
integer, intent(in) :: n_basis
real(kind=long), intent(in) :: lc(:,:)
real(kind=long), intent(out) :: c_teint(size(lc,2),n_basis,n_basis*(n_basis+1)/2)
real(kind=long), optional :: accuracy

public subroutine calculate_dipole_gaussian_basis_set(env, atm, bas, first_bf_of_shell, dipole)

calculates dipole matrix elements

Arguments

Type IntentOptional Attributes Name
real(kind=long), intent(in) :: env(:)
integer, intent(in) :: atm(:,:)
integer, intent(in) :: bas(:,:)
integer, intent(in) :: first_bf_of_shell(:)
real(kind=long), intent(out) :: dipole(:,:,:)

public subroutine construct_fock_matrix_direct(env, atm, bas, first_bf_of_shell, dmat, h0, fout, accuracy)

This procedure constructs the RHF Fock Matrix assumes that density matrix is symmetric

Arguments

Type IntentOptional Attributes Name
real(kind=long), intent(in) :: env(:)
integer, intent(in) :: atm(:,:)
integer, intent(in) :: bas(:,:)
integer, intent(in) :: first_bf_of_shell(:)
real(kind=long), intent(in) :: dmat(:,:)
real(kind=long), intent(in) :: h0(:,:)
real(kind=long), intent(out) :: fout(:,:)
real(kind=long), optional :: accuracy

public subroutine construct_uhf_fock_matrix_direct(env, atm, bas, first_bf_of_shell, dmata, dmatb, h0, fouta, foutb, accuracy)

construct the UHF Fock Matrix assumes that density matrix is symmetric

Arguments

Type IntentOptional Attributes Name
real(kind=long), intent(in) :: env(:)
integer, intent(in) :: atm(:,:)
integer, intent(in) :: bas(:,:)
integer, intent(in) :: first_bf_of_shell(:)
real(kind=long), intent(in) :: dmata(:,:)
real(kind=long), intent(in) :: dmatb(:,:)
real(kind=long), intent(in) :: h0(:,:)
real(kind=long), intent(out) :: fouta(:,:)
real(kind=long), intent(out) :: foutb(:,:)
real(kind=long), optional :: accuracy

public subroutine add_electrostatic_interaction(env, atm, bas, first_bf_of_shell, d, h)

this is for the direct SCF GTO calculation adds all electrostatic interaction matrix elements to the one-electron Hamiltonian H

Arguments

Type IntentOptional Attributes Name
real(kind=long), intent(in) :: env(:)
integer, intent(in) :: atm(:,:)
integer, intent(in) :: bas(:,:)
integer, intent(in) :: first_bf_of_shell(:)
real(kind=long), intent(in) :: d(:,:)
real(kind=long), intent(inout) :: h(:,:)

public subroutine overlap_gradient(env, atm, bas, first_bf_of_shell, ograd)

compute gradient of overlap matrix

Arguments

Type IntentOptional Attributes Name
real(kind=long), intent(in) :: env(:)
integer, intent(in) :: atm(:,:)
integer, intent(in) :: bas(:,:)
integer, intent(in) :: first_bf_of_shell(:)
real(kind=long), intent(out) :: ograd(:,:,:)

public subroutine overlap_energy_gradient_rohf(env, atm, bas, first_bf_of_shell, lc, fa, fb, n_occ, gradient)

computes the gradient due to basis overlap

Arguments

Type IntentOptional Attributes Name
real(kind=long), intent(in) :: env(:)
integer, intent(in) :: atm(:,:)
integer, intent(in) :: bas(:,:)
integer, intent(in) :: first_bf_of_shell(:)
real(kind=long), intent(in) :: lc(:,:)
real(kind=long), intent(in) :: fa(:,:)
real(kind=long), intent(in) :: fb(:,:)
integer, intent(in) :: n_occ(:)
real(kind=long), intent(out) :: gradient(:)

public subroutine overlap_energy_gradient_uhf(env, atm, bas, first_bf_of_shell, lca, lcb, fa, fb, n_occ, gradient)

computes the gradient due to basis overlap

Arguments

Type IntentOptional Attributes Name
real(kind=long), intent(in) :: env(:)
integer, intent(in) :: atm(:,:)
integer, intent(in) :: bas(:,:)
integer, intent(in) :: first_bf_of_shell(:)
real(kind=long), intent(in) :: lca(:,:)
real(kind=long), intent(in) :: lcb(:,:)
real(kind=long), intent(in) :: fa(:,:)
real(kind=long), intent(in) :: fb(:,:)
integer, intent(in) :: n_occ(:)
real(kind=long), intent(out) :: gradient(:)

public subroutine overlap_energy_gradient(env, atm, bas, first_bf_of_shell, lc, ce, f_occ, gradient)

computes derivative of overlap matrix with respect to geometry

Arguments

Type IntentOptional Attributes Name
real(kind=long), intent(in) :: env(:)
integer, intent(in) :: atm(:,:)
integer, intent(in) :: bas(:,:)
integer, intent(in) :: first_bf_of_shell(:)
real(kind=long), intent(in) :: lc(:,:)
real(kind=long), intent(in) :: ce(:)
real(kind=long), intent(in) :: f_occ(:)
real(kind=long), intent(out) :: gradient(:)

public subroutine one_electron_gradient(env, atm, bas, first_bf_of_shell, d, gradient)

returns the gradient resulting from one-electron contributions (h0)

Arguments

Type IntentOptional Attributes Name
real(kind=long), intent(in) :: env(:)
integer, intent(in) :: atm(:,:)
integer, intent(in) :: bas(:,:)
integer, intent(in) :: first_bf_of_shell(:)
real(kind=long), intent(in) :: d(:,:)
real(kind=long), intent(out) :: gradient(:)

public subroutine rinv_gradient(env, atm, bas, first_bf_of_shell, pos, at, gradient)

calculates where R is given by pos(1:3) and has atom index 'at'

Arguments

Type IntentOptional Attributes Name
real(kind=long), intent(in) :: env(:)
integer, intent(in) :: atm(:,:)
integer, intent(in) :: bas(:,:)
integer, intent(in) :: first_bf_of_shell(:)
real(kind=long), intent(in) :: pos(:)
integer, intent(in) :: at
real(kind=long), intent(out) :: gradient(:,:,:)

public subroutine rinv_gradient_int(env, atm, bas, first_bf_of_shell, pos, gradient)

calculates where R_ is given by pos(1:3)

Arguments

Type IntentOptional Attributes Name
real(kind=long), intent(in) :: env(:)
integer, intent(in) :: atm(:,:)
integer, intent(in) :: bas(:,:)
integer, intent(in) :: first_bf_of_shell(:)
real(kind=long), intent(in) :: pos(:)
real(kind=long), intent(out) :: gradient(:,:,:)

public subroutine two_electron_gradient(env, atm, bas, first_bf_of_shell, d, gradient, fock, accuracy)

returns the gradient resulting from two-electron contributions

Arguments

Type IntentOptional Attributes Name
real(kind=long), intent(in) :: env(:)
integer, intent(in) :: atm(:,:)
integer, intent(in) :: bas(:,:)
integer, intent(in) :: first_bf_of_shell(:)
real(kind=long), intent(in) :: d(:,:)
real(kind=long), intent(out) :: gradient(:)
logical, intent(in), optional :: fock
real(kind=long), intent(in), optional :: accuracy

public subroutine two_electron_gradient_rohf(env, atm, bas, first_bf_of_shell, da, db, gradient)

returns the gradient resulting from two-electron contributions for ROHF

Arguments

Type IntentOptional Attributes Name
real(kind=long), intent(in) :: env(:)
integer, intent(in) :: atm(:,:)
integer, intent(in) :: bas(:,:)
integer, intent(in) :: first_bf_of_shell(:)
real(kind=long), intent(in) :: da(:,:)
real(kind=long), intent(in) :: db(:,:)
real(kind=long), intent(out) :: gradient(:)

public subroutine jk_gradient(env, atm, bas, first_bf_of_shell, d, jd, kd, symmetricd, accuracy)

computes derivatives of J and K

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=long), intent(in) :: env(:)
integer, intent(in) :: atm(:,:)
integer, intent(in) :: bas(:,:)
integer, intent(in) :: first_bf_of_shell(:)
real(kind=long), intent(in) :: d(:,:)
real(kind=long), intent(out) :: jd(:,:,:)
real(kind=long), intent(out) :: kd(:,:,:)
logical, intent(in), optional :: symmetricd
real(kind=long), intent(in), optional :: accuracy

public subroutine h0_gradient(env, atm, bas, first_bf_of_shell, h_grad)

computes gradient of h0 $omp parallel shared(atm, bas, env, first_bf_of_shell) private(sh1,sh2,atom,atom1,atom2,shls,ret1,n1,n2,data1,m1,m2,i,j,factor,my_env,mu) default(none) reduction(+:v,t) $omp do $omp end do $omp end parallel

Arguments

Type IntentOptional Attributes Name
real(kind=long), intent(in) :: env(:)
integer, intent(in) :: atm(:,:)
integer, intent(in) :: bas(:,:)
integer, intent(in) :: first_bf_of_shell(:)
real(kind=long), intent(out) :: h_grad(:,:,:)