TABLE OF CONTENTS


ABINIT/m_paw_dmft [ Modules ]

[ Top ] [ Modules ]

NAME

  m_paw_dmft

FUNCTION

COPYRIGHT

 Copyright (C) 2006-2024 ABINIT group (BAmadon)
 This file is distributed under the terms of the
 GNU General Public License, see ~abinit/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .

INPUTS

OUTPUT

SOURCE

19 #if defined HAVE_CONFIG_H
20 #include "config.h"
21 #endif
22 
23 #include "abi_common.h"
24 
25 MODULE m_paw_dmft
26 
27  use defs_basis
28  use m_CtqmcInterface
29  use m_errors
30  use m_abicore
31  use m_xmpi
32  use m_data4entropyDMFT
33  use m_dtset
34 
35  use defs_datatypes, only : pseudopotential_type
36  use defs_abitypes, only : MPI_type
37  use m_io_tools,  only : open_file
38  use m_pawtab,    only : pawtab_type
39 
40  implicit none
41 
42  private
43 
44  public :: init_dmft
45  public :: init_sc_dmft
46  public :: construct_nwli_dmft
47  public :: destroy_dmft
48  public :: destroy_sc_dmft
49  public :: print_dmft
50  public :: print_sc_dmft
51  public :: saveocc_dmft
52  public :: readocc_dmft
53 
54  private :: init_sc_dmft_paralkgb, destroy_sc_dmft_paralkgb

m_paw_dmft/construct_nwli_dmft [ Functions ]

[ Top ] [ m_paw_dmft ] [ Functions ]

NAME

 construct_nwli_dmft

FUNCTION

  Compute linear frequencies

INPUTS

  paw_dmft=structure for dmft
  nwli=number of linear frequencies

 OUTPUTS
  omegali(1:nwli)=computed frequencies

SOURCE

921 subroutine construct_nwli_dmft(paw_dmft,nwli,omega_li)
922 
923  type(paw_dmft_type), intent(in) :: paw_dmft
924  integer, intent(in) :: nwli
925  real(dp), intent(out) :: omega_li(:)
926  !fortran2003 ?
927  !real(dp), allocatable, intent(inout) :: omega_li(:)
928  integer :: ifreq
929  real(dp) :: factor
930  character(len=100) :: message
931 
932 ! if (allocated(omega_li)) then
933    if (size(omega_li) .ne. nwli) then
934      write(message,'(2a,i8,a,i8)') ch10, "Number of linear frequencies asked is", &
935        &    nwli, "whereas dimension of array omega_li is", size(omega_li)
936      ABI_BUG(message)
937 !     ABI_FREE(omega_li)
938 !     ABI_MALLOC(omega_li,(nwli))
939 !     write(*,*) "RESIZE"
940 !     call flush(6)
941    endif
942 !     write(*,*) "NOTHING"
943 !     call flush(6)
944 ! else
945 !     write(*,*) "ALLOCATE"
946 !     call flush(6)
947 !   ABI_MALLOC(omega_li,(nwli))
948 ! endif
949 
950 ! Set up linear frequencies
951  factor = pi*paw_dmft%temp
952  do ifreq=1,nwli
953    omega_li(ifreq)=factor*(real(2*ifreq-1,kind=dp))
954    ! (2(ifreq-1)+1 = 2ifreq-1
955  enddo
956 end subroutine construct_nwli_dmft

m_paw_dmft/construct_nwlo_dmft [ Functions ]

[ Top ] [ m_paw_dmft ] [ Functions ]

NAME

 construct_nwlo_dmft

FUNCTION

  Allocate log frequencies if used and compute them as well as their weight

INPUTS

  paw_dmft=structure for dmft calculation

SOURCE

 972 !! NOTE
 973 !! The part of the code which deals
 974 !! with the use of logarithmic frequencies
 975 !! is a modification of the GNU GPL
 976 !! code available on http://dmft.rutgers.edu/ and
 977 !! described in the  RMP paper written by
 978 !! G.Kotliar,  S.Y.Savrasov, K.Haule, V.S.Oudovenko, O.Parcollet, C.A.Marianetti.
 979 !!
 980 
 981 subroutine construct_nwlo_dmft(paw_dmft)
 982  use m_splines
 983 
 984  type(paw_dmft_type), intent(inout) :: paw_dmft
 985  integer :: cubic_freq
 986  integer :: ifreq,ifreq2
 987  ! for parallel
 988  integer :: myproc, nproc, spacecomm
 989  integer :: deltaw, residu, omegaBegin, omegaEnd
 990  ! end
 991  character(len=500) :: message
 992  real(dp) :: deltaomega,expfac,omegamaxmin,prefacexp,AA,BB,CC,nlin,nlog,t1
 993  real(dp) :: wl
 994  complex(dpc):: ybcbeg,ybcend
 995  integer, allocatable :: select_log(:)
 996  real(dp), allocatable :: omega_lo_tmp(:)
 997  real(dp), allocatable :: omega_li(:)
 998  real(dp), allocatable :: wgt_wlo(:)
 999  complex(dpc), allocatable :: tospline_lo(:), splined_li(:),ysplin2_lo(:)
1000 
1001  ABI_MALLOC(omega_lo_tmp,(paw_dmft%dmft_nwlo))
1002  ABI_MALLOC(wgt_wlo,(paw_dmft%dmft_nwlo))
1003 
1004 !==  Variables for DMFT related to frequencies
1005 ! the part of the code which deals
1006 ! with the use of logarithmic frequencies
1007 ! is a modification of the GNU GPL
1008 ! code available on http://dmft.rutgers.edu/ and
1009 ! described in the  RMP paper written by
1010 ! G.Kotliar, S.Y.Savrasov, K.Haule, V.S.Oudovenko, O.Parcollet, C.A.Marianetti
1011 
1012 !========================================
1013 !== construct log. freq.
1014  if(paw_dmft%dmft_log_freq==1) then
1015 !=======================================
1016    cubic_freq=0
1017    !omegamaxmin=paw_dmft%omega_li(paw_dmft%dmft_nwli)-paw_dmft%omega_li(paw_dmft%dmftqmc_l+1)
1018    omegamaxmin=pi*paw_dmft%temp*two*real(paw_dmft%dmft_nwli-paw_dmft%dmftqmc_l-1,kind=dp)
1019 
1020    if(cubic_freq==1) then
1021 
1022      if (paw_dmft%dmft_solv .eq. 5 ) then
1023        write(message, '(2a)') ch10, "Warning : Cubish Mesh not tested with CT-QMC"
1024        ABI_WARNING(message)
1025      end if
1026 !  ------------  CUBIC MESH MESH
1027 !    useless
1028      nlin=dble(paw_dmft%dmft_nwli)
1029      nlog=dble(paw_dmft%dmft_nwlo)
1030      AA=((nlin-one)/nlin/(nlog**2-one)-one/(three*nlin))/((nlog**3-one)/(nlog**2-one)-seven/three)
1031      BB=(one/nlin - seven*AA)/three
1032      CC=-AA-BB
1033 !    AA=((nlin-one)/nlin/(nlog-one)-one/(nlin))/((nlog**2-one)/(nlog-one)-three)
1034 !    BB=(one/nlin - three*AA)
1035 !    CC=-AA-BB
1036      write(message, '(a,16x,2(2x,a))') ch10,"  Cubic Mesh Parameters are"
1037      call wrtout(std_out,message,'COLL')
1038      write(message, '(3x,a,3(2x,e13.5))') "AA,BB,CC",AA,BB,CC
1039      call wrtout(std_out,message,'COLL')
1040      do ifreq=1,paw_dmft%dmft_nwlo
1041        t1=dble(ifreq)
1042        !omega_lo_tmp(ifreq)=(AA*t1**3+BB*t1**2+CC)*omegamaxmin+paw_dmft%omega_li(1)
1043        omega_lo_tmp(ifreq)=(AA*t1**3+BB*t1**2+CC)*omegamaxmin+paw_dmft%temp*pi
1044 !       paw_dmft%omega_lo(ifreq)=(AA*t1**2+BB*t1+CC)*omegamaxmin+paw_dmft%omega_li(1)
1045 !     write(69,*) paw_dmft%omega_lo(ifreq),0.5
1046      enddo
1047    else
1048      if(paw_dmft%dmft_solv<4) paw_dmft%dmftqmc_l=0
1049 
1050 !   ------------  LOGARITHMIC MESH
1051      deltaomega=0.5_dp
1052      expfac=log(omegamaxmin/deltaomega)/(float(paw_dmft%dmft_nwlo-paw_dmft%dmftqmc_l-1)/two)
1053      prefacexp=omegamaxmin/(exp(expfac*float(paw_dmft%dmft_nwlo-paw_dmft%dmftqmc_l-1))-one)
1054      ABI_MALLOC(select_log,(paw_dmft%dmft_nwlo))
1055      select_log=0
1056 
1057 !   ------------ IMPOSE LINEAR MESH for w < 2*w_n=(2*l-1)pi/beta
1058 !         Check variables (Already done in chkinp if dmft_solv==5)
1059      if (paw_dmft%dmftqmc_l .gt. paw_dmft%dmft_nwlo) then
1060        write(message, '(a,a,i6)' )ch10,&
1061 &       ' ERROR: dmft_nwlo has to be at least equal to 2xdmftqmc_l :',2*paw_dmft%dmftqmc_l
1062        ABI_ERROR(message)
1063      end if
1064 !         End Check
1065 
1066      call construct_nwli_dmft(paw_dmft,paw_dmft%dmftqmc_l,omega_lo_tmp(1:paw_dmft%dmftqmc_l))
1067      select_log(1:paw_dmft%dmftqmc_l) = (/ (ifreq,ifreq=1,paw_dmft%dmftqmc_l) /)
1068 
1069      !do ifreq=1,paw_dmft%dmftqmc_l
1070      !  omega_lo_tmp(ifreq)=(two*DBLE(ifreq-1)+one)*pi*paw_dmft%temp
1071      !  select_log(ifreq)=ifreq
1072      !enddo
1073 
1074 !   ------------ COMPLETE FREQUENCIES WITH LOG MESH
1075      wl = paw_dmft%temp*pi*real(2*paw_dmft%dmftqmc_l+1,kind=dp)
1076      do ifreq=1,paw_dmft%dmft_nwlo-paw_dmft%dmftqmc_l
1077        !omega_lo_tmp(ifreq+paw_dmft%dmftqmc_l)=prefacexp*(exp(expfac*float(ifreq-1))-one)+paw_dmft%omega_li(paw_dmft%dmftqmc_l+1)
1078        omega_lo_tmp(ifreq+paw_dmft%dmftqmc_l)=prefacexp*(exp(expfac*real(ifreq-1,kind=dp))-one) + wl
1079 !       -------- Impose that the each frequency of the logarithmic mesh is on a Matsubara frequency
1080 ! FIXME : This may be done for all solver, not only for QMCs
1081        if(paw_dmft%dmft_solv>=4) then
1082          ! Compute the integer "n" of iwn
1083          ifreq2 = nint((omega_lo_tmp(ifreq+paw_dmft%dmftqmc_l)/(paw_dmft%temp*pi)-one)*half)
1084          ! compute freq
1085          omega_lo_tmp(ifreq+paw_dmft%dmftqmc_l)= (dble(ifreq2)*two+one)*pi*paw_dmft%temp
1086 
1087          if((ifreq2+1)>paw_dmft%dmft_nwli) then
1088            write(message, '(a,a,i8)' )ch10,&
1089 &          ' BUG: init_dmft,   dimension  of array select_log is about to be overflown',&
1090 &          (ifreq2+1)
1091            ABI_BUG(message)
1092          endif
1093          select_log(paw_dmft%dmftqmc_l+ifreq)=ifreq2+1
1094        endif
1095      enddo
1096 
1097 !       -------- Suppress duplicate frequencies
1098 ! FIXME : So this also should be done for all solver and remove useless
1099 ! frequencies
1100      if(paw_dmft%dmft_solv>=4) then
1101        ifreq2=1
1102        do ifreq=2,paw_dmft%dmft_nwlo-1
1103          if(select_log(ifreq2).ne.select_log(ifreq)) then
1104            ifreq2=ifreq2+1
1105            omega_lo_tmp(ifreq2)=omega_lo_tmp(ifreq)
1106            select_log(ifreq2)=select_log(ifreq)
1107          endif
1108        enddo
1109        paw_dmft%dmft_nwlo=ifreq2+1
1110      endif
1111 
1112      omega_lo_tmp(1)=paw_dmft%temp*pi
1113      omega_lo_tmp(paw_dmft%dmft_nwlo)=paw_dmft%temp*pi*real(2*paw_dmft%dmft_nwli-1,kind=dp)
1114    endif
1115 
1116 !=======================
1117 !== construct weight for log. freq.
1118 !=======================
1119    ABI_MALLOC(tospline_lo,(paw_dmft%dmft_nwlo))
1120    ABI_MALLOC(splined_li,(paw_dmft%dmft_nwli))
1121    ABI_MALLOC(ysplin2_lo,(paw_dmft%dmft_nwlo))
1122    if (allocated(omega_li)) then
1123      ABI_FREE(omega_li)
1124    endif
1125    ABI_MALLOC(omega_li,(1:paw_dmft%dmft_nwli))
1126    call construct_nwli_dmft(paw_dmft,paw_dmft%dmft_nwli,omega_li)
1127   
1128    !Parallelisation over frequencies!
1129    ! ============= Set up =============
1130    myproc = paw_dmft%myproc
1131    nproc = paw_dmft%nproc
1132    spacecomm = paw_dmft%spacecomm
1133    deltaw = paw_dmft%dmft_nwlo / nproc
1134    residu = paw_dmft%dmft_nwlo - nproc*deltaw
1135    if ( myproc .LT. nproc - residu ) then
1136      omegaBegin = 1 + myproc*deltaw
1137      omegaEnd   = (myproc + 1)*deltaw
1138    else
1139      omegaBegin = 1 + myproc*(deltaw + 1) -nproc + residu
1140      omegaEnd = omegaBegin + deltaw
1141    end if
1142    wgt_wlo(:)=zero
1143    ! ============= END Set up =============
1144    do ifreq=omegaBegin,omegaEnd
1145      tospline_lo=cmplx(0_dp,0_dp,kind=dp)
1146 !    do ifreq1=1,paw_dmft%dmft_nwlo
1147      tospline_lo(ifreq)=cmplx(1_dp,0_dp,kind=dp)
1148 !    tospline_lo(ifreq1)=ifreq1**2-ifreq1
1149 !    enddo
1150      splined_li=cmplx(0_dp,0_dp,kind=dp)
1151 !    ybcbeg=cmplx(one/tol16**2,zero)
1152 !    ybcend=cmplx(one/tol16**2,zero)
1153      ybcbeg=czero
1154      ybcend=czero
1155 
1156 
1157 !==         spline delta function
1158      call spline_complex( omega_lo_tmp, tospline_lo, paw_dmft%dmft_nwlo, &
1159      & ybcbeg, ybcend, ysplin2_lo)
1160 !   do ifreq1=1,paw_dmft%dmft_nwlo
1161 !    write(6588,*) paw_dmft%omega_lo(ifreq1),ysplin2_lo(ifreq1)
1162 !   enddo
1163 
1164      call splint_complex( paw_dmft%dmft_nwlo, omega_lo_tmp, tospline_lo,&
1165      & ysplin2_lo, paw_dmft%dmft_nwli, omega_li, splined_li)
1166 
1167 !==         accumulate weights
1168      wgt_wlo(ifreq)=zero
1169      do ifreq2=1,paw_dmft%dmft_nwli
1170        wgt_wlo(ifreq)=wgt_wlo(ifreq)+real(splined_li(ifreq2),kind=dp)
1171      enddo
1172 ! do ifreq1=1,paw_dmft%dmft_nwlo
1173 !  write(6688,*) paw_dmft%omega_lo(ifreq1),tospline_lo(ifreq1)
1174 ! enddo
1175 ! do ifreq1=1,paw_dmft%dmft_nwli
1176 !  write(6788,*) paw_dmft%omega_li(ifreq1),splined_li(ifreq1)
1177    enddo
1178   ! ============= Gatherall  =============
1179    call xmpi_sum(wgt_wlo, spacecomm, residu)
1180   ! ============= END Gatherall ==========
1181   ! end parallelisation over frequencies
1182 
1183    ABI_FREE(tospline_lo)
1184    ABI_FREE(splined_li)
1185    ABI_FREE(ysplin2_lo)
1186 ! if(abs(dtset%pawprtvol)>=3) then
1187    write(message, '(a,18x,2(2x,a21))') ch10,"Log. Freq","weight"
1188    call wrtout(std_out,message,'COLL')
1189    do ifreq=1,paw_dmft%dmft_nwlo
1190      write(message, '(3x,a9,i6,2(2x,e21.14))') "--ifreq--",ifreq,omega_lo_tmp(ifreq),wgt_wlo(ifreq)
1191      call wrtout(std_out,message,'COLL')
1192    enddo
1193    write(message, '(3x,a,i6)') "  Total number of log frequencies is", paw_dmft%dmft_nwlo
1194    call wrtout(std_out,message,'COLL')
1195    ifreq2 = 1
1196    do ifreq=1,min(30,paw_dmft%dmft_nwlo)
1197      write(message, '(3x,a9,i6,2(2x,e21.14))') "--ifreq--",ifreq,omega_li(ifreq)
1198      call wrtout(std_out,message,'COLL')
1199      if ( select_log(ifreq2) .eq. ifreq ) then
1200        write(message, '(3x,a,i4,2(2x,i5))') "--sel_log",1
1201        ifreq2 = ifreq+1
1202      else
1203        write(message, '(3x,a,i4,2(2x,i5))') "--sel_log",0
1204      end if
1205      call wrtout(std_out,message,'COLL')
1206    enddo
1207    write(message, '(3x,2a)') "--ifreq--","..."
1208    call wrtout(std_out,message,'COLL')
1209    write(message, '(3x,a,i6,2(2x,e13.5))') "--ifreq--",paw_dmft%dmft_nwli,omega_li(paw_dmft%dmft_nwli)
1210    call wrtout(std_out,message,'COLL')
1211 !   endif
1212    ABI_FREE(select_log)
1213    ABI_FREE(omega_li)
1214 
1215 !=========================================================
1216 !== do not construct log. freq. and use linear frequencies
1217  else
1218 !=========================================================
1219    write(message, '(a,10x,2(2x,a))') ch10,"   Use of linear frequencies for DMFT calculation"
1220    call wrtout(std_out,message,'COLL')
1221    call construct_nwli_dmft(paw_dmft,paw_dmft%dmft_nwli,omega_lo_tmp)
1222    wgt_wlo=one
1223  endif
1224 
1225  ! Should be check but since type definition does not initialize pointer with
1226  ! =>null() (fortran95 and later) it produces conditional jump in valgrind
1227  !if ( associated(paw_dmft%omega_lo) ) then
1228  !  ABI_FREE(paw_dmft%omega_lo)
1229  !endif
1230  !if ( associated(paw_dmft%wgt_wlo) ) then
1231  !  ABI_FREE(paw_dmft%wgt_wlo)
1232  !endif
1233  ABI_MALLOC(paw_dmft%omega_lo,(paw_dmft%dmft_nwlo))
1234  ABI_MALLOC(paw_dmft%wgt_wlo,(paw_dmft%dmft_nwlo))
1235  paw_dmft%omega_lo(1:paw_dmft%dmft_nwlo) = omega_lo_tmp(1:paw_dmft%dmft_nwlo)
1236  paw_dmft%wgt_wlo(1:paw_dmft%dmft_nwlo) = wgt_wlo(1:paw_dmft%dmft_nwlo)
1237  ABI_FREE(omega_lo_tmp)
1238  ABI_FREE(wgt_wlo)
1239 end subroutine construct_nwlo_dmft

m_paw_dmft/destroy_dmft [ Functions ]

[ Top ] [ m_paw_dmft ] [ Functions ]

NAME

 destroy_dmft

FUNCTION

  deallocate some variables related to paw_dmft

INPUTS

  paw_dmft

OUTPUT

SOURCE

1256 subroutine destroy_dmft(paw_dmft)
1257 
1258 !Arguments ------------------------------------
1259 !scalars
1260  type(paw_dmft_type),intent(inout) :: paw_dmft
1261 
1262 !Local variables-------------------------------
1263  integer :: iatom
1264 
1265 ! *********************************************************************
1266 
1267    if (paw_dmft%dmft_solv == 5 .and. allocated(paw_dmft%hybrid)) then
1268      do iatom=1, size(paw_dmft%hybrid) !paw_dmft%natom
1269        !if(paw_dmft%lpawu(iatom)/=-1) then
1270          call ctqmcinterface_finalize(paw_dmft%hybrid(iatom))
1271        !endif
1272      enddo
1273      ABI_FREE(paw_dmft%hybrid)
1274    endif
1275    if (allocated(paw_dmft%psichi))  then
1276      ABI_FREE(paw_dmft%psichi)
1277    end if
1278 !   paw_dmft%wtk is only an explicit pointer =>dtset%wtk
1279 !   if (associated(paw_dmft%wtk)) deallocate(paw_dmft%wtk)
1280    paw_dmft%wtk => null()
1281    paw_dmft%fixed_self => null()
1282    if (allocated(paw_dmft%eigen_dft))  then
1283      ABI_FREE(paw_dmft%eigen_dft)
1284    endif
1285    if (associated(paw_dmft%omega_lo))  then
1286      ABI_FREE(paw_dmft%omega_lo)
1287    end if
1288    if (associated(paw_dmft%omega_r))  then
1289      ABI_FREE(paw_dmft%omega_r)
1290    end if
1291    if (associated(paw_dmft%wgt_wlo))  then
1292      ABI_FREE(paw_dmft%wgt_wlo)
1293    end if
1294    if (allocated(paw_dmft%lpawu))  then
1295      ABI_FREE(paw_dmft%lpawu)
1296    end if
1297 
1298 end subroutine destroy_dmft

m_paw_dmft/destroy_sc_dmft [ Functions ]

[ Top ] [ m_paw_dmft ] [ Functions ]

NAME

 destroy_sc_dmft

FUNCTION

  deallocate paw_dmft

INPUTS

  paw_dmft

OUTPUT

SOURCE

1315 subroutine destroy_sc_dmft(paw_dmft)
1316 
1317 !Arguments ------------------------------------
1318 !scalars
1319  type(paw_dmft_type),intent(inout) :: paw_dmft
1320 
1321 !Local variables-------------------------------
1322  character(len=500) :: message
1323 
1324 ! *********************************************************************
1325 
1326  if (( .not. allocated(paw_dmft%occnd) .or. .not. allocated(paw_dmft%band_in) &
1327 &  .or. .not. allocated(paw_dmft%include_bands) .or. .not. allocated(paw_dmft%exclude_bands)) &
1328 &  .and. paw_dmft%use_dmft == 1 )  then
1329   write(message, '(a,a,a)' )&
1330 &  '  an array is not allocated and is not deallocated with use_dmft==1 ',ch10, &
1331 &  '  Action : check the code'
1332   ABI_WARNING(message)
1333  endif
1334  if ( allocated(paw_dmft%occnd) )          then
1335    ABI_FREE(paw_dmft%occnd)
1336  end if
1337  if ( allocated(paw_dmft%band_in) )        then
1338    ABI_FREE(paw_dmft%band_in)
1339  end if
1340  if ( allocated(paw_dmft%include_bands) )  then
1341    ABI_FREE(paw_dmft%include_bands)
1342  end if
1343  if ( allocated(paw_dmft%exclude_bands) )  then
1344    ABI_FREE(paw_dmft%exclude_bands)
1345  end if
1346 
1347  call destroy_sc_dmft_paralkgb(paw_dmft)
1348 
1349 end subroutine destroy_sc_dmft

m_paw_dmft/destroy_sc_dmft_paralkgb [ Functions ]

[ Top ] [ m_paw_dmft ] [ Functions ]

NAME

 destroy_sc_dmft_paralkgb

FUNCTION

   deallocate bandc_proc and use_bandc

INPUTS

  paw_dmft   = data structure

OUTPUT

SOURCE

1687 subroutine destroy_sc_dmft_paralkgb(paw_dmft)
1688 
1689 !Arguments ------------------------------------
1690  type(paw_dmft_type),intent(inout) :: paw_dmft
1691 ! *********************************************************************
1692 
1693  if ( allocated(paw_dmft%bandc_proc) )  then
1694    ABI_FREE(paw_dmft%bandc_proc)
1695  end if
1696 
1697  if ( allocated(paw_dmft%use_bandc) )  then
1698    ABI_FREE(paw_dmft%use_bandc)
1699  end if
1700 
1701 end subroutine destroy_sc_dmft_paralkgb

m_paw_dmft/init_dmft [ Functions ]

[ Top ] [ m_paw_dmft ] [ Functions ]

NAME

 init_dmft

FUNCTION

  Allocate variables and setup DFT hamiltonian and related data
  (init_sc_dmft has to been called before)

INPUTS

  dmatpawu   = fixed occupation matrix of correlated orbitals
  eigen      = DFT eigenvalues
  fermie_dft = DFT Fermi level
  psichi     = <chi|Psi> projection of KS states over atomic !wavefunction
  nkpt       = number of k-points
  nsppol     = number of spin polarisation
  nspinor    = number of spinorial component

SOURCE

554 !! NOTE
555 !! The part of the code which deals
556 !! with the use of logarithmic frequencies
557 !! is a modification of the GNU GPL
558 !! code available on http://dmft.rutgers.edu/ and
559 !! described in the  RMP paper written by
560 !! G.Kotliar,  S.Y.Savrasov, K.Haule, V.S.Oudovenko, O.Parcollet, C.A.Marianetti.
561 !!
562 
563 subroutine init_dmft(dmatpawu, dtset, fermie_dft, fnametmp_app, fnamei, nspinor, paw_dmft, pawtab, psps, typat)
564 
565  use m_splines
566  !use m_CtqmcInterface
567 
568 !Arguments ------------------------------------
569 !scalars
570  integer, intent(in)  :: nspinor
571  real(dp), intent(in) :: fermie_dft
572 !type
573  type(pseudopotential_type), intent(in) :: psps
574  type(dataset_type),target,intent(in) :: dtset
575  type(pawtab_type),intent(in)  :: pawtab(psps%ntypat*psps%usepaw)
576  type(paw_dmft_type),intent(inout) :: paw_dmft
577  character(len=fnlen), intent(in) :: fnametmp_app
578  character(len=fnlen), intent(in) :: fnamei
579 !arrays
580  integer,intent(in) :: typat(dtset%natom)
581  real(dp),intent(in),target :: dmatpawu(:,:,:,:)
582 !Local variables ------------------------------------
583  integer :: grid_unt,ikpt,isym,itypat,nsppol,mbandc,maxlpawu, iatom, ifreq, ioerr
584  integer :: nflavor,ngrid,iexist2
585  character(len=fnlen) :: tmpfil
586  real(dp) :: sumwtk
587  character(len=500) :: message
588  real(dp) :: unit_e,step
589  logical :: lexist
590 ! *********************************************************************
591 
592  nsppol = dtset%nsppol
593  if(dtset%ucrpa==0) then
594  write(message,'(6a)') ch10,' ====================================', &
595 &                      ch10,' =====  Start of DMFT calculation', &
596 &                      ch10,' ===================================='
597  else if(dtset%ucrpa>0) then
598  write(message,'(6a)') ch10,' ============================================================', &
599 &                      ch10,' =====  Initialize construction of Wannier in DMFT routines',&
600 &                      ch10,' ============================================================'
601  endif
602  call wrtout(std_out,message,'COLL')
603 
604  unit_e=2_dp
605 !=======================
606 !==  Check sym
607 !=======================
608  do isym=1,dtset%nsym
609    if(dtset%symafm(isym)<0) then
610      message = 'symafm negative is not implemented in DMFT '
611      ABI_ERROR(message)
612    endif
613  enddo
614 
615 ! paw_dmft%dmft_mag=0
616 ! do iatom=1,dtset%natom
617 !   do  ii=1,3
618 !     if ( dtset(ii,iatom) > 0.001 ) paw_dmft%dmft_mag=1
619 !   enddo
620 ! enddo
621 
622 !=======================
623 !==  Define integers and reals
624 !=======================
625  paw_dmft%fermie_dft=fermie_dft ! in Ha
626  paw_dmft%fermie= fermie_dft
627  if(nspinor==1) then
628    if(paw_dmft%nsppol==2) then
629      paw_dmft%nelectval= dtset%nelect-float(paw_dmft%dmftbandi-1)*paw_dmft%nsppol
630    else if(paw_dmft%nsppol==1) then
631      paw_dmft%nelectval= dtset%nelect-float(paw_dmft%dmftbandi-1)*paw_dmft%nsppol*2
632    endif
633  else if (nspinor==2) then
634    paw_dmft%nelectval= dtset%nelect-float(paw_dmft%dmftbandi-1)*paw_dmft%nsppol
635  endif
636  paw_dmft%filapp= fnametmp_app
637  paw_dmft%filnamei= fnamei
638  paw_dmft%natpawu=dtset%natpawu
639  paw_dmft%natom=dtset%natom
640  paw_dmft%temp=dtset%tsmear!*unit_e
641  paw_dmft%dmft_iter=dtset%dmft_iter
642  paw_dmft%dmft_entropy=dtset%dmft_entropy
643  paw_dmft%dmft_kspectralfunc=dtset%dmft_kspectralfunc
644  paw_dmft%dmft_dc=dtset%dmft_dc
645  paw_dmft%dmft_wanorthnorm=dtset%dmft_wanorthnorm
646  !paw_dmft%idmftloop=0
647  paw_dmft%prtvol = dtset%prtvol
648  paw_dmft%prtdos = dtset%prtdos
649  paw_dmft%dmft_tolfreq = dtset%dmft_tolfreq
650  paw_dmft%dmft_lcpr = dtset%dmft_tollc
651  paw_dmft%dmft_charge_prec = dtset%dmft_charge_prec
652 
653 ! for entropy (alternate external calculation)
654  paw_dmft%ientropy  =  0
655  paw_dmft%u_for_s   =  4.1_dp
656  paw_dmft%j_for_s   =  0.5_dp
657 
658 !=======================
659 !==  Fixed self for input
660 !=======================
661  paw_dmft%use_fixed_self=dtset%usedmatpu
662  paw_dmft%fixed_self=>dmatpawu
663 
664 !=======================
665 !==  Choose solver
666 !=======================
667 
668  paw_dmft%dmft_solv=dtset%dmft_solv
669  paw_dmft%dmft_blockdiag=0
670  if(paw_dmft%dmft_solv==-2) then
671    paw_dmft%dmft_solv=2
672    paw_dmft%dmft_blockdiag=1
673  endif
674 !  0: DFT, no solver
675 !  1: DFT+U
676 ! -1: DFT+U but DFT values are not renormalized !
677 ! if((paw_dmft%dmft_solv==0.and.paw_dmft%prtvol>4).or.&
678 !&   (paw_dmft%dmft_solv>=-1.and.paw_dmft%dmft_solv<=2)) then
679 !   call wrtout(std_out,message,'COLL')
680 !   call wrtout(ab_out,message,'COLL')
681 ! endif
682 
683  if(paw_dmft%dmft_solv==0) then
684    do itypat=1,psps%ntypat
685      if(pawtab(itypat)%lpawu/=-1) then
686        if((pawtab(itypat)%upawu)>tol5.or.(pawtab(itypat)%jpawu)>tol5) then
687           write(message, '(2a,i5,2a,2e15.6)' )ch10,&
688 &          ' option dmft_solv=0 requires upaw=jpaw=0 for species',itypat,ch10,&
689 &          ' Value of upawu and jpawu are here',pawtab(itypat)%upawu,pawtab(itypat)%jpawu
690           ABI_ERROR(message)
691         endif
692      endif
693    enddo
694  endif
695 
696 ! todo_ab: why upaw and jpawu are not zero (on bigmac) if lpawu==-1 ?
697 ! if(paw_dmft%dmft_solv==0.and.&
698 !& (maxval(abs(pawtab(:)%upawu))>tol5.or.maxval(abs(pawtab(:)%jpawu))>tol5)) then
699 !   write(message, '(a,a,2f12.3)' )ch10,&
700 !&   ' option dmft_solv=0 requires upaw=jpaw=0',maxval(abs(pawtab(:)%upawu)),maxval(abs(pawtab(:)%jpawu))
701 !    ABI_WARNING(message)
702 ! endif
703 
704  paw_dmft%dmftcheck=dtset%dmftcheck
705 
706  if(paw_dmft%dmftcheck==-1) then
707    message = ' init_dmft: dmftcheck=-1 should not happend here'
708    ABI_BUG(message)
709  endif
710  paw_dmft%dmft_log_freq=1 ! use logarithmic frequencies.
711  if(paw_dmft%dmft_solv==6.or.paw_dmft%dmft_solv==7.or.paw_dmft%dmft_solv==9) then
712    paw_dmft%dmft_log_freq=0 ! do not use logarithmic frequencies.
713  endif
714  paw_dmft%dmft_nwli=dtset%dmft_nwli
715  if(paw_dmft%dmft_log_freq==1) then
716    paw_dmft%dmft_nwlo=dtset%dmft_nwlo
717  else
718    paw_dmft%dmft_nwlo=dtset%dmft_nwli
719  endif
720  paw_dmft%dmft_nwr=800
721  paw_dmft%dmft_rslf=dtset%dmft_rslf
722  paw_dmft%dmft_mxsf=dtset%dmft_mxsf
723  paw_dmft%dmftqmc_l=dtset%dmftqmc_l
724  paw_dmft%dmftqmc_n=dtset%dmftqmc_n
725  paw_dmft%dmftqmc_seed=dtset%dmftqmc_seed
726  paw_dmft%dmftqmc_therm=dtset%dmftqmc_therm
727  paw_dmft%dmftqmc_t2g=dtset%dmft_t2g
728 
729 !paw_dmft%dmftqmc_x2my2d=dtset%dmft_x2my2d
730  paw_dmft%dmftqmc_x2my2d=0
731 
732  paw_dmft%dmftctqmc_basis =dtset%dmftctqmc_basis
733  paw_dmft%dmftctqmc_check =dtset%dmftctqmc_check
734  paw_dmft%dmftctqmc_correl=dtset%dmftctqmc_correl
735  paw_dmft%dmftctqmc_gmove =dtset%dmftctqmc_gmove
736  paw_dmft%dmftctqmc_grnns =dtset%dmftctqmc_grnns
737  paw_dmft%dmftctqmc_meas  =dtset%dmftctqmc_meas
738  paw_dmft%dmftctqmc_mrka  =dtset%dmftctqmc_mrka
739  paw_dmft%dmftctqmc_mov   =dtset%dmftctqmc_mov
740  paw_dmft%dmftctqmc_order =dtset%dmftctqmc_order
741  paw_dmft%dmftctqmc_config =dtset%dmftctqmc_config
742  paw_dmft%dmftctqmc_triqs_nleg =dtset%dmftctqmc_triqs_nleg
743 
744  if ( paw_dmft%dmft_solv >= 4 ) then
745  write(message, '(a,a,i6)' )ch10,&
746 &   '=> Seed for QMC inside DMFT is dmftqmc_seed=',paw_dmft%dmftqmc_seed
747    call wrtout(std_out,message,'COLL')
748  endif
749 
750 
751 !=======================
752 !==  Variables for DMFT itself
753 !=======================
754 
755  mbandc = paw_dmft%mbandc
756 
757  ABI_MALLOC(paw_dmft%eigen_dft,(paw_dmft%nsppol,paw_dmft%nkpt,paw_dmft%mbandc))
758  paw_dmft%eigen_dft=zero
759 
760 ! allocate(paw_dmft%wtk(paw_dmft%nkpt))
761  paw_dmft%wtk=>dtset%wtk
762  if(dtset%iscf<0) then
763    paw_dmft%wtk=one/float(dtset%nkpt)
764  endif
765  sumwtk=0
766  do ikpt=1,paw_dmft%nkpt
767    sumwtk=sumwtk+paw_dmft%wtk(ikpt)
768  enddo
769  if(abs(sumwtk-1_dp)>tol11.and.dtset%iscf>=0) then
770    write(message, '(a,f15.11)' )' sum of k-point is incorrect',sumwtk
771    ABI_ERROR(message)
772  endif
773  ABI_MALLOC(paw_dmft%lpawu,(paw_dmft%natom))
774  do iatom=1,paw_dmft%natom
775    paw_dmft%lpawu(iatom)=pawtab(typat(iatom))%lpawu
776    if(paw_dmft%dmftqmc_t2g==1.and.paw_dmft%lpawu(iatom)==2) paw_dmft%lpawu(iatom)=1
777    if(paw_dmft%dmftqmc_x2my2d==1.and.paw_dmft%lpawu(iatom)==2) paw_dmft%lpawu(iatom)=0
778  enddo
779  paw_dmft%maxlpawu=maxval(paw_dmft%lpawu(:))
780  maxlpawu = paw_dmft%maxlpawu
781 
782  ABI_MALLOC(paw_dmft%psichi,(nsppol,dtset%nkpt,mbandc,nspinor,dtset%natom,(2*maxlpawu+1)))
783 
784  paw_dmft%psichi=cmplx(zero,zero,kind=dp)
785  paw_dmft%lpsichiortho=0
786 
787 
788 !=======================
789 ! Real      frequencies
790 !=======================
791 
792  iexist2=1
793  if(dtset%iscf<0.and.(paw_dmft%dmft_solv==5.or.paw_dmft%dmft_solv==8)) then
794      tmpfil = trim(paw_dmft%filapp)//'_spectralfunction_realfrequencygrid'
795      inquire(file=trim(tmpfil),exist=lexist)!,recl=nrecl)
796    !  write(6,*) "inquire",lexist
797    grid_unt=2000
798      if((.not.lexist)) then
799        iexist2=0
800        write(message,'(4x,a,i5,3a)') "File number",grid_unt,&
801 &       " called ",trim(tmpfil)," does not exist"
802        call wrtout(std_out,message,'COLL')
803        message = "Cannot continue: the missing file coming from Maxent code is needed"
804        ABI_WARNING(message)
805      endif
806 
807      if(iexist2==1) then
808 #ifdef FC_NAG
809        open (unit=grid_unt,file=trim(tmpfil),status='unknown',form='formatted',recl=ABI_RECL)
810 #else
811        open (unit=grid_unt,file=trim(tmpfil),status='unknown',form='formatted')
812 #endif
813        rewind(grid_unt)
814       ! if (open_file(tmpfil, message, newunit=grid_unt, status='unknown', form='formatted') /= 0) then
815       !   ABI_ERROR(message)
816       ! end if
817        write(message,'(3a)') ch10,"  == Read  grid frequency in file ",trim(tmpfil)
818        call wrtout(std_out,message,'COLL')
819        write(message,'(a,a,a,i4)') 'opened file : ', trim(tmpfil), ' unit', grid_unt
820        call wrtout(std_out,message,'COLL')
821        read(grid_unt,*) ngrid
822        ABI_MALLOC(paw_dmft%omega_r,(ngrid))
823        if(ioerr<0) then
824          message = "Error reading grid file"
825          ABI_ERROR(message)
826        endif
827        do ifreq=1,ngrid
828          read(grid_unt,*) paw_dmft%omega_r(ifreq)
829          paw_dmft%omega_r(ifreq)=paw_dmft%omega_r(ifreq)
830        enddo
831        if(ioerr<0) then
832          message = "Error reading grid file"
833          ABI_ERROR(message)
834        endif
835      endif
836  else
837   ABI_MALLOC(paw_dmft%omega_r,(2*paw_dmft%dmft_nwr))
838   ! Set up real frequencies for spectral function in Hubbard one.
839    step=0.00005_dp
840    paw_dmft%omega_r(2*paw_dmft%dmft_nwr)=pi*step*(two*float(paw_dmft%dmft_nwr-1)+one)
841    do ifreq=1,2*paw_dmft%dmft_nwr-1
842     paw_dmft%omega_r(ifreq)=pi*step*(two*float(ifreq-1)+one)-paw_dmft%omega_r(2*paw_dmft%dmft_nwr)
843   !  write(std_out,*) ifreq,paw_dmft%omega_r(ifreq)
844    enddo
845 
846  endif
847 !=======================
848 ! Imaginary frequencies
849 !=======================
850 ! Set up log frequencies
851  if(dtset%ucrpa==0) call construct_nwlo_dmft(paw_dmft)
852 
853 !=========================================================
854 !== if we use ctqmc impurity solver
855 !=========================================================
856 ! IMPORTANT : paw_dmft%hybrid is corrupted somewhere in DMFT routines on
857 ! tikal_psc and max2_open64. Use a local hybrid in qmc_prep even if not optimal.
858 ! Anyway Initializing ctqmc here is not good and produce the same result for
859 ! dmft_iter =1 which speed up the convergency ...
860 ! FIXME : Move this to init_sc_dmft and find bug
861  if(paw_dmft%dmft_solv==5) then ! CTQMC initialisation
862    write(message,'(a,2x,a,f13.5)') ch10,&
863 &  " == Initializing CTQMC"
864 !   call wrtout(std_out,message,'COLL')
865 
866    ABI_MALLOC(paw_dmft%hybrid,(paw_dmft%natom))
867    do iatom=1,paw_dmft%natom
868      if(paw_dmft%lpawu(iatom)/=-1) then
869        nflavor=2*(2*paw_dmft%lpawu(iatom)+1)
870 #ifdef HAVE_MPI
871        call CtqmcInterface_init(paw_dmft%hybrid(iatom),paw_dmft%dmftqmc_seed,paw_dmft%dmftqmc_n, &
872 &       paw_dmft%dmftqmc_therm, paw_dmft%dmftctqmc_meas,nflavor,paw_dmft%dmftqmc_l,one/paw_dmft%temp,zero,&
873 &       std_out,paw_dmft%spacecomm,nspinor=paw_dmft%nspinor)
874 #else
875        call CtqmcInterface_init(paw_dmft%hybrid(iatom),paw_dmft%dmftqmc_seed,paw_dmft%dmftqmc_n, &
876 &       paw_dmft%dmftqmc_therm, paw_dmft%dmftctqmc_meas,nflavor,paw_dmft%dmftqmc_l,one/paw_dmft%temp,zero,&
877 &       std_out,nspinor=paw_dmft%nspinor)
878 #endif
879        call CtqmcInterface_setOpts(paw_dmft%hybrid(iatom),&
880                                    opt_Fk      =1,&
881 &                                  opt_order   =paw_dmft%dmftctqmc_order ,&
882 &                                  opt_histo   =paw_dmft%dmftctqmc_config ,&
883 &                                  opt_movie   =paw_dmft%dmftctqmc_mov   ,&
884 &                                  opt_analysis=paw_dmft%dmftctqmc_correl,&
885 &                                  opt_check   =paw_dmft%dmftctqmc_check ,&
886 &                                  opt_noise   =paw_dmft%dmftctqmc_grnns ,&
887 &                                  opt_spectra =paw_dmft%dmftctqmc_mrka  ,&
888 &                                  opt_gmove   =paw_dmft%dmftctqmc_gmove )
889      end if
890    enddo
891    write(message,'(a,2x,a,f13.5)') ch10,&
892 &  " == Initialization CTQMC done"
893    !call wrtout(std_out,message,'COLL')
894  endif
895 
896  if(paw_dmft%dmftcheck==1.and.paw_dmft%dmft_solv<4) then
897    paw_dmft%dmftqmc_l=64
898  endif
899 
900 !************************************************************************
901 end subroutine init_dmft

m_paw_dmft/init_sc_dmft [ Functions ]

[ Top ] [ m_paw_dmft ] [ Functions ]

NAME

 init_sc_dmft

FUNCTION

  Allocate variables used in type paw_dmft_type.

INPUTS

 dmftbandi = lower bound for band states included in DMFT calculation
 dmftbandf = upper bound for band states included in DMFT calculation
 mband     = max number of bands
 nband     = number of bands for each k-point
 nkpt      = number of k-points
 nsppol    = number of spin polarisation
 occ       =  occupations
 usedmft  = ==1 if dmft is activated
 use_sc_dmft = for self-consistency in dmft

 OUTPUTS
 paw_dmft  = structure of data for dmft

SOURCE

362 subroutine init_sc_dmft(bandkss,dmftbandi,dmftbandf,dmft_read_occnd,mband,nband,nkpt,nspden,&
363 &nspinor,nsppol,occ,usedmft,paw_dmft,use_sc_dmft,dmft_solv,mpi_enreg)
364 
365 !Arguments ------------------------------------
366 !scalars
367  integer, intent(in) :: bandkss,dmft_read_occnd,dmftbandi,dmftbandf,mband,nkpt,nspden,&
368 & nspinor,nsppol,usedmft,use_sc_dmft,dmft_solv
369 !type
370  type(paw_dmft_type),intent(out) :: paw_dmft
371  type(MPI_type), intent(in) :: mpi_enreg
372 ! arrays
373  integer,intent(in) :: nband(nkpt*nsppol)
374  real(dp),intent(in) :: occ(mband*nkpt*nsppol)
375 !Local variables ------------------------------------
376  integer :: iband,icb,ikpt,isppol,nband_k,bdtot_index
377  integer :: myproc,nproc,spacecomm,use_dmft
378 ! integer :: ie,nb_procs
379  character(len=500) :: message
380 
381 !************************************************************************
382  use_dmft=abs(usedmft)
383 
384 
385 ! Check processors for DMFT
386 ! Initialise spaceComm, myproc, and nproc
387 
388  spacecomm=mpi_enreg%comm_cell
389  myproc=mpi_enreg%me_cell
390  nproc=mpi_enreg%nproc_cell
391  spacecomm=mpi_enreg%comm_world
392  myproc=mpi_enreg%me
393  nproc=mpi_enreg%nproc
394  !print *, " spacecomm,myproc,nproc",spacecomm,myproc,nproc
395  paw_dmft%spacecomm=spacecomm
396  paw_dmft%myproc=myproc
397  paw_dmft%nproc=nproc
398 
399  ! Do not comment these lines: it guarantees the parallelism in DMFT/HI or QMC will work.
400  if ((use_dmft/=0).and.(xmpi_comm_size(xmpi_world) /= xmpi_comm_size(mpi_enreg%comm_world))) then
401    ABI_ERROR("Someone changed the k-point parallelism again")
402  end if
403 
404  if(use_dmft/=0) then
405    write(message,'(2a,i3)') ch10,'-       ( number of procs used in dmft ) =', nproc
406    call wrtout(std_out,message,'COLL')
407    call wrtout(ab_out,message,'COLL')
408    write(std_out_default,'(2a,i3)') ch10,'       ( current proc is        ) =', myproc
409   ! write(ab_out_default,'(2a,i3)') ch10,'       ( current proc is        ) =', myproc
410    if(myproc==nproc-1) then
411      write(std_out_default,'(2a,i3)') ch10,'      ( last proc            ) =', myproc
412   !   write(ab_out_default,'(2a,i3)') ch10,'       ( last proc            ) =', myproc
413    endif
414  endif
415 !#ifdef HAVE_MPI
416 ! call MPI_COMM_SIZE(MPI_COMM_WORLD,nb_procs,ie)
417 ! write(6,*) "nprocs,nb_procs",nproc,nb_procs
418 ! if(nb_procs/=nproc)  then
419 !   message = ' Number of procs used in DMFT is erroneously computed '
420 !   ABI_ERROR(message)
421 ! endif
422 !#endif
423 
424  paw_dmft%mband       = mband
425  paw_dmft%dmftbandf   = dmftbandf
426  paw_dmft%dmftbandi   = dmftbandi
427  paw_dmft%nkpt        = nkpt
428 
429 !  Spin related variables and check
430  paw_dmft%nsppol      = nsppol
431  paw_dmft%nspinor     = nspinor
432  paw_dmft%nspden      = nspden
433  if(nspinor==2.and.nspden==1.and.use_dmft/=0) then
434    message = ' nspinor==2 and nspden =1 and usedmft=1 is not implemented'
435    ABI_ERROR(message)
436  endif
437 
438 ! if(nspinor==1.and.nspden==1.and.use_dmft/=0) then
439 !   message = ' nspinor==1 and nspden =1 and usedmft=1 is not implemented'
440 !   ABI_ERROR(message)
441 ! endif
442 
443  paw_dmft%use_dmft    = use_dmft
444  if (bandkss/=0) then
445    paw_dmft%use_sc_dmft = 0
446  else
447    paw_dmft%use_sc_dmft = use_sc_dmft
448  endif
449  paw_dmft%dmft_read_occnd = dmft_read_occnd
450  paw_dmft%idmftloop=0
451  paw_dmft%mbandc  = 0
452  ABI_MALLOC(paw_dmft%occnd,(2,mband,mband,nkpt,nsppol*use_dmft))
453  ABI_MALLOC(paw_dmft%band_in,(mband*use_dmft))
454  ABI_MALLOC(paw_dmft%include_bands,((dmftbandf-dmftbandi+1)*use_dmft))
455  ABI_MALLOC(paw_dmft%exclude_bands,(mband*use_dmft))
456 ! allocate(paw_dmft%ph0phiiint()
457  paw_dmft%band_in(:)=.false.
458  paw_dmft%occnd=zero
459  icb=0
460  if(use_dmft==1) then
461   do iband=1, mband
462    if(iband>=paw_dmft%dmftbandi.and.iband<=paw_dmft%dmftbandf) then
463     paw_dmft%band_in(iband)=.true.
464     paw_dmft%mbandc = paw_dmft%mbandc+1
465     paw_dmft%include_bands(paw_dmft%mbandc) = iband
466    else
467     icb=icb+1
468     paw_dmft%exclude_bands(icb)=iband
469    endif
470   enddo
471   bdtot_index=1
472   do isppol=1,nsppol
473    do ikpt=1,nkpt
474     nband_k=nband(ikpt+(isppol-1)*nkpt)
475     do iband=1,nband_k
476      paw_dmft%occnd(1,iband,iband,ikpt,isppol)=occ(bdtot_index)
477      bdtot_index=bdtot_index+1
478     end do
479    end do
480   end do
481  else
482   paw_dmft%mbandc = 0
483  endif
484 
485  if(paw_dmft%use_sc_dmft /= 0 .and. mpi_enreg%paral_kgb/=0) then
486    call init_sc_dmft_paralkgb(paw_dmft, mpi_enreg)
487  end if
488 
489  if(paw_dmft%use_dmft > 0 .and. paw_dmft%mbandc /= dmftbandf-dmftbandi+1) then
490   write(message, '(3a)' )&
491 &  ' WARNING init_sc_dmft',ch10,&
492 &  '  number of bands in dmft is not correctly computed ',ch10, &
493 &  '  Action : check the code'
494   ABI_WARNING(message)
495  endif
496  if(use_dmft>=1) then
497    write(message, '(7a)' ) ch10,ch10," ******************************************", &
498 &   ch10," DFT+DMFT Method is used", &
499 &   ch10," ******************************************"
500    call wrtout(std_out,  message,'COLL')
501    call wrtout(ab_out,  message,'COLL')
502  endif
503 
504  if(use_dmft>=1) then
505    if(dmft_solv==0) then
506      write(message, '(a,a)') ch10,' DMFT check: no solver and U=J=0'
507    else if(dmft_solv==1) then
508      write(message, '(a,a)') ch10,' DMFT check: static solver'
509    else if(dmft_solv==-1) then
510      write(message, '(a,a)') ch10,' DMFT check: static solver without renormalization of projectors: should recover DFT+U'
511    else if(dmft_solv==2) then
512      write(message, '(a,a)') ch10,' DMFT uses the Hubbard one solver'
513    else if(dmft_solv==4) then
514      write(message, '(a,a)') ch10,' DMFT uses the Hirsch Fye solver'
515    else if(dmft_solv==5) then
516      write(message, '(a,a)') ch10,' DMFT uses the Continuous Time Quantum Monte Carlo solver of ABINIT'
517    else if(dmft_solv==6) then
518      write(message, '(a,a)') ch10,' DMFT uses the Continuous Time Quantum Monte Carlo solver of TRIQS&
519      & (with density density interactions)'
520    else if(dmft_solv==7) then
521      write(message, '(a,a)') ch10,' DMFT uses the Continuous Time Quantum Monte Carlo solver of TRIQS&
522      & (with rotationaly invariant interaction)'
523    else if(dmft_solv==9) then
524      write(message, '(a,a)') ch10,' DMFT uses the python invocation of TRIQS, for which you need to &
525      &give your personal script'
526   endif
527   call wrtout(std_out,message,'COLL')
528   call wrtout(ab_out,message,'COLL')
529  endif
530 
531 end subroutine init_sc_dmft

m_paw_dmft/init_sc_dmft_paralkgb [ Functions ]

[ Top ] [ m_paw_dmft ] [ Functions ]

NAME

 init_sc_dmft_paralkgb

FUNCTION

  Init some values used with KGB parallelism in self consistent DMFT
  calculation.

INPUTS

  paw_dmft   = data structure

OUTPUT

  paw_dmft: bandc_proc, use_bandc

SOURCE

1643 subroutine init_sc_dmft_paralkgb(paw_dmft,mpi_enreg)
1644 
1645 !Arguments ------------------------------------
1646  type(paw_dmft_type),intent(inout) :: paw_dmft
1647  type(MPI_type), intent(in) :: mpi_enreg
1648 
1649 !Local variables-------------------------------
1650 !scalars
1651  integer :: nproc, ib, ibc, proc
1652 
1653 ! *********************************************************************
1654  nproc = mpi_enreg%nproc_band
1655 
1656  ABI_MALLOC(paw_dmft%bandc_proc,(paw_dmft%mbandc))
1657  ABI_MALLOC(paw_dmft%use_bandc,(nproc))
1658  paw_dmft%bandc_proc = 0
1659  paw_dmft%use_bandc = .false.
1660 
1661  do ibc=1,paw_dmft%mbandc
1662    ib = paw_dmft%include_bands(ibc)
1663    proc = mod((ib-1)/mpi_enreg%bandpp,nproc)
1664 
1665    paw_dmft%bandc_proc(ibc) = proc
1666 
1667    paw_dmft%use_bandc(proc+1) = .true.
1668  end do
1669 end subroutine init_sc_dmft_paralkgb

m_paw_dmft/paw_dmft_type [ Types ]

[ Top ] [ m_paw_dmft ] [ Types ]

NAME

  paw_dmft_type

FUNCTION

  This structured datatype contains the necessary data for the link
  between dmft and paw.
  occnd(non-diagonal band occupations for self-consistency), band_in
  (say which band are taken into account in the calculation), and the

SOURCE

 72  type, public :: paw_dmft_type
 73 
 74   integer :: dmft_dc
 75   ! Type of double counting used in DMFT
 76 
 77   integer :: dmft_iter
 78   ! Nb of iterations for dmft
 79 
 80   integer :: dmft_kspectralfunc
 81   ! =0 Default
 82   ! =1 Activate calculation of k-resolved spectral function
 83 
 84   integer :: dmft_solv
 85   ! choice of solver for DMFT
 86 
 87   integer :: dmftcheck
 88   ! Check various part of the implementation
 89 
 90   integer :: dmft_entropy
 91   ! = 0: do not compute entropy
 92   ! = 1: compute entropy
 93 
 94   integer :: dmft_log_freq
 95   ! = 0: do not use log frequencies
 96   ! = 1: use log frequencies
 97 
 98 !  integer :: dmft_mag
 99 !  ! 0 if non magnetic calculation, 1 if magnetic calculation
100 
101   integer :: dmft_nwlo
102   ! dmft frequencies
103 
104   integer :: dmft_nwr
105   ! dmft frequencies
106 
107   integer :: dmft_nwli
108   ! dmft frequencies
109 
110   integer :: dmftqmc_l
111   ! qmc related input
112 
113   integer :: dmftqmc_seed
114   ! qmc seed
115 
116   integer :: dmftqmc_therm
117   ! qmc thermalization
118 
119   integer :: dmftctqmc_basis
120   ! CTQMC: Basis to perform the CTQMC calculation
121   ! for historical reasons and tests
122   ! 0 : Slm basis, 1 : diagonalise local Hamiltonian, 2: diago density matrix
123   !
124   integer :: dmft_blockdiag
125   !  Block diagonalize Hamiltonian in the local basis
126 
127   integer :: dmftctqmc_check
128   ! CTQMC: perform a check on the impurity and/or bath operator
129   ! only for debug
130   ! 0 : nothing, 1 : impurity, 2 : bath, 3 : both
131 
132   integer :: dmftctqmc_correl
133   ! CTQMC: Gives analysis for CTQMC
134   ! 0 : nothing, 1 : activated Correlations.dat
135 
136   integer :: dmftctqmc_gmove
137   ! CTQMC: add global move every dmftctqmc_gmove sweeps
138   ! >= 0 ; warning inside CT-QMC with warning
139   ! == 0 ; no global moves
140 
141   integer :: dmftctqmc_grnns
142   ! CTQMC: compute green function noise for each imaginary time
143   ! 0 : nothing, 1 : activated
144 
145   integer :: dmftctqmc_config
146   ! CTQMC: Enables histogram of occupations.
147   ! 0 : nothing, 1 enabled
148 
149   integer :: dmftctqmc_meas
150   ! CTQMC: every each dmftctqmc_meas step energy is measured
151   ! Speed up caltion about 10% for dmftctqmc_meas = 2
152 
153   integer :: dmftctqmc_mrka
154   ! CTQMC: Write a temporary file Spectra_RANK.dat with the sweep evolution of
155   ! the number of electron for each flavor
156   ! The measurement is done every dmftctqmc_meas*dmftctqmc_mrka sweep
157   ! e.g. : meas=2 mrka=10 -> every 20 sweeps sum_i c+(ti)c(t'i) is mesured
158 
159   integer :: dmftctqmc_mov
160   ! CTQMC: Gives movie for CTQMC
161   ! 0 : nothing, 1 : 1 file Movie_RANK.tex for each cpu
162 
163   integer :: dmftctqmc_order
164   ! CTQMC: Gives order in perturbation for CTQMC solver
165   ! 0 : nothing, >=1 max order evaluated in Perturbation.dat
166 
167   integer :: dmftctqmc_triqs_nleg
168   ! CTQMC of TRIQS: Nb of Legendre polynomial used to compute the
169   ! Green's function (Phys. Rev. B 84, 075145) [[cite:Boehnke2011]]. Default is 30.
170 
171   ! 0 : nothing, >=1 max order evaluated in Perturbation.dat
172 
173   real(dp) :: dmftqmc_n
174   ! qmc number of sweeps
175 
176   integer :: dmftqmc_x2my2d
177   ! for doing qmc with x2my2d only (for testing purposes)
178 
179   integer :: dmftqmc_t2g
180   ! for doing qmc with t2g only (for testing purposes)
181 
182   integer :: dmftbandi
183   ! Number of bands
184 
185   integer :: dmftbandf
186   ! Number of bands
187 
188   integer :: dmft_read_occnd
189   ! Number of bands
190 
191   integer :: dmft_rslf
192   ! Number of bands
193 
194   integer :: dmft_prgn
195   ! Precise the way of printing the green function.
196   !  =1   print green
197   !  =2   print self
198 
199   integer :: dmft_wanorthnorm
200   !  1 orthonormalisation of Wannier functions k-point after k-point
201   !  2 orthonormalisation over the sum over k-points.
202 
203   integer :: idmftloop
204   ! current iteration in the dmft loop
205 
206   integer :: maxlpawu         ! Number of correlated atoms
207 
208 
209   integer :: mband
210   ! Number of bands
211 
212   integer :: mbandc
213   ! Total number of bands in the Kohn-Sham Basis for PAW+DMFT
214 
215   integer :: natom
216   ! Number of atom
217 
218   integer :: natpawu         ! Number of correlated atoms
219 
220   integer :: nkpt
221   ! Number of k-point in the IBZ.
222 
223   integer :: nspden
224 
225   integer :: nspinor
226 
227   integer :: nsppol
228 
229   integer :: prtdos
230 
231   integer :: prtvol
232 
233   integer  :: lpsichiortho
234 
235   integer  :: use_fixed_self
236 
237   integer :: ientropy
238   ! activate evaluation of terms for alternative calculation of entropy in DMFT
239 
240   real(dp) :: edmft
241 
242   real(dp) :: dmft_charge_prec
243   ! Precision on charge required for determination of fermi level (fermi_green) with newton method
244 
245   real(dp) :: dmft_fermi_prec
246   ! Required precision on Fermi level (fermi_green) during the DMFT SCF cycle, (=> ifermie_cv)
247   ! used also for self (new_self)  (=> iself_cv).
248 
249   real(dp) :: dmft_mxsf
250   ! Mixing coefficient for Self-Energy during the SCF DMFT cycle.
251 
252   real(dp) :: dmft_tolfreq
253   ! Required precision on local correlated density matrix  (depends on
254   ! frequency mesh), used in m_dmft/dmft_solve
255 
256   real(dp) :: dmft_lcpr
257   ! Required precision on local correlated charge  in order to stop SCF
258   ! DMFT cycle (integrate_green) => ichargeloc_cv
259 
260   real(dp) :: fermie
261 
262   real(dp) :: u_for_s
263   ! Variable for evaluation of correlation energy for U=0 in the entropic
264   ! calculation
265 
266   real(dp) :: j_for_s
267   ! Variable for evaluation of correlation energy for U=0 in the entropic
268   ! calculation
269 
270 
271   real(dp) :: fermie_dft
272 
273   real(dp) :: nelectval
274 
275   character(len=fnlen) :: filapp
276 
277   character(len=fnlen) :: filnamei
278 
279   real(dp) :: temp
280 
281   integer, allocatable :: lpawu(:)
282 
283   integer, allocatable :: include_bands(:)
284   ! for each bands included in the calculation (1..mbandc), include_bands
285   ! gives the index in the full band index  (1...mband)
286 
287   integer, allocatable :: exclude_bands(:)
288   ! gives the bands than are not in the DMFT calculations.
289 
290   real(dp), allocatable :: occnd(:,:,:,:,:)
291   ! non diagonal band-occupation for each k-point, polarisation.
292 
293 !  real(dp), allocatable :: phi0phiiint(:)
294 !  ! non diagonal band-occupation for each k-point, polarisation.
295 
296   logical, allocatable :: band_in(:)
297   ! true for each band included in the calculation.
298 
299   integer, allocatable :: bandc_proc(:)
300   ! proc index (on comm_band) for each correlated band in DMFT
301 
302   logical, allocatable :: use_bandc(:)
303   ! true for each proc wich has at least one band involved in DMFT non diagonal
304   ! occupations on band parallelism
305 
306   integer :: use_dmft
307   ! 1 if non diagonal occupations are used, else 0
308 
309   integer :: use_sc_dmft
310   ! 1 if calculations have to be carried out self-consistently in the
311   ! electronic density.
312 
313   complex(dpc), allocatable :: psichi(:,:,:,:,:,:)
314 
315   real(dp), allocatable :: eigen_dft(:,:,:)
316 
317   real(dp), pointer :: wtk(:) => null()
318   real(dp), pointer :: fixed_self(:,:,:,:) => null()
319   real(dp), pointer :: omega_lo(:) => null()
320   real(dp), pointer :: omega_r(:) => null()
321   real(dp), pointer :: wgt_wlo(:) => null()
322 
323   type(CtqmcInterface), allocatable :: hybrid(:)
324   type(data4entropyDMFT_t) :: forentropyDMFT
325 
326   ! MPI realated variables
327   integer :: myproc
328   integer :: nproc
329   integer :: spacecomm
330 
331  end type paw_dmft_type

m_paw_dmft/print_dmft [ Functions ]

[ Top ] [ m_paw_dmft ] [ Functions ]

NAME

 print_dmft

FUNCTION

INPUTS

OUTPUT

SOURCE

1364 subroutine print_dmft(paw_dmft,pawprtvol)
1365 
1366 !Arguments ------------------------------------
1367 !type
1368  type(paw_dmft_type),intent(in) :: paw_dmft
1369  integer :: pawprtvol
1370 
1371 !Local variables-------------------------------
1372  integer :: ikpt,iband,ifreq,isppol
1373  character(len=500) :: message
1374 ! *********************************************************************
1375 
1376  if( abs(pawprtvol) >= 3 )  then
1377   write(message,'(4a,3(a,2x,e21.14,a))') &
1378 &   "  -------------------------------------------------",ch10,&
1379 &   "  --- Data for DMFT ",ch10,&
1380 &   "  --- paw_dmft%fermie     = ",paw_dmft%fermie    ,ch10,&
1381 &   "  --- paw_dmft%fermie_dft = ",paw_dmft%fermie_dft,ch10,&
1382 &   "  --- paw_dmft%temp       = ",paw_dmft%temp      ,ch10
1383   call wrtout(std_out,message,'COLL')
1384   write(message,'(7(a,15x,i8,a),a,2x,e21.14,2a)') &
1385 &   "  --- paw_dmft%natpawu    = ",paw_dmft%natpawu   ,ch10,&
1386 &   "  --- paw_dmft%dmft_iter  = ",paw_dmft%dmft_iter ,ch10,&
1387 &   "  --- paw_dmft%dmft_solv  = ",paw_dmft%dmft_solv ,ch10,&
1388 &   "  --- paw_dmft%dmft_nwlo  = ",paw_dmft%dmft_nwlo ,ch10,&
1389 &   "  --- paw_dmft%dmft_nwli  = ",paw_dmft%dmft_nwli ,ch10,&
1390 &   "  --- paw_dmft%dmft_dc    = ",paw_dmft%dmft_dc   ,ch10,&
1391 &   "  --- paw_dmft%dmftqmc_l  = ",paw_dmft%dmftqmc_l ,ch10,&
1392 &   "  --- paw_dmft%dmftqmc_n  = ",paw_dmft%dmftqmc_n ,ch10,&
1393 &   "  -------------------------------------------------"
1394   call wrtout(std_out,message,'COLL')
1395 
1396 !  write(message,'(4a,3(a,2x,f8.3,a),8(a,2x,i8,a),a)') "-----------------------------------------------",ch10,&
1397 !&   "--- Data for DMFT ",ch10,&
1398 !&   "--- paw_dmft%fermie     = ",paw_dmft%fermie    ,ch10,&
1399 !&   "--- paw_dmft%fermie_dft = ",paw_dmft%fermie_dft,ch10,&
1400 !&   "--- paw_dmft%temp       = ",paw_dmft%temp      ,ch10,&
1401 !&   "--- paw_dmft%natpawu    = ",paw_dmft%natpawu   ,ch10,&
1402 !&   "--- paw_dmft%dmft_iter  = ",paw_dmft%dmft_iter ,ch10,&
1403 !&   "--- paw_dmft%dmft_solv  = ",paw_dmft%dmft_solv ,ch10,&
1404 !&   "--- paw_dmft%dmft_nwlo  = ",paw_dmft%dmft_nwlo ,ch10,&
1405 !&   "--- paw_dmft%dmft_nwli  = ",paw_dmft%dmft_nwli ,ch10,&
1406 !&   "--- paw_dmft%dmft_dc    = ",paw_dmft%dmft_dc   ,ch10,&
1407 !&   "--- paw_dmft%dmftqmc_l  = ",paw_dmft%dmftqmc_l ,ch10,&
1408 !&   "--- paw_dmft%dmftqmc_n  = ",paw_dmft%dmftqmc_n ,ch10,&
1409 !&   "-----------------------------------------------"
1410   if(abs(pawprtvol)>10) then
1411    call wrtout(std_out,message,'COLL')
1412    write(message, '(a)') " DFT Eigenvalues "
1413    do isppol=1,paw_dmft%nsppol
1414     write(message, '(a,i4)') "--isppol--",isppol
1415     call wrtout(std_out,message,'COLL')
1416     do ikpt=1,paw_dmft%nkpt
1417      write(message, '(a,i4,2x,f14.5,a)') "  -k-pt--",ikpt,paw_dmft%wtk(ikpt),"(<-weight(k-pt))"
1418 
1419      call wrtout(std_out,message,'COLL')
1420      do iband=1,paw_dmft%mbandc
1421       write(message, '(a,i4,f10.5)') "   -iband--",iband,paw_dmft%eigen_dft(isppol,ikpt,iband)
1422       call wrtout(std_out,message,'COLL')
1423      enddo
1424     enddo
1425    enddo
1426    write(message, '(3x,a)') "Log. Freq"
1427    call wrtout(std_out,message,'COLL')
1428    do ifreq=1,paw_dmft%dmft_nwlo
1429     write(message, '(3x,a,i4,2(2x,e13.5))') "--ifreq--",ifreq,paw_dmft%omega_lo(ifreq),paw_dmft%wgt_wlo(ifreq)
1430     call wrtout(std_out,message,'COLL')
1431    enddo
1432   endif
1433  endif
1434 
1435 end subroutine print_dmft

m_paw_dmft/print_sc_dmft [ Functions ]

[ Top ] [ m_paw_dmft ] [ Functions ]

NAME

 print_sc_dmft

FUNCTION

INPUTS

OUTPUT

SOURCE

1450 subroutine print_sc_dmft(paw_dmft,pawprtvol)
1451 
1452 !Arguments ------------------------------------
1453 !type
1454  type(paw_dmft_type),intent(in) :: paw_dmft
1455  integer :: pawprtvol
1456 
1457 !Local variables-------------------------------
1458  integer :: iband
1459  character(len=500) :: message
1460 ! *********************************************************************
1461 
1462  if( abs(pawprtvol) >= 3 )  then
1463    write(message,'(5a,8(a,2x,i5,a),a)')ch10,"-----------------------------------------------",ch10,&
1464 &    "--- Data for SC DMFT ",ch10,&
1465 &    "--- paw_dmft%mband       = ",paw_dmft%mband,ch10,&
1466 &    "--- paw_dmft%dmftbandf   = ",paw_dmft%dmftbandf,ch10,&
1467 &    "--- paw_dmft%dmftbandi   = ",paw_dmft%dmftbandi,ch10,&
1468 &    "--- paw_dmft%nkpt        = ",paw_dmft%nkpt,ch10,&
1469 &    "--- paw_dmft%nsppol      = ",paw_dmft%nsppol,ch10,&
1470 &    "--- paw_dmft%use_dmft    = ",paw_dmft%use_dmft,ch10,&
1471 &    "--- paw_dmft%use_sc_dmft = ",paw_dmft%use_sc_dmft,ch10,&
1472 &    "--- paw_dmft%mbandc      = ",paw_dmft%mbandc,ch10,&
1473 &    "-----------------------------------------------"
1474    call wrtout(std_out,message,'COLL')
1475    write(message, '(a)') " paw_dmft%band_in"
1476    call wrtout(std_out,message,'COLL')
1477    write(message, '(100i5)') (iband,iband=1,min(paw_dmft%mband,100))
1478    call wrtout(std_out,message,'COLL')
1479    write(message, '(100L3)') (paw_dmft%band_in(iband),iband=1,min(paw_dmft%mband,100))
1480    call wrtout(std_out,message,'COLL')
1481    do iband=1,paw_dmft%mbandc
1482      write(message,*) "include_bands",iband,paw_dmft%include_bands(iband)
1483      call wrtout(std_out,message,'COLL')
1484    enddo
1485    write(message, '(a,a,i4,a)' )ch10,&
1486 &    'The',paw_dmft%mband-paw_dmft%dmftbandf+paw_dmft%dmftbandi-1,&
1487 &    '  Following bands are excluded from the DMFT calculation  '
1488    call wrtout(std_out,message,'COLL')
1489    write(message,'(100i5)') (paw_dmft%exclude_bands(iband),iband=1,min(paw_dmft%mband-paw_dmft%dmftbandf+paw_dmft%dmftbandi-1,100))
1490    call wrtout(std_out,message,'COLL')
1491  endif
1492 
1493 end subroutine print_sc_dmft

m_paw_dmft/readocc_dmft [ Functions ]

[ Top ] [ m_paw_dmft ] [ Functions ]

NAME

 readocc_dmft

FUNCTION

  read occnd on disk

INPUTS

  paw_dmft   = data structure
  filnam_ds3 = root for filname to read (input)
  filnam_ds4 = root for filname to read (output)

OUTPUT

  paw_dmft: occnd

SOURCE

1567 subroutine readocc_dmft(paw_dmft,filnam_ds3,filnam_ds4)
1568 
1569 !Arguments ------------------------------------
1570  type(paw_dmft_type),intent(inout) :: paw_dmft
1571  character(len=fnlen) :: filnam_ds3
1572  character(len=fnlen) :: filnam_ds4
1573 
1574 !Local variables-------------------------------
1575 !scalars
1576  character(len=fnlen) :: tmpfil
1577  integer :: ib,ib1,ikpt,is,unitsaveocc,dum1,dum2,dum3,dum4,ioerr
1578  logical :: lexist
1579  character(len=500) :: message
1580  character(len=4) :: chtemp
1581 ! *********************************************************************
1582  if(paw_dmft%dmft_read_occnd==0) return
1583  if(paw_dmft%dmft_read_occnd==1) tmpfil=trim(filnam_ds3)//'_DMFTOCCND'
1584  if(paw_dmft%dmft_read_occnd==2) tmpfil=trim(filnam_ds4)//'_DMFTOCCND'
1585  inquire(file=trim(tmpfil),exist=lexist)!,recl=nrecl)
1586  unitsaveocc=679
1587  if (lexist) then
1588    if (open_file(tmpfil,message,unit=unitsaveocc,status='unknown',form='formatted') /= 0) then
1589      ABI_ERROR(message)
1590    end if
1591    rewind(unitsaveocc)
1592    write(message,'(3a)') ch10,"  == Read DMFT non diagonal occupations on disk"
1593    call wrtout(std_out,message,'COLL')
1594    read(unitsaveocc,*)
1595    read(unitsaveocc,*,iostat=ioerr)&
1596 &              chtemp,dum1,dum2,dum3,dum4
1597    if(ioerr<0) then
1598      write(std_out,*) "read",dum1,dum2,dum3,dum4
1599    endif
1600    write(message,'(2a,4i4)') ch10,"  == natom, nsppol, nbandc, nkpt  read are",dum1,dum2,dum3,dum4
1601    call wrtout(std_out,message,'COLL')
1602    do is = 1 , paw_dmft%nsppol
1603      do ikpt = 1, paw_dmft%nkpt
1604        do ib = 1, paw_dmft%mbandc
1605          do ib1 = 1, paw_dmft%mbandc
1606            read(unitsaveocc,*) dum1,dum2,dum3,dum4,&
1607 &           paw_dmft%occnd(1,paw_dmft%include_bands(ib),paw_dmft%include_bands(ib1),ikpt,is),&
1608 &           paw_dmft%occnd(2,paw_dmft%include_bands(ib),paw_dmft%include_bands(ib1),ikpt,is)
1609          enddo
1610        enddo
1611      enddo
1612    enddo
1613 !   write(read,'(3a)') "# end of record",ch10&
1614 !&                ,"####  1234 "
1615 !   call wrtout(unitsaveocc,message,'COLL')
1616  else
1617    write(message,'(2a,2x,2a)') ch10,"   File",trim(tmpfil),"is not available"
1618    call wrtout(std_out,message,'COLL')
1619    write(message,'(4a)') ch10,"  ==> DMFT Occupations not available for restart", ch10, &
1620 &   "      -> The calculation is started with Fermi Dirac scheme for occupations"
1621    call wrtout(std_out,message,'COLL')
1622  endif
1623 
1624 end subroutine readocc_dmft

m_paw_dmft/saveocc_dmft [ Functions ]

[ Top ] [ m_paw_dmft ] [ Functions ]

NAME

 saveocc_dmft

FUNCTION

  save occnd on disk

INPUTS

  paw_dmft

OUTPUT

SOURCE

1510 subroutine saveocc_dmft(paw_dmft)
1511 
1512 !Arguments ------------------------------------
1513  type(paw_dmft_type),intent(inout) :: paw_dmft
1514 
1515 !Local variables-------------------------------
1516 !scalars
1517  character(len=fnlen) :: tmpfil
1518  integer :: ib,ib1,ikpt,is,unitsaveocc
1519  character(len=500) :: message
1520 ! *********************************************************************
1521  tmpfil = trim(paw_dmft%filapp)//'_DMFTOCCND'
1522  if (open_file(tmpfil,message,newunit=unitsaveocc,status='unknown',form='formatted') /= 0) then
1523    ABI_ERROR(message)
1524  end if
1525 
1526  rewind(unitsaveocc)
1527  write(message,'(2a)') ch10,"  == Print DMFT non diagonal occupations on disk"
1528  call wrtout(std_out,message,'COLL')
1529  write(message,'(3a,2x,4i5)') "# natom,nsppol,mbandc,nkpt",ch10&
1530 &              ,"####",paw_dmft%natom,paw_dmft%nsppol,paw_dmft%mbandc,paw_dmft%nkpt
1531  call wrtout(unitsaveocc,message,'COLL')
1532  do is = 1 , paw_dmft%nsppol
1533    do ikpt = 1, paw_dmft%nkpt
1534      do ib = 1, paw_dmft%mbandc
1535        do ib1 = 1, paw_dmft%mbandc
1536          write(unitsaveocc,*) is,ikpt,ib,ib1,paw_dmft%occnd(1,paw_dmft%include_bands(ib),paw_dmft%include_bands(ib1),ikpt,is),&
1537 &         paw_dmft%occnd(2,paw_dmft%include_bands(ib),paw_dmft%include_bands(ib1),ikpt,is)
1538        enddo
1539      enddo
1540    enddo
1541  enddo
1542  write(message,'(3a)') "# end of record",ch10&
1543 &              ,"####  1234 "
1544  call wrtout(unitsaveocc,message,'COLL')
1545  close(unitsaveocc)
1546 
1547 end subroutine saveocc_dmft