mirror of
https://github.com/LCPQ/quantum_package
synced 2025-01-03 18:16:12 +01:00
ASSERT Fixed
This commit is contained in:
parent
c06a4c2ecd
commit
3d51dde718
@ -63,13 +63,14 @@ subroutine davidson_slave_work(zmq_to_qp_run_socket, zmq_socket_push, N_st, sze,
|
||||
character*(512) :: msg
|
||||
integer :: imin, imax, ishift, istep
|
||||
|
||||
integer, allocatable :: psi_det_read(:,:,:)
|
||||
double precision, allocatable :: v_0(:,:), s_0(:,:), u_t(:,:)
|
||||
!DIR$ ATTRIBUTES ALIGN : $IRP_ALIGN :: u_t, v_0, s_0
|
||||
|
||||
! Get wave function (u_t)
|
||||
! -----------------------
|
||||
|
||||
integer :: rc
|
||||
integer :: rc
|
||||
integer :: N_states_read, N_det_read, psi_det_size_read
|
||||
integer :: N_det_selectors_read, N_det_generators_read
|
||||
double precision :: energy(N_st)
|
||||
@ -107,12 +108,11 @@ subroutine davidson_slave_work(zmq_to_qp_run_socket, zmq_socket_push, N_st, sze,
|
||||
TOUCH N_det
|
||||
endif
|
||||
|
||||
allocate(v_0(sze,N_st), s_0(sze,N_st),u_t(N_st,N_det_read))
|
||||
|
||||
allocate(v_0(sze,N_st), s_0(sze,N_st),u_t(N_st,N_det))
|
||||
|
||||
rc = f77_zmq_recv(zmq_to_qp_run_socket,psi_det,N_int*2*N_det*bit_kind,0)
|
||||
if (rc /= N_int*2*N_det*bit_kind) then
|
||||
print *, 'f77_zmq_recv(zmq_to_qp_run_socket,psi_det,N_int*2*N_det*bit_kind,0)'
|
||||
rc = f77_zmq_recv(zmq_to_qp_run_socket,psi_det,N_int*2*N_det_read*bit_kind,0)
|
||||
if (rc /= N_int*2*N_det_read*bit_kind) then
|
||||
print *, 'f77_zmq_recv(zmq_to_qp_run_socket,psi_det,N_int*2*N_det_read*bit_kind,0)'
|
||||
stop 'error'
|
||||
endif
|
||||
|
||||
@ -129,11 +129,10 @@ subroutine davidson_slave_work(zmq_to_qp_run_socket, zmq_socket_push, N_st, sze,
|
||||
stop 'error'
|
||||
endif
|
||||
|
||||
PROVIDE psi_bilinear_matrix_columns psi_bilinear_matrix_transp_rows_loc
|
||||
|
||||
! Run tasks
|
||||
! ---------
|
||||
|
||||
|
||||
do
|
||||
v_0 = 0.d0
|
||||
s_0 = 0.d0
|
||||
@ -295,10 +294,6 @@ subroutine H_S2_u_0_nstates_zmq(v_0,s_0,u_0,N_st,sze)
|
||||
|
||||
if(N_st /= N_states_diag .or. sze < N_det) stop "assert fail in H_S2_u_0_nstates"
|
||||
|
||||
ASSERT (Nint > 0)
|
||||
ASSERT (Nint == N_int)
|
||||
ASSERT (n>0)
|
||||
|
||||
call new_parallel_job(zmq_to_qp_run_socket,'davidson')
|
||||
|
||||
character*(512) :: task
|
||||
|
@ -9,12 +9,13 @@ subroutine u_0_H_u_0(e_0,u_0,n,keys_tmp,Nint,N_st,sze)
|
||||
END_DOC
|
||||
integer, intent(in) :: n,Nint, N_st, sze
|
||||
double precision, intent(out) :: e_0(N_st)
|
||||
double precision, intent(inout):: u_0(sze,N_st)
|
||||
double precision, intent(inout) :: u_0(sze,N_st)
|
||||
integer(bit_kind),intent(in) :: keys_tmp(Nint,2,n)
|
||||
|
||||
double precision, allocatable :: v_0(:,:), s_0(:,:)
|
||||
double precision :: u_dot_u,u_dot_v,diag_H_mat_elem
|
||||
integer :: i,j
|
||||
integer :: i,j
|
||||
|
||||
allocate (v_0(sze,N_st),s_0(sze,N_st))
|
||||
call H_S2_u_0_nstates_openmp(v_0,s_0,u_0,N_st,sze)
|
||||
do i=1,N_st
|
||||
@ -158,11 +159,11 @@ subroutine H_S2_u_0_nstates_openmp_work_$N_int(v_0,s_0,u_t,N_st,sze,istart,iend,
|
||||
!$OMP psi_bilinear_matrix_transp_order, N_st, &
|
||||
!$OMP psi_bilinear_matrix_order_transp_reverse, &
|
||||
!$OMP psi_bilinear_matrix_columns_loc, &
|
||||
!$OMP istart, iend, istep, &
|
||||
!$OMP istart, iend, istep, irp_here, &
|
||||
!$OMP ishift, idx0, u_t, maxab, v_0, s_0) &
|
||||
!$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 buffer, doubles, n_doubles, &
|
||||
!$OMP tmp_det2, hij, sij, idx, l, kcol_prev, v_t, &
|
||||
!$OMP singles_a, n_singles_a, singles_b, &
|
||||
!$OMP n_singles_b, s_t, k8)
|
||||
@ -181,12 +182,18 @@ subroutine H_S2_u_0_nstates_openmp_work_$N_int(v_0,s_0,u_t,N_st,sze,istart,iend,
|
||||
v_t = 0.d0
|
||||
s_t = 0.d0
|
||||
|
||||
ASSERT (iend <= N_det)
|
||||
ASSERT (istart > 0)
|
||||
ASSERT (istep > 0)
|
||||
|
||||
!$OMP DO SCHEDULE(dynamic,64)
|
||||
do k_a=istart+ishift,iend,istep
|
||||
|
||||
krow = psi_bilinear_matrix_rows(k_a)
|
||||
ASSERT (krow <= N_det_alpha_unique)
|
||||
|
||||
kcol = psi_bilinear_matrix_columns(k_a)
|
||||
ASSERT (kcol <= N_det_beta_unique)
|
||||
|
||||
tmp_det(1:$N_int,1) = psi_det_alpha_unique(1:$N_int, krow)
|
||||
tmp_det(1:$N_int,2) = psi_det_beta_unique (1:$N_int, kcol)
|
||||
@ -208,10 +215,15 @@ subroutine H_S2_u_0_nstates_openmp_work_$N_int(v_0,s_0,u_t,N_st,sze,istart,iend,
|
||||
tmp_det2(1:$N_int,2) = psi_det_beta_unique(1:$N_int, lcol)
|
||||
|
||||
l_a = psi_bilinear_matrix_columns_loc(lcol)
|
||||
ASSERT (l_a <= N_det)
|
||||
|
||||
do j=1,psi_bilinear_matrix_columns_loc(lcol+1) - l_a
|
||||
lrow = psi_bilinear_matrix_rows(l_a)
|
||||
ASSERT (lrow <= N_det_alpha_unique)
|
||||
|
||||
buffer(1:$N_int,j) = psi_det_alpha_unique(1:$N_int, lrow)
|
||||
|
||||
ASSERT (l_a <= N_det)
|
||||
idx(j) = l_a
|
||||
l_a = l_a+1
|
||||
enddo
|
||||
@ -226,7 +238,11 @@ subroutine H_S2_u_0_nstates_openmp_work_$N_int(v_0,s_0,u_t,N_st,sze,istart,iend,
|
||||
|
||||
do k = 1,n_singles_a
|
||||
l_a = singles_a(k)
|
||||
ASSERT (l_a <= N_det)
|
||||
|
||||
lrow = psi_bilinear_matrix_rows(l_a)
|
||||
ASSERT (lrow <= N_det_alpha_unique)
|
||||
|
||||
tmp_det2(1:$N_int,1) = psi_det_alpha_unique(1:$N_int, lrow)
|
||||
call i_H_j_double_alpha_beta(tmp_det,tmp_det2,$N_int,hij)
|
||||
call get_s2(tmp_det,tmp_det2,$N_int,sij)
|
||||
@ -255,7 +271,10 @@ subroutine H_S2_u_0_nstates_openmp_work_$N_int(v_0,s_0,u_t,N_st,sze,istart,iend,
|
||||
! -----------------------------------------------------------------------
|
||||
|
||||
krow = psi_bilinear_matrix_rows(k_a)
|
||||
ASSERT (krow <= N_det_alpha_unique)
|
||||
|
||||
kcol = psi_bilinear_matrix_columns(k_a)
|
||||
ASSERT (kcol <= N_det_beta_unique)
|
||||
|
||||
tmp_det(1:$N_int,1) = psi_det_alpha_unique(1:$N_int, krow)
|
||||
tmp_det(1:$N_int,2) = psi_det_beta_unique (1:$N_int, kcol)
|
||||
@ -264,19 +283,22 @@ subroutine H_S2_u_0_nstates_openmp_work_$N_int(v_0,s_0,u_t,N_st,sze,istart,iend,
|
||||
! ----------------------------------------------------------------------
|
||||
|
||||
k_b = psi_bilinear_matrix_order_transp_reverse(k_a)
|
||||
ASSERT (k_b <= N_det)
|
||||
|
||||
spindet(1:$N_int) = tmp_det(1:$N_int,1)
|
||||
|
||||
! Loop inside the beta column to gather all the connected alphas
|
||||
l_a = k_a+1
|
||||
do i=1,N_det_alpha_unique
|
||||
if (l_a > N_det) exit
|
||||
lcol = psi_bilinear_matrix_columns(l_a)
|
||||
if (lcol /= kcol) exit
|
||||
lrow = psi_bilinear_matrix_rows(l_a)
|
||||
ASSERT (lrow <= N_det_alpha_unique)
|
||||
|
||||
buffer(1:$N_int,i) = psi_det_alpha_unique(1:$N_int, lrow)
|
||||
idx(i) = l_a
|
||||
l_a = l_a+1
|
||||
if (l_a > N_det) exit
|
||||
enddo
|
||||
i = i-1
|
||||
|
||||
@ -290,9 +312,14 @@ subroutine H_S2_u_0_nstates_openmp_work_$N_int(v_0,s_0,u_t,N_st,sze,istart,iend,
|
||||
tmp_det2(1:$N_int,2) = psi_det_beta_unique (1:$N_int, kcol)
|
||||
do i=1,n_singles_a
|
||||
l_a = singles_a(i)
|
||||
ASSERT (l_a <= N_det)
|
||||
|
||||
lrow = psi_bilinear_matrix_rows(l_a)
|
||||
ASSERT (lrow <= N_det_alpha_unique)
|
||||
|
||||
tmp_det2(1:$N_int,1) = psi_det_alpha_unique(1:$N_int, lrow)
|
||||
call i_H_j_mono_spin( tmp_det, tmp_det2, $N_int, 1, hij)
|
||||
|
||||
do l=1,N_st
|
||||
v_t(l,l_a) = v_t(l,l_a) + hij * u_t(l,k_a)
|
||||
v_t(l,k_a) = v_t(l,k_a) + hij * u_t(l,l_a)
|
||||
@ -306,7 +333,11 @@ subroutine H_S2_u_0_nstates_openmp_work_$N_int(v_0,s_0,u_t,N_st,sze,istart,iend,
|
||||
|
||||
do i=1,n_doubles
|
||||
l_a = doubles(i)
|
||||
ASSERT (l_a <= N_det)
|
||||
|
||||
lrow = psi_bilinear_matrix_rows(l_a)
|
||||
ASSERT (lrow <= N_det_alpha_unique)
|
||||
|
||||
call i_H_j_double_spin( tmp_det(1,1), psi_det_alpha_unique(1, lrow), $N_int, hij)
|
||||
do l=1,N_st
|
||||
v_t(l,l_a) = v_t(l,l_a) + hij * u_t(l,k_a)
|
||||
@ -335,17 +366,20 @@ subroutine H_S2_u_0_nstates_openmp_work_$N_int(v_0,s_0,u_t,N_st,sze,istart,iend,
|
||||
! -----------------------------------------------------------------------
|
||||
|
||||
k_b = psi_bilinear_matrix_order_transp_reverse(k_a)
|
||||
ASSERT (k_b <= N_det)
|
||||
|
||||
! Loop inside the alpha row to gather all the connected betas
|
||||
l_b = k_b+1
|
||||
do i=1,N_det_beta_unique
|
||||
if (l_b > N_det) exit
|
||||
lrow = psi_bilinear_matrix_transp_rows(l_b)
|
||||
if (lrow /= krow) exit
|
||||
lcol = psi_bilinear_matrix_transp_columns(l_b)
|
||||
ASSERT (lcol <= N_det_beta_unique)
|
||||
|
||||
buffer(1:$N_int,i) = psi_det_beta_unique(1:$N_int, lcol)
|
||||
idx(i) = l_b
|
||||
l_b = l_b+1
|
||||
if (l_b > N_det) exit
|
||||
enddo
|
||||
i = i-1
|
||||
|
||||
@ -359,10 +393,15 @@ subroutine H_S2_u_0_nstates_openmp_work_$N_int(v_0,s_0,u_t,N_st,sze,istart,iend,
|
||||
tmp_det2(1:$N_int,1) = psi_det_alpha_unique(1:$N_int, krow)
|
||||
do i=1,n_singles_b
|
||||
l_b = singles_b(i)
|
||||
ASSERT (l_b <= N_det)
|
||||
|
||||
lcol = psi_bilinear_matrix_transp_columns(l_b)
|
||||
ASSERT (lcol <= N_det_beta_unique)
|
||||
|
||||
tmp_det2(1:$N_int,2) = psi_det_beta_unique (1:$N_int, lcol)
|
||||
call i_H_j_mono_spin( tmp_det, tmp_det2, $N_int, 2, hij)
|
||||
l_a = psi_bilinear_matrix_transp_order(l_b)
|
||||
ASSERT (l_a <= N_det)
|
||||
do l=1,N_st
|
||||
v_t(l,l_a) = v_t(l,l_a) + hij * u_t(l,k_a)
|
||||
v_t(l,k_a) = v_t(l,k_a) + hij * u_t(l,l_a)
|
||||
@ -375,9 +414,15 @@ subroutine H_S2_u_0_nstates_openmp_work_$N_int(v_0,s_0,u_t,N_st,sze,istart,iend,
|
||||
|
||||
do i=1,n_doubles
|
||||
l_b = doubles(i)
|
||||
ASSERT (l_b <= N_det)
|
||||
|
||||
lcol = psi_bilinear_matrix_transp_columns(l_b)
|
||||
ASSERT (lcol <= N_det_beta_unique)
|
||||
|
||||
call i_H_j_double_spin( tmp_det(1,2), psi_det_beta_unique(1, lcol), $N_int, hij)
|
||||
l_a = psi_bilinear_matrix_transp_order(l_b)
|
||||
ASSERT (l_a <= N_det)
|
||||
|
||||
do l=1,N_st
|
||||
v_t(l,l_a) = v_t(l,l_a) + hij * u_t(l,k_a)
|
||||
v_t(l,k_a) = v_t(l,k_a) + hij * u_t(l,l_a)
|
||||
@ -394,7 +439,10 @@ subroutine H_S2_u_0_nstates_openmp_work_$N_int(v_0,s_0,u_t,N_st,sze,istart,iend,
|
||||
! -----------------------------------------------------------------------
|
||||
|
||||
krow = psi_bilinear_matrix_rows(k_a)
|
||||
ASSERT (krow <= N_det_alpha_unique)
|
||||
|
||||
kcol = psi_bilinear_matrix_columns(k_a)
|
||||
ASSERT (kcol <= N_det_beta_unique)
|
||||
|
||||
tmp_det(1:$N_int,1) = psi_det_alpha_unique(1:$N_int, krow)
|
||||
tmp_det(1:$N_int,2) = psi_det_beta_unique (1:$N_int, kcol)
|
||||
|
@ -1929,7 +1929,7 @@ subroutine a_operator(iorb,ispin,key,hjj,Nint,na,nb)
|
||||
ASSERT (Nint > 0)
|
||||
|
||||
k = ishft(iorb-1,-bit_kind_shift)+1
|
||||
ASSERT (k > 0)
|
||||
ASSERT (k>0)
|
||||
l = iorb - ishft(k-1,bit_kind_shift)-1
|
||||
key(k,ispin) = ibclr(key(k,ispin),l)
|
||||
other_spin = iand(ispin,1)+1
|
||||
@ -1977,11 +1977,12 @@ subroutine ac_operator(iorb,ispin,key,hjj,Nint,na,nb)
|
||||
!DIR$ FORCEINLINE
|
||||
call bitstring_to_list_ab(key, occ, tmp, Nint)
|
||||
ASSERT (tmp(1) == elec_alpha_num)
|
||||
ASSERT (tmp(2) == elec_beta_num)
|
||||
ASSERT (tmp(2) == elec_alpha_num)
|
||||
|
||||
k = ishft(iorb-1,-bit_kind_shift)+1
|
||||
ASSERT (k > 0)
|
||||
ASSERT (k >0)
|
||||
l = iorb - ishft(k-1,bit_kind_shift)-1
|
||||
ASSERT (l >= 0)
|
||||
key(k,ispin) = ibset(key(k,ispin),l)
|
||||
other_spin = iand(ispin,1)+1
|
||||
|
||||
|
@ -189,9 +189,7 @@ integer function get_index_in_psi_det_alpha_unique(key,Nint)
|
||||
enddo
|
||||
i += 1
|
||||
|
||||
if (i > N_det_alpha_unique) then
|
||||
return
|
||||
endif
|
||||
ASSERT (i <= N_det_alpha_unique)
|
||||
|
||||
!DIR$ FORCEINLINE
|
||||
do while (spin_det_search_key(psi_det_alpha_unique(1,i),Nint) == det_ref)
|
||||
@ -213,6 +211,7 @@ integer function get_index_in_psi_det_alpha_unique(key,Nint)
|
||||
endif
|
||||
i += 1
|
||||
if (i > N_det_alpha_unique) then
|
||||
ASSERT (get_index_in_psi_det_alpha_unique > 0)
|
||||
return
|
||||
endif
|
||||
|
||||
@ -270,9 +269,7 @@ integer function get_index_in_psi_det_beta_unique(key,Nint)
|
||||
enddo
|
||||
i += 1
|
||||
|
||||
if (i > N_det_beta_unique) then
|
||||
return
|
||||
endif
|
||||
ASSERT (i <= N_det_beta_unique)
|
||||
|
||||
!DIR$ FORCEINLINE
|
||||
do while (spin_det_search_key(psi_det_beta_unique(1,i),Nint) == det_ref)
|
||||
@ -294,6 +291,7 @@ integer function get_index_in_psi_det_beta_unique(key,Nint)
|
||||
endif
|
||||
i += 1
|
||||
if (i > N_det_beta_unique) then
|
||||
ASSERT (get_index_in_psi_det_beta_unique > 0)
|
||||
return
|
||||
endif
|
||||
|
||||
@ -413,15 +411,20 @@ BEGIN_PROVIDER [ double precision, psi_bilinear_matrix_values, (N_det,N_states)
|
||||
!$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(i,j,k,l)
|
||||
do k=1,N_det
|
||||
i = get_index_in_psi_det_alpha_unique(psi_det(1,1,k),N_int)
|
||||
ASSERT (i>0)
|
||||
ASSERT (i<=N_det_alpha_unique)
|
||||
|
||||
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'
|
||||
ASSERT (j>0)
|
||||
ASSERT (j<=N_det_alpha_unique)
|
||||
|
||||
do l=1,N_states
|
||||
psi_bilinear_matrix_values(k,l) = psi_coef(k,l)
|
||||
enddo
|
||||
psi_bilinear_matrix_rows(k) = i
|
||||
psi_bilinear_matrix_columns(k) = j
|
||||
to_sort(k) = int(N_det_alpha_unique,8) * int(j-1,8) + int(i,8)
|
||||
ASSERT (to_sort(k) > 0_8)
|
||||
psi_bilinear_matrix_order(k) = k
|
||||
enddo
|
||||
!$OMP END PARALLEL DO
|
||||
@ -432,6 +435,12 @@ BEGIN_PROVIDER [ double precision, psi_bilinear_matrix_values, (N_det,N_states)
|
||||
call dset_order(psi_bilinear_matrix_values(1,l),psi_bilinear_matrix_order,N_det)
|
||||
enddo
|
||||
deallocate(to_sort)
|
||||
ASSERT (minval(psi_bilinear_matrix_rows) == 1)
|
||||
ASSERT (minval(psi_bilinear_matrix_columns) == 1)
|
||||
ASSERT (minval(psi_bilinear_matrix_order) == 1)
|
||||
ASSERT (maxval(psi_bilinear_matrix_rows) == N_det_alpha_unique)
|
||||
ASSERT (maxval(psi_bilinear_matrix_columns) == N_det_beta_unique)
|
||||
ASSERT (maxval(psi_bilinear_matrix_order) == N_det)
|
||||
END_PROVIDER
|
||||
|
||||
|
||||
@ -447,6 +456,8 @@ BEGIN_PROVIDER [ integer, psi_bilinear_matrix_order_reverse , (N_det) ]
|
||||
psi_bilinear_matrix_order_reverse(psi_bilinear_matrix_order(k)) = k
|
||||
enddo
|
||||
!$OMP END PARALLEL DO
|
||||
ASSERT (minval(psi_bilinear_matrix_order) == 1)
|
||||
ASSERT (maxval(psi_bilinear_matrix_order) == N_det)
|
||||
END_PROVIDER
|
||||
|
||||
|
||||
@ -477,6 +488,8 @@ BEGIN_PROVIDER [ integer, psi_bilinear_matrix_columns_loc, (N_det_beta_unique+1)
|
||||
endif
|
||||
enddo
|
||||
psi_bilinear_matrix_columns_loc(N_det_beta_unique+1) = N_det+1
|
||||
ASSERT (minval(psi_bilinear_matrix_columns_loc) == 1)
|
||||
ASSERT (maxval(psi_bilinear_matrix_columns_loc) == N_det+1)
|
||||
END_PROVIDER
|
||||
|
||||
BEGIN_PROVIDER [ double precision, psi_bilinear_matrix_transp_values, (N_det,N_states) ]
|
||||
@ -510,16 +523,17 @@ BEGIN_PROVIDER [ double precision, psi_bilinear_matrix_transp_values, (N_det,N_
|
||||
!$OMP DO
|
||||
do k=1,N_det
|
||||
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
|
||||
ASSERT (psi_bilinear_matrix_transp_columns(k) > 0)
|
||||
ASSERT (psi_bilinear_matrix_transp_columns(k) <= N_det)
|
||||
|
||||
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
|
||||
ASSERT (psi_bilinear_matrix_transp_rows(k) > 0)
|
||||
ASSERT (psi_bilinear_matrix_transp_rows(k) <= N_det)
|
||||
|
||||
i = psi_bilinear_matrix_transp_columns(k)
|
||||
j = psi_bilinear_matrix_transp_rows (k)
|
||||
to_sort(k) = int(N_det_beta_unique,8) * int(j-1,8) + int(i,8)
|
||||
ASSERT (to_sort(k) > 0)
|
||||
psi_bilinear_matrix_transp_order(k) = k
|
||||
enddo
|
||||
!$OMP ENDDO
|
||||
@ -531,6 +545,12 @@ BEGIN_PROVIDER [ double precision, psi_bilinear_matrix_transp_values, (N_det,N_
|
||||
call dset_order(psi_bilinear_matrix_transp_values(1,l),psi_bilinear_matrix_transp_order,N_det)
|
||||
enddo
|
||||
deallocate(to_sort)
|
||||
ASSERT (minval(psi_bilinear_matrix_transp_columns) == 1)
|
||||
ASSERT (minval(psi_bilinear_matrix_transp_rows) == 1)
|
||||
ASSERT (minval(psi_bilinear_matrix_transp_order) == 1)
|
||||
ASSERT (maxval(psi_bilinear_matrix_transp_columns) == N_det_beta_unique)
|
||||
ASSERT (maxval(psi_bilinear_matrix_transp_rows) == N_det_alpha_unique)
|
||||
ASSERT (maxval(psi_bilinear_matrix_transp_order) == N_det)
|
||||
END_PROVIDER
|
||||
|
||||
BEGIN_PROVIDER [ integer, psi_bilinear_matrix_transp_rows_loc, (N_det_alpha_unique+1) ]
|
||||
@ -552,6 +572,8 @@ BEGIN_PROVIDER [ integer, psi_bilinear_matrix_transp_rows_loc, (N_det_alpha_uniq
|
||||
endif
|
||||
enddo
|
||||
psi_bilinear_matrix_transp_rows_loc(N_det_alpha_unique+1) = N_det+1
|
||||
ASSERT (minval(psi_bilinear_matrix_transp_rows_loc) == 1)
|
||||
ASSERT (maxval(psi_bilinear_matrix_transp_rows_loc) == N_det+1)
|
||||
END_PROVIDER
|
||||
|
||||
BEGIN_PROVIDER [ integer, psi_bilinear_matrix_order_transp_reverse , (N_det) ]
|
||||
@ -562,11 +584,14 @@ BEGIN_PROVIDER [ integer, psi_bilinear_matrix_order_transp_reverse , (N_det) ]
|
||||
END_DOC
|
||||
integer :: k
|
||||
|
||||
psi_bilinear_matrix_order_transp_reverse = -1
|
||||
!$OMP PARALLEL DO DEFAULT(SHARED) PRIVATE(k)
|
||||
do k=1,N_det
|
||||
psi_bilinear_matrix_order_transp_reverse(psi_bilinear_matrix_transp_order(k)) = k
|
||||
enddo
|
||||
!$OMP END PARALLEL DO
|
||||
ASSERT (minval(psi_bilinear_matrix_order_transp_reverse) == 1)
|
||||
ASSERT (maxval(psi_bilinear_matrix_order_transp_reverse) == N_det)
|
||||
END_PROVIDER
|
||||
|
||||
|
||||
|
Loading…
Reference in New Issue
Block a user