Skip to content

Commit

Permalink
extracted runtime critical loop in get_AA_R to new subroutine; see po…
Browse files Browse the repository at this point in the history
  • Loading branch information
Thomas Ponweiser committed Feb 24, 2017
1 parent 569dd8a commit 57098a8
Showing 1 changed file with 44 additions and 18 deletions.
62 changes: 44 additions & 18 deletions src/postw90/get_oper.F90
Original file line number Diff line number Diff line change
Expand Up @@ -249,7 +249,7 @@ subroutine get_AA_R
use w90_parameters, only : num_kpts,nntot,num_wann,wb,bk,timing_level,&
num_bands,ndimwin,nnlist,have_disentangled,&
transl_inv,nncell,effective_model
use w90_postw90_common, only : nrpts,v_matrix
use w90_postw90_common, only : nrpts
use w90_io, only : stdout,io_file_unit,io_error,io_stopwatch,&
seedname
use w90_comms, only : on_root,comms_bcast
Expand All @@ -258,7 +258,7 @@ subroutine get_AA_R
complex(kind=dp), allocatable :: AA_q_diag(:,:)
complex(kind=dp), allocatable :: S_o(:,:)
complex(kind=dp), allocatable :: S(:,:)
integer :: n,m,i,ii,j,jj,winmin_q,winmin_qb,&
integer :: n,m,i,j,&
ik,ik2,ik_prev,nn,inn,nnl,nnm,nnn,&
idir,ncount,nn_count,mmn_in,&
nb_tmp,nkp_tmp,nntot_tmp,file_unit,&
Expand Down Expand Up @@ -419,22 +419,7 @@ subroutine get_AA_R

! Wannier-gauge overlap matrix S in the projected subspace
!
call get_win_min(ik,winmin_q)
call get_win_min(nnlist(ik,nn),winmin_qb)
S=cmplx_0
do m=1,num_wann
do n=1,num_wann
do i=1,num_states(ik)
ii=winmin_q+i-1
do j=1,num_states(nnlist(ik,nn))
jj=winmin_qb+j-1
S(n,m)=S(n,m)&
+conjg(v_matrix(i,n,ik))*S_o(ii,jj)&
*v_matrix(j,m,nnlist(ik,nn))
enddo
enddo
enddo
enddo
call get_gauge_overlap_matrix(ik,nn,num_states,S_o,S)

! Berry connection matrix
! Assuming all neighbors of a given point are read in sequence!
Expand Down Expand Up @@ -1215,4 +1200,45 @@ subroutine get_win_min(ik,win_min)

end subroutine get_win_min

!==========================================================!
subroutine get_gauge_overlap_matrix(ik,nn,num_states,S_o,S)
!==========================================================!
! !
! Wannier-gauge overlap matrix S in the projected subspace !
! !
!==========================================================!

use w90_constants, only : dp, cmplx_0
use w90_postw90_common, only : v_matrix
use w90_parameters, only : nnlist, num_wann

integer, intent(in) :: ik
integer, intent(in) :: nn

integer, dimension(:), intent(in) :: num_states
complex(kind=dp), dimension(:,:), intent(in) :: S_o
complex(kind=dp), dimension(:,:), intent(out) :: S

integer :: winmin_q, winmin_qb, &
m, n, i, ii, j, jj


call get_win_min(ik,winmin_q)
call get_win_min(nnlist(ik,nn),winmin_qb)
S=cmplx_0
do m=1,num_wann
do n=1,num_wann
do i=1,num_states(ik)
ii=winmin_q+i-1
do j=1,num_states(nnlist(ik,nn))
jj=winmin_qb+j-1
S(n,m)=S(n,m)&
+conjg(v_matrix(i,n,ik))*S_o(ii,jj)&
*v_matrix(j,m,nnlist(ik,nn))
end do
end do
end do
end do
end subroutine get_gauge_overlap_matrix

end module w90_get_oper

0 comments on commit 57098a8

Please sign in to comment.