diff --git a/plugins/local/dmc_dress/EZFIO.cfg b/plugins/local/dmc_dress/EZFIO.cfg new file mode 100644 index 00000000..88bc4ff3 --- /dev/null +++ b/plugins/local/dmc_dress/EZFIO.cfg @@ -0,0 +1,6 @@ +[dmc_delta_h] +type: double precision +doc: Dressing matrix obtained from DMC +size: (determinants.n_det) +interface: ezfio, provider + diff --git a/plugins/local/dmc_dress/NEED b/plugins/local/dmc_dress/NEED new file mode 100644 index 00000000..657d6cfd --- /dev/null +++ b/plugins/local/dmc_dress/NEED @@ -0,0 +1,3 @@ +selectors_full +generators_full +davidson_dressed diff --git a/plugins/local/dmc_dress/dmc_dress.irp.f b/plugins/local/dmc_dress/dmc_dress.irp.f new file mode 100644 index 00000000..f79b98bf --- /dev/null +++ b/plugins/local/dmc_dress/dmc_dress.irp.f @@ -0,0 +1,19 @@ +program diagonalize_h + implicit none + BEGIN_DOC +! Program that extracts the lowest states of the Hamiltonian dressed by the QMC +! dressing vector stored in :option:`dmc_dressing dmc_delta_h` +! + END_DOC + read_wf = .True. + touch read_wf + call routine + call save_wavefunction_general(N_det,N_states,psi_det_sorted,size(psi_coef_sorted,1),psi_coef_sorted) +end + +subroutine routine + implicit none + psi_coef(1:N_det,1) = ci_eigenvectors_dressed(1:N_det,1) + print*,'N_det = ',N_det + SOFT_TOUCH psi_coef +end diff --git a/plugins/local/dmc_dress/dressing_vector.irp.f b/plugins/local/dmc_dress/dressing_vector.irp.f new file mode 100644 index 00000000..9f799967 --- /dev/null +++ b/plugins/local/dmc_dress/dressing_vector.irp.f @@ -0,0 +1,24 @@ + BEGIN_PROVIDER [ double precision, dressing_column_h, (N_det,N_states) ] +&BEGIN_PROVIDER [ double precision, dressing_column_s, (N_det,N_states) ] + implicit none + BEGIN_DOC + ! \Delta_{state-specific}. \Psi + END_DOC + + integer :: i,ii,k,j, l + double precision :: f, tmp + double precision, external :: u_dot_v + logical, external :: detEq + + dressing_column_h(:,:) = 0.d0 + dressing_column_s(:,:) = 0.d0 + + l = dressed_column_idx(1) + do j = 1, n_det + dressing_column_h(j,1) = 0.5d0*dmc_delta_h(j) + dressing_column_h(l,1) -= 0.5d0 * psi_coef(j,1) * dmc_delta_h(j) /psi_coef(l,1) + enddo +END_PROVIDER + + +