From eac877dd44705f0a5cfad315cc984f8a5c79ef01 Mon Sep 17 00:00:00 2001 From: Anthony Scemama Date: Wed, 26 Feb 2020 23:30:38 +0100 Subject: [PATCH] Removed dead code --- scripts/generate_h_apply.py | 7 ------- src/tools/print_energy.irp.f | 2 +- src/tools/print_hamiltonian.irp.f | 29 +++++++++++++++++++++++++++++ 3 files changed, 30 insertions(+), 8 deletions(-) create mode 100644 src/tools/print_hamiltonian.irp.f diff --git a/scripts/generate_h_apply.py b/scripts/generate_h_apply.py index 0f63308b..8c399e55 100644 --- a/scripts/generate_h_apply.py +++ b/scripts/generate_h_apply.py @@ -11,7 +11,6 @@ declarations decls_main deinit_thread init_main -filter_integrals filter2p filter2h2p_double filter2h2p_single @@ -106,12 +105,6 @@ class H_apply(object): s["do_double_excitations"] = d[do_double_exc] s["keys_work"] += "call fill_H_apply_buffer_no_selection(key_idx,keys_out,N_int,iproc)" - s["filter_integrals"] = "array_pairs = .True." - if SingleRef: - s["filter_integrals"] = """ - call get_mo_bielec_integrals_existing_ik(i_a,j_a,mo_num,array_pairs,mo_integrals_map) - """ - s["generate_psi_guess"] = """ ! Sort H_jj to find the N_states lowest states integer :: i diff --git a/src/tools/print_energy.irp.f b/src/tools/print_energy.irp.f index be7075c1..4703e7d4 100644 --- a/src/tools/print_energy.irp.f +++ b/src/tools/print_energy.irp.f @@ -1,4 +1,4 @@ -program print_wf +program print_energy implicit none BEGIN_DOC ! Prints the energy of the wave function stored in the |EZFIO| directory. diff --git a/src/tools/print_hamiltonian.irp.f b/src/tools/print_hamiltonian.irp.f new file mode 100644 index 00000000..207161dd --- /dev/null +++ b/src/tools/print_hamiltonian.irp.f @@ -0,0 +1,29 @@ +program print_hamiltonian + implicit none + BEGIN_DOC + ! Prints the Hamiltonian matrix defined in the space of determinants + ! present in the |EZFIO| directory. + END_DOC + + ! this has to be done in order to be sure that N_det, psi_det and + ! psi_coef_sorted are the wave function stored in the |EZFIO| directory. + read_wf = .True. + touch read_wf + call run +end + +subroutine run + implicit none + integer :: i, j + double precision :: hij + + do j=1,N_det + do i=1,N_det + call i_H_j(psi_det(1,1,i), psi_det(1,1,j), N_int, hij) + if (dabs(hij) > 1.d-20) then + print *, i, j, hij + endif + enddo + enddo + +end