mirror of
https://github.com/LCPQ/quantum_package
synced 2025-01-11 05:28:29 +01:00
Completed SVD ezfio for QMC=Chem
This commit is contained in:
parent
3284b9d645
commit
be8818538a
@ -42,8 +42,10 @@ subroutine save_dets_qmcchem
|
||||
enddo
|
||||
close(31)
|
||||
call system('gzip -f '//trim(ezfio_filename)//'/mo_basis/mo_classif')
|
||||
|
||||
end
|
||||
|
||||
program save_for_qmc
|
||||
call save_dets_qmcchem
|
||||
call write_spindeterminants
|
||||
end
|
||||
|
@ -7,4 +7,8 @@ spindeterminants
|
||||
psi_det_alpha integer*8 (spindeterminants_n_int*spindeterminants_bit_kind/8,spindeterminants_n_det_alpha)
|
||||
psi_det_beta integer*8 (spindeterminants_n_int*spindeterminants_bit_kind/8,spindeterminants_n_det_beta)
|
||||
psi_coef_matrix double precision (spindeterminants_n_det_alpha,spindeterminants_n_det_beta,spindeterminants_n_states)
|
||||
n_svd_coefs integer
|
||||
psi_svd_alpha double precision (spindeterminants_n_det_alpha,spindeterminants_n_svd_coefs,spindeterminants_n_states)
|
||||
psi_svd_beta double precision (spindeterminants_n_det_beta,spindeterminants_n_svd_coefs,spindeterminants_n_states)
|
||||
psi_svd_coefs double precision (spindeterminants_n_svd_coefs,spindeterminants_n_states)
|
||||
|
||||
|
@ -26,7 +26,7 @@ subroutine write_spindeterminants
|
||||
enddo
|
||||
call ezfio_set_spindeterminants_psi_det_alpha(psi_det_alpha_unique)
|
||||
deallocate(tmpdet)
|
||||
|
||||
|
||||
allocate(tmpdet(N_int2,N_det_beta_unique))
|
||||
do i=1,N_det_beta_unique
|
||||
do k=1,N_int
|
||||
@ -38,7 +38,54 @@ subroutine write_spindeterminants
|
||||
enddo
|
||||
call ezfio_set_spindeterminants_psi_det_beta(psi_det_beta_unique)
|
||||
deallocate(tmpdet)
|
||||
|
||||
|
||||
call ezfio_set_spindeterminants_psi_coef_matrix(psi_svd_matrix)
|
||||
|
||||
integer :: n_svd_coefs
|
||||
double precision :: norm, f
|
||||
f = 1.d0/dble(N_states)
|
||||
norm = 1.d0
|
||||
do n_svd_coefs=1,N_det_alpha_unique
|
||||
do k=1,N_states
|
||||
norm -= psi_svd_coefs(n_svd_coefs,k)*psi_svd_coefs(n_svd_coefs,k)
|
||||
enddo
|
||||
if (norm < 1.d-6) then
|
||||
exit
|
||||
endif
|
||||
enddo
|
||||
n_svd_coefs -= 1
|
||||
call ezfio_set_spindeterminants_n_svd_coefs(n_svd_coefs)
|
||||
|
||||
double precision, allocatable :: dtmp(:,:,:)
|
||||
allocate(dtmp(N_det_alpha_unique,n_svd_coefs,N_states))
|
||||
do k=1,N_states
|
||||
do j=1,n_svd_coefs
|
||||
do i=1,N_det_alpha_unique
|
||||
dtmp(i,j,k) = psi_svd_alpha(i,j,k)
|
||||
enddo
|
||||
enddo
|
||||
enddo
|
||||
call ezfio_set_spindeterminants_psi_svd_alpha(dtmp)
|
||||
deallocate(dtmp)
|
||||
|
||||
allocate(dtmp(N_det_beta_unique,n_svd_coefs,N_states))
|
||||
do k=1,N_states
|
||||
do j=1,n_svd_coefs
|
||||
do i=1,N_det_beta_unique
|
||||
dtmp(i,j,k) = psi_svd_beta(i,j,k)
|
||||
enddo
|
||||
enddo
|
||||
enddo
|
||||
call ezfio_set_spindeterminants_psi_svd_beta(dtmp)
|
||||
deallocate(dtmp)
|
||||
|
||||
allocate(dtmp(n_svd_coefs,N_states,1))
|
||||
do k=1,N_states
|
||||
do j=1,n_svd_coefs
|
||||
dtmp(j,k,1) = psi_svd_coefs(j,k)
|
||||
enddo
|
||||
enddo
|
||||
call ezfio_set_spindeterminants_psi_svd_coefs(dtmp)
|
||||
deallocate(dtmp)
|
||||
|
||||
end
|
||||
|
Loading…
Reference in New Issue
Block a user