init_hop_trancorr_fac_cached_vec Subroutine

private subroutine init_hop_trancorr_fac_cached_vec(J_fac, in_lat)

Arguments

Type IntentOptional Attributes Name
real(kind=dp), intent(in) :: J_fac
class(lattice), intent(in) :: in_lat

Contents


Source Code

    subroutine init_hop_trancorr_fac_cached_vec(J_fac, in_lat)
        ! maybe it is better to cache it to be accessible by the r-vecs,
        ! since this is the main use of this functionality and then we would
        ! not need to map r-vectors to orbital indices..
        real(dp), intent(in) :: J_fac
        class(lattice), intent(in) :: in_lat
#ifdef DEBUG_
        character(*), parameter :: this_routine = "init_hop_trancorr_fac_cached_vec"
#endif
        integer :: n_sites, i, j, k, ri(3), rj(3), r_vec(3), r_min(3), r_max(3)
        integer :: k_vec(3)
        complex(dp) :: temp, temp_m

        ! similar to Kais way of initializing the BZ zone i need to find the
        ! lower and upper bounds of our arrays..
        ! if the bounds are not .and. allocated(hop_transcorr_factor_cached_vec
        ! then init them
        call in_lat%init_hop_cache_bounds(r_min, r_max)

        if (.not. (t_recalc_umat .or. t_recalc_tmat) .and. allocated(hop_transcorr_factor_cached_vec)) then
            ! then we have already done that..
            return
        end if

        if (allocated(hop_transcorr_factor_cached_vec)) deallocate(hop_transcorr_factor_cached_vec)
        if (allocated(hop_transcorr_factor_cached_vec_m)) deallocate(hop_transcorr_factor_cached_vec_m)

        allocate (hop_transcorr_factor_cached_vec( &
                  r_min(1):r_max(1), r_min(2):r_max(2), r_min(3):r_max(3)))

        ! i also need one for -J!
        allocate (hop_transcorr_factor_cached_vec_m( &
                  r_min(1):r_max(1), r_min(2):r_max(2), r_min(3):r_max(3)))

        hop_transcorr_factor_cached_vec = 0.0_dp
        hop_transcorr_factor_cached_vec_m = 0.0_dp

        n_sites = in_lat%get_nsites()

        do i = 1, n_sites
            ri = in_lat%get_r_vec(i)
            do j = 1, n_sites
                rj = in_lat%get_r_vec(j)
                r_vec = ri - rj

                temp = 0.0_dp
                temp_m = 0.0_dp

                do k = 1, n_sites
                    k_vec = in_lat%get_k_vec(k)

                    temp = temp + exp(imag_unit * lat%dot_prod(k_vec, r_vec)) * &
                           exp(-J_fac * epsilon_kvec(k))

                    temp_m = temp_m + exp(imag_unit * lat%dot_prod(k_vec, r_vec)) * &
                             exp(J_fac * epsilon_kvec(k))
                end do

                ASSERT(abs(aimag(temp)) < EPS)
                ASSERT(abs(aimag(temp_m)) < EPS)

                hop_transcorr_factor_cached_vec(r_vec(1), r_vec(2), r_vec(3)) = &
                    real(temp) / real(n_sites, dp)

                hop_transcorr_factor_cached_vec_m(r_vec(1), r_vec(2), r_vec(3)) = &
                    real(temp_m) / real(n_sites, dp)
            end do
        end do

    end subroutine init_hop_trancorr_fac_cached_vec