BASIS_FUNC module to define basis-set functions from atomic orbitals represented on multicenter molecular grids written by Son, Sang-Kil and Hao, Yajiang in Mar. 2014
USAGE: type(Basis) :: B
Ref) Becke, J. Chem. Phys. 88, 2547 (1988) Becke & Dickson, J. Chem. Phys. 89, 2993 (1988) Becke & Dickson, J. Chem. Phys. 92, 3610 (1990) Dickson & Becke, J. Chem. Phys. 99, 3898 (1993)
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
integer, | public, | parameter | :: | n_symop | = | 8 | |
integer, | public, | parameter | :: | symop_e | = | 1 |
1 E |
integer, | public, | parameter | :: | symop_c2x | = | 2 |
2 C2X |
integer, | public, | parameter | :: | symop_c2y | = | 3 |
3 C2Y |
integer, | public, | parameter | :: | symop_c2z | = | 4 |
4 C2Z |
integer, | public, | parameter | :: | symop_i | = | 5 |
5 i |
integer, | public, | parameter | :: | symop_sigxy | = | 6 |
6 SIGXY |
integer, | public, | parameter | :: | symop_sigxz | = | 7 |
7 SIGXZ |
integer, | public, | parameter | :: | symop_sigyz | = | 8 |
8 SIGYZ |
character(len=5), | public, | parameter | :: | sym_op_name(8) | = | (/'E ', 'C2X ', 'C2Y ', 'C2Z ', 'i ', 'SIGXY', 'SIGXZ', 'SIGYZ'/) |
basis information
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
logical, | public | :: | yes_initialized | = | .false. | ||
integer, | public | :: | n_basis |
total # of basis functions |
|||
integer, | public | :: | nr_basis |
total # of basis functions projected out |
|||
real(kind=long), | public, | allocatable | :: | func(:,:) |
(1:N_grid,1:N_basis) |
||
real(kind=long), | public, | allocatable | :: | dfdr(:,:,:) |
(1:N_grid,1:N_basis,3) |
||
real(kind=long), | public, | allocatable | :: | d2f(:,:) |
(1:N_grid,1:N_basis) |
||
logical, | public, | allocatable | :: | yes_skip(:,:) |
(1:N_nuc,1:N_basis) skip if all basis functions (and its deriv.) are zero for a given atom A it helps to speed up because we know that B%func are zero if r > A%rmax and we don't need to evaluate anything beyond A%rmax --> initialized in place_AO_on_Grid3D() |
||
integer, | public, | allocatable | :: | info(:,:) |
basis function info B%info(i,1) : atomic index for the i-th basis B%info(i,2) : nuclear charge for the i-th basis B%info(i,3) : n of AO for the i-th basis B%info(i,4) : l of AO for the i-th basis B%info(i,5) : m of AO for the i-th basis |
||
character(len=10), | public, | allocatable | :: | label(:) | |||
real(kind=long), | public, | allocatable | :: | s(:,:) |
overlap matrix |
||
real(kind=long), | public, | allocatable | :: | ssqrt(:,:) |
sqrt overlap matrix |
||
real(kind=long), | public, | allocatable | :: | sinvsqrt(:,:) |
sqrt of inverse overlap matrix |
||
real(kind=long), | public, | allocatable | :: | sinv(:,:) |
inverse overlap matrix |
||
real(kind=long), | public, | allocatable | :: | x(:,:) | |||
real(kind=long), | public, | allocatable | :: | xxts(:,:) | |||
real(kind=long), | public, | allocatable | :: | hkin(:,:) |
Kinetic operator matrix, (1:N_basis,1:N_basis) |
||
real(kind=long), | public, | allocatable | :: | vext(:,:) |
Nuclear potential (1:N_basis,1:N_basis) |
||
real(kind=long), | public, | allocatable | :: | h0(:,:) |
Kinetic + nuclear potential |
||
real(kind=long), | public, | allocatable | :: | teint(:) |
Two-Electron integrals |
||
logical, | public, | allocatable | :: | yes_overlap(:,:) |
yes_overlap(k_nuc,ij): whether a multiplication of two basis functions (i and j) contributes to the k_nuc-th atomic sphere [k_nuc=1...N_nuc] If this variable is yes, then the k_nuc-th atom will be included in multicenter integration involving the i-th and j-th basis functions. If no, then it is not included, which means "save the CPU time!" |
||
real(kind=long), | public, | allocatable | :: | overlap(:,:) | |||
real(kind=long), | public | :: | eps_overlap | ||||
real(kind=long), | public | :: | eps_overlap_single | ||||
real(kind=long), | public | :: | r_max_overlap | ||||
real(kind=long), | public | :: | error_single |
numerical accuracy for single-center and multicenter integration and integration optimization by truncation |
|||
real(kind=long), | public | :: | error_multi |
numerical accuracy for single-center and multicenter integration and integration optimization by truncation |
|||
real(kind=long), | public | :: | error_trunc |
numerical accuracy for single-center and multicenter integration and integration optimization by truncation |
|||
integer, | public, | allocatable | :: | first_bf_on_atom(:) | |||
integer, | public, | allocatable | :: | num_bf_on_atom(:) | |||
logical, | public | :: | yes_symmetry(n_symop) | ||||
integer, | public, | allocatable | :: | irrep_prod(:,:) | |||
real(kind=long), | public, | allocatable | :: | sym_bf(:,:) |
conversion matrix from symmetric bf to atomic orbitals |
||
integer, | public, | allocatable | :: | irrep_sym_bf(:) |
index of irreducible representation per sym bf |
||
integer, | public, | allocatable | :: | irrep_sym_bfr(:) |
index of irreducible representation per sym bf |
||
character(len=5), | public | :: | sym_label |
Schoenflies symbol |
|||
integer, | public | :: | nirrep |
number of irrep |
|||
character(len=4), | public, | allocatable | :: | irrep_name(:) |
name of irrep |
||
real(kind=long), | public | :: | minimal_s_eigenvalue |
minimal overlap matrix eigenvalue |
|||
real(kind=long), | public | :: | minimal_s_remaining_eigenvalue |
minimal overlap matrix eigenval that is not projected out |
|||
real(kind=long), | public, | allocatable | :: | env(:) | |||
integer, | public, | allocatable | :: | atm(:,:) | |||
integer, | public, | allocatable | :: | bas(:,:) | |||
integer, | public, | allocatable | :: | first_bf_of_shell(:) |
braket(): where i, j indicate i-th and j-th basis function braket_Laplacian():
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(grid3d), | intent(in) | :: | g | |||
type(basis), | intent(in) | :: | b | |||
integer, | intent(in) | :: | i | |||
integer, | intent(in) | :: | j | |||
real(kind=long), | intent(in) | :: | f(:) |
construct (or reconstruct) B from A and G all necessary informations are obtained from - A : atomic orbitals of individual atoms - G : molecular grids
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(atom_set), | intent(in) | :: | a(:) |
atom set |
||
type(grid3d), | intent(in) | :: | g |
grid |
||
type(basis), | intent(out) | :: | b |
basis |
||
real(kind=long), | intent(in), | optional | :: | eps_overlap |
truncation parameters |
|
real(kind=long), | intent(in), | optional | :: | eps_overlap_single |
truncation parameters |
|
logical, | intent(in), | optional | :: | yes_sym |
use symmetry? |
|
character(len=*), | intent(in), | optional | :: | gto_file_path |
path to file with gto informations |
|
logical, | intent(in), | optional | :: | yes_gradient |
prepare for calculating gradients? |
|
real(kind=long), | intent(in), | optional | :: | basissetoverlaplim |
truncation parameters |
deallocate basis structure
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(basis), | intent(inout) | :: | b |
print grid parameters and basis function (numerical atomic orbital) info it is supposed to be in read_param.f90, but it requires A(:) and G wait a minute... even it doesn't contain F!!! why should it be here? it has beend moved to basis_func
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(param), | intent(in) | :: | p | |||
type(atom_set), | intent(in) | :: | a(:) | |||
type(grid3d), | intent(in) | :: | g | |||
type(basis), | intent(in) | :: | b |
calculates the transition dipole matrix elements between basisfunctions and stores them in tdipole_bf( mu, nu, p)
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(param), | intent(in) | :: | p | |||
type(basis), | intent(in) | :: | b | |||
type(grid3d), | intent(in) | :: | g | |||
real(kind=long), | intent(out), | dimension(p%n_basis, p%n_basis, 3) | :: | tdipole_bf |
Do a Gram-Schmidt orthogonalization for the coefficients given by C according to the metric given by the overlap matrix B%S, i.e. such that < C(:,i) S(:,:) C(:,j) > = delta_i,j B are the basis set parameters
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(basis), | intent(in) | :: | b | |||
real(kind=long), | intent(inout), | dimension(:, :) | :: | c | ||
integer, | intent(in), | optional | :: | irrep_mo(:) |