10
0
mirror of https://github.com/LCPQ/quantum_package synced 2024-12-22 12:23:48 +01:00

Integrated in main program

This commit is contained in:
Anthony Scemama 2017-10-02 11:40:36 +02:00
parent 827e6933d4
commit 1c58249f53
4 changed files with 15 additions and 57 deletions

View File

@ -1,53 +0,0 @@
program FourIdx
use map_module
implicit none
BEGIN_DOC
! Performs a four index transformation of the two-electron integrals
END_DOC
type(map_type) :: test_map
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
call map_init(test_map,sze)
! call four_index_transform(ao_integrals_map,test_map, &
! mo_coef, size(mo_coef,1), &
! 1, 1, 1, 1, ao_num, ao_num, ao_num, ao_num, &
! 1, 1, 1, 1, mo_tot_num, mo_tot_num, mo_tot_num, mo_tot_num)
double precision :: t0,t1
call wall_time(t0)
call four_index_transform_sym(ao_integrals_map,test_map, &
mo_coef, size(mo_coef,1), &
1, 1, 1, 1, ao_num, ao_num, ao_num, ao_num, &
1, 1, 1, 1, mo_tot_num, mo_tot_num, mo_tot_num, mo_tot_num)
call wall_time(t1)
print *, 'Time: ', t1-t0, 's'
integer :: i,j,k,l
real(integral_kind) :: integral1, integral2
provide mo_bielec_integrals_in_map
do i=1,mo_tot_num
do j=1,mo_tot_num
do k=1,mo_tot_num
do l=1,mo_tot_num
call bielec_integrals_index(i,j,k,l,key_max)
call map_get(test_map,key_max,integral1)
call map_get(mo_integrals_map,key_max,integral2)
if (dabs(integral2) >=1.d-10 ) then
if (dabs(integral1 / integral2 -1.d0) > .001d0) then
print *, i,j,k,l
print *, integral1, integral2
print *, ''
endif
endif
enddo
enddo
enddo
enddo
end

View File

@ -67,7 +67,7 @@ subroutine four_index_transform_sym(map_a,map_c,matrix_B,LDB, &
allocate(l_pointer(l_start:l_end+1), value((i_max*k_max)) )
ii = 1_8
!$OMP PARALLEL DEFAULT(SHARED) PRIVATE(i,j,k,l,ik,idx,ii)
!$OMP PARALLEL DEFAULT(SHARED) PRIVATE(i,j,k,l,ik,idx)
do l=l_start,l_end
!$OMP SINGLE
l_pointer(l) = ii
@ -102,8 +102,10 @@ subroutine four_index_transform_sym(map_a,map_c,matrix_B,LDB, &
!$OMP END SINGLE
enddo
enddo
!$OMP END PARALLEL
!$OMP SINGLE
l_pointer(l_end+1) = ii
!$OMP END SINGLE
!$OMP END PARALLEL
deallocate(value)
!INPUT DATA

View File

@ -1 +1 @@
Perturbation Selectors_full Generators_full ZMQ
Perturbation Selectors_full Generators_full ZMQ FourIdx

View File

@ -117,7 +117,16 @@ BEGIN_PROVIDER [ logical, mo_bielec_integrals_in_map ]
endif
else
call add_integrals_to_map(full_ijkl_bitmask_4)
! call add_integrals_to_map(full_ijkl_bitmask_4)
call four_index_transform_sym(ao_integrals_map,mo_integrals_map, &
mo_coef, size(mo_coef,1), &
1, 1, 1, 1, ao_num, ao_num, ao_num, ao_num, &
1, 1, 1, 1, mo_tot_num, mo_tot_num, mo_tot_num, mo_tot_num)
integer*8 :: get_mo_map_size, mo_map_size
mo_map_size = get_mo_map_size()
print*,'Molecular integrals provided'
endif
if (write_mo_integrals) then
call ezfio_set_work_empty(.False.)