mirror of
https://github.com/LCPQ/quantum_package
synced 2024-12-23 04:43:50 +01:00
Added runtime checks
This commit is contained in:
parent
1c3d8f6a09
commit
767074ab95
@ -63,7 +63,9 @@ subroutine run_wf
|
|||||||
! --------
|
! --------
|
||||||
|
|
||||||
print *, 'Davidson'
|
print *, 'Davidson'
|
||||||
call davidson_slave_tcp(i)
|
call omp_set_nested(.True.)
|
||||||
|
call davidson_slave_tcp(0)
|
||||||
|
call omp_set_nested(.False.)
|
||||||
print *, 'Davidson done'
|
print *, 'Davidson done'
|
||||||
|
|
||||||
else if (trim(zmq_state) == 'pt2') then
|
else if (trim(zmq_state) == 'pt2') then
|
||||||
|
@ -129,6 +129,8 @@ subroutine davidson_slave_work(zmq_to_qp_run_socket, zmq_socket_push, N_st, sze,
|
|||||||
stop 'error'
|
stop 'error'
|
||||||
endif
|
endif
|
||||||
|
|
||||||
|
PROVIDE psi_bilinear_matrix_columns psi_bilinear_matrix_transp_rows_loc
|
||||||
|
|
||||||
! Run tasks
|
! Run tasks
|
||||||
! ---------
|
! ---------
|
||||||
|
|
||||||
|
@ -10,6 +10,7 @@ program davidson_slave
|
|||||||
|
|
||||||
call provide_everything
|
call provide_everything
|
||||||
call switch_qp_run_to_master
|
call switch_qp_run_to_master
|
||||||
|
call omp_set_nested(.True.)
|
||||||
|
|
||||||
zmq_context = f77_zmq_ctx_new ()
|
zmq_context = f77_zmq_ctx_new ()
|
||||||
zmq_state = 'davidson'
|
zmq_state = 'davidson'
|
||||||
|
@ -20,7 +20,7 @@ integer*8 function spin_det_search_key(det,Nint)
|
|||||||
do i=2,Nint
|
do i=2,Nint
|
||||||
spin_det_search_key = ieor(spin_det_search_key,det(i))
|
spin_det_search_key = ieor(spin_det_search_key,det(i))
|
||||||
enddo
|
enddo
|
||||||
spin_det_search_key = spin_det_search_key+1_bit_kind-unsigned_shift
|
spin_det_search_key = spin_det_search_key+unsigned_shift
|
||||||
end
|
end
|
||||||
|
|
||||||
|
|
||||||
@ -414,7 +414,8 @@ BEGIN_PROVIDER [ double precision, psi_bilinear_matrix_values, (N_det,N_states)
|
|||||||
do k=1,N_det
|
do k=1,N_det
|
||||||
i = get_index_in_psi_det_alpha_unique(psi_det(1,1,k),N_int)
|
i = get_index_in_psi_det_alpha_unique(psi_det(1,1,k),N_int)
|
||||||
j = get_index_in_psi_det_beta_unique (psi_det(1,2,k),N_int)
|
j = get_index_in_psi_det_beta_unique (psi_det(1,2,k),N_int)
|
||||||
|
if (i < 1) stop 'i<1'
|
||||||
|
if (j < 1) stop 'j<1'
|
||||||
do l=1,N_states
|
do l=1,N_states
|
||||||
psi_bilinear_matrix_values(k,l) = psi_coef(k,l)
|
psi_bilinear_matrix_values(k,l) = psi_coef(k,l)
|
||||||
enddo
|
enddo
|
||||||
@ -471,6 +472,9 @@ BEGIN_PROVIDER [ integer, psi_bilinear_matrix_columns_loc, (N_det_beta_unique+1)
|
|||||||
l = psi_bilinear_matrix_columns(k)
|
l = psi_bilinear_matrix_columns(k)
|
||||||
psi_bilinear_matrix_columns_loc(l) = k
|
psi_bilinear_matrix_columns_loc(l) = k
|
||||||
endif
|
endif
|
||||||
|
if (psi_bilinear_matrix_columns(k) < 1) then
|
||||||
|
stop '(psi_bilinear_matrix_columns(k) < 1)'
|
||||||
|
endif
|
||||||
enddo
|
enddo
|
||||||
psi_bilinear_matrix_columns_loc(N_det_beta_unique+1) = N_det+1
|
psi_bilinear_matrix_columns_loc(N_det_beta_unique+1) = N_det+1
|
||||||
END_PROVIDER
|
END_PROVIDER
|
||||||
@ -496,17 +500,23 @@ BEGIN_PROVIDER [ double precision, psi_bilinear_matrix_transp_values, (N_det,N_
|
|||||||
integer*8, allocatable :: to_sort(:)
|
integer*8, allocatable :: to_sort(:)
|
||||||
allocate(to_sort(N_det))
|
allocate(to_sort(N_det))
|
||||||
!$OMP PARALLEL DEFAULT(SHARED) PRIVATE(i,j,k,l)
|
!$OMP PARALLEL DEFAULT(SHARED) PRIVATE(i,j,k,l)
|
||||||
!$OMP DO COLLAPSE(2)
|
|
||||||
do l=1,N_states
|
do l=1,N_states
|
||||||
|
!$OMP DO
|
||||||
do k=1,N_det
|
do k=1,N_det
|
||||||
psi_bilinear_matrix_transp_values (k,l) = psi_bilinear_matrix_values (k,l)
|
psi_bilinear_matrix_transp_values (k,l) = psi_bilinear_matrix_values (k,l)
|
||||||
enddo
|
enddo
|
||||||
|
!$OMP ENDDO
|
||||||
enddo
|
enddo
|
||||||
!$OMP ENDDO
|
|
||||||
!$OMP DO
|
!$OMP DO
|
||||||
do k=1,N_det
|
do k=1,N_det
|
||||||
psi_bilinear_matrix_transp_columns(k) = psi_bilinear_matrix_columns(k)
|
psi_bilinear_matrix_transp_columns(k) = psi_bilinear_matrix_columns(k)
|
||||||
|
if (psi_bilinear_matrix_transp_columns(k) < 1) then
|
||||||
|
stop '(psi_bilinear_matrix_transp_columns(k) < 1)'
|
||||||
|
endif
|
||||||
psi_bilinear_matrix_transp_rows (k) = psi_bilinear_matrix_rows (k)
|
psi_bilinear_matrix_transp_rows (k) = psi_bilinear_matrix_rows (k)
|
||||||
|
if (psi_bilinear_matrix_transp_rows(k) < 1) then
|
||||||
|
stop '(psi_bilinear_matrix_transp_rows(k) < 1)'
|
||||||
|
endif
|
||||||
i = psi_bilinear_matrix_transp_columns(k)
|
i = psi_bilinear_matrix_transp_columns(k)
|
||||||
j = psi_bilinear_matrix_transp_rows (k)
|
j = psi_bilinear_matrix_transp_rows (k)
|
||||||
to_sort(k) = int(N_det_beta_unique,8) * int(j-1,8) + int(i,8)
|
to_sort(k) = int(N_det_beta_unique,8) * int(j-1,8) + int(i,8)
|
||||||
|
Loading…
Reference in New Issue
Block a user