2020-03-04 14:43:22 +01:00
|
|
|
program import_integrals_mo
|
2020-03-04 16:17:22 +01:00
|
|
|
PROVIDE mo_num
|
2020-03-04 14:43:22 +01:00
|
|
|
call run
|
|
|
|
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:
|
|
|
|
!
|
|
|
|
! E.qp : Contains the nuclear repulsion energy
|
|
|
|
!
|
|
|
|
! Tmo.qp : kinetic energy integrals
|
|
|
|
!
|
|
|
|
! Smo.qp : overlap matrix
|
|
|
|
!
|
|
|
|
! Pmo.qp : pseudopotential integrals
|
|
|
|
!
|
|
|
|
! Vmo.qp : electron-nucleus potential
|
|
|
|
!
|
|
|
|
! Wmo.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), allocatable :: buffer_i(:)
|
|
|
|
real(integral_kind), allocatable :: buffer_values(:)
|
|
|
|
|
2020-03-04 16:17:22 +01:00
|
|
|
allocate(buffer_i(mo_num**3), buffer_values(mo_num**3))
|
2020-03-04 14:43:22 +01:00
|
|
|
allocate (A(mo_num,mo_num))
|
|
|
|
|
|
|
|
A(1,1) = huge(1.d0)
|
|
|
|
iunit = getunitandopen('E.qp','r')
|
|
|
|
read (iunit,*,end=9) A(1,1)
|
|
|
|
9 continue
|
|
|
|
close(iunit)
|
|
|
|
if (A(1,1) /= huge(1.d0)) then
|
|
|
|
call ezfio_set_nuclei_nuclear_repulsion(A(1,1))
|
2020-03-04 16:17:22 +01:00
|
|
|
call ezfio_set_nuclei_io_nuclear_repulsion('Read')
|
2020-03-04 14:43:22 +01:00
|
|
|
endif
|
|
|
|
|
|
|
|
A = 0.d0
|
|
|
|
iunit = getunitandopen('Tmo.qp','r')
|
|
|
|
do
|
|
|
|
read (iunit,*,end=10) i,j, integral
|
2020-03-04 16:17:22 +01:00
|
|
|
if (i<0 .or. i>mo_num) then
|
|
|
|
print *, i
|
|
|
|
stop 'i out of bounds in Tmo.qp'
|
|
|
|
endif
|
|
|
|
if (j<0 .or. j>mo_num) then
|
|
|
|
print *, j
|
|
|
|
stop 'j out of bounds in Tmo.qp'
|
|
|
|
endif
|
2020-03-04 14:43:22 +01:00
|
|
|
A(i,j) = integral
|
|
|
|
enddo
|
|
|
|
10 continue
|
|
|
|
close(iunit)
|
2020-03-04 16:17:22 +01:00
|
|
|
call ezfio_set_mo_one_e_ints_mo_integrals_kinetic(A)
|
|
|
|
call ezfio_set_mo_one_e_ints_io_mo_integrals_kinetic('Read')
|
2020-03-04 14:43:22 +01:00
|
|
|
|
|
|
|
A = 0.d0
|
|
|
|
iunit = getunitandopen('Pmo.qp','r')
|
|
|
|
do
|
|
|
|
read (iunit,*,end=14) i,j, integral
|
2020-03-04 16:17:22 +01:00
|
|
|
if (i<0 .or. i>mo_num) then
|
|
|
|
print *, i
|
|
|
|
stop 'i out of bounds in Pmo.qp'
|
|
|
|
endif
|
|
|
|
if (j<0 .or. j>mo_num) then
|
|
|
|
print *, j
|
|
|
|
stop 'j out of bounds in Pmo.qp'
|
|
|
|
endif
|
2020-03-04 14:43:22 +01:00
|
|
|
A(i,j) = integral
|
|
|
|
enddo
|
|
|
|
14 continue
|
|
|
|
close(iunit)
|
2020-03-04 16:17:22 +01:00
|
|
|
call ezfio_set_mo_one_e_ints_mo_integrals_pseudo(A)
|
|
|
|
call ezfio_set_mo_one_e_ints_io_mo_integrals_pseudo('Read')
|
2020-03-04 14:43:22 +01:00
|
|
|
|
|
|
|
A = 0.d0
|
|
|
|
iunit = getunitandopen('Vmo.qp','r')
|
|
|
|
do
|
|
|
|
read (iunit,*,end=12) i,j, integral
|
2020-03-04 16:17:22 +01:00
|
|
|
if (i<0 .or. i>mo_num) then
|
|
|
|
print *, i
|
|
|
|
stop 'i out of bounds in Vmo.qp'
|
|
|
|
endif
|
|
|
|
if (j<0 .or. j>mo_num) then
|
|
|
|
print *, j
|
|
|
|
stop 'j out of bounds in Vmo.qp'
|
|
|
|
endif
|
2020-03-04 14:43:22 +01:00
|
|
|
A(i,j) = integral
|
|
|
|
enddo
|
|
|
|
12 continue
|
|
|
|
close(iunit)
|
2020-03-04 16:17:22 +01:00
|
|
|
call ezfio_set_mo_one_e_ints_mo_integrals_e_n(A)
|
|
|
|
call ezfio_set_mo_one_e_ints_io_mo_integrals_e_n('Read')
|
2020-03-04 14:43:22 +01:00
|
|
|
|
|
|
|
iunit = getunitandopen('Wmo.qp','r')
|
|
|
|
n_integrals=0
|
|
|
|
buffer_values = 0.d0
|
|
|
|
do
|
|
|
|
read (iunit,*,end=13) i,j,k,l, integral
|
2020-03-04 16:17:22 +01:00
|
|
|
if (i<0 .or. i>mo_num) then
|
|
|
|
print *, i
|
|
|
|
stop 'i out of bounds in Wmo.qp'
|
|
|
|
endif
|
|
|
|
if (j<0 .or. j>mo_num) then
|
|
|
|
print *, j
|
|
|
|
stop 'j out of bounds in Wmo.qp'
|
|
|
|
endif
|
|
|
|
if (k<0 .or. k>mo_num) then
|
|
|
|
print *, k
|
|
|
|
stop 'k out of bounds in Wmo.qp'
|
|
|
|
endif
|
|
|
|
if (l<0 .or. l>mo_num) then
|
|
|
|
print *, l
|
|
|
|
stop 'l out of bounds in Wmo.qp'
|
|
|
|
endif
|
2020-03-04 14:43:22 +01:00
|
|
|
n_integrals += 1
|
2020-03-04 16:17:22 +01:00
|
|
|
call mo_two_e_integrals_index(i, j, k, l, buffer_i(n_integrals) )
|
2020-03-04 14:43:22 +01:00
|
|
|
buffer_values(n_integrals) = integral
|
|
|
|
if (n_integrals == size(buffer_i)) then
|
2020-03-04 16:17:22 +01:00
|
|
|
call insert_into_mo_integrals_map(n_integrals, buffer_i, buffer_values, &
|
|
|
|
real(mo_integrals_threshold,integral_kind))
|
2020-03-04 14:43:22 +01:00
|
|
|
n_integrals = 0
|
|
|
|
endif
|
|
|
|
enddo
|
|
|
|
13 continue
|
|
|
|
close(iunit)
|
|
|
|
|
|
|
|
if (n_integrals > 0) then
|
2020-03-04 16:17:22 +01:00
|
|
|
call insert_into_mo_integrals_map(n_integrals, buffer_i, buffer_values, &
|
|
|
|
real(mo_integrals_threshold,integral_kind))
|
2020-03-04 14:43:22 +01:00
|
|
|
endif
|
|
|
|
|
2020-03-04 16:17:22 +01:00
|
|
|
call map_merge(mo_integrals_map)
|
2020-03-04 14:43:22 +01:00
|
|
|
|
|
|
|
call map_save_to_disk(trim(ezfio_filename)//'/work/mo_ints',mo_integrals_map)
|
|
|
|
call ezfio_set_mo_two_e_ints_io_mo_two_e_integrals('Read')
|
|
|
|
|
|
|
|
end
|