wham_get_occ_mat_list Subroutine

public subroutine wham_get_occ_mat_list(UU, f_list, g_list, eig, occ)

Uses

  • proc~~wham_get_occ_mat_list~~UsesGraph proc~wham_get_occ_mat_list wham_get_occ_mat_list module~w90_postw90_common w90_postw90_common proc~wham_get_occ_mat_list->module~w90_postw90_common module~w90_constants w90_constants proc~wham_get_occ_mat_list->module~w90_constants module~w90_parameters w90_parameters proc~wham_get_occ_mat_list->module~w90_parameters module~w90_io w90_io proc~wham_get_occ_mat_list->module~w90_io module~w90_postw90_common->module~w90_constants module~w90_comms w90_comms module~w90_postw90_common->module~w90_comms module~w90_parameters->module~w90_constants module~w90_parameters->module~w90_io module~w90_io->module~w90_constants module~w90_comms->module~w90_constants module~w90_comms->module~w90_io

Occupation matrix f, and g=1-f for a list of Fermi energies

Arguments

Type IntentOptional AttributesName
complex(kind=dp), intent(in), dimension(:, :):: UU
complex(kind=dp), intent(out), dimension(:, :, :):: f_list
complex(kind=dp), intent(out), dimension(:, :, :):: g_list
real(kind=dp), intent(in), optional dimension(:):: eig
real(kind=dp), intent(in), optional dimension(:):: occ

Calls

proc~~wham_get_occ_mat_list~~CallsGraph proc~wham_get_occ_mat_list wham_get_occ_mat_list proc~io_error io_error proc~wham_get_occ_mat_list->proc~io_error

Called by

proc~~wham_get_occ_mat_list~~CalledByGraph proc~wham_get_occ_mat_list wham_get_occ_mat_list proc~berry_get_imfgh_klist berry_get_imfgh_klist proc~berry_get_imfgh_klist->proc~wham_get_occ_mat_list proc~berry_get_imf_klist berry_get_imf_klist proc~berry_get_imf_klist->proc~berry_get_imfgh_klist proc~k_slice k_slice proc~k_slice->proc~berry_get_imfgh_klist proc~k_slice->proc~berry_get_imf_klist proc~k_path k_path proc~k_path->proc~berry_get_imfgh_klist proc~k_path->proc~berry_get_imf_klist proc~gyrotropic_get_k_list gyrotropic_get_k_list proc~gyrotropic_get_k_list->proc~berry_get_imfgh_klist proc~gyrotropic_get_k_list->proc~berry_get_imf_klist proc~berry_main berry_main proc~berry_main->proc~berry_get_imfgh_klist proc~berry_main->proc~berry_get_imf_klist proc~gyrotropic_main gyrotropic_main proc~gyrotropic_main->proc~gyrotropic_get_k_list

Contents

Source Code


Source Code

  subroutine wham_get_occ_mat_list(UU, f_list, g_list, eig, occ)
!  subroutine wham_get_occ_mat_list(eig,UU,f_list,g_list)
    !================================!
    !                                !
    !! Occupation matrix f, and g=1-f
    !! for a list of Fermi energies
    ! Tsirkin: !now optionally either eig or occ parameters may be supplied  !
    !    (Changed consistently the calls from the Berry module)        !
    !================================!

    use w90_constants, only: dp, cmplx_0, cmplx_1
    use w90_parameters, only: num_wann, nfermi, fermi_energy_list
    use w90_postw90_common, only: pw90common_get_occ
    use w90_io, only: io_error

    ! Arguments
    !
    complex(kind=dp), dimension(:, :), intent(in)  :: UU
    complex(kind=dp), dimension(:, :, :), intent(out) :: f_list
    complex(kind=dp), dimension(:, :, :), intent(out) :: g_list
    real(kind=dp), intent(in), optional, dimension(:) :: eig
    real(kind=dp), intent(in), optional, dimension(:) :: occ

    integer       :: n, m, i, if, nfermi_loc
    real(kind=dp), allocatable :: occ_list(:, :)

    if (present(occ)) then
      nfermi_loc = 1
    else
      nfermi_loc = nfermi
    endif
    allocate (occ_list(num_wann, nfermi_loc))

    if (present(occ) .and. present(eig)) then
      call io_error( &
        'occ_list and eig cannot be both arguments in get_occ_mat_list')
    elseif (.not. present(occ) .and. .not. present(eig)) then
      call io_error( &
        'either occ_list or eig must be passed as arguments to get_occ_mat_list')
    endif

    if (present(occ)) then
      occ_list(:, 1) = occ(:)
    else
      do if = 1, nfermi_loc
        call pw90common_get_occ(eig, occ_list(:, if), fermi_energy_list(if))
      enddo
    endif

    f_list = cmplx_0
    do if = 1, nfermi_loc
      do n = 1, num_wann
        do m = 1, num_wann
          do i = 1, num_wann
            f_list(n, m, if) = f_list(n, m, if) &
                               + UU(n, i)*occ_list(i, if)*conjg(UU(m, i))
          enddo
          g_list(n, m, if) = -f_list(n, m, if)
          if (m == n) g_list(n, n, if) = g_list(n, n, if) + cmplx_1
        enddo
      enddo
    enddo

  end subroutine wham_get_occ_mat_list