@brief Calculate AO overlap between two different geometries @param[in,out] infos Information structure containing all necessary data
This subroutine calculates the overlap between basis sets (Atomic Orbitals) of two different geometries. It computes both AO and derived MO overlaps and stores the results in the infos structure.
Note
Output AO overlap is written to mol.data["OQP::overlap_ao_non_orthogonal"] Output MO overlap is written to mol.data["OQP::overlap_mo_non_orthogonal"] Current geometry (xyz) is taken from mol.data["OQP::xyz"] Current MO coefficients (mo_a) are taken from mol.data["OQP::VEC_MO_A"] Current MO energies (e_a) are taken from mol.data["OQP::E_MO_A"] Old geometry (xyz_old) is taken from mol.data["OQP::xyz_old"] Old MO coefficients (mo_a_old) are taken from mol.data["OQP::VEC_MO_A_old"] Old MO energies (e_a_old) are taken from mol.data["OQP::E_MO_A_old"]
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(information), | intent(inout), | target | :: | infos |
Information structure containing all necessary data |
subroutine get_structures_ao_overlap(infos) use precision, only: dp use io_constants, only: iw use oqp_tagarray_driver use types, only: information use strings, only: Cstring, fstring use basis_tools, only: basis_set use atomic_structure_m, only: atomic_structure use messages, only: show_message, with_abort use int1, only: basis_overlap use constants, only: tol_int use util, only: measure_time implicit none character(len=*), parameter :: subroutine_name = "get_structures_ao_overlap" !> Information structure containing all necessary data type(information), target, intent(inout) :: infos type(basis_set), pointer :: basis type(basis_set), allocatable :: basis_old type(atomic_structure), allocatable, target :: atoms_old integer :: i, nbf ! Tagarray definitions and data pointers character(len=*), parameter :: tags_general(*) = (/ character(len=80) :: & OQP_XYZ_old, OQP_VEC_MO_A, OQP_E_MO_A, OQP_VEC_MO_A_old, OQP_E_MO_A_old /) character(len=*), parameter :: tags_alloc(*) = (/ character(len=80) :: & OQP_overlap_mo, OQP_overlap_ao/) real(kind=dp), pointer :: xyz_old(:,:), overlap_ao_out(:,:), overlap_mo_out(:,:), & mo_a(:,:), mo_a_old(:,:), e_a(:), e_a_old(:) ! Open log file open (unit=IW, file=infos%log_filename, position="append") ! Initialize basis sets and atomic structures allocate(basis_old, source=infos%basis) allocate(atoms_old, source=infos%atoms) basis => infos%basis basis%atoms => infos%atoms nbf = basis%nbf ! Allocate and prepare data for output call infos%dat%remove_records(tags_alloc) call infos%dat%reserve_data(OQP_overlap_mo, TA_TYPE_REAL64, & nbf*nbf, (/ nbf, nbf /), comment=OQP_overlap_mo_comment) call infos%dat%reserve_data(OQP_overlap_ao, TA_TYPE_REAL64, & nbf*nbf, (/ nbf, nbf /), comment=OQP_overlap_ao_comment) call data_has_tags(infos%dat, tags_alloc, module_name, subroutine_name, with_abort) call tagarray_get_data(infos%dat, OQP_overlap_mo, overlap_mo_out) call tagarray_get_data(infos%dat, OQP_overlap_ao, overlap_ao_out) ! Load data from python level call data_has_tags(infos%dat, tags_general, module_name, subroutine_name, with_abort) call tagarray_get_data(infos%dat, OQP_xyz_old, xyz_old) call tagarray_get_data(infos%dat, OQP_VEC_MO_A, mo_a) call tagarray_get_data(infos%dat, OQP_E_MO_A, e_a) call tagarray_get_data(infos%dat, OQP_VEC_MO_A_old, mo_a_old) call tagarray_get_data(infos%dat, OQP_E_MO_A_old, e_a_old) ! Apply old geometry to old basis atoms_old%xyz = xyz_old basis_old%atoms => atoms_old call print_geo(basis_old, "Previous geometry") call print_geo(basis, " Current geometry") overlap_ao_out = 0.0_dp overlap_mo_out = 0.0_dp ! Calculate AO overlap between old and new basis sets call basis_overlap(overlap_ao_out, basis, basis_old, tol=log(10.0_dp)*tol_int) ! Normalize AO overlap do i = 1, nbf overlap_ao_out(:,i) = overlap_ao_out(:,i) * basis%bfnrm(i) * basis_old%bfnrm(:) end do ! Calculate MO overlap: <old(I)|new(J')> = transpose[Cold(AI)] . <old(A)|new(A')> . Cnew(A'J') call mo_overlap(overlap_mo_out, mo_a, mo_a_old, overlap_ao_out, nbf) ! Output results: overlap between old and current MOs call print_results(overlap_mo_out, e_a, e_a_old, nbf, & infos%mol_prop%nelec_a, iw) ! Print timings call measure_time(print_total=1, log_unit=iw) call flush(iw) close(iw) end subroutine get_structures_ao_overlap