10
0
mirror of https://github.com/LCPQ/quantum_package synced 2024-11-07 22:53:57 +01:00
quantum_package/src/Bielec_integrals/map_integrals.irp.f

548 lines
14 KiB
Fortran
Raw Normal View History

2014-04-17 23:50:51 +02:00
use map_module
!! AO Map
!! ======
BEGIN_PROVIDER [ type(map_type), ao_integrals_map ]
implicit none
BEGIN_DOC
! AO integrals
END_DOC
2015-03-22 20:48:32 +01:00
integer(key_kind) :: key_max
integer(map_size_kind) :: sze
call bielec_integrals_index(ao_num,ao_num,ao_num,ao_num,key_max)
sze = key_max
2014-04-17 23:50:51 +02:00
call map_init(ao_integrals_map,sze)
2015-03-31 09:58:25 +02:00
write(output_bielec_integrals,*) 'AO map initialized'
2014-04-17 23:50:51 +02:00
END_PROVIDER
subroutine bielec_integrals_index(i,j,k,l,i1)
2015-03-22 20:48:32 +01:00
use map_module
2014-04-17 23:50:51 +02:00
implicit none
integer, intent(in) :: i,j,k,l
2015-03-22 20:48:32 +01:00
integer(key_kind), intent(out) :: i1
integer(key_kind) :: p,q,r,s,i2
2014-04-17 23:50:51 +02:00
p = min(i,k)
r = max(i,k)
p = p+ishft(r*r-r,-1)
q = min(j,l)
s = max(j,l)
q = q+ishft(s*s-s,-1)
i1 = min(p,q)
i2 = max(p,q)
i1 = i1+ishft(i2*i2-i2,-1)
end
2014-09-19 11:32:45 +02:00
subroutine bielec_integrals_index_reverse(i,j,k,l,i1)
2015-03-22 20:48:32 +01:00
use map_module
2014-09-19 11:32:45 +02:00
implicit none
2014-09-20 23:45:18 +02:00
integer, intent(out) :: i(8),j(8),k(8),l(8)
2015-03-22 20:48:32 +01:00
integer(key_kind), intent(in) :: i1
integer(key_kind) :: i2,i3
2014-09-20 23:45:18 +02:00
i = 0
2014-10-03 16:40:18 +02:00
i2 = ceiling(0.5d0*(dsqrt(8.d0*dble(i1)+1.d0)-1.d0))
l(1) = ceiling(0.5d0*(dsqrt(8.d0*dble(i2)+1.d0)-1.d0))
2014-09-22 23:38:49 +02:00
i3 = i1 - ishft(i2*i2-i2,-1)
2014-10-03 16:40:18 +02:00
k(1) = ceiling(0.5d0*(dsqrt(8.d0*dble(i3)+1.d0)-1.d0))
2014-09-20 23:45:18 +02:00
j(1) = i2 - ishft(l(1)*l(1)-l(1),-1)
i(1) = i3 - ishft(k(1)*k(1)-k(1),-1)
!ijkl
2014-09-22 23:38:49 +02:00
i(2) = i(1) !ilkj
j(2) = l(1)
k(2) = k(1)
l(2) = j(1)
2014-09-20 23:45:18 +02:00
2014-09-22 23:38:49 +02:00
i(3) = k(1) !kjil
j(3) = j(1)
k(3) = i(1)
l(3) = l(1)
2014-09-20 23:45:18 +02:00
i(4) = k(1) !klij
j(4) = l(1)
k(4) = i(1)
l(4) = j(1)
i(5) = j(1) !jilk
j(5) = i(1)
k(5) = l(1)
l(5) = k(1)
2014-09-22 23:38:49 +02:00
i(6) = j(1) !jkli
j(6) = k(1)
k(6) = l(1)
l(6) = i(1)
2014-09-20 23:45:18 +02:00
2014-09-22 23:38:49 +02:00
i(7) = l(1) !lijk
j(7) = i(1)
k(7) = j(1)
l(7) = k(1)
2014-09-20 23:45:18 +02:00
i(8) = l(1) !lkji
j(8) = k(1)
k(8) = j(1)
l(8) = i(1)
integer :: ii, jj
do ii=2,8
do jj=1,ii-1
if ( (i(ii) == i(jj)).and. &
(j(ii) == j(jj)).and. &
(k(ii) == k(jj)).and. &
(l(ii) == l(jj)) ) then
2014-09-22 20:19:56 +02:00
i(ii) = 0
2014-09-20 23:45:18 +02:00
exit
endif
enddo
enddo
2014-10-27 12:40:30 +01:00
do ii=1,8
if (i(ii) /= 0) then
call bielec_integrals_index(i(ii),j(ii),k(ii),l(ii),i2)
if (i1 /= i2) then
print *, i1, i2
print *, i(ii), j(jj), k(jj), l(jj)
stop 'bielec_integrals_index_reverse failed'
endif
endif
enddo
2014-09-20 23:45:18 +02:00
2014-09-19 11:32:45 +02:00
end
2014-04-17 23:50:51 +02:00
double precision function get_ao_bielec_integral(i,j,k,l,map)
use map_module
implicit none
BEGIN_DOC
! Gets one AO bi-electronic integral from the AO map
END_DOC
integer, intent(in) :: i,j,k,l
2015-03-22 20:48:32 +01:00
integer(key_kind) :: idx
2014-04-17 23:50:51 +02:00
type(map_type), intent(inout) :: map
real(integral_kind) :: tmp
PROVIDE ao_bielec_integrals_in_map
!DIR$ FORCEINLINE
2014-09-19 11:35:53 +02:00
if (ao_overlap_abs(i,k)*ao_overlap_abs(j,l) < ao_integrals_threshold ) then
tmp = 0.d0
else if (ao_bielec_integral_schwartz(i,k)*ao_bielec_integral_schwartz(j,l) < ao_integrals_threshold) then
tmp = 0.d0
else
call bielec_integrals_index(i,j,k,l,idx)
call map_get(map,idx,tmp)
endif
2014-04-17 23:50:51 +02:00
get_ao_bielec_integral = tmp
end
subroutine get_ao_bielec_integrals(j,k,l,sze,out_val)
use map_module
BEGIN_DOC
! Gets multiple AO bi-electronic integral from the AO map .
! All i are retrieved for j,k,l fixed.
END_DOC
implicit none
integer, intent(in) :: j,k,l, sze
real(integral_kind), intent(out) :: out_val(sze)
integer :: i
2015-03-22 20:48:32 +01:00
integer(key_kind) :: hash
2014-04-17 23:50:51 +02:00
double precision :: thresh
2014-10-15 15:19:34 +02:00
PROVIDE ao_bielec_integrals_in_map ao_integrals_map
2014-04-17 23:50:51 +02:00
thresh = ao_integrals_threshold
if (ao_overlap_abs(j,l) < thresh) then
out_val = 0.d0
return
endif
do i=1,sze
if (ao_overlap_abs(i,k)*ao_overlap_abs(j,l) < thresh ) then
out_val(i) = 0.d0
2014-09-19 11:35:53 +02:00
else if (ao_bielec_integral_schwartz(i,k)*ao_bielec_integral_schwartz(j,l) < thresh) then
out_val(i)=0.d0
2014-04-17 23:50:51 +02:00
else
!DIR$ FORCEINLINE
call bielec_integrals_index(i,j,k,l,hash)
call map_get(ao_integrals_map, hash, out_val(i))
endif
enddo
end
subroutine get_ao_bielec_integrals_non_zero(j,k,l,sze,out_val,out_val_index,non_zero_int)
use map_module
implicit none
BEGIN_DOC
! Gets multiple AO bi-electronic integral from the AO map .
! All non-zero i are retrieved for j,k,l fixed.
END_DOC
integer, intent(in) :: j,k,l, sze
real(integral_kind), intent(out) :: out_val(sze)
integer, intent(out) :: out_val_index(sze),non_zero_int
integer :: i
2015-03-22 20:48:32 +01:00
integer(key_kind) :: hash
2014-04-17 23:50:51 +02:00
double precision :: thresh,tmp
PROVIDE ao_bielec_integrals_in_map
thresh = ao_integrals_threshold
2014-08-21 11:14:30 +02:00
non_zero_int = 0
2014-04-17 23:50:51 +02:00
if (ao_overlap_abs(j,l) < thresh) then
out_val = 0.d0
return
endif
2014-08-21 11:14:30 +02:00
2014-04-17 23:50:51 +02:00
non_zero_int = 0
do i=1,sze
2014-09-19 11:35:53 +02:00
integer, external :: ao_l4
double precision, external :: ao_bielec_integral
2014-04-17 23:50:51 +02:00
!DIR$ FORCEINLINE
2014-09-19 11:35:53 +02:00
if (ao_bielec_integral_schwartz(i,k)*ao_bielec_integral_schwartz(j,l) < thresh) then
cycle
endif
2014-04-17 23:50:51 +02:00
call bielec_integrals_index(i,j,k,l,hash)
call map_get(ao_integrals_map, hash,tmp)
if (dabs(tmp) < thresh ) cycle
non_zero_int = non_zero_int+1
out_val_index(non_zero_int) = i
out_val(non_zero_int) = tmp
enddo
end
2015-03-22 20:48:32 +01:00
function get_ao_map_size()
2014-04-17 23:50:51 +02:00
implicit none
2015-03-22 20:48:32 +01:00
integer (map_size_kind) :: get_ao_map_size
2014-04-17 23:50:51 +02:00
BEGIN_DOC
! Returns the number of elements in the AO map
END_DOC
get_ao_map_size = ao_integrals_map % n_elements
end
subroutine clear_ao_map
implicit none
BEGIN_DOC
! Frees the memory of the AO map
END_DOC
call map_deinit(ao_integrals_map)
FREE ao_integrals_map
end
!! MO Map
!! ======
BEGIN_PROVIDER [ type(map_type), mo_integrals_map ]
implicit none
BEGIN_DOC
! MO integrals
END_DOC
2015-03-22 20:48:32 +01:00
integer(key_kind) :: key_max
integer(map_size_kind) :: sze
call bielec_integrals_index(mo_tot_num,mo_tot_num,mo_tot_num,mo_tot_num,key_max)
sze = key_max
2014-04-17 23:50:51 +02:00
call map_init(mo_integrals_map,sze)
2015-03-31 09:58:25 +02:00
write(output_bielec_integrals,*) 'MO map initialized'
2014-04-17 23:50:51 +02:00
END_PROVIDER
subroutine insert_into_ao_integrals_map(n_integrals, &
buffer_i, buffer_values)
use map_module
implicit none
BEGIN_DOC
! Create new entry into AO map
END_DOC
2015-03-22 20:48:32 +01:00
integer, intent(in) :: n_integrals
integer(key_kind), intent(inout) :: buffer_i(n_integrals)
2014-04-17 23:50:51 +02:00
real(integral_kind), intent(inout) :: buffer_values(n_integrals)
call map_append(ao_integrals_map, buffer_i, buffer_values, n_integrals)
end
subroutine insert_into_mo_integrals_map(n_integrals, &
buffer_i, buffer_values, thr)
use map_module
implicit none
BEGIN_DOC
! Create new entry into MO map, or accumulate in an existing entry
END_DOC
2015-03-22 20:48:32 +01:00
integer, intent(in) :: n_integrals
integer(key_kind), intent(inout) :: buffer_i(n_integrals)
2014-04-17 23:50:51 +02:00
real(integral_kind), intent(inout) :: buffer_values(n_integrals)
2015-03-22 20:48:32 +01:00
real(integral_kind), intent(in) :: thr
2014-04-17 23:50:51 +02:00
call map_update(mo_integrals_map, buffer_i, buffer_values, n_integrals, thr)
end
double precision function get_mo_bielec_integral(i,j,k,l,map)
use map_module
implicit none
BEGIN_DOC
! Returns one integral <ij|kl> in the MO basis
END_DOC
integer, intent(in) :: i,j,k,l
2015-03-22 20:48:32 +01:00
integer(key_kind) :: idx
2014-04-17 23:50:51 +02:00
type(map_type), intent(inout) :: map
real(integral_kind) :: tmp
PROVIDE mo_bielec_integrals_in_map
!DIR$ FORCEINLINE
call bielec_integrals_index(i,j,k,l,idx)
call map_get(map,idx,tmp)
get_mo_bielec_integral = dble(tmp)
end
double precision function mo_bielec_integral(i,j,k,l)
implicit none
BEGIN_DOC
! Returns one integral <ij|kl> in the MO basis
END_DOC
integer, intent(in) :: i,j,k,l
double precision :: get_mo_bielec_integral
PROVIDE mo_bielec_integrals_in_map
mo_bielec_integral = get_mo_bielec_integral(i,j,k,l,mo_integrals_map)
return
end
subroutine get_mo_bielec_integrals(j,k,l,sze,out_val,map)
use map_module
implicit none
BEGIN_DOC
! Returns multiple integrals <ij|kl> in the MO basis, all
! i for j,k,l fixed.
END_DOC
integer, intent(in) :: j,k,l, sze
2015-03-22 20:48:32 +01:00
double precision, intent(out) :: out_val(sze)
2014-04-17 23:50:51 +02:00
type(map_type), intent(inout) :: map
integer :: i
2015-03-22 20:48:32 +01:00
integer(key_kind) :: hash(sze)
real(integral_kind) :: tmp_val(sze)
2014-04-17 23:50:51 +02:00
PROVIDE mo_bielec_integrals_in_map
do i=1,sze
!DIR$ FORCEINLINE
call bielec_integrals_index(i,j,k,l,hash(i))
enddo
2015-03-22 20:48:32 +01:00
if (key_kind == 8) then
call map_get_many(map, hash, out_val, sze)
else
call map_get_many(map, hash, tmp_val, sze)
! Conversion to double precision
do i=1,sze
out_val(i) = tmp_val(i)
enddo
endif
2014-04-17 23:50:51 +02:00
end
2014-06-06 16:19:14 +02:00
subroutine get_mo_bielec_integrals_existing_ik(j,l,sze,out_array,map)
use map_module
implicit none
BEGIN_DOC
! Returns multiple integrals <ij|kl> in the MO basis, all
! i(1)j(1) 1/r12 k(2)l(2)
2014-06-06 16:19:14 +02:00
! i for j,k,l fixed.
END_DOC
integer, intent(in) :: j,l, sze
logical, intent(out) :: out_array(sze,sze)
type(map_type), intent(inout) :: map
integer :: i,k,kk,ll,m
2015-03-22 20:48:32 +01:00
integer(key_kind),allocatable :: hash(:)
2014-06-06 16:19:14 +02:00
integer ,allocatable :: pairs(:,:), iorder(:)
PROVIDE mo_bielec_integrals_in_map
allocate (hash(sze*sze), pairs(2,sze*sze),iorder(sze*sze))
kk=0
do k=1,sze
do i=1,sze
kk += 1
!DIR$ FORCEINLINE
call bielec_integrals_index(i,j,k,l,hash(kk))
pairs(1,kk) = i
pairs(2,kk) = k
iorder(kk) = kk
enddo
enddo
logical :: integral_is_in_map
2015-04-09 21:46:28 +02:00
if (key_kind == 8) then
2015-03-22 20:48:32 +01:00
call i8radix_sort(hash,iorder,kk,-1)
2015-04-09 21:46:28 +02:00
else if (key_kind == 4) then
2015-03-22 20:48:32 +01:00
call iradix_sort(hash,iorder,kk,-1)
2015-04-09 21:46:28 +02:00
else if (key_kind == 2) then
2015-03-22 20:48:32 +01:00
call i2radix_sort(hash,iorder,kk,-1)
endif
2014-06-06 16:19:14 +02:00
call map_exists_many(mo_integrals_map, hash, kk)
do ll=1,kk
m = iorder(ll)
i=pairs(1,m)
k=pairs(2,m)
out_array(i,k) = (hash(ll) /= 0_8)
enddo
deallocate(pairs,hash,iorder)
end
2014-04-17 23:50:51 +02:00
integer*8 function get_mo_map_size()
implicit none
BEGIN_DOC
! Return the number of elements in the MO map
END_DOC
get_mo_map_size = mo_integrals_map % n_elements
end
subroutine clear_mo_map
implicit none
BEGIN_DOC
! Frees the memory of the MO map
END_DOC
call map_deinit(mo_integrals_map)
FREE mo_integrals_map
end
BEGIN_TEMPLATE
subroutine dump_$ao_integrals(filename)
use map_module
implicit none
BEGIN_DOC
! Save to disk the $ao integrals
END_DOC
character*(*), intent(in) :: filename
integer(cache_key_kind), pointer :: key(:)
real(integral_kind), pointer :: val(:)
integer*8 :: i,j, n
call ezfio_set_work_empty(.False.)
2014-04-17 23:50:51 +02:00
open(unit=66,file=filename,FORM='unformatted')
write(66) integral_kind, key_kind
write(66) $ao_integrals_map%sorted, $ao_integrals_map%map_size, &
$ao_integrals_map%n_elements
do i=0_8,$ao_integrals_map%map_size
write(66) $ao_integrals_map%map(i)%sorted, $ao_integrals_map%map(i)%map_size,&
$ao_integrals_map%map(i)%n_elements
enddo
do i=0_8,$ao_integrals_map%map_size
key => $ao_integrals_map%map(i)%key
val => $ao_integrals_map%map(i)%value
n = $ao_integrals_map%map(i)%n_elements
write(66) (key(j), j=1,n), (val(j), j=1,n)
enddo
close(66)
end
IRP_IF COARRAY
subroutine communicate_$ao_integrals()
use map_module
implicit none
BEGIN_DOC
! Communicate the $ao integrals with co-array
END_DOC
integer(cache_key_kind), pointer :: key(:)
real(integral_kind), pointer :: val(:)
integer*8 :: i,j, k, nmax
integer*8, save :: n[*]
integer :: copy_n
real(integral_kind), allocatable :: buffer_val(:)[:]
integer(cache_key_kind), allocatable :: buffer_key(:)[:]
real(integral_kind), allocatable :: copy_val(:)
2015-03-22 20:48:32 +01:00
integer(key_kind), allocatable :: copy_key(:)
n = 0_8
do i=0_8,$ao_integrals_map%map_size
n = max(n,$ao_integrals_map%map(i)%n_elements)
enddo
sync all
nmax = 0_8
do j=1,num_images()
nmax = max(nmax,n[j])
enddo
allocate( buffer_key(nmax)[*], buffer_val(nmax)[*])
allocate( copy_key(nmax), copy_val(nmax))
do i=0_8,$ao_integrals_map%map_size
key => $ao_integrals_map%map(i)%key
val => $ao_integrals_map%map(i)%value
n = $ao_integrals_map%map(i)%n_elements
do j=1,n
buffer_key(j) = key(j)
buffer_val(j) = val(j)
enddo
sync all
do j=1,num_images()
if (j /= this_image()) then
copy_n = n[j]
do k=1,copy_n
copy_val(k) = buffer_val(k)[j]
copy_key(k) = buffer_key(k)[j]
copy_key(k) = copy_key(k)+ishft(i,-map_shift)
enddo
call map_append($ao_integrals_map, copy_key, copy_val, copy_n )
endif
enddo
sync all
enddo
deallocate( buffer_key, buffer_val, copy_val, copy_key)
end
IRP_ENDIF
2014-04-17 23:50:51 +02:00
integer function load_$ao_integrals(filename)
implicit none
BEGIN_DOC
! Read from disk the $ao integrals
END_DOC
character*(*), intent(in) :: filename
integer*8 :: i
integer(cache_key_kind), pointer :: key(:)
real(integral_kind), pointer :: val(:)
integer :: iknd, kknd, n, j
load_$ao_integrals = 1
open(unit=66,file=filename,FORM='unformatted',STATUS='UNKNOWN')
read(66,err=98,end=98) iknd, kknd
if (iknd /= integral_kind) then
print *, 'Wrong integrals kind in file :', iknd
stop 1
endif
if (kknd /= key_kind) then
print *, 'Wrong key kind in file :', kknd
stop 1
endif
read(66,err=98,end=98) $ao_integrals_map%sorted, $ao_integrals_map%map_size,&
$ao_integrals_map%n_elements
do i=0_8, $ao_integrals_map%map_size
read(66,err=99,end=99) $ao_integrals_map%map(i)%sorted, &
$ao_integrals_map%map(i)%map_size, $ao_integrals_map%map(i)%n_elements
call cache_map_reallocate($ao_integrals_map%map(i),$ao_integrals_map%map(i)%map_size)
enddo
do i=0_8, $ao_integrals_map%map_size
key => $ao_integrals_map%map(i)%key
val => $ao_integrals_map%map(i)%value
n = $ao_integrals_map%map(i)%n_elements
read(66,err=99,end=99) (key(j), j=1,n), (val(j), j=1,n)
enddo
call map_sort($ao_integrals_map)
load_$ao_integrals = 0
return
99 continue
call map_deinit($ao_integrals_map)
FREE $ao_integrals_map
if (.True.) then
PROVIDE $ao_integrals_map
endif
stop 'Problem reading $ao_integrals_map file in work/'
98 continue
end
SUBST [ ao_integrals_map, ao_integrals, ao_num , get_ao_bielec_integral ]
ao_integrals_map ; ao_integrals ; ao_num ; get_ao_bielec_integral ;;
2014-05-13 13:57:58 +02:00
mo_integrals_map ; mo_integrals ; mo_tot_num ; get_mo_bielec_integral ;;
2014-04-17 23:50:51 +02:00
END_TEMPLATE