2020-03-23 17:34:43 +01:00
|
|
|
program export_integrals_mo
|
2021-03-17 14:00:56 +01:00
|
|
|
PROVIDE mo_two_e_integrals_in_map
|
2020-03-23 17:34:43 +01:00
|
|
|
call run
|
|
|
|
end
|
|
|
|
|
|
|
|
subroutine write_2d(name,A)
|
|
|
|
implicit none
|
|
|
|
character*(*), intent(in) :: name
|
|
|
|
double precision, intent(in) :: A(mo_num,mo_num)
|
|
|
|
integer :: i, j
|
|
|
|
integer :: iunit
|
|
|
|
integer :: getunitandopen
|
|
|
|
|
|
|
|
iunit = getunitandopen(name,'w')
|
|
|
|
do j=1,mo_num
|
|
|
|
do i=1,mo_num
|
|
|
|
if (A(i,j) /= 0.d0) then
|
|
|
|
write (iunit,*) i,j, A(i,j)
|
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
close(iunit)
|
|
|
|
end
|
|
|
|
|
|
|
|
subroutine run
|
|
|
|
use map_module
|
|
|
|
implicit none
|
|
|
|
BEGIN_DOC
|
|
|
|
! Program to import integrals in the MO basis.
|
|
|
|
!
|
|
|
|
! one-electron integrals : format is : i j value
|
|
|
|
! two-electron integrals : format is : i j k l value
|
|
|
|
! Dirac's notation is used : <ij|kl> is <r1 r2|r1 r2>
|
|
|
|
!
|
|
|
|
! These files are read:
|
|
|
|
!
|
|
|
|
! T_mo.qp : kinetic energy integrals
|
|
|
|
!
|
|
|
|
! P_mo.qp : pseudopotential integrals
|
|
|
|
!
|
|
|
|
! V_mo.qp : electron-nucleus potential
|
|
|
|
!
|
|
|
|
! W_mo.qp : electron repulsion integrals
|
|
|
|
!
|
|
|
|
END_DOC
|
|
|
|
|
|
|
|
integer :: iunit
|
|
|
|
integer :: getunitandopen
|
|
|
|
|
|
|
|
integer :: i,j,k,l
|
|
|
|
double precision :: integral
|
|
|
|
double precision, allocatable :: A(:,:)
|
|
|
|
|
|
|
|
integer :: n_integrals
|
|
|
|
integer(key_kind) :: idx
|
|
|
|
|
2021-03-17 14:08:24 +01:00
|
|
|
double precision, external :: mo_two_e_integral
|
2020-03-23 17:34:43 +01:00
|
|
|
|
2021-03-17 14:00:56 +01:00
|
|
|
|
2020-03-23 17:34:43 +01:00
|
|
|
allocate (A(mo_num,mo_num))
|
|
|
|
call ezfio_set_nuclei_nuclear_repulsion(nuclear_repulsion)
|
|
|
|
|
|
|
|
call write_2d('T_mo.qp',mo_kinetic_integrals)
|
|
|
|
call write_2d('P_mo.qp',mo_pseudo_integrals)
|
|
|
|
call write_2d('V_mo.qp',mo_integrals_n_e)
|
2021-03-12 13:13:49 +01:00
|
|
|
call write_2d('H_mo.qp',mo_one_e_integrals)
|
2020-03-23 17:34:43 +01:00
|
|
|
|
|
|
|
iunit = getunitandopen('W_mo.qp','w')
|
|
|
|
do l=1,mo_num
|
|
|
|
do k=1,mo_num
|
2021-03-17 14:08:24 +01:00
|
|
|
do j=l,mo_num
|
|
|
|
do i=k,mo_num
|
|
|
|
if (i==j .and. k<l) cycle
|
|
|
|
integral = mo_two_e_integral(i,j,k,l)
|
2020-03-23 17:34:43 +01:00
|
|
|
if (integral /= 0.d0) then
|
|
|
|
write (iunit,'(4(I5,2X),E22.15)') i,j,k,l, integral
|
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
enddo
|
|
|
|
close(iunit)
|
|
|
|
|
|
|
|
end
|