TABLE OF CONTENTS


ABINIT/complete_gamma [ Functions ]

[ Top ] [ Functions ]

NAME

 complete_gamma

FUNCTION

 Use the set of special q points calculated by the Monkhorst & Pack Technique.
 Check if all the information for the q points are present in the input gamma matrices.
 Generate the gamma matrices (already summed over the FS) of the set of q points which
 samples homogeneously the entire Brillouin zone.

INPUTS

 qpttoqpt = qpoint index mapping under symops

OUTPUT

 gamma_qpt = in/out: set of gamma matrix elements completed and symmetrized
    gamma_qpt(2,nbranch**2,nsppol,nqpt_full)

SOURCE

642 subroutine complete_gamma(Cryst,nbranch,nsppol,nqptirred,nqpt_full,ep_scalprod,qirredtofull,qpttoqpt,gamma_qpt)
643 
644 !Arguments ------------------------------------
645 !scalars
646  integer,intent(in) :: nsppol,nbranch,nqptirred,nqpt_full,ep_scalprod
647  type(crystal_t),intent(in) :: Cryst
648 !arrays
649  integer,intent(in) :: qirredtofull(nqptirred)
650  integer,intent(in) :: qpttoqpt(2,Cryst%nsym,nqpt_full)
651  real(dp), intent(inout) :: gamma_qpt(2,nbranch**2,nsppol,nqpt_full)
652 
653 !Local variables-------------------------------
654 !scalars
655  integer :: ibranch,ieqqpt,ii,natom,nsym,iqpt,isppol,isym
656  integer :: itim,jbranch,jj,kk,ll,neqqpt,iatom,ancestor_iatom,iqpt_fullbz
657 !arrays
658  integer :: symmetrized_qpt(nqpt_full)
659  integer :: gkk_flag(nbranch,nbranch,nsppol,nqpt_full)
660  real(dp) :: ss(3,3)
661  real(dp) :: tmp_mat(2,nbranch,nbranch)
662  real(dp) :: tmp_mat2(2,nbranch,nbranch)
663  real(dp) :: ss_allatoms(2,nbranch,nbranch)
664  complex(dpc) :: c_one, c_zero
665  real(dp),allocatable :: gkk_qpt_new(:,:,:),gkk_qpt_tmp(:,:,:)
666 
667 ! *********************************************************************
668 
669  c_one = dcmplx(one,zero)
670  c_zero = dcmplx(zero,zero)
671 
672  natom = Cryst%natom
673  nsym  = Cryst%nsym
674 
675 !Generation of the gkk matrices relative to the q points
676 !of the set which samples the entire Brillouin zone
677 
678 !set up flags for gamma_qpt matrices we have
679  gkk_flag = -1
680  do iqpt=1,nqptirred
681    iqpt_fullbz = qirredtofull(iqpt)
682    gkk_flag(:,:,:,iqpt_fullbz) = 1
683  end do
684 
685  symmetrized_qpt(:) = -1
686 
687  ABI_MALLOC(gkk_qpt_new,(2,nbranch**2,nsppol))
688  ABI_MALLOC(gkk_qpt_tmp,(2,nbranch**2,nsppol))
689 
690  do iqpt=1,nqpt_full
691 !
692 !  Already symmetrized?
693    if (symmetrized_qpt(iqpt) == 1) cycle
694 
695    gkk_qpt_new(:,:,:) = zero
696 
697 !  loop over qpoints equivalent to iqpt
698    neqqpt=0
699 !  do not use time reversal symmetry to complete the qpoints:
700 !  do not know what happens to the gamma matrices
701 !  11/2011: MJV: time reversal is needed here if inversion is absent
702 !  - used in read_gkk and all reductions of q-points by symmetry.
703 
704    do itim=1,2
705      do isym=1,nsym
706 !      ieqqpt is sent onto iqpt by itim/isym
707        ieqqpt = qpttoqpt(itim,isym,iqpt)
708 
709 
710        if (gkk_flag(1,1,1,ieqqpt) == -1) cycle
711 !      if we have information on this qpt
712 !      iqpt is equivalent to ieqqpt: get it from file or memory
713        gkk_qpt_tmp(:,:,:) = gamma_qpt(:,:,:,ieqqpt)
714 
715        neqqpt=neqqpt+1
716 
717 !
718 !      MJV note 02/2010:
719 !      the correspondence of symrel and symrec in the different cases, symmetrizing there
720 !      and back, has been fixed in the cases with and without scalprod (ie cartesian
721 !      and reduced real space coordinates) with respect to a calculation with no symmetries
722 !      I believe everything is settled, but still do not know why the 2 versions of the ss
723 !      matrices here use different rel/rec, instead of just being multiplied by the rprim gprim...
724 !
725        if (ep_scalprod==1) then
726          do ii=1,3
727            do jj=1,3
728              ss(ii,jj)=zero
729              do kk=1,3
730                do ll=1,3
731                  ss(ii,jj)=ss(ii,jj)+Cryst%rprimd(ii,kk)*Cryst%symrel(kk,ll,isym)*Cryst%gprimd(ll,jj)
732                end do
733              end do
734            end do
735          end do
736        else
737          do ii=1,3
738            do jj=1,3
739              ss(ii,jj) = Cryst%symrec(ii,jj,isym)
740            end do
741          end do
742        end if
743 
744        ss_allatoms(:,:,:) = zero
745        do iatom=1,natom
746          ancestor_iatom = Cryst%indsym(4,isym,iatom)
747          ss_allatoms(1, (ancestor_iatom-1)*3+1:(ancestor_iatom-1)*3+3, (iatom-1)*3+1:(iatom-1)*3+3) = ss(1:3,1:3)
748        end do
749 
750 
751 !      NOTE   ssinv(ii,jj)=ssinv(ii,jj)+Cryst%gprimd(ii,kk)*rprimd(jj,ll)*Cryst%symrec(ll,kk,isym)
752 
753        do isppol=1,nsppol
754 !        multiply by the ss matrices
755          tmp_mat2(:,:,:) = zero
756          tmp_mat(:,:,:) = reshape(gkk_qpt_tmp(:,:,isppol),(/2,nbranch,nbranch/))
757 
758          call ZGEMM ('N','N',nbranch,nbranch,nbranch,&
759 &         c_one,ss_allatoms,nbranch,tmp_mat,nbranch,c_zero,tmp_mat2,nbranch)
760 
761          call ZGEMM ('N','T',nbranch,nbranch,nbranch,&
762 &         c_one,tmp_mat2,nbranch,ss_allatoms,nbranch,c_zero,tmp_mat,nbranch)
763 
764 !        add to gkk_qpt_new
765          do ibranch =1,nbranch
766            do jbranch =1,nbranch
767              gkk_qpt_new(:,(jbranch-1)*nbranch+ibranch,isppol) = &
768 &             gkk_qpt_new(:,(jbranch-1)*nbranch+ibranch,isppol) + tmp_mat(:,jbranch,ibranch)
769            end do
770          end do
771        end do ! isppol
772 !
773      end do ! isym
774    end do ! itim
775 !
776    ABI_CHECK(neqqpt>0,'no q-points found equivalent to iqpt ')
777 !  Divide by number of equivalent qpts found.
778    gkk_qpt_new(:,:,:) = gkk_qpt_new(:,:,:)/neqqpt
779 
780 !  copy the symmetrized version into all the equivalent qpoints, appropriately transformed
781    do itim=1,2
782      do isym=1,nsym
783 !      ieqqpt is sent onto iqpt by itim/isym
784        ieqqpt = qpttoqpt(itim,isym,iqpt)
785 
786        if (symmetrized_qpt(ieqqpt) /= -1) cycle
787        gkk_qpt_tmp(:,:,:) = zero
788 
789 !      use symrec matrices to get inverse transform from isym^{-1}
790        if (ep_scalprod==1) then
791          do ii=1,3
792            do jj=1,3
793              ss(ii,jj)=zero
794              do kk=1,3
795                do ll=1,3
796 !                Use inverse of symop matrix here to get back to ieqqpt (inv+transpose is in symrec and in gprimd)
797                  ss(ii,jj)=ss(ii,jj)+Cryst%rprimd(ii,kk)*Cryst%symrec(ll,kk,isym)*Cryst%gprimd(ll,jj)
798                end do
799              end do
800            end do
801          end do
802        else
803          do ii=1,3
804            do jj=1,3
805              ss(ii,jj) = Cryst%symrel(jj,ii,isym)
806            end do
807          end do
808        end if
809 
810        ss_allatoms(:,:,:) = zero
811        do iatom=1,natom
812          ancestor_iatom = Cryst%indsym(4,isym,iatom)
813          ss_allatoms(1, (ancestor_iatom-1)*3+1:(ancestor_iatom-1)*3+3, (iatom-1)*3+1:(iatom-1)*3+3) = ss(1:3,1:3)
814        end do
815 
816 !      ! Use inverse of symop matrix here to get back to ieqqpt
817 !      ssinv(ii,jj)=ssinv(ii,jj)+gprimd(ii,kk)*rprimd(jj,ll)*Cryst%symrel(kk,ll,isym)
818 
819        do isppol=1,nsppol
820 !        multiply by the ss^{-1} matrices
821          tmp_mat2(:,:,:) = zero
822          tmp_mat(:,:,:) = reshape(gkk_qpt_new(:,:,isppol),(/2,nbranch,nbranch/))
823 
824          call ZGEMM ('N','N',nbranch,nbranch,nbranch,&
825 &         c_one,ss_allatoms,nbranch,tmp_mat,nbranch,c_zero,tmp_mat2,nbranch)
826 
827          call ZGEMM ('N','T',nbranch,nbranch,nbranch,&
828 &         c_one,tmp_mat2,nbranch,ss_allatoms,nbranch,c_zero,tmp_mat,nbranch)
829 
830 !        FIXME: the following could just be a reshape
831          do ibranch =1,nbranch
832            do jbranch =1,nbranch
833              gkk_qpt_tmp(:,(jbranch-1)*nbranch+ibranch,isppol) =&
834 &             tmp_mat(:,jbranch,ibranch)
835            end do
836          end do
837          if (gkk_flag (1,1,isppol,ieqqpt) == -1) gkk_flag (:,:,isppol,ieqqpt) = 0
838        end do ! end isppol do
839 
840 !      save symmetrized matrices for qpt ieqqpt
841        gamma_qpt(:,:,:,ieqqpt) = gkk_qpt_tmp(:,:,:)
842 
843        symmetrized_qpt(ieqqpt) = 1
844 
845      end do !isym
846    end do !itim
847  end do !iqpt
848 
849  ABI_FREE(gkk_qpt_new)
850  ABI_FREE(gkk_qpt_tmp)
851 
852 end subroutine complete_gamma

ABINIT/complete_gamma_tr [ Functions ]

[ Top ] [ Functions ]

NAME

 complete_gamma_tr

FUNCTION

 Use the set of special q points calculated by the Monkhorst & Pack Technique.
 Check if all the information for the q points are present in
 the input gamma transport matrices.
 Generate the gamma transport matrices (already summed over the FS) of the set of q points which
 samples homogeneously the entire Brillouin zone.

INPUTS

 crystal<crystal_t>=data type gathering info on the crystalline structure.
 ep_scalprod= flag for scalar product of gkk with phonon displacement vectors
 nbranch=number of phonon branches = 3*natom
 nqptirred=nqpt irred BZ
 nqpt_full=nqpt full BZ
 nsppol=number of spins
 qirredtofull= mapping irred to full qpoints
 qpttoqpt = qpoint index mapping under symops

OUTPUT

 gamma_qpt_tr = in/out: set of gamma matrix elements completed and symmetrized
    gamma_qpt_tr(2,9,nbranch*nbranch,nsppol,nqpt_full)

SOURCE

 883 subroutine complete_gamma_tr(crystal,ep_scalprod,nbranch,nqptirred,nqpt_full,nsppol,gamma_qpt_tr,qirredtofull,qpttoqpt)
 884 
 885  use m_linalg_interfaces
 886 
 887 !Arguments ------------------------------------
 888 !scalars
 889  integer, intent(in) :: nbranch,nqptirred,nqpt_full,nsppol, ep_scalprod
 890  type(crystal_t),intent(in) :: crystal
 891 !arrays
 892  integer,intent(in) :: qpttoqpt(2,crystal%nsym,nqpt_full)
 893  integer,intent(in) :: qirredtofull(nqptirred)
 894  real(dp), intent(inout) :: gamma_qpt_tr(2,9,nbranch*nbranch,nsppol,nqpt_full)
 895 
 896 !Local variables-------------------------------
 897 !scalars
 898  integer :: ieqqpt,ii,iqpt,isppol,isym
 899  integer :: itim,jj,kk,ll,neqqpt
 900  integer :: iatom,ancestor_iatom
 901  integer :: iqpt_fullbz,imode, itensor
 902  real(dp),parameter :: tol=2.d-8
 903 !arrays
 904  integer :: symrel(3,3,crystal%nsym),symrec(3,3,crystal%nsym)
 905  integer :: symmetrized_qpt(nqpt_full)
 906  integer :: gkk_flag(nbranch,nbranch,nsppol,nqpt_full)
 907  real(dp) :: gprimd(3,3),rprimd(3,3)
 908  real(dp) :: ss(3,3), sscart(3,3)
 909  real(dp) :: tmp_mat(2,nbranch,nbranch)
 910  real(dp) :: tmp_mat2(2,nbranch,nbranch)
 911  real(dp) :: tmp_tensor(2,3,3)
 912  real(dp) :: tmp_tensor2(2,3,3)
 913  real(dp) :: ss_allatoms(nbranch,nbranch)
 914  real(dp),allocatable :: gkk_qpt_new(:,:,:,:),gkk_qpt_tmp(:,:,:,:)
 915 
 916 ! *********************************************************************
 917 
 918  gprimd = crystal%gprimd
 919  rprimd = crystal%rprimd
 920 
 921  symrec =  crystal%symrec
 922  symrel =  crystal%symrel
 923 
 924 !Generation of the gkk matrices relative to the q points
 925 !of the set which samples the entire Brillouin zone
 926 
 927 !set up flags for gamma_qpt matrices we have
 928  gkk_flag = -1
 929  do iqpt=1,nqptirred
 930    iqpt_fullbz = qirredtofull(iqpt)
 931    gkk_flag(:,:,:,iqpt_fullbz) = 1
 932  end do
 933 
 934  symmetrized_qpt(:) = -1
 935 ! isppol=1
 936 
 937  ABI_MALLOC(gkk_qpt_new,(2,9,nbranch*nbranch, nsppol))
 938  ABI_MALLOC(gkk_qpt_tmp,(2,9,nbranch*nbranch, nsppol))
 939 
 940  do iqpt=1,nqpt_full
 941 
 942 !  Already symmetrized?
 943    if (symmetrized_qpt(iqpt) == 1) cycle
 944 
 945    gkk_qpt_new(:,:,:,:) = zero
 946 
 947 !  loop over qpoints equivalent to iqpt
 948    neqqpt=0
 949 !  do not use time reversal symmetry to complete the qpoints:
 950 !  do not know what happens to the gamma matrices
 951 
 952    do itim=1,2
 953      do isym=1,crystal%nsym
 954 !      ieqqpt is sent onto iqpt by itim/isym
 955        ieqqpt = qpttoqpt(itim,isym,iqpt)
 956 
 957        if (gkk_flag(1,1,1,ieqqpt) == -1) cycle
 958 !      if we have information on this qpt
 959 !      iqpt is equivalent to ieqqpt: get it from file or memory
 960        gkk_qpt_tmp(:,:,:,:) = gamma_qpt_tr(:,:,:,:,ieqqpt)
 961 
 962        neqqpt=neqqpt+1
 963 
 964 !
 965 !      MJV note 02/2010:
 966 !      the correspondence of symrel and symrec in the different cases, symmetrizing there
 967 !      and back, has been fixed in the cases with and without scalprod (ie cartesian
 968 !      and reduced real space coordinates) with respect to a calculation with no symmetries
 969 !      I believe everything is settled, but still do not know why the 2 versions of the ss
 970 !      matrices here use different rel/rec, instead of just being multiplied by the rprim gprim...
 971 !
 972        do ii=1,3
 973          do jj=1,3
 974            sscart(ii,jj)=0.0_dp
 975            do kk=1,3
 976              do ll=1,3
 977                sscart(ii,jj)=sscart(ii,jj)+rprimd(ii,kk)*symrel(kk,ll,isym)*gprimd(ll,jj)
 978 !              sscart(ii,jj)=sscart(ii,jj)+rprimd(ii,kk)*symrel(kk,ll,isym)*gprimd(ll,jj)
 979              end do
 980            end do
 981          end do
 982        end do
 983        if (ep_scalprod==1) then
 984          do ii=1,3
 985            do jj=1,3
 986              ss(ii,jj)=0.0_dp
 987              do kk=1,3
 988                do ll=1,3
 989                  ss(ii,jj)=ss(ii,jj)+rprimd(ii,kk)*symrel(kk,ll,isym)*gprimd(ll,jj)
 990                end do
 991              end do
 992            end do
 993          end do
 994        else
 995          do ii=1,3
 996            do jj=1,3
 997              ss(ii,jj) = symrec(ii,jj,isym)
 998            end do
 999          end do
1000        end if
1001 
1002        ss_allatoms(:,:) = zero
1003        do iatom=1,crystal%natom
1004          ancestor_iatom = crystal%indsym(4,isym,iatom)
1005          ss_allatoms((ancestor_iatom-1)*3+1:(ancestor_iatom-1)*3+3,&
1006 &         (iatom-1)*3+1:         (iatom-1)*3+3) = ss(1:3,1:3)
1007        end do
1008 
1009 
1010 !      NOTE   ssinv(ii,jj)=ssinv(ii,jj)+gprimd(ii,kk)*rprimd(jj,ll)*symrec(ll,kk,isym)
1011 
1012        do isppol=1,nsppol
1013 
1014 !        for each tensor component, rotate the cartesian directions of phonon modes
1015          do itensor = 1, 9
1016 !          multiply by the ss matrices
1017            tmp_mat2(:,:,:) = zero
1018            tmp_mat(:,:,:) = reshape(gkk_qpt_tmp(:,itensor,:,isppol),&
1019 &           (/2,nbranch,nbranch/))
1020            call DGEMM ('N','N',nbranch,nbranch,nbranch,&
1021 &           one,ss_allatoms,nbranch,tmp_mat(1,:,:),nbranch,zero,&
1022 &           tmp_mat2(1,:,:),nbranch)
1023            call DGEMM ('N','N',nbranch,nbranch,nbranch,&
1024 &           one,ss_allatoms,nbranch,tmp_mat(2,:,:),nbranch,zero,&
1025 &           tmp_mat2(2,:,:),nbranch)
1026 
1027            call DGEMM ('N','T',nbranch,nbranch,nbranch,&
1028 &           one,tmp_mat2(1,:,:),nbranch,ss_allatoms,nbranch,zero,&
1029 &           tmp_mat(1,:,:),nbranch)
1030            call DGEMM ('N','T',nbranch,nbranch,nbranch,&
1031 &           one,tmp_mat2(2,:,:),nbranch,ss_allatoms,nbranch,zero,&
1032 &           tmp_mat(2,:,:),nbranch)
1033 
1034            gkk_qpt_tmp(:,itensor,:,isppol) = reshape (tmp_mat, (/2,nbranch*nbranch/))
1035          end do ! itensor
1036 
1037 !        for each cartesian direction/phonon mode, rotate the tensor components
1038          do imode = 1, nbranch*nbranch
1039            tmp_tensor2(:,:,:) = zero
1040            tmp_tensor(:,:,:) = reshape(gkk_qpt_tmp(:,:,imode,isppol),&
1041 &           (/2,3,3/))
1042            call DGEMM ('N','N',3,3,3,&
1043 &           one,sscart,3,tmp_tensor(1,:,:),3,zero,&
1044 &           tmp_tensor2(1,:,:),3)
1045            call DGEMM ('N','T',3,3,3,&
1046 &           one,tmp_tensor2(1,:,:),3,sscart,3,zero,&
1047 &           tmp_tensor(1,:,:),3)
1048 
1049            call DGEMM ('N','N',3,3,3,&
1050 &           one,sscart,3,tmp_tensor(2,:,:),3,zero,&
1051 &           tmp_tensor2(2,:,:),3)
1052            call DGEMM ('N','T',3,3,3,&
1053 &           one,tmp_tensor2(2,:,:),3,sscart,3,zero,&
1054 &           tmp_tensor(2,:,:),3)
1055 
1056            gkk_qpt_tmp(:,:,imode,isppol) = reshape (tmp_tensor, (/2,9/)) ! modified by BX
1057          end do ! imode
1058 
1059 !        add to gkk_qpt_new
1060          gkk_qpt_new(:,:,:,isppol) = gkk_qpt_new(:,:,:,isppol) + gkk_qpt_tmp(:,:,:,isppol)
1061 
1062        end do ! end isppol do
1063 
1064      end do ! end isym do
1065    end do ! end itim do
1066 
1067    ABI_CHECK(neqqpt>0,'no q-points found equivalent to iqpt ')
1068 
1069 !  divide by number of equivalent qpts found
1070    gkk_qpt_new = gkk_qpt_new/neqqpt
1071 
1072 
1073 !  copy the symmetrized version into all the equivalent qpoints, appropriately transformed
1074    do itim=1,2
1075      do isym=1,crystal%nsym
1076 !      ieqqpt is sent onto iqpt by itim/isym
1077        ieqqpt = qpttoqpt(itim,isym,iqpt)
1078 
1079        if (symmetrized_qpt(ieqqpt) /= -1) cycle
1080        gkk_qpt_tmp = zero
1081 
1082 !      use symrec matrices to get inverse transform from isym^{-1}
1083        do ii=1,3
1084          do jj=1,3
1085            sscart(ii,jj)=0.0_dp
1086            do kk=1,3
1087              do ll=1,3
1088 !              Use inverse of symop matrix here to get back to ieqqpt (inv+transpose is in symrec and in gprimd)
1089 !              sscart(ii,jj)=sscart(ii,jj)+rprimd(ii,kk)*symrec(ll,kk,isym)*gprimd(ll,jj)
1090                sscart(ii,jj)=sscart(ii,jj)+rprimd(ii,kk)*symrec(ll,kk,isym)*gprimd(ll,jj)
1091              end do
1092            end do
1093          end do
1094        end do
1095        if (ep_scalprod==1) then
1096          do ii=1,3
1097            do jj=1,3
1098              ss(ii,jj)=0.0_dp
1099              do kk=1,3
1100                do ll=1,3
1101 !                Use inverse of symop matrix here to get back to ieqqpt (inv+transpose is in symrec and in gprimd)
1102                  ss(ii,jj)=ss(ii,jj)+rprimd(ii,kk)*symrec(ll,kk,isym)*gprimd(ll,jj)
1103                end do
1104              end do
1105            end do
1106          end do
1107        else
1108          do ii=1,3
1109            do jj=1,3
1110              ss(ii,jj) = symrel(jj,ii,isym)
1111            end do
1112          end do
1113        end if
1114 
1115        ss_allatoms(:,:) = zero
1116        do iatom=1,crystal%natom
1117          ancestor_iatom = crystal%indsym(4,isym,iatom)
1118          ss_allatoms((ancestor_iatom-1)*3+1:(ancestor_iatom-1)*3+3,&
1119 &         (iatom-1)*3+1:          (iatom-1)*3+3) = ss(1:3,1:3)
1120        end do
1121 
1122 !      ! Use inverse of symop matrix here to get back to ieqqpt
1123 !      ssinv(ii,jj)=ssinv(ii,jj)+gprimd(ii,kk)*rprimd(jj,ll)*symrel(kk,ll,isym)
1124 
1125        do isppol=1,nsppol
1126          do itensor = 1, 9
1127 !          multiply by the ss^{-1} matrices
1128            tmp_mat2(:,:,:) = zero
1129            tmp_mat(:,:,:) = reshape(gkk_qpt_new(:,itensor,:,isppol),&
1130 &           (/2,nbranch,nbranch/))
1131 
1132 
1133            call DGEMM ('N','N',nbranch,nbranch,nbranch,&
1134 &           one,ss_allatoms,nbranch,tmp_mat(1,:,:),nbranch,zero,&
1135 &           tmp_mat2(1,:,:),nbranch)
1136            call DGEMM ('N','N',nbranch,nbranch,nbranch,&
1137 &           one,ss_allatoms,nbranch,tmp_mat(2,:,:),nbranch,zero,&
1138 &           tmp_mat2(2,:,:),nbranch)
1139 
1140            call DGEMM ('N','T',nbranch,nbranch,nbranch,&
1141 &           one,tmp_mat2(1,:,:),nbranch,ss_allatoms,nbranch,zero,&
1142 &           tmp_mat(1,:,:),nbranch)
1143            call DGEMM ('N','T',nbranch,nbranch,nbranch,&
1144 &           one,tmp_mat2(2,:,:),nbranch,ss_allatoms,nbranch,zero,&
1145 &           tmp_mat(2,:,:),nbranch)
1146 
1147 
1148            gkk_qpt_tmp(:,itensor,:,isppol) = reshape (tmp_mat, (/2,nbranch*nbranch/))
1149          end do ! itensor
1150 
1151 !        for each cartesian direction/phonon mode, rotate the tensor components
1152          do imode = 1, nbranch*nbranch
1153            tmp_tensor2(:,:,:) = zero
1154            tmp_tensor(:,:,:) = reshape(gkk_qpt_tmp(:,:,imode,isppol),&
1155 &           (/2,3,3/))
1156            call DGEMM ('N','N',3,3,3,&
1157 &           one,sscart,3,tmp_tensor(1,:,:),3,zero,&
1158 &           tmp_tensor2(1,:,:),3)
1159            call DGEMM ('N','T',3,3,3,&
1160 &           one,tmp_tensor2(1,:,:),3,sscart,3,zero,&
1161 &           tmp_tensor(1,:,:),3)
1162 
1163            call DGEMM ('N','N',3,3,3,&
1164 &           one,sscart,3,tmp_tensor(2,:,:),3,zero,&
1165 &           tmp_tensor2(2,:,:),3)
1166            call DGEMM ('N','T',3,3,3,&
1167 &           one,tmp_tensor2(2,:,:),3,sscart,3,zero,&
1168 &           tmp_tensor(2,:,:),3)
1169 
1170 !          gkk_qpt_new(:,:,imode,isppol) = reshape (tmp_tensor, (/2,9/)) ! Modified by BX
1171            gkk_qpt_tmp(:,:,imode,isppol) = reshape (tmp_tensor, (/2,9/)) ! Modified by BX
1172          end do ! imode
1173 
1174          if (gkk_flag (1,1,isppol,ieqqpt) == -1) then
1175            gkk_flag (:,:,isppol,ieqqpt) = 0
1176          end if
1177 
1178        end do ! end isppol do
1179 
1180 
1181 !      save symmetrized matrices for qpt ieqqpt
1182        gamma_qpt_tr(:,:,:,:,ieqqpt) = gkk_qpt_tmp(:,:,:,:)
1183 
1184        symmetrized_qpt(ieqqpt) = 1
1185 
1186      end do ! end isym do
1187    end do ! end itim do
1188 
1189  end do
1190 !end iqpt do
1191 
1192  ABI_FREE(gkk_qpt_new)
1193  ABI_FREE(gkk_qpt_tmp)
1194 
1195 end subroutine complete_gamma_tr

ABINIT/defs_elphon [ Modules ]

[ Top ] [ Modules ]

NAME

 defs_elphon

FUNCTION

 This module contains the datastructures for elphon
  the different (huge) matrices will either be allocated and
  used, or be written to disk. All combinations should be feasible.

COPYRIGHT

 Copyright (C) 2004-2024 ABINIT group (MVer, MG)
 This file is distributed under the terms of the
 GNU General Public Licence, see ~abinit/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .
 For the initials of contributors, see ~abinit/doc/developers/contributors.txt .

NOTES

  Contains the following datastructures:
   1) elph_type contains data and dimensions for the kpoints near the
      fermi surface and the $g_{k k+q}$ matrix elements

SOURCE

25 #if defined HAVE_CONFIG_H
26 #include "config.h"
27 #endif
28 
29 #include "abi_common.h"
30 
31 module defs_elphon
32 
33  use defs_basis
34  use m_abicore
35  use m_errors
36  use m_xmpi
37  use m_krank
38  use m_crystal
39 
40  implicit none
41 
42  private

defs_elphon/elph_ds_clean [ Functions ]

[ Top ] [ defs_elphon ] [ Functions ]

NAME

   elph_ds_clean

FUNCTION

   deallocate remaining arrays in the elph_ds datastructure

INPUTS

  elph_ds = elphon datastructure

NOTES

SOURCE

313 subroutine elph_ds_clean(elph_ds)
314 
315 !Arguments ------------------------------------
316 !scalars
317  type(elph_type), intent(inout) :: elph_ds
318 
319 ! *************************************************************************
320 
321  !@elph_type
322  ABI_SFREE(elph_ds%qirredtofull)
323  ABI_SFREE(elph_ds%wtq)
324  ABI_SFREE(elph_ds%n0)
325  ABI_SFREE(elph_ds%qpt_full)
326  ABI_SFREE(elph_ds%gkk_intweight)
327  ABI_SFREE(elph_ds%gkk_qpt)
328  ABI_SFREE(elph_ds%gkk_rpt)
329  ABI_SFREE(elph_ds%gkk2)
330  ABI_SFREE(elph_ds%gamma_qpt)
331  ABI_SFREE(elph_ds%gamma_rpt)
332  ABI_SFREE(elph_ds%phfrq)
333  ABI_SFREE(elph_ds%a2f)
334  ABI_SFREE(elph_ds%qgrid_data)
335 
336  call elph_k_destroy (elph_ds%k_phon)
337  call elph_k_destroy (elph_ds%k_fine)
338 
339  call elph_ds%k_fine%krank%free()
340 
341 end subroutine elph_ds_clean

defs_elphon/elph_k_copy [ Functions ]

[ Top ] [ defs_elphon ] [ Functions ]

NAME

   elph_k_copy

FUNCTION

   allocate and copy arrays in the elph_k datastructure

INPUTS

  elph_k = elphon k-points datastructure

NOTES

SOURCE

417 subroutine elph_k_copy(elph_k_in, elph_k_out)
418 
419 !Arguments ------------------------------------
420 !scalars
421  type(elph_kgrid_type), intent(in) :: elph_k_in
422  type(elph_kgrid_type), intent(out) :: elph_k_out
423 
424 ! *************************************************************************
425 
426  !@elph_kgrid_type
427  elph_k_out%nband = elph_k_in%nband
428  elph_k_out%nsppol = elph_k_in%nsppol
429  elph_k_out%nsym = elph_k_in%nsym
430 
431  elph_k_out%nkpt = elph_k_in%nkpt
432  elph_k_out%nkptirr = elph_k_in%nkptirr
433 
434  elph_k_out%my_nkpt = elph_k_in%my_nkpt
435 
436  ABI_MALLOC(elph_k_out%my_kpt,(elph_k_out%nkpt))
437  elph_k_out%my_kpt = elph_k_in%my_kpt
438 
439  ABI_MALLOC(elph_k_out%my_ikpt,(elph_k_out%my_nkpt))
440  elph_k_out%my_ikpt = elph_k_in%my_ikpt
441 
442  ABI_MALLOC(elph_k_out%kptirr,(3,elph_k_out%nkptirr))
443  elph_k_out%kptirr = elph_k_in%kptirr
444  ABI_MALLOC(elph_k_out%wtkirr,(elph_k_out%nkptirr))
445  elph_k_out%wtkirr = elph_k_in%wtkirr
446 
447  ABI_MALLOC(elph_k_out%wtk,(elph_k_out%nband,elph_k_out%nkpt,elph_k_out%nsppol))
448  elph_k_out%wtk = elph_k_in%wtk
449  ABI_MALLOC(elph_k_out%kpt,(3,elph_k_out%nkpt))
450  elph_k_out%kpt = elph_k_in%kpt
451 
452  elph_k_out%krank = elph_k_in%krank%copy()
453 
454  ABI_MALLOC(elph_k_out%irr2full,(elph_k_out%nkptirr))
455  elph_k_out%irr2full = elph_k_in%irr2full
456  ABI_MALLOC(elph_k_out%full2irr,(3,elph_k_out%nkpt))
457  elph_k_out%full2irr = elph_k_in%full2irr
458  ABI_MALLOC(elph_k_out%full2full,(2,elph_k_out%nsym,elph_k_out%nkpt))
459  elph_k_out%full2full = elph_k_in%full2full
460 
461  ABI_MALLOC(elph_k_out%irredtoGS,(elph_k_out%nkptirr))
462  elph_k_out%irredtoGS = elph_k_in%irredtoGS
463 
464 end subroutine elph_k_copy

defs_elphon/elph_k_destroy [ Functions ]

[ Top ] [ defs_elphon ] [ Functions ]

NAME

   elph_k_destroy

FUNCTION

   deallocate arrays in the elph_k datastructure

INPUTS

  elph_k = elphon k-points datastructure

NOTES

SOURCE

484 subroutine elph_k_destroy(elph_k)
485 
486 !Arguments ------------------------------------
487 !scalars
488  type(elph_kgrid_type), intent(inout) :: elph_k
489 
490 ! *************************************************************************
491 
492  !@elph_kgrid_type
493  ABI_SFREE(elph_k%irr2full)
494  ABI_SFREE(elph_k%full2irr)
495  ABI_SFREE(elph_k%full2full)
496  ABI_SFREE(elph_k%irredtoGS)
497  ABI_SFREE(elph_k%new_irredtoGS)
498  ABI_SFREE(elph_k%kpt)
499  ABI_SFREE(elph_k%kptirr)
500  ABI_SFREE(elph_k%new_kptirr)
501  ABI_SFREE(elph_k%my_kpt)
502  ABI_SFREE(elph_k%my_ikpt)
503  ABI_SFREE(elph_k%wtk)
504  ABI_SFREE(elph_k%wtq)
505  ABI_SFREE(elph_k%wtkirr)
506  ABI_SFREE(elph_k%new_wtkirr)
507  ABI_SFREE(elph_k%velocwtk)
508  ABI_SFREE(elph_k%vvelocwtk)
509 
510  call elph_k%krank%free()
511 
512 end subroutine elph_k_destroy

defs_elphon/elph_k_procs [ Functions ]

[ Top ] [ defs_elphon ] [ Functions ]

NAME

   elph_k_procs

FUNCTION

   allocate kpt to processors, in the elph_k datastructure

INPUTS

  nproc = number of k-parallel processors
  elph_k = elphon k-points datastructure

NOTES

SOURCE

533 subroutine elph_k_procs(nproc, elph_k)
534 
535 !Arguments ------------------------------------
536 !scalars
537  integer, intent(in) :: nproc
538  type(elph_kgrid_type), intent(inout) :: elph_k
539 
540  integer :: ikpt, me, ik_this_proc
541 ! *************************************************************************
542 
543  if (allocated(elph_k%my_kpt)) then
544    ABI_FREE (elph_k%my_kpt)
545  end if
546  ABI_MALLOC (elph_k%my_kpt, (elph_k%nkpt))
547 
548  elph_k%my_kpt = 0
549  elph_k%my_nkpt = 0
550  me = xmpi_comm_rank(xmpi_world)
551  do ikpt = 1, elph_k%nkpt
552    elph_k%my_kpt(ikpt) = MOD(ikpt-1, nproc)
553    if (elph_k%my_kpt(ikpt) == me) elph_k%my_nkpt = elph_k%my_nkpt + 1
554  end do
555 
556 ! create inverse mapping from ik_this_proc to ikpt
557  if (allocated(elph_k%my_ikpt)) then
558    ABI_FREE (elph_k%my_ikpt)
559  end if
560  ABI_MALLOC (elph_k%my_ikpt, (elph_k%my_nkpt))
561  elph_k%my_ikpt = 0
562 
563  ik_this_proc = 0
564  do ikpt = 1, elph_k%nkpt
565    if (elph_k%my_kpt(ikpt) == me) then
566      ik_this_proc = ik_this_proc + 1
567      elph_k%my_ikpt(ik_this_proc) = ikpt
568    end if
569  end do
570  ABI_CHECK(ik_this_proc == elph_k%my_nkpt, 'found inconsistent k distribution in processors')
571 
572  write (std_out,*) 'elph_k_procs : nkpt, distrib = ', elph_k%my_nkpt
573  write (std_out,*) elph_k%my_kpt
574  write (std_out,*) elph_k%my_ikpt
575 
576 end subroutine elph_k_procs

defs_elphon/elph_kgrid_type [ Types ]

[ Top ] [ defs_elphon ] [ Types ]

NAME

 elph_kgrid_type

FUNCTION

 elph_kgrid_type contains k-point grid data and dimensions
  this is a sub object of elph_type

SOURCE

61   type,public :: elph_kgrid_type
62 
63    integer :: nband                           ! number of bands for weights
64    integer :: nsppol                          ! number of spin pol for weights
65    integer :: nsym                            ! number of symmetry operations
66    integer :: nkpt                            ! number of k-points in full grid
67    integer :: nkptirr                         ! number of k-points in irreducible grid
68    integer :: new_nkptirr                     ! number of k-points in irreducible grid
69    integer :: my_nkpt                         ! number of k-points on present processor
70 
71    type(krank_t) :: krank            ! ranking of all kpoints on phonon calculation grid, and inverse rank
72 
73    integer, allocatable :: irr2full(:)            ! correspondence of irred kpoints to a full one
74    integer, allocatable :: full2irr(:,:)          ! correspondence of full k to one irred kpoints through sym and timrev
75    integer, allocatable :: full2full(:,:,:)       ! symmetry mapping of kpoints
76    integer, allocatable :: my_kpt(:)              ! flag for k-points belonging to present proc (= me index of proc for each k-point)
77    integer, allocatable :: my_ikpt(:)             ! flag for k-points belonging to present proc (= me index of proc for each k-point)
78    integer, allocatable :: irredtoGS(:)           ! (nkptirr)
79    integer, allocatable :: new_irredtoGS(:)       ! (new_nkptirr)
80 
81    real(dp),allocatable :: kpt(:,:)               ! coordinates of the full kpoints from phonon calculation
82    real(dp),allocatable :: kptirr(:,:)            ! irreducible k-points, for preliminary set up
83    real(dp),allocatable :: new_kptirr(:,:)        ! irreducible k-points, for preliminary set up
84    real(dp),allocatable :: wtk(:,:,:)             ! integration weights (see also gkk_intweight)
85    real(dp),allocatable :: wtq(:,:,:)             ! integration weights (see also gkk_intweight)
86    real(dp),allocatable :: wtkirr(:)              ! weights for irreducible kpoints, to sum over _whole_ BZ (not just Fermi Surface)
87    real(dp),allocatable :: new_wtkirr(:)          ! weights for irreducible kpoints, to sum over _whole_ BZ (not just Fermi Surface)
88    real(dp),allocatable :: velocwtk(:,:,:,:)      ! (nFSband,nkpt_fine,3,nsppol), v(k)*wtk
89    real(dp),allocatable :: vvelocwtk(:,:,:,:,:)   ! (nFSband,nkpt_fine,3,3,nsppol), v(k)*v(k)*wtk
90 
91   end type elph_kgrid_type
92 
93   public :: elph_k_copy
94   public :: elph_k_procs
95   public :: elph_k_destroy

defs_elphon/elph_tr_ds_clean [ Functions ]

[ Top ] [ defs_elphon ] [ Functions ]

NAME

   elph_tr_ds_clean

FUNCTION

   deallocate remaining arrays in the elph_tr_ds datastructure

INPUTS

  elph_tr_ds = elphon transport datastructure

NOTES

SOURCE

361 subroutine elph_tr_ds_clean(elph_tr_ds)
362 
363 !Arguments ------------------------------------
364 !scalars
365  type(elph_tr_type), intent(inout) :: elph_tr_ds
366 
367 ! *************************************************************************
368 
369  !@elph_tr_type
370  ABI_SFREE(elph_tr_ds%el_veloc)
371  ABI_SFREE(elph_tr_ds%FSelecveloc_sq)
372  ABI_SFREE(elph_tr_ds%veloc_sq0)
373  ABI_SFREE(elph_tr_ds%veloc_sq)
374  ABI_SFREE(elph_tr_ds%dos_n0)
375  ABI_SFREE(elph_tr_ds%dos_n)
376  ABI_SFREE(elph_tr_ds%en_all)
377  ABI_SFREE(elph_tr_ds%de_all)
378  ABI_SFREE(elph_tr_ds%gamma_qpt_tr)
379  ABI_SFREE(elph_tr_ds%gamma_qpt_trin)
380  ABI_SFREE(elph_tr_ds%gamma_qpt_trout)
381  ABI_SFREE(elph_tr_ds%gamma_rpt_tr)
382  ABI_SFREE(elph_tr_ds%gamma_rpt_trin)
383  ABI_SFREE(elph_tr_ds%gamma_rpt_trout)
384  ABI_SFREE(elph_tr_ds%a2f_1d_tr)
385  ABI_SFREE(elph_tr_ds%a2f_1d_trin)
386  ABI_SFREE(elph_tr_ds%a2f_1d_trout)
387  ABI_SFREE(elph_tr_ds%tmp_gkk_intweight)
388  ABI_SFREE(elph_tr_ds%tmp_gkk_intweight1)
389  ABI_SFREE(elph_tr_ds%tmp_gkk_intweight2)
390  ABI_SFREE(elph_tr_ds%tmp_velocwtk)
391  ABI_SFREE(elph_tr_ds%tmp_velocwtk1)
392  ABI_SFREE(elph_tr_ds%tmp_velocwtk2)
393  ABI_SFREE(elph_tr_ds%tmp_vvelocwtk)
394  ABI_SFREE(elph_tr_ds%tmp_vvelocwtk1)
395  ABI_SFREE(elph_tr_ds%tmp_vvelocwtk2)
396 
397 end subroutine elph_tr_ds_clean

defs_elphon/elph_tr_type [ Types ]

[ Top ] [ defs_elphon ] [ Types ]

NAME

 elph_tr_type

FUNCTION

 elph_tr_ds contains the necessary data for the transport properties

SOURCE

242   type,public :: elph_tr_type
243 
244      integer :: ifltransport
245      integer :: unitgkq_trin,unitgkq_trout
246      integer :: gkqwrite,gkqexist
247      integer :: onegkksize
248 
249      character(len=fnlen) :: ddkfilename
250 
251      real(dp),allocatable :: dos_n0(:,:)                  ! (nT,nsppol) DOS at the Fermi level (states/Ha/spin) at input temperatures
252      real(dp),allocatable :: dos_n(:,:)                   ! (nE,nsppol) DOS at the selected energies (states/Ha/spin)
253      real(dp),allocatable :: en_all(:,:)                  ! (nE,nsppol) selected energies
254      real(dp),allocatable :: de_all(:,:)                  ! (nE,nsppol) differences between selected energies
255      real(dp),allocatable :: veloc_sq0(:,:,:)             ! (3,nsppol,nT)
256      real(dp),allocatable :: veloc_sq(:,:,:)              ! (3,nsppol,nE)
257 
258      real(dp),allocatable :: el_veloc(:,:,:,:)        ! nkpt nband 3 nsppol
259 ! the 9 = 3x3 is for the full tensorial transport coefficients
260      real(dp),allocatable :: gamma_qpt_tr(:,:,:,:,:)    ! 2 9 branches**2 nsppol qpt
261      real(dp),allocatable :: gamma_qpt_trin(:,:,:,:,:)  !idem
262      real(dp),allocatable :: gamma_qpt_trout(:,:,:,:,:) !idem
263 
264      real(dp),allocatable :: gamma_rpt_tr(:,:,:,:,:,:,:)    !idem
265      real(dp),allocatable :: gamma_rpt_trin(:,:,:,:,:)  !idem
266      real(dp),allocatable :: gamma_rpt_trout(:,:,:,:,:) !idem
267 
268      real(dp),allocatable :: a2f_1d_tr(:,:,:,:,:,:)           ! nfreq 9 nsppol 4 n_pair ntemp
269      real(dp),allocatable :: a2f_1d_trin(:,:,:)
270      real(dp),allocatable :: a2f_1d_trout(:,:,:)
271 
272      real(dp),allocatable :: FSelecveloc_sq(:,:)       ! 3 nsppol
273 
274      real(dp),allocatable :: tmp_gkk_intweight(:,:,:,:)
275      real(dp),allocatable :: tmp_gkk_intweight1(:,:,:)
276      real(dp),allocatable :: tmp_gkk_intweight2(:,:,:)
277 
278      real(dp),allocatable :: tmp_velocwtk(:,:,:,:,:)
279      real(dp),allocatable :: tmp_velocwtk1(:,:,:,:)
280      real(dp),allocatable :: tmp_velocwtk2(:,:,:,:)
281 
282      real(dp),allocatable :: tmp_vvelocwtk(:,:,:,:,:,:)
283      real(dp),allocatable :: tmp_vvelocwtk1(:,:,:,:,:)
284      real(dp),allocatable :: tmp_vvelocwtk2(:,:,:,:,:)
285 
286   end type elph_tr_type
287 
288  public :: elph_tr_ds_clean

defs_elphon/elph_type [ Types ]

[ Top ] [ defs_elphon ] [ Types ]

NAME

 elph_type

FUNCTION

 elph_type contains data and dimensions for the kpoints near the
 fermi surface and the $g_{k k+q}$ matrix elements

SOURCE

110   type,public :: elph_type
111 
112    type(elph_kgrid_type) :: k_phon               ! object for k-grid of phonon calculation
113    type(elph_kgrid_type) :: k_fine               ! object for fine k-grid for FS integration
114 
115    integer :: natom,nbranch,nFSband,nband
116    integer :: minFSband,maxFSband                !Index of lower and upper bands used for the FS integration
117 
118    integer :: ngkkband                           !Number of bands kept in final gkk matrix elements:
119                                                  !either 1 if sum is performed immediately
120                                                  !or = nband if all elements are kept based on flag ep_keepbands
121 
122    integer :: nenergy                            ! number of points in energy
123                                                  !   space for the electron band energies
124                                                  !   energies from
125                                                  !   ef-nenergy*delta_e to
126                                                  !   ef+nenergy*delta_e
127 
128    integer :: n_pair                             ! number of pairs considered
129 
130    integer :: nqpt_full                          !number of q in full BZ
131    integer :: nqptirred                          !number of irred q-points
132 
133 
134    integer :: unita2f,unit_gkk2,unit_gkk_rpt
135    integer :: unitgkq                            !units for file output
136 
137    integer :: gkqwrite
138    integer :: gkk2write
139    integer :: gkk_rptwrite
140 
141    integer :: ep_scalprod                        !flag to perform the scalar product
142    integer :: symgkq                             !flag to symmetrize gkq matrix elements
143    integer :: ep_keepbands                       !flag to sum over bands or not
144    integer :: ep_lova                            ! 1 for lova, and 0 for general
145    integer :: tuniformgrid                       !flag to expect uniform grid of q or not
146    integer :: prtbltztrp                         !flag to output BoltzTraP input files
147 
148    integer :: na2f                               !dimensions and increments for a2F function
149    integer :: nsppol                             ! number of spin polarization channels
150    integer :: nspinor                            ! number of spinorial components
151    integer :: telphint                           ! flag for integration over the FS with 0=tetrahedra 1=gaussians
152    integer :: ep_nspline                         ! scale factor for spline interpolation in RTA
153    integer :: ntemper                            ! number of temperature points
154    integer :: use_k_fine                         ! flag for using fine k-grids for eigenvalues and velocities. 0=no 1=yes
155    integer :: ep_int_gkk                         ! flag for interpolate gkk(1) or gamma (0)
156    integer :: ep_b_min                           ! first band taken into account in FS integration (if telphint==2)
157    integer :: ep_b_max                           ! last band taken into account in FS integration (if telphint==2)
158    integer :: kptrlatt(3,3)                      ! kpoint grid generating vectors, as in abinit
159    integer :: kptrlatt_fine(3,3)                 ! kpoint grid generating vectors, for fine grid used in FS integration
160 
161    real(dp) :: delta_e                           ! step in electronic energies, around Fermi level
162    real(dp) :: omega_min,omega_max
163    real(dp) :: a2fsmear,domega
164    real(dp) :: nelect                            ! number of electrons per unit cell, eventually with extra charges for carriers in semiconductors.
165    real(dp) :: occ_factor                        ! normalization for integrals over FS, for num of spins, spinors, etc...
166 
167    real(dp) :: mustar                            ! mustar parameter
168    real(dp) :: fermie                            ! Fermi energy (Ha), either comes from wfk file or from anaddb input file
169    real(dp) :: elphsmear                         ! smearing width for gaussian integration or buffer in energy for
170                                                  ! calculations with tetrahedra (telphint=0)
171    real(dp) :: tempermin                         ! minimum temperature at which resistivity etc are calculated (in K)
172    real(dp) :: temperinc                         ! interval temperature grid on which resistivity etc are calculated (in K)
173 
174    character(len=fnlen) :: elph_base_name        !base name for output files
175 
176    integer,allocatable :: qirredtofull(:)            !mapping between the qpoints found in the GGK file
177                                                  !and the array of qpoints generated by the code
178 
179    real(dp),allocatable :: wtq(:)                    !weight for each qpoint in the full grid spqt
180                                                  !if a point is not in the IBZ ==>  wtq=0
181                                                  !MG we can also use indqpt
182 
183    real(dp),allocatable :: n0(:)                     !DOS at the Fermi level (states/Ha/spin)
184    real(dp),allocatable :: qpt_full(:,:)             !special q points obtained by the Monkhorst & Pack method,
185                                                  !in reduced coordinates
186 
187 
188    real(dp),allocatable :: gkk_intweight(:,:,:)      ! (nFSband,nkpt_fine,nsppol)
189                                                  !integration weights for gkk matrix elements on FS:
190                                                  !if ep_keepbands == 0 all are 1
191                                                  !if ep_keepbands == 1 then = to wtk_phon in elphon
192                                                  !DOES NOT INCLUDE FACTOR OF 1/nkpt_phon
193 
194    real(dp),allocatable :: gkk_velocwtk(:,:,:)      ! (nFSband,nkpt_fine,nsppol)
195 
196    real(dp),allocatable :: gkk_vvelocwtk(:,:,:)      ! (nFSband,nkpt_fine,nsppol)
197 
198    real(dp),allocatable :: gkk_qpt(:,:,:,:,:,:)      ! (2, ngkkband*ngkkband, nbranch*nbranch, nkpt_phon, nsppol, nqptirred)
199                                                  !Now gkq contains gkk2 matrices on basic qpts,
200                                                  !summed over bands if ngkkband==1
201 
202 
203    real(dp),allocatable :: gkk_rpt(:,:,:,:,:,:)      ! (2, ngkkband**2, nbranch**2, nkpt_phon, nsppol, nrpt)
204                                                  !For the moment, gkk_rpt in memory is out of the question
205    real(dp),allocatable :: gkk2(:,:,:,:,:,:)         ! (nbranch, ngkkband,ngkkband, nkpt_phon, nkpt_phon, nsppol)
206 
207    real(dp),allocatable :: gamma_qpt(:,:,:,:)        !gamma matrices integrated over kpoint coeff
208                                                  !  and bands: still depends on qpt
209                                                  ! dims= 2, nbranch**2, nsppol, nqpt
210    real(dp),allocatable :: gamma_rpt(:,:,:,:)
211                                                  ! dims= 2, nbranch**2, nsppol, nrpt
212 !NOTE: choice to put nsppol before or after nqpt is a bit arbitrary
213 !   abinit uses nband,nkpt,nsppol, but here for convenience nkpt_phon,nsppol,nqpt
214 !   as interpolation is on qpt
215 
216    real(dp),allocatable :: phfrq(:,:)                !phonon frequencies
217    real(dp),allocatable :: a2f(:,:,:)                !a2f function
218 
219    real(dp),allocatable :: qgrid_data(:,:,:,:)       !e-ph values calculated over the irreducible part of the q-grid:
220                                                  !first entry  =  index of the q-point,
221                                                  !second index =  branch index
222                                                  !the third slice contains the frequency, the linewidth and lambda(q,nu)
223                                                  !for that particular phonon mode
224                                                  ! dims= nqptirred,elph_ds%nbranch,nsppol,3
225 
226  end type elph_type
227 
228  public :: elph_ds_clean

defs_elphon/gam_mult_displ [ Functions ]

[ Top ] [ defs_elphon ] [ Functions ]

NAME

 gam_mult_displ

FUNCTION

 This routine takes the bare gamma matrices and multiplies them
  by the displ_red matrices (related to the scalprod variable)

INPUTS

   nbranch = number of phonon branches (3*natom)
   displ_red = phonon mode displacement vectors in reduced coordinates.
   gam_bare = bare gamma matrices before multiplication

OUTPUT

   gam_now = output gamma matrices multiplied by displacement matrices

SOURCE

600 subroutine gam_mult_displ(nbranch, displ_red, gam_bare, gam_now)
601 
602 !Arguments -------------------------------
603  integer, intent(in)  :: nbranch
604  real(dp), intent(in)  :: displ_red(2,nbranch,nbranch)
605  real(dp), intent(in)  :: gam_bare(2,nbranch,nbranch)
606  real(dp), intent(out) :: gam_now(2,nbranch,nbranch)
607 
608 !Local variables -------------------------
609  real(dp) :: zgemm_tmp_mat(2,nbranch,nbranch)
610 
611 ! *********************************************************************
612 
613  gam_now = zero
614 
615  call zgemm('c','n',nbranch,nbranch,nbranch,cone,displ_red,nbranch,gam_bare,nbranch,czero,zgemm_tmp_mat,nbranch)
616 
617  call zgemm('n','n',nbranch,nbranch,nbranch,cone,zgemm_tmp_mat,nbranch,displ_red,nbranch,czero,gam_now,nbranch)
618 
619 end subroutine gam_mult_displ

m_fstab/mkqptequiv [ Functions ]

[ Top ] [ m_fstab ] [ Functions ]

NAME

 mkqptequiv

FUNCTION

 This routine determines the equivalence between
   1) qpoints and fermi surface kpoints
   2) qpoints under symmetry operations

INPUTS

   Cryst<crystal_t>=Info on unit cell and symmetries.
   kpt_phon = fermi surface kpoints
   nkpt_phon = number of kpoints in the full FS set
   nqpt = number of qpoints
   qpt_full = qpoint coordinates

OUTPUT

   FSfullpqtofull = mapping of k + q onto k' for k and k' in full BZ
   qpttoqpt(itim,isym,iqpt) = qpoint index which transforms to iqpt under isym and with time reversal itim.

NOTES

   REMOVED 3/6/2008: much too large matrix, and not used at present
       FStoqpt = mapping of kpoint pairs (1 irreducible and 1 full) to qpoints

SOURCE

1226 subroutine mkqptequiv(FSfullpqtofull,Cryst,kpt_phon,nkpt_phon,nqpt,qpttoqpt,qpt_full,mqtofull)
1227 
1228 !Arguments ------------------------------------
1229 !scalars
1230  integer,intent(in) :: nkpt_phon,nqpt
1231  type(crystal_t),intent(in) :: Cryst
1232 !arrays
1233  integer,intent(out) :: FSfullpqtofull(nkpt_phon,nqpt),qpttoqpt(2,Cryst%nsym,nqpt)
1234  integer,intent(out),optional :: mqtofull(nqpt)
1235  real(dp),intent(in) :: kpt_phon(3,nkpt_phon),qpt_full(3,nqpt)
1236 
1237 !Local variables-------------------------------
1238 !scalars
1239  integer :: ikpt_phon,iFSqpt,iqpt,isym,symrankkpt_phon
1240  !character(len=500) :: message
1241  type(krank_t) :: krank
1242 !arrays
1243  real(dp) :: tmpkpt(3),gamma_kpt(3)
1244 
1245 ! *************************************************************************
1246 
1247  call wrtout(std_out,' mkqptequiv : making rankkpt_phon and invrankkpt_phon',"COLL")
1248 
1249  krank = krank_new(nkpt_phon, kpt_phon)
1250 
1251  FSfullpqtofull = -999
1252  gamma_kpt(:) = zero
1253 
1254  do ikpt_phon=1,nkpt_phon
1255    do iqpt=1,nqpt
1256      ! tmpkpt = jkpt = ikpt + qpt
1257      tmpkpt(:) = kpt_phon(:,ikpt_phon) + qpt_full(:,iqpt)
1258 
1259      ! which kpt is it among the full FS kpts?
1260      symrankkpt_phon = krank%get_rank(tmpkpt)
1261 
1262      FSfullpqtofull(ikpt_phon,iqpt) = krank%invrank(symrankkpt_phon)
1263      if (FSfullpqtofull(ikpt_phon, iqpt) == -1) then
1264        ABI_ERROR("looks like no kpoint equiv to k+q !!!")
1265      end if
1266 
1267    end do
1268  end do
1269 
1270  if (present(mqtofull)) then
1271    do iqpt=1,nqpt
1272      tmpkpt(:) = gamma_kpt(:) - qpt_full(:,iqpt)
1273 
1274      ! which kpt is it among the full FS kpts?
1275      symrankkpt_phon = krank%get_rank(tmpkpt)
1276 
1277      mqtofull(iqpt) = krank%invrank(symrankkpt_phon)
1278      if (mqtofull(iqpt) == -1) then
1279        ABI_ERROR("looks like no kpoint equiv to -q !!!")
1280      end if
1281    end do
1282  end if
1283 
1284  call krank%free()
1285 
1286  ! start over with q grid
1287  call wrtout(std_out,' mkqptequiv : FSfullpqtofull made. Do qpttoqpt',"COLL")
1288 
1289  krank = krank_new(nqpt, qpt_full)
1290 
1291  qpttoqpt(:,:,:) = -1
1292  do iFSqpt=1,nqpt
1293    do isym=1,Cryst%nsym
1294      tmpkpt(:) =  Cryst%symrec(:,1,isym)*qpt_full(1,iFSqpt) &
1295                 + Cryst%symrec(:,2,isym)*qpt_full(2,iFSqpt) &
1296                 + Cryst%symrec(:,3,isym)*qpt_full(3,iFSqpt)
1297 
1298      symrankkpt_phon = krank%get_rank(tmpkpt)
1299      if (krank%invrank(symrankkpt_phon) == -1) then
1300        ABI_ERROR("looks like no kpoint equiv to q by symmetry without time reversal!!!")
1301      end if
1302      qpttoqpt(1,isym,krank%invrank(symrankkpt_phon)) = iFSqpt
1303 
1304      tmpkpt = -tmpkpt
1305      symrankkpt_phon = krank%get_rank(tmpkpt)
1306      if (krank%invrank(symrankkpt_phon) == -1) then
1307        ABI_ERROR('looks like no kpoint equiv to q by symmetry with time reversal!!!')
1308      end if
1309      qpttoqpt(2,isym,krank%invrank(symrankkpt_phon)) = iFSqpt
1310    end do
1311  end do
1312 
1313  call krank%free()
1314 
1315 end subroutine mkqptequiv