10
0
mirror of https://github.com/LCPQ/quantum_package synced 2024-11-12 17:13:54 +01:00

Corrected bug when reading basis set

This commit is contained in:
Anthony Scemama 2014-10-09 12:26:17 +02:00
parent 6eaadc9b0c
commit c1d8e6d29a
3 changed files with 64 additions and 10 deletions

View File

@ -1,3 +1,5 @@
open Core.Std;;
exception ElementError of string
type t =
@ -8,7 +10,8 @@ type t =
|K |Ca|Sc|Ti|V |Cr|Mn|Fe|Co|Ni|Cu|Zn|Ga|Ge|As|Se|Br|Kr
;;
let of_string = function
let of_string x =
match (String.capitalize (String.lowercase x)) with
| "X" | "Dummy" -> X
| "H" | "Hydrogen" -> H
| "He" | "Helium" -> He

View File

@ -310,12 +310,28 @@ subroutine save_wavefunction_general(ndet,nstates,psidet,psicoef)
integer :: i,k
PROVIDE progress_bar
call start_progress(7,'Saving wfunction',0.d0)
progress_bar(1) = 1
progress_value = dble(progress_bar(1))
call ezfio_set_determinants_N_int(N_int)
progress_bar(1) = 2
progress_value = dble(progress_bar(1))
call ezfio_set_determinants_bit_kind(bit_kind)
progress_bar(1) = 3
progress_value = dble(progress_bar(1))
call ezfio_set_determinants_N_det(ndet)
progress_bar(1) = 4
progress_value = dble(progress_bar(1))
call ezfio_set_determinants_n_states(nstates)
progress_bar(1) = 5
progress_value = dble(progress_bar(1))
call ezfio_set_determinants_mo_label(mo_label)
progress_bar(1) = 6
progress_value = dble(progress_bar(1))
N_int2 = (N_int*bit_kind)/8
allocate (psi_det_save(N_int2,2,ndet))
do i=1,ndet
@ -336,6 +352,8 @@ subroutine save_wavefunction_general(ndet,nstates,psidet,psicoef)
call ezfio_set_determinants_psi_det(psi_det_save)
deallocate (psi_det_save)
progress_bar(7) = 7
progress_value = dble(progress_bar(1))
allocate (psi_coef_save(ndet,nstates))
do k=1,nstates
do i=1,ndet
@ -344,5 +362,6 @@ subroutine save_wavefunction_general(ndet,nstates,psidet,psicoef)
enddo
call ezfio_set_determinants_psi_coef(psi_coef_save)
call write_int(output_dets,ndet,'Saved determinants')
call stop_progress
deallocate (psi_coef_save)
end

View File

@ -953,23 +953,36 @@ END_PROVIDER
BEGIN_PROVIDER [ integer*8, det_connections, (N_con_int,N_det) ]
implicit none
BEGIN_DOC
!
! Build connection proxy between determinants
END_DOC
integer :: i,j
integer :: degree
integer :: j_int, j_k, j_l
integer, allocatable :: idx(:)
integer :: thread_num
!$ integer :: omp_get_thread_num
PROVIDE progress_bar
call start_progress(N_det,'Det connections',0.d0)
select case(N_int)
case(1)
!$OMP PARALLEL DEFAULT (NONE) &
!$OMP SHARED(N_det, N_con_int, psi_det,N_int, det_connections)&
!$OMP PRIVATE(i,j_int,j_k,j_l,j,degree,idx)
!$OMP SHARED(N_det, N_con_int, psi_det,N_int, det_connections, &
!$OMP progress_bar,progress_value)&
!$OMP PRIVATE(i,j_int,j_k,j_l,j,degree,idx,thread_num)
!$ thread_num = omp_get_thread_num()
allocate (idx(0:N_det))
!$OMP DO SCHEDULE(guided)
do i=1,N_det
if (thread_num == 0) then
progress_bar(1) = i
progress_value = dble(i)
endif
do j_int=1,N_con_int
det_connections(j_int,i) = 0_8
j_k = ishft(j_int-1,11)
@ -992,11 +1005,17 @@ BEGIN_PROVIDER [ integer*8, det_connections, (N_con_int,N_det) ]
case(2)
!$OMP PARALLEL DEFAULT (NONE) &
!$OMP SHARED(N_det, N_con_int, psi_det,N_int, det_connections)&
!$OMP PRIVATE(i,j_int,j_k,j_l,j,degree,idx)
!$OMP SHARED(N_det, N_con_int, psi_det,N_int, det_connections,&
!$OMP progress_bar,progress_value)&
!$OMP PRIVATE(i,j_int,j_k,j_l,j,degree,idx,thread_num)
!$ thread_num = omp_get_thread_num()
allocate (idx(0:N_det))
!$OMP DO SCHEDULE(guided)
do i=1,N_det
if (thread_num == 0) then
progress_bar(1) = i
progress_value = dble(i)
endif
do j_int=1,N_con_int
det_connections(j_int,i) = 0_8
j_k = ishft(j_int-1,11)
@ -1021,11 +1040,17 @@ BEGIN_PROVIDER [ integer*8, det_connections, (N_con_int,N_det) ]
case(3)
!$OMP PARALLEL DEFAULT (NONE) &
!$OMP SHARED(N_det, N_con_int, psi_det,N_int, det_connections)&
!$OMP PRIVATE(i,j_int,j_k,j_l,j,degree,idx)
!$OMP SHARED(N_det, N_con_int, psi_det,N_int, det_connections,&
!$OMP progress_bar,progress_value)&
!$OMP PRIVATE(i,j_int,j_k,j_l,j,degree,idx,thread_num)
!$ thread_num = omp_get_thread_num()
allocate (idx(0:N_det))
!$OMP DO SCHEDULE(guided)
do i=1,N_det
if (thread_num == 0) then
progress_bar(1) = i
progress_value = dble(i)
endif
do j_int=1,N_con_int
det_connections(j_int,i) = 0_8
j_k = ishft(j_int-1,11)
@ -1053,11 +1078,17 @@ BEGIN_PROVIDER [ integer*8, det_connections, (N_con_int,N_det) ]
!$OMP PARALLEL DEFAULT (NONE) &
!$OMP SHARED(N_det, N_con_int, psi_det,N_int, det_connections)&
!$OMP PRIVATE(i,j_int,j_k,j_l,j,degree,idx)
!$OMP SHARED(N_det, N_con_int, psi_det,N_int, det_connections,&
!$OMP progress_bar,progress_value)&
!$OMP PRIVATE(i,j_int,j_k,j_l,j,degree,idx,thread_num)
!$ thread_num = omp_get_thread_num()
allocate (idx(0:N_det))
!$OMP DO SCHEDULE(guided)
do i=1,N_det
if (thread_num == 0) then
progress_bar(1) = i
progress_value = dble(i)
endif
do j_int=1,N_con_int
det_connections(j_int,i) = 0_8
j_k = ishft(j_int-1,11)
@ -1078,6 +1109,7 @@ BEGIN_PROVIDER [ integer*8, det_connections, (N_con_int,N_det) ]
!$OMP END PARALLEL
end select
call stop_progress
END_PROVIDER