1
0
mirror of https://gitlab.com/scemama/qp_plugins_scemama.git synced 2024-06-02 03:15:25 +02:00
qp_plugins_scemama/devel/import_integrals/import_integrals_ao.irp.f
2020-03-24 09:33:51 +01:00

191 lines
4.4 KiB
Fortran

program import_integrals_ao
print *, 'Number of AOs?'
read(*,*) ao_num
TOUCH ao_num
call run
end
subroutine run
use map_module
implicit none
BEGIN_DOC
! Program to import integrals in the AO 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
!
! T.qp : kinetic energy integrals
!
! S.qp : overlap matrix
!
! P.qp : pseudopotential integrals
!
! V.qp : electron-nucleus potential
!
! W.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(:)
call ezfio_set_ao_basis_ao_num(ao_num)
allocate (A(ao_num,ao_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))
call ezfio_set_nuclei_io_nuclear_repulsion('Read')
endif
A = 0.d0
i = 1
j = 1
iunit = getunitandopen('T.qp','r')
do
read (iunit,*,end=10) i,j, integral
if (i<0 .or. i>ao_num) then
print *, i
stop 'i out of bounds in T.qp'
endif
if (j<0 .or. j>ao_num) then
print *, j
stop 'j out of bounds in T.qp'
endif
A(i,j) = integral
enddo
10 continue
close(iunit)
call ezfio_set_ao_one_e_ints_ao_integrals_kinetic(A(1:ao_num, 1:ao_num))
call ezfio_set_ao_one_e_ints_io_ao_integrals_kinetic('Read')
A = 0.d0
i = 1
j = 1
iunit = getunitandopen('S.qp','r')
do
read (iunit,*,end=11) i,j, integral
if (i<0 .or. i>ao_num) then
print *, i
stop 'i out of bounds in S.qp'
endif
if (j<0 .or. j>ao_num) then
print *, j
stop 'j out of bounds in S.qp'
endif
A(i,j) = integral
enddo
11 continue
close(iunit)
call ezfio_set_ao_one_e_ints_ao_integrals_overlap(A(1:ao_num, 1:ao_num))
call ezfio_set_ao_one_e_ints_io_ao_integrals_overlap('Read')
A = 0.d0
i = 1
j = 1
iunit = getunitandopen('P.qp','r')
do
if (i<0 .or. i>ao_num) then
print *, i
stop 'i out of bounds in P.qp'
endif
if (j<0 .or. j>ao_num) then
print *, j
stop 'j out of bounds in P.qp'
endif
read (iunit,*,end=14) i,j, integral
A(i,j) = integral
enddo
14 continue
close(iunit)
call ezfio_set_ao_one_e_ints_ao_integrals_pseudo(A(1:ao_num,1:ao_num))
call ezfio_set_ao_one_e_ints_io_ao_integrals_pseudo('Read')
A = 0.d0
i = 1
j = 1
iunit = getunitandopen('V.qp','r')
do
read (iunit,*,end=12) i,j, integral
if (i<0 .or. i>ao_num) then
print *, i
stop 'i out of bounds in V.qp'
endif
if (j<0 .or. j>ao_num) then
print *, j
stop 'j out of bounds in V.qp'
endif
A(i,j) = integral
enddo
12 continue
close(iunit)
call ezfio_set_ao_one_e_ints_ao_integrals_n_e(A(1:ao_num, 1:ao_num))
call ezfio_set_ao_one_e_ints_io_ao_integrals_n_e('Read')
allocate(buffer_i(ao_num**3), buffer_values(ao_num**3))
iunit = getunitandopen('W.qp','r')
n_integrals=0
i = 1
j = 1
k = 1
l = 1
buffer_values = 0.d0
do
read (iunit,*,end=13) i,j,k,l, integral
if (i<0 .or. i>ao_num) then
print *, i
stop 'i out of bounds in W.qp'
endif
if (j<0 .or. j>ao_num) then
print *, j
stop 'j out of bounds in W.qp'
endif
if (k<0 .or. k>ao_num) then
print *, k
stop 'k out of bounds in W.qp'
endif
if (l<0 .or. l>ao_num) then
print *, l
stop 'l out of bounds in W.qp'
endif
n_integrals += 1
call two_e_integrals_index(i, j, k, l, buffer_i(n_integrals) )
buffer_values(n_integrals) = integral
if (n_integrals == size(buffer_i)) then
call insert_into_ao_integrals_map(n_integrals,buffer_i,buffer_values)
n_integrals = 0
endif
enddo
13 continue
close(iunit)
if (n_integrals > 0) then
call insert_into_ao_integrals_map(n_integrals,buffer_i,buffer_values)
endif
call map_sort(ao_integrals_map)
call map_unique(ao_integrals_map)
call map_save_to_disk(trim(ezfio_filename)//'/work/ao_ints',ao_integrals_map)
call ezfio_set_ao_two_e_ints_io_ao_two_e_integrals('Read')
end