mirror of
https://gitlab.com/scemama/qp_plugins_scemama.git
synced 2024-08-30 07:53:39 +02:00
212 lines
5.1 KiB
Fortran
212 lines
5.1 KiB
Fortran
program import_integrals_mo
|
|
PROVIDE mo_num
|
|
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
|
|
!
|
|
! T_mo.qp : kinetic energy integrals
|
|
!
|
|
! P_mo.qp : pseudopotential integrals
|
|
!
|
|
! V_mo.qp : electron-nucleus potential
|
|
!
|
|
! H_mo.qp : core hamiltonian. If hmo exists, the other 1-e files are not read
|
|
!
|
|
! 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
|
|
logical :: exists
|
|
integer(key_kind), allocatable :: buffer_i(:)
|
|
real(integral_kind), allocatable :: buffer_values(:)
|
|
|
|
allocate(buffer_i(mo_num**3), buffer_values(mo_num**3))
|
|
allocate (A(mo_num,mo_num))
|
|
PROVIDE mo_integrals_map
|
|
PROVIDE mo_integrals_threshold
|
|
|
|
! Nuclear repulsion
|
|
|
|
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
|
|
|
|
! One-electron integrals
|
|
|
|
exists = .False.
|
|
inquire(file='H_mo.qp',exist=exists)
|
|
if (exists) then
|
|
|
|
A = 0.d0
|
|
i = 1
|
|
j = 1
|
|
iunit = getunitandopen('H_mo.qp','r')
|
|
do
|
|
read (iunit,*,end=8) i,j, integral
|
|
if (i<0 .or. i>mo_num) then
|
|
print *, i
|
|
stop 'i out of bounds in H_mo.qp'
|
|
endif
|
|
if (j<0 .or. j>mo_num) then
|
|
print *, j
|
|
stop 'j out of bounds in H_mo.qp'
|
|
endif
|
|
A(i,j) = integral
|
|
enddo
|
|
8 continue
|
|
close(iunit)
|
|
call ezfio_set_mo_one_e_ints_mo_one_e_integrals(A)
|
|
call ezfio_set_mo_one_e_ints_io_mo_one_e_integrals('Read')
|
|
call ezfio_set_mo_one_e_ints_io_mo_integrals_kinetic('None')
|
|
call ezfio_set_mo_one_e_ints_io_mo_integrals_pseudo('None')
|
|
call ezfio_set_mo_one_e_ints_io_mo_integrals_n_e('None')
|
|
|
|
else
|
|
|
|
call ezfio_set_mo_one_e_ints_io_mo_one_e_integrals('None')
|
|
|
|
A = 0.d0
|
|
i = 1
|
|
j = 1
|
|
iunit = getunitandopen('T_mo.qp','r')
|
|
do
|
|
read (iunit,*,end=10) i,j, integral
|
|
if (i<0 .or. i>mo_num) then
|
|
print *, i
|
|
stop 'i out of bounds in T_mo.qp'
|
|
endif
|
|
if (j<0 .or. j>mo_num) then
|
|
print *, j
|
|
stop 'j out of bounds in T_mo.qp'
|
|
endif
|
|
A(i,j) = integral
|
|
enddo
|
|
10 continue
|
|
close(iunit)
|
|
call ezfio_set_mo_one_e_ints_mo_integrals_kinetic(A)
|
|
call ezfio_set_mo_one_e_ints_io_mo_integrals_kinetic('Read')
|
|
|
|
A = 0.d0
|
|
i = 1
|
|
j = 1
|
|
iunit = getunitandopen('P_mo.qp','r')
|
|
do
|
|
read (iunit,*,end=14) i,j, integral
|
|
if (i<0 .or. i>mo_num) then
|
|
print *, i
|
|
stop 'i out of bounds in P_mo.qp'
|
|
endif
|
|
if (j<0 .or. j>mo_num) then
|
|
print *, j
|
|
stop 'j out of bounds in P_mo.qp'
|
|
endif
|
|
A(i,j) = integral
|
|
enddo
|
|
14 continue
|
|
close(iunit)
|
|
call ezfio_set_mo_one_e_ints_mo_integrals_pseudo(A)
|
|
call ezfio_set_mo_one_e_ints_io_mo_integrals_pseudo('Read')
|
|
|
|
A = 0.d0
|
|
i = 1
|
|
j = 1
|
|
iunit = getunitandopen('V_mo.qp','r')
|
|
do
|
|
read (iunit,*,end=12) i,j, integral
|
|
if (i<0 .or. i>mo_num) then
|
|
print *, i
|
|
stop 'i out of bounds in V_mo.qp'
|
|
endif
|
|
if (j<0 .or. j>mo_num) then
|
|
print *, j
|
|
stop 'j out of bounds in V_mo.qp'
|
|
endif
|
|
A(i,j) = integral
|
|
enddo
|
|
12 continue
|
|
close(iunit)
|
|
call ezfio_set_mo_one_e_ints_mo_integrals_n_e(A)
|
|
call ezfio_set_mo_one_e_ints_io_mo_integrals_n_e('Read')
|
|
|
|
end if
|
|
|
|
|
|
! Two-electron integrals
|
|
|
|
iunit = getunitandopen('W_mo.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>mo_num) then
|
|
print *, i
|
|
stop 'i out of bounds in W_mo.qp'
|
|
endif
|
|
if (j<0 .or. j>mo_num) then
|
|
print *, j
|
|
stop 'j out of bounds in W_mo.qp'
|
|
endif
|
|
if (k<0 .or. k>mo_num) then
|
|
print *, k
|
|
stop 'k out of bounds in W_mo.qp'
|
|
endif
|
|
if (l<0 .or. l>mo_num) then
|
|
print *, l
|
|
stop 'l out of bounds in W_mo.qp'
|
|
endif
|
|
n_integrals += 1
|
|
call mo_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 map_append(mo_integrals_map, buffer_i, buffer_values, n_integrals)
|
|
n_integrals = 0
|
|
endif
|
|
enddo
|
|
13 continue
|
|
close(iunit)
|
|
|
|
if (n_integrals > 0) then
|
|
call map_append(mo_integrals_map, buffer_i, buffer_values, n_integrals)
|
|
endif
|
|
|
|
call map_sort(mo_integrals_map)
|
|
call map_unique(mo_integrals_map)
|
|
|
|
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
|
|
|
|
|