fixed qp_convert

This commit is contained in:
Anthony Scemama 2019-01-12 02:41:16 +01:00
parent cf2fa908e2
commit 43a9f7f302
6 changed files with 154 additions and 41 deletions

1
TODO
View File

@ -52,7 +52,6 @@ Refaire les benchmarks
# Documentation de /etc
# Toto
singles_alpha_csc_idx, singles_alpha_size | git diff
Environment variables dans qp_run (prefix,etc)
Si un provider est un programme, generer une page a lui tout seul avec le man

View File

@ -1,6 +1,6 @@
#!/bin/bash
export QP_ROOT=$(dirname $0)
export QP_ROOT=$(dirname $0)/..
exec bash --init-file <(cat << EOF

4
configure vendored
View File

@ -3,7 +3,7 @@
# Quantum Package configuration script
#
export QP_ROOT="$( cd "$(dirname "$0")" ; pwd --physical )"
export QP_ROOT="$( cd "$(dirname "$0")" ; pwd -P )"
echo "QP_ROOT="$QP_ROOT
@ -234,7 +234,7 @@ EOF
--yes --comp=4.07.0
eval $(${QP_ROOT}/bin/opam env)
opam install -y ${OCAML_PACKAGES}
opam install -y ${OCAML_PACKAGES} || exit 1
elif [[ ${PACKAGE} = ezfio ]] ; then

View File

@ -3,12 +3,11 @@
convert output of gamess/GAU$$IAN to ezfio
Usage:
qp_convert_output_to_ezfio.py <file.out> [-o|--output <EZFIO_DIRECTORY>]
qp_convert_output_to_ezfio [-o EZFIO_FILE|--output=EZFIO_FILE] <file.out>
Option:
file.out is the file to check (like gamess.out)
EZFIO_DIRECTORY is the name of the produced directory
(by default is file.out.ezfio)
Options:
-o EZFIO_FILE --output=EZFIO_FILE Produced directory
by default is file.out.ezfio
"""
@ -358,8 +357,10 @@ if __name__ == '__main__':
file_ = get_full_path(arguments['<file.out>'])
if arguments["-o"]:
ezfio_file = get_full_path(arguments["<ezfio_directory>"])
if arguments["--output"]:
ezfio_file = get_full_path(arguments["--output"])
elif arguments["-o"]:
ezfio_file = get_full_path(arguments["-o"])
else:
ezfio_file = "{0}.ezfio".format(file_)

View File

@ -138,9 +138,19 @@ subroutine H_S2_u_0_nstates_openmp_work_$N_int(v_t,s_t,u_t,N_st,sze,istart,iend,
integer :: maxab, n_singles_a, n_singles_b, kcol_prev
integer*8 :: k8
logical :: compute_singles
integer*8 :: last_found, left, right
!TODO
compute_singles = .True.
double precision :: rss, mem
call resident_memory(rss)
mem = dble(singles_beta_csc_size) / 1024.d0**3
compute_singles = (mem+rss > qp_max_mem)
! compute_singles = (iend-istart < 100000).and.(mem+rss < qp_max_mem)
if (.not.compute_singles) then
! provide singles_alpha_csc
provide singles_beta_csc
endif
maxab = max(N_det_alpha_unique, N_det_beta_unique)+1
@ -154,7 +164,7 @@ subroutine H_S2_u_0_nstates_openmp_work_$N_int(v_t,s_t,u_t,N_st,sze,istart,iend,
! -------------------------------------------------
PROVIDE N_int nthreads_davidson
!$OMP PARALLEL DEFAULT(NONE) NUM_THREADS(nthreads_davidson) &
!$OMP PARALLEL DEFAULT(SHARED) NUM_THREADS(nthreads_davidson) &
!$OMP SHARED(psi_bilinear_matrix_rows, N_det, &
!$OMP psi_bilinear_matrix_columns, &
!$OMP psi_det_alpha_unique, psi_det_beta_unique, &
@ -166,13 +176,15 @@ subroutine H_S2_u_0_nstates_openmp_work_$N_int(v_t,s_t,u_t,N_st,sze,istart,iend,
!$OMP psi_bilinear_matrix_columns_loc, &
!$OMP psi_bilinear_matrix_transp_rows_loc, &
!$OMP istart, iend, istep, irp_here, v_t, s_t, &
!$OMP ishift, idx0, u_t, maxab, compute_singles) &
!$OMP ishift, idx0, u_t, maxab, compute_singles, &
!$OMP singles_alpha_csc,singles_alpha_csc_idx, &
!$OMP singles_beta_csc,singles_beta_csc_idx) &
!$OMP PRIVATE(krow, kcol, tmp_det, spindet, k_a, k_b, i, &
!$OMP lcol, lrow, l_a, l_b, &
!$OMP buffer, doubles, n_doubles, &
!$OMP tmp_det2, hij, sij, idx, l, kcol_prev, &
!$OMP singles_a, n_singles_a, singles_b, &
!$OMP n_singles_b, k8)
!$OMP n_singles_b, k8, last_found,left,right)
! Alpha/Beta double excitations
! =============================
@ -201,12 +213,18 @@ subroutine H_S2_u_0_nstates_openmp_work_$N_int(v_t,s_t,u_t,N_st,sze,istart,iend,
tmp_det(1:$N_int,1) = psi_det_alpha_unique(1:$N_int, krow)
if (kcol /= kcol_prev) then
tmp_det(1:$N_int,2) = psi_det_beta_unique (1:$N_int, kcol)
if (compute_singles) then
tmp_det(1:$N_int,2) = psi_det_beta_unique (1:$N_int, kcol)
call get_all_spin_singles_$N_int( &
psi_det_beta_unique, idx0, &
tmp_det(1,2), N_det_beta_unique, &
singles_b, n_singles_b)
else
n_singles_b = 0
do k8=singles_beta_csc_idx(kcol),singles_beta_csc_idx(kcol+1)-1
n_singles_b = n_singles_b+1
singles_b(n_singles_b) = singles_beta_csc(k8)
enddo
endif
endif
kcol_prev = kcol
@ -222,7 +240,7 @@ subroutine H_S2_u_0_nstates_openmp_work_$N_int(v_t,s_t,u_t,N_st,sze,istart,iend,
l_a = psi_bilinear_matrix_columns_loc(lcol)
ASSERT (l_a <= N_det)
if (compute_singles) then
! if (compute_singles) then
do j=1,psi_bilinear_matrix_columns_loc(lcol+1) - psi_bilinear_matrix_columns_loc(lcol)
lrow = psi_bilinear_matrix_rows(l_a)
ASSERT (lrow <= N_det_alpha_unique)
@ -238,7 +256,34 @@ subroutine H_S2_u_0_nstates_openmp_work_$N_int(v_t,s_t,u_t,N_st,sze,istart,iend,
call get_all_spin_singles_$N_int( &
buffer, idx, tmp_det(1,1), j, &
singles_a, n_singles_a )
endif
! else
!
! n_singles_a = 0
! last_found = singles_alpha_csc_idx(krow)
! do j=1,psi_bilinear_matrix_columns_loc(lcol+1) - psi_bilinear_matrix_columns_loc(lcol)
! lrow = psi_bilinear_matrix_rows(l_a)
! ASSERT (lrow <= N_det_alpha_unique)
!
! left = last_found
! right = singles_alpha_csc_idx(krow+1)
! do while (right-left>0_8)
! k8 = shiftr(right+left,1)
! if (singles_alpha_csc(k8) > lrow) then
! right = k8
! else if (singles_alpha_csc(k8) < lrow) then
! left = k8 + 1_8
! else
! n_singles_a += 1
! singles_a(n_singles_a) = l_a
! last_found = k8+1_8
! exit
! endif
! enddo
! l_a = l_a+1
! enddo
! j = j-1
!
! endif
! Loop over alpha singles
! -----------------------

View File

@ -149,7 +149,7 @@ END_TEMPLATE
integer function get_index_in_psi_det_alpha_unique(key,Nint)
use bitmasks
BEGIN_DOC
! Returns the index of the determinant in the ``psi_det_alpha_unique`` array
! Returns the index of the determinant in the :c:data:`psi_det_alpha_unique` array
END_DOC
implicit none
@ -230,7 +230,7 @@ end
integer function get_index_in_psi_det_beta_unique(key,Nint)
use bitmasks
BEGIN_DOC
! Returns the index of the determinant in the ``psi_det_beta_unique`` array
! Returns the index of the determinant in the :c:data:`psi_det_beta_unique` array
END_DOC
implicit none
@ -466,7 +466,7 @@ BEGIN_PROVIDER [ integer, psi_bilinear_matrix_order_reverse , (N_det) ]
use bitmasks
implicit none
BEGIN_DOC
! Order which allows to go from ``psi_bilinear_matrix`` to ``psi_det``
! Order which allows to go from :c:data:`psi_bilinear_matrix` to :c:data:`psi_det`
END_DOC
integer :: k
!$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(k)
@ -489,7 +489,7 @@ BEGIN_PROVIDER [ integer, psi_bilinear_matrix_columns_loc, (N_det_beta_unique+1)
!
! Rows are $\alpha$ determinants and columns are $\beta$.
!
! Order refers to ``psi_det``
! Order refers to :c:data:`psi_det`
END_DOC
integer :: i,j,k, l
@ -521,7 +521,7 @@ END_PROVIDER
use bitmasks
implicit none
BEGIN_DOC
! Transpose of ``psi_bilinear_matrix``
! Transpose of :c:data:`psi_bilinear_matrix`
!
! $D_\beta^\dagger.C^\dagger.D_\alpha$
!
@ -582,7 +582,7 @@ BEGIN_PROVIDER [ integer, psi_bilinear_matrix_transp_rows_loc, (N_det_alpha_uniq
use bitmasks
implicit none
BEGIN_DOC
! Location of the columns in the ``psi_bilinear_matrix``
! Location of the columns in the :c:data:`psi_bilinear_matrix`
END_DOC
integer :: i,j,k, l
@ -608,8 +608,8 @@ BEGIN_PROVIDER [ integer, psi_bilinear_matrix_order_transp_reverse , (N_det) ]
use bitmasks
implicit none
BEGIN_DOC
! Order which allows to go from ``psi_bilinear_matrix_order_transp`` to
! ``psi_bilinear_matrix``
! Order which allows to go from :c:data:`psi_bilinear_matrix_order_transp` to
! :c:data:`psi_bilinear_matrix`
END_DOC
integer :: k
psi_bilinear_matrix_order_transp_reverse = -1
@ -730,7 +730,7 @@ subroutine generate_all_alpha_beta_det_products
!$OMP PRIVATE(i,j,k,l,tmp_det,iproc)
!$ iproc = omp_get_thread_num()
allocate (tmp_det(N_int,2,N_det_alpha_unique))
!$OMP DO SCHEDULE(static,1)
!$OMP DO SCHEDULE(static,8)
do j=1,N_det_beta_unique
l = 1
do i=1,N_det_alpha_unique
@ -858,7 +858,8 @@ end
subroutine copy_psi_bilinear_to_psi(psi, isize)
implicit none
BEGIN_DOC
! Overwrites ``psi_det`` and ``psi_coef`` with the wave function in bilinear order
! Overwrites :c:data:`psi_det` and :c:data:`psi_coef` with the wave function
! in bilinear order
END_DOC
integer, intent(in) :: isize
integer(bit_kind), intent(out) :: psi(N_int,2,isize)
@ -871,19 +872,14 @@ subroutine copy_psi_bilinear_to_psi(psi, isize)
enddo
end
BEGIN_PROVIDER [ integer, singles_alpha_size ]
implicit none
BEGIN_DOC
! Dimension of the ``singles_alpha`` array
END_DOC
singles_alpha_size = elec_alpha_num * (mo_num - elec_alpha_num)
END_PROVIDER
BEGIN_PROVIDER [ integer*8, singles_alpha_csc_idx, (N_det_alpha_unique+1) ]
&BEGIN_PROVIDER [ integer*8, singles_alpha_csc_size ]
implicit none
BEGIN_DOC
! Dimension of the ``singles_alpha`` array
! singles_alpha_csc_size : Dimension of the :c:data:`singles_alpha_csc` array
!
! singles_alpha_csc_idx : Index where the single excitations of determinant i start
END_DOC
integer :: i,j
integer, allocatable :: idx0(:), s(:)
@ -895,10 +891,10 @@ END_PROVIDER
!$OMP PARALLEL DEFAULT(NONE) &
!$OMP SHARED(N_det_alpha_unique, psi_det_alpha_unique, &
!$OMP idx0, N_int, singles_alpha_csc, &
!$OMP singles_alpha_size, singles_alpha_csc_idx) &
!$OMP elec_alpha_num, mo_num, singles_alpha_csc_idx) &
!$OMP PRIVATE(i,s,j)
allocate (s(singles_alpha_size))
!$OMP DO SCHEDULE(static,1)
allocate (s(elec_alpha_num * (mo_num-elec_alpha_num) ))
!$OMP DO SCHEDULE(static,64)
do i=1, N_det_alpha_unique
call get_all_spin_singles( &
psi_det_alpha_unique, idx0, psi_det_alpha_unique(1,i), N_int,&
@ -921,7 +917,7 @@ END_PROVIDER
BEGIN_PROVIDER [ integer, singles_alpha_csc, (singles_alpha_csc_size) ]
implicit none
BEGIN_DOC
! Dimension of the singles_alpha array
! Indices of all single excitations
END_DOC
integer :: i, k
integer, allocatable :: idx0(:)
@ -948,6 +944,78 @@ END_PROVIDER
BEGIN_PROVIDER [ integer*8, singles_beta_csc_idx, (N_det_beta_unique+1) ]
&BEGIN_PROVIDER [ integer*8, singles_beta_csc_size ]
implicit none
BEGIN_DOC
! singles_beta_csc_size : Dimension of the :c:data:`singles_beta_csc` array
!
! singles_beta_csc_idx : Index where the single excitations of determinant i start
END_DOC
integer :: i,j
integer, allocatable :: idx0(:), s(:)
allocate (idx0(N_det_beta_unique))
do i=1, N_det_beta_unique
idx0(i) = i
enddo
!$OMP PARALLEL DEFAULT(NONE) &
!$OMP SHARED(N_det_beta_unique, psi_det_beta_unique, &
!$OMP idx0, N_int, singles_beta_csc, &
!$OMP elec_beta_num, mo_num, singles_beta_csc_idx) &
!$OMP PRIVATE(i,s,j)
allocate (s(elec_beta_num*(mo_num-elec_beta_num)))
!$OMP DO SCHEDULE(static,1)
do i=1, N_det_beta_unique
call get_all_spin_singles( &
psi_det_beta_unique, idx0, psi_det_beta_unique(1,i), N_int,&
N_det_beta_unique, s, j)
singles_beta_csc_idx(i+1) = int(j,8)
enddo
!$OMP END DO
deallocate(s)
!$OMP END PARALLEL
deallocate(idx0)
singles_beta_csc_idx(1) = 1_8
do i=2, N_det_beta_unique+1
singles_beta_csc_idx(i) = singles_beta_csc_idx(i) + singles_beta_csc_idx(i-1)
enddo
singles_beta_csc_size = singles_beta_csc_idx(N_det_beta_unique+1)
END_PROVIDER
BEGIN_PROVIDER [ integer, singles_beta_csc, (singles_beta_csc_size) ]
implicit none
BEGIN_DOC
! Indices of all single excitations
END_DOC
integer :: i, k
integer, allocatable :: idx0(:)
allocate (idx0(N_det_beta_unique))
do i=1, N_det_beta_unique
idx0(i) = i
enddo
!$OMP PARALLEL DO DEFAULT(NONE) &
!$OMP SHARED(N_det_beta_unique, psi_det_beta_unique, &
!$OMP idx0, N_int, singles_beta_csc, singles_beta_csc_idx)&
!$OMP PRIVATE(i,k) SCHEDULE(static,64)
do i=1, N_det_beta_unique
call get_all_spin_singles( &
psi_det_beta_unique, idx0, psi_det_beta_unique(1,i), N_int,&
N_det_beta_unique, singles_beta_csc(singles_beta_csc_idx(i)),&
k)
enddo
!$OMP END PARALLEL DO
deallocate(idx0)
END_PROVIDER
subroutine get_all_spin_singles_and_doubles_1(buffer, idx, spindet, size_buffer, singles, doubles, n_singles, n_doubles)
use bitmasks
implicit none