gyrotropic_get_NOA_Bnl_orb Subroutine

private subroutine gyrotropic_get_NOA_Bnl_orb(eig, del_eig, AA, num_occ, occ_list, num_unocc, unocc_list, Bnl)

Uses

  • proc~~gyrotropic_get_noa_bnl_orb~~UsesGraph proc~gyrotropic_get_noa_bnl_orb gyrotropic_get_NOA_Bnl_orb module~w90_parameters w90_parameters proc~gyrotropic_get_noa_bnl_orb->module~w90_parameters module~w90_constants w90_constants proc~gyrotropic_get_noa_bnl_orb->module~w90_constants module~w90_parameters->module~w90_constants module~w90_io w90_io module~w90_parameters->module~w90_io module~w90_io->module~w90_constants

Arguments

Type IntentOptional AttributesName
real(kind=dp), intent(in), dimension(:):: eig
real(kind=dp), intent(in), dimension(:, :):: del_eig
complex(kind=dp), intent(in), dimension(:, :, :):: AA
integer, intent(in) :: num_occ
integer, intent(in), dimension(:):: occ_list
integer, intent(in) :: num_unocc
integer, intent(in), dimension(:):: unocc_list
complex(kind=dp), intent(out), dimension(:, :, :, :):: Bnl

Called by

proc~~gyrotropic_get_noa_bnl_orb~~CalledByGraph proc~gyrotropic_get_noa_bnl_orb gyrotropic_get_NOA_Bnl_orb proc~gyrotropic_get_noa_k gyrotropic_get_NOA_k proc~gyrotropic_get_noa_k->proc~gyrotropic_get_noa_bnl_orb proc~gyrotropic_get_k_list gyrotropic_get_k_list proc~gyrotropic_get_k_list->proc~gyrotropic_get_noa_k proc~gyrotropic_main gyrotropic_main proc~gyrotropic_main->proc~gyrotropic_get_k_list

Contents


Source Code

  subroutine gyrotropic_get_NOA_Bnl_orb(eig, del_eig, AA, &
                                        num_occ, occ_list, num_unocc, unocc_list, Bnl)
    !====================================================================!
    !                                                                    !
    ! Calculating the matrix                                             !
    ! B_{nl,ac}(num_occ,num_unocc,3,3)=                    !
    !      -sum_m(  (E_m-E_n)A_nma*Amlc +(E_l-E_m)A_nmc*A_mla -          !
    !      -i( del_a (E_n+E_l) A_nlc                                     !
    !   in units eV*Ang^2                                                !
    !====================================================================!

    use w90_constants, only: dp, cmplx_i, cmplx_0
    use w90_parameters, only: gyrotropic_num_bands, gyrotropic_band_list

    ! Arguments
    !
    integer, intent(in) ::num_occ, num_unocc
    integer, dimension(:), intent(in) ::occ_list, unocc_list
    real(kind=dp), dimension(:), intent(in) ::eig     !  n
    real(kind=dp), dimension(:, :), intent(in) ::del_eig !  n
    complex(kind=dp), dimension(:, :, :), intent(in) ::AA      !  n,l,a
    complex(kind=dp), dimension(:, :, :, :), intent(out)::Bnl     !   n,l,a,c
    integer n, m, l, a, c, n1, m1, l1

    Bnl(:, :, :, :) = cmplx_0

    do a = 1, 3
      do c = 1, 3
        do n1 = 1, num_occ
          n = occ_list(n1)
          do l1 = 1, num_unocc
            l = unocc_list(l1)
            Bnl(n1, l1, a, c) = -cmplx_i*(del_eig(n, a) + del_eig(l, a))*AA(n, l, c)
            do m1 = 1, gyrotropic_num_bands
              m = gyrotropic_band_list(m1)
              Bnl(n1, l1, a, c) = Bnl(n1, l1, a, c) + &
                                  (eig(n) - eig(m))*AA(n, m, a)*AA(m, l, c) - &
                                  (eig(l) - eig(m))*AA(n, m, c)*AA(m, l, a)
            enddo ! m1
          enddo !l1
        enddo !n1
      enddo !c
    enddo !a

  end subroutine gyrotropic_get_NOA_Bnl_orb