calc_contribs_spawn Subroutine

public subroutine calc_contribs_spawn(nvecs, krylov_array, krylov_ht, nspawns_this_proc, fac, est_matrix)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: nvecs
integer(kind=n_int), intent(in) :: krylov_array(0:,:)
type(ll_node), intent(in), pointer :: krylov_ht(:)
integer, intent(in) :: nspawns_this_proc
real(kind=dp), intent(in) :: fac
real(kind=dp), intent(inout) :: est_matrix(:,:)

Contents

Source Code


Source Code

    subroutine calc_contribs_spawn(nvecs, krylov_array, krylov_ht, nspawns_this_proc, fac, est_matrix)

        use bit_rep_data, only: IlutBits
        use bit_reps, only: decode_bit_det
        use FciMCData, only: SpawnedParts, ll_node
        use hash, only: FindWalkerHash
        use SystemData, only: nel

        integer, intent(in) :: nvecs
        integer(n_int), intent(in) :: krylov_array(0:, :)
        type(ll_node), pointer, intent(in) :: krylov_ht(:)
        integer, intent(in) :: nspawns_this_proc
        real(dp), intent(in) :: fac
        real(dp), intent(inout) :: est_matrix(:, :)

        integer :: idet, det_hash, det_ind, i, j
        integer(n_int) :: ilut_spawn(0:NIfTot)
        integer :: nI_spawn(nel)
        integer(n_int) :: int_sign(lenof_all_signs)
        real(dp) :: real_sign_1(lenof_all_signs), real_sign_2(lenof_all_signs)
        type(ll_node), pointer :: temp_node
        logical :: tDetFound

        ilut_spawn = 0_n_int

        do idet = 1, nspawns_this_proc

            ilut_spawn(0:nifd) = SpawnedParts(0:nifd, idet)
            int_sign = SpawnedParts(IlutBits%ind_pop:IlutBits%ind_pop + lenof_all_signs - 1, idet)
            real_sign_1 = transfer(int_sign, real_sign_1)
            call decode_bit_det(nI_spawn, ilut_spawn)
            det_hash = FindWalkerHash(nI_spawn, size(krylov_ht))
            ! Point to the first node with this hash value in krylov_array.
            temp_node => krylov_ht(det_hash)

            if (temp_node%ind == 0) then
                ! If there are no determinants at all with this hash value in krylov_array.
                cycle
            else
                tDetFound = .false.
                do while (associated(temp_node))
                    if (all(ilut_spawn(0:nifd) == krylov_array(0:nifd, temp_node%ind))) then
                        ! If this CurrentDets determinant has been found in krylov_array.
                        det_ind = temp_node%ind
                        tDetFound = .true.
                        exit
                    end if
                    ! Move on to the next determinant with this hash value.
                    temp_node => temp_node%next
                end do
                if (tDetFound) then
                    int_sign = &
                        krylov_array(IlutBits%ind_pop:IlutBits%ind_pop + lenof_all_signs - 1, &
                                     det_ind)
                    real_sign_2 = transfer(int_sign, real_sign_1)
                    if (IsUnoccDet(real_sign_2)) cycle

                    ! Finally, add in the contribution to the projected Hamiltonian for each pair of Krylov vectors.
                    do i = 1, nvecs
                        do j = i, nvecs
                            est_matrix(i, j) = est_matrix(i, j) + &
                                               fac * (real_sign_1(kp_ind_1(i)) * real_sign_2(kp_ind_2(j)) + &
                                                      real_sign_1(kp_ind_2(i)) * real_sign_2(kp_ind_1(j))) / 2.0_dp
                        end do
                    end do
                end if
            end if
        end do

    end subroutine calc_contribs_spawn