gen_guga_heisenberg_cum_list Subroutine

private subroutine gen_guga_heisenberg_cum_list(csf_i, id, cum_arr, tgt, tgt_pgen)

Arguments

Type IntentOptional Attributes Name
type(CSF_Info_t), intent(in) :: csf_i
integer, intent(in) :: id
real(kind=dp), intent(out), allocatable :: cum_arr(:)
integer, intent(in), optional :: tgt
real(kind=dp), intent(out), optional :: tgt_pgen

Contents


Source Code

    subroutine gen_guga_heisenberg_cum_list(csf_i, id, cum_arr, tgt, tgt_pgen)
        ! make a routine for this, for easy pgen recalculation
        type(CSF_Info_t), intent(in) :: csf_i
        integer, intent(in) :: id
        real(dp), allocatable, intent(out) :: cum_arr(:)
        integer, intent(in), optional :: tgt
        real(dp), intent(out), optional :: tgt_pgen

        integer :: step, i, n
        integer, allocatable :: neighbors(:)
        real(dp) :: cum_sum, tmp

        neighbors = lat%get_neighbors(lat%get_site_index(id))

        ! then check if the neighbors are available for exchange
        allocate(cum_arr(size(neighbors)), source=0.0_dp)

        cum_sum = 0.0_dp
        tmp = 0.0_dp

        step = csf_i%stepvector(id)

        if (present(tgt)) then

            tgt_pgen = 0.0_dp
            ! if target is not in the neighbors list we can exit
            if (.not. any(tgt == neighbors)) return

            do i = 1, size(neighbors)
                n = neighbors(i)

                if (csf_i%stepvector(n) == 0) then
                    cycle

                else if (csf_i%stepvector(n) == step) then
                    if (abs(id - n) == 1) then
                        cycle
                    else

                        if (step == 1 .and. count_alpha_orbs_ij(csf_i, min(id, n), max(id, n)) == 0) cycle

                        if (step == 2 .and. count_beta_orbs_ij(csf_i, min(id, n), max(id, n)) == 0) cycle

                        tmp = 1.0_dp

                        cum_sum = cum_sum + tmp

                        if (tgt == n) tgt_pgen = tmp
                    end if

                else if (csf_i%stepvector(n) /= 0 &
                         .and. csf_i%stepvector(n) /= step) then

                    if (id - n == -1 &
                            .and. csf_i%stepvector(id) == 1 &
                            .and.  csf_i%B_int(id) == 1) then
                        cycle
                    end if

                    if (id - n == 1 &
                            .and. csf_i%stepvector(n) == 1 &
                            .and. csf_i%B_int(n) == 1) then
                        cycle
                    end if

                    tmp = 1.0_dp

                    cum_sum = cum_sum + tmp

                    if (tgt == n) tgt_pgen = tmp
                end if
            end do

            if (.not. near_zero(cum_sum)) then
                tgt_pgen = tgt_pgen / cum_sum
            end if

        else
            do i = 1, size(neighbors)
                n = neighbors(i)
                if (csf_i%stepvector(n) == 0) then
                    ! t-J case
                    cum_arr(i) = cum_sum
                    cycle

                else if (csf_i%stepvector(n) == step) then
                    ! this can only be if we have space
                    ! between the step-vectors
                    if (abs(id - n) == 1) then
                        cum_arr(i) = cum_sum
                        cycle
                    else
                        if (step == 1 &
                                .and. count_alpha_orbs_ij(csf_i, min(id, n), max(id, n)) == 0) then
                            cum_arr(i) = cum_sum
                            cycle
                        end if

                        if (step == 2 &
                                .and. count_beta_orbs_ij(csf_i, min(id, n), max(id, n)) == 0) then
                            cum_arr(i) = cum_sum
                            cycle
                        end if

                        cum_arr(i) = cum_sum + 1.0_dp
                        cum_sum = cum_sum + 1.0_dp
                    end if

                else if (csf_i%stepvector(n) /= 0 .and. csf_i%stepvector(n) /= step) then

                    if (id - n == -1 &
                            .and. csf_i%stepvector(id) == 1 &
                            .and.  csf_i%B_int(id) == 1) then
                        cum_arr(i) = cum_sum
                        cycle
                    end if

                    if (id - n == 1 &
                            .and. csf_i%stepvector(n) == 1 &
                            .and. csf_i%B_int(n) == 1) then
                        cum_arr(i) = cum_sum
                        cycle
                    end if

                    cum_arr(i) = cum_sum + 1.0_dp
                    cum_sum = cum_sum + 1.0_dp
                end if
            end do
        end if

    end subroutine gen_guga_heisenberg_cum_list