init_exact_spectrum Subroutine

public subroutine init_exact_spectrum(ndets)

Arguments

Type IntentOptional Attributes Name
integer, intent(out) :: ndets

Contents

Source Code


Source Code

    subroutine init_exact_spectrum(ndets)

        use bit_rep_data, only: NIfTot, nifd
        use DetBitOps, only: EncodeBitDet, ilut_lt, ilut_gt
        use exact_diag, only: calculate_full_hamiltonian
        use gndts_mod, only: gndts
        use sort_mod, only: sort
        use SystemData, only: nbasis, nel, BRR, nBasisMax, G1, tSpn, lms, tParity, SymRestrict

        integer, intent(out) :: ndets
        integer, allocatable :: nI_list(:, :)
        integer(n_int), allocatable :: ilut_list(:, :)
        integer(n_int) :: ilut(0:NIfTot)
        integer :: i, hf_ind, temp(1, 1), ierr
        character(len=*), parameter :: t_r = 'init_exact_spectrum'

        write(stdout, '(/,1x,a48,/)') "Beginning calculation of exact spectral density."
        call neci_flush(stdout)

        write(stdout, '(1x,a56)', advance='yes') "Enumerating and storing all determinants in the space..."
        call neci_flush(stdout)

        ! Calculate the number of determinants.
        call gndts(nel, nbasis, BRR, nBasisMax, temp, .true., G1, tSpn, lms, tParity, SymRestrict, ndets, hf_ind)
        allocate(nI_list(nel, ndets))
        ! Now generate them again and store them this time.
        call gndts(nel, nbasis, BRR, nBasisMax, nI_list, .false., G1, tSpn, lms, tParity, SymRestrict, ndets, hf_ind)

        write(stdout, '(1x,a33,i7)') "Number of determinants in space: ", ndets

        allocate(ilut_list(0:NIfTot, ndets))
        ilut_list = 0_n_int

        do i = 1, ndets
            call EncodeBitDet(nI_list(:, i), ilut)
            ilut_list(0:nifd, i) = ilut(0:nifd)
        end do

        call sort(ilut_list, ilut_lt, ilut_gt)

        allocate(eigv_es(ndets), stat=ierr)
        if (ierr /= 0) then
            write(stdout, '(1x,a11,1x,i5)') "Error code:", ierr
            call stop_all(t_r, "Error allocating eigenvalue array.")
        end if
        eigv_es = 0.0_dp

        allocate(trans_amps_left(ndets), stat=ierr)
        if (ierr /= 0) then
            write(stdout, '(1x,a11,1x,i5)') "Error code:", ierr
            call stop_all(t_r, "Error allocating array to hold left transition amplitudes.")
        end if
        trans_amps_left = 0.0_dp

        allocate(trans_amps_right(ndets), stat=ierr)
        if (ierr /= 0) then
            write(stdout, '(1x,a11,1x,i5)') "Error code:", ierr
            call stop_all(t_r, "Error allocating array to hold right transition amplitudes.")
        end if
        trans_amps_right = 0.0_dp

        allocate(pert_ground_left(ndets), stat=ierr)
        if (ierr /= 0) then
            write(stdout, '(1x,a11,1x,i5)') "Error code:", ierr
            call stop_all(t_r, "Error allocating array to hold left perturbed ground state components.")
        end if

        allocate(pert_ground_right(ndets), stat=ierr)
        if (ierr /= 0) then
            write(stdout, '(1x,a11,1x,i5)') "Error code:", ierr
            call stop_all(t_r, "Error allocating array to hold right perturbed ground state components.")
        end if

        write(stdout, '(1x,a48)') "Allocating and calculating Hamiltonian matrix..."
        call neci_flush(stdout)
        allocate(hamiltonian(ndets, ndets), stat=ierr)
        if (ierr /= 0) then
            write(stdout, '(1x,a11,1x,i5)') "Error code:", ierr
            call stop_all(t_r, "Error allocating Hamiltonian array.")
        end if
        write(stdout, '(1x,a46)') "Hamiltonian allocation completed successfully."
        call neci_flush(stdout)
        call calculate_full_hamiltonian(ilut_list, hamiltonian)
        write(stdout, '(1x,a33)') "Hamiltonian calculation complete."
        call neci_flush(stdout)

        call return_perturbed_ground_spec(left_perturb_spectral, ilut_list, pert_ground_left, left_pert_norm)
        call return_perturbed_ground_spec(right_perturb_spectral, ilut_list, pert_ground_right, right_pert_norm)

        deallocate(ilut_list)
        deallocate(nI_list)

    end subroutine init_exact_spectrum