get_gauge_overlap_matrix Subroutine

public subroutine get_gauge_overlap_matrix(ik_a, ns_a, ik_b, ns_b, S_o, S, H)

Uses

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

Arguments

Type IntentOptional AttributesName
integer, intent(in) :: ik_a
integer, intent(in) :: ns_a
integer, intent(in) :: ik_b
integer, intent(in) :: ns_b
complex(kind=dp), intent(in), dimension(:, :):: S_o
complex(kind=dp), intent(out), optional dimension(:, :):: S
complex(kind=dp), intent(out), optional dimension(:, :):: H

Calls

proc~~get_gauge_overlap_matrix~~CallsGraph proc~get_gauge_overlap_matrix get_gauge_overlap_matrix proc~get_win_min get_win_min proc~get_gauge_overlap_matrix->proc~get_win_min

Contents


Source Code

  subroutine get_gauge_overlap_matrix(ik_a, ns_a, ik_b, ns_b, S_o, S, H)
    !==========================================================
    !
    ! Wannier-gauge overlap matrix S in the projected subspace
    !
    ! TODO: Update this documentation of this routine and
    ! possibliy give it a better name. The routine has been
    ! generalized multiple times.
    !
    !==========================================================

    use w90_constants, only: dp, cmplx_0
    use w90_postw90_common, only: v_matrix
    use w90_parameters, only: num_wann, eigval
    use w90_utility, only: utility_zgemmm

    integer, intent(in) :: ik_a, ns_a, ik_b, ns_b

    complex(kind=dp), dimension(:, :), intent(in)            :: S_o
    complex(kind=dp), dimension(:, :), intent(out), optional :: S, H

    integer :: wm_a, wm_b

    call get_win_min(ik_a, wm_a)
    call get_win_min(ik_b, wm_b)

    call utility_zgemmm(v_matrix(1:ns_a, 1:num_wann, ik_a), 'C', &
                        S_o(wm_a:wm_a + ns_a - 1, wm_b:wm_b + ns_b - 1), 'N', &
                        v_matrix(1:ns_b, 1:num_wann, ik_b), 'N', &
                        S, eigval(wm_a:wm_a + ns_a - 1, ik_a), H)
  end subroutine get_gauge_overlap_matrix