get_proj_e_for_preconditioner Subroutine

public subroutine get_proj_e_for_preconditioner(ValidSpawned, proj_energy)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: ValidSpawned
real(kind=dp), intent(out) :: proj_energy(lenof_sign)

Contents


Source Code

    subroutine get_proj_e_for_preconditioner(ValidSpawned, proj_energy)

        integer, intent(in) :: ValidSpawned
        real(dp), intent(out) :: proj_energy(lenof_sign)

        integer :: i, run, ierr
        real(dp) :: SignTemp(lenof_sign), ref_pop(lenof_sign)
        logical :: ref_found(lenof_sign), tSuccess
        integer :: PartInd, DetHash
        character(*), parameter :: t_r = "get_proj_e_for_preconditioner"

        proj_energy = 0.0_dp
        ref_found = .false.

        ! Find the weight spawned on the Hartree--Fock determinant.
        if (tSemiStochastic) then
            if (.not. cs_replicas(core_run)%t_global) then
                call stop_all(t_r, "Invalid core space for preconditioning. Global core space required")
            end if

            do run = 1, lenof_sign
                associate(rep => cs_replicas(core_run))
                do i = 1, rep%determ_sizes(iProcIndex)
                    if (DetBitEQ(rep%core_space(0:nifd, rep%determ_displs(iProcIndex) + i), iLutRef(:, run), nifd)) then
                        ! This might need adjustment for the complex case
                        proj_energy(run) = -rep%partial_determ_vecs(run, i)
                        ref_found(run) = .true.
                    end if
                end do
                end associate
            end do
        end if

        do i = 1, ValidSpawned
            do run = 1, lenof_sign
                if (DetBitEQ(SpawnedParts(:, i), iLutRef(:, run), nifd)) then
                    call extract_sign(SpawnedParts(:, i), SignTemp)
                    proj_energy(run) = proj_energy(run) - SignTemp(run)
                    ref_found(run) = .true.
                end if
            end do
        end do

        do run = 1, lenof_sign
            if (iProcIndex == iRefProc(run)) then
                if ((.not. ref_found(run))) then
                    proj_energy(run) = -0.01_dp
                    write(stdout, '("WARNING: The reference determinant was not spawned to in the last iteration.")')
                else if (abs(proj_energy(run)) < 1.e-12_dp) then
                    proj_energy(run) = -0.01_dp
                    write(stdout, '("WARNING: The projected energy from the last iteration was zero. Setting to -0.1.")')
                end if

                call hash_table_lookup(ProjEDet(:, run), ilutRef(:, run), nifd, HashIndex, &
                                       CurrentDets, PartInd, DetHash, tSuccess)

                if (tSuccess) then
                    call extract_sign(CurrentDets(:, PartInd), ref_pop)
                    proj_energy(run) = proj_energy(run) / ref_pop(run)
                else
                    write(stdout, '("WARNING: Reference determinant not found in main walker list.")')
                end if
            end if

            call MPIBarrier(ierr)
            call MPIBCast(proj_energy(run), 1, iRefProc(run))
        end do

        ! Remove time step
        proj_energy = proj_energy / tau

    end subroutine get_proj_e_for_preconditioner