TABLE OF CONTENTS
ABINIT/cchi0 [ Functions ]
NAME
cchi0
FUNCTION
Main calculation of the independent-particle susceptibility chi0 for qpoint != 0
INPUTS
use_tr=If .TRUE. valence states are stored in Wfs_val and only resonant transitions are calculated (time reversal is assumed) Dtset <type(dataset_type)>=all input variables in this dataset Cryst<crystal_t>= data type gathering info on symmetries and unit cell %natom=number of atoms %nsym=number of symmetries %xred(3,natom)=reduced coordinated of atoms %typat(natom)=type of each atom %rprimd(3,3)=dimensional primitive translations in real space (bohr) %timrev= 2 if time reversal can be used, 1 otherwise qpoint(3)=reciprocal space coordinates of the q wavevector Ep<type(em1params_t_type)>= Parameters related to the calculation of the inverse dielectric matrix. %nbnds=number of bands summed over %npwe=number of planewaves for the irreducible polarizability X^0_GGp %npwvec=maximum number of G vectors used to define the dimension of some arrays e.g igfft %nsppol=1 for unpolarized, 2 for spin-polarized %nomega=total number of frequencies in X^0 (both real and imaginary) %nomegasf=number of real frequencies used to sample the imaginary part of X^0 (spectral method) %spmeth=1 if we use the spectral method, 0 for standard Adler-Wiser expression %spsmear=gaussian broadening used to approximate the delta distribution %zcut=small imaginary shift to avoid poles in X^0 Psps <type(pseudopotential_type)>=variables related to pseudopotentials Kmesh <kmesh_t>= datatype gathering parameters related to the k-point sampling %nibz=number of k-points in the IBZ %nbz=number of k-points in the BZ %bz(3,nbz)=reduced coordinates for k-points in the full Brillouin zone %ibz(3,nibz)=reduced coordinates for k-points in the irreducible wedge %tab(nbz)=mapping between a kpt in the BZ (array bz) and the irred point in the array ibz %tabi(nbz)= -1 if inversion is needed to obtain this particular kpt in the BZ, 1 means identity %tabo(nbz)= for each point in the BZ, the index of the symmetry operation S in reciprocal space which rotates k_IBZ onto \pm k_BZ (depending on tabi) %tabp(nbz)= For each k_BZ, it gives the phase factors associated to non-symmorphic operations, i.e e^{-i 2 \pi k_IBZ \cdot R{^-1}t} == e{-i 2\pi k_BZ cdot t} where : \transpose R{-1}=S and (S k_IBZ) = \pm k_BZ (depending on ktabi) %tabr(nfftot,nbz) For each point r on the real mesh and for each k-point in the BZ, tabr gives the index of (R^-1 (r-t)) in the FFT array where R=\transpose S^{-1} and k_BZ=S k_IBZ. t is the fractional translation associated to R Gsph_epsG0<gsphere_t data type> The G-sphere used to describe chi0/eps. (including umklapp G0 vectors) %ng=number of G vectors for chi0 %rottbm1(ng,2,nsym)=contains the index (IS^{-1}) G %phmGt(Ep%npwe,nsym)=phase factors e^{-iG \cdot t} needed to symmetrize oscillator matrix elements and epsilon %gprimd(3,3)=dimensional reciprocal space primitive translations (b^-1) %gmet(3,3)=reciprocal space metric ($\textrm{bohr}^{-2}$). nbvw=number of bands in the arrays wfrv ngfft_gw(18)= array containing all the information for 3D FFT for the oscillator strengths (see input variable) nfftot_gw=Total number of points in the GW FFT grid Ltg_q<Little group>=Data type gathering information on the little group of the q-points. Pawtab(Psps%ntypat) <type(pawtab_type)>=paw tabulated starting data Pawang<pawang_type> angular mesh discretization and related data: qp_ebands<ebands_t>=Quasiparticle energies and occupations (for the moment real quantities) %mband=MAX number of bands over k-points and spin (==Ep%nbnds) %occ(mband,nkpt,nsppol)=QP occupation numbers, for each k point in IBZ, and each band %eig(mband,nkpt,nsppol)=GW energies, for self-consistency purposes Paw_pwff<pawpwff_t>=Form factor used to calculate the onsite mat. elements of a plane wave. Wfd<wfdgw_t>=Object used to access the wavefunctions
OUTPUT
chi0(Ep%npwe,Ep%npwe,Ep%nomega)=independent-particle susceptibility matrix at wavevector qpoint and each frequeny defined by Ep%omega and Ep%nomega
SOURCE
1136 subroutine cchi0(use_tr,Dtset,Cryst,qpoint,Ep,Psps,Kmesh,qp_ebands,Gsph_epsG0,& 1137 Pawtab,Pawang,Paw_pwff,Pawfgrtab,Paw_onsite,nbvw,ngfft_gw,nfftot_gw,ngfftf,nfftf_tot,& 1138 chi0,ktabr,ktabrf,Ltg_q,chi0_sumrule,Wfd,Wfdf,wan) 1139 1140 !Arguments ------------------------------------ 1141 !scalars 1142 integer,intent(in) :: nbvw,nfftot_gw,nfftf_tot 1143 logical,intent(in) :: use_tr 1144 type(ebands_t),target,intent(in) :: qp_ebands 1145 type(kmesh_t),intent(in) :: Kmesh 1146 type(crystal_t),intent(in) :: Cryst 1147 type(Dataset_type),intent(in) :: Dtset 1148 type(em1params_t),intent(in) :: Ep 1149 type(gsphere_t),intent(in) :: Gsph_epsG0 1150 type(littlegroup_t),intent(in) :: Ltg_q 1151 type(Pawang_type),intent(in) :: Pawang 1152 type(Pseudopotential_type),intent(in) :: Psps 1153 type(wfdgw_t),target,intent(inout) :: Wfd,Wfdf 1154 !arrays 1155 integer,intent(in) :: ktabr(nfftot_gw,Kmesh%nbz),ktabrf(nfftf_tot*Dtset%pawcross,Kmesh%nbz) 1156 integer,intent(in) :: ngfft_gw(18),ngfftf(18) 1157 real(dp),intent(in) :: qpoint(3) 1158 real(dp),intent(out) :: chi0_sumrule(Ep%npwe) 1159 complex(gwpc),intent(out) :: chi0(Ep%npwe*Ep%nI,Ep%npwe*Ep%nJ,Ep%nomega) 1160 type(Pawtab_type),intent(in) :: Pawtab(Psps%ntypat*Psps%usepaw) 1161 type(pawpwff_t),intent(in) :: Paw_pwff(Psps%ntypat*Psps%usepaw) 1162 type(pawfgrtab_type),intent(inout) :: Pawfgrtab(Cryst%natom*Psps%usepaw) 1163 type(plowannier_type),intent(inout) :: wan 1164 type(paw_pwaves_lmn_t),intent(in) :: Paw_onsite(Cryst%natom) 1165 1166 !Local variables ------------------------------ 1167 !scalars 1168 integer,parameter :: tim_fourdp1 = 1, two_poles = 2, one_pole = 1, ndat1 = 1 1169 integer :: bandinf,bandsup,dim_rtwg,band1,band2,ierr,band1c,band2c 1170 integer :: ig1,ig2,iat1,iat2,iat,ik_bz,ik_ibz,ikmq_bz,ikmq_ibz 1171 integer :: io,iomegal,iomegar,ispinor1,ispinor2,isym_k,itypatcor,nfft,il1,il2 1172 integer :: isym_kmq,itim_k,itim_kmq,m1,m2,my_wl,my_wr,size_chi0 1173 integer :: nfound,nkpt_summed,nspinor,nsppol,mband 1174 integer :: comm,gw_mgfft,use_padfft,gw_fftalga,lcor,mgfftf,use_padfftf 1175 integer :: my_nbbp,my_nbbpks,spin,nbmax,dummy 1176 real(dp) :: cpu_time,wall_time,gflops 1177 real(dp) :: deltaeGW_b1kmq_b2k,deltaeGW_enhigh_b2k,deltaf_b1kmq_b2k 1178 real(dp) :: e_b1_kmq,en_high,fac,fac2,fac3,f_b1_kmq,factor,max_rest,min_rest,my_max_rest 1179 real(dp) :: my_min_rest,numerator,spin_fact,weight,wl,wr 1180 real(dp) :: gw_gsq,memreq 1181 complex(dpc) :: ph_mkmqt,ph_mkt 1182 complex(gwpc) :: local_czero_gw 1183 logical :: qzero,isirred_k,isirred_kmq,luwindow,is_metallic 1184 character(len=500) :: msg,allup 1185 type(gsphere_t) :: Gsph_FFT 1186 !arrays 1187 integer :: G0(3),umklp_k(3),umklp_kmq(3) 1188 integer :: ucrpa_bands(2) 1189 integer :: wtk_ltg(Kmesh%nbz),got(Wfd%nproc) 1190 integer,allocatable :: tabr_k(:),tabr_kmq(:),tabrf_k(:),tabrf_kmq(:) 1191 integer,allocatable :: igfftepsG0(:),gspfft_igfft(:),igfftepsG0f(:) 1192 integer,allocatable :: gw_gfft(:,:),gw_gbound(:,:),dummy_gbound(:,:),gboundf(:,:) 1193 integer,allocatable :: bbp_ks_distrb(:,:,:,:) 1194 real(dp) :: kbz(3),kmq_bz(3),spinrot_k(4),spinrot_kmq(4),q0(3),tsec(2) 1195 real(dp),ABI_CONTIGUOUS pointer :: qp_eig(:,:,:),qp_occ(:,:,:) 1196 real(dp),allocatable :: omegasf(:) 1197 complex(dpc),allocatable :: green_enhigh_w(:),green_w(:),kkweight(:,:) 1198 complex(gwpc),allocatable :: sf_chi0(:,:,:),rhotwg(:) 1199 complex(gwpc),allocatable :: ur1_kmq_ibz(:),ur2_k_ibz(:),wfwfg(:) 1200 complex(gwpc),allocatable :: usr1_kmq(:),ur2_k(:) 1201 complex(gwpc),allocatable :: ur_ae1(:),ur_ae_onsite1(:),ur_ps_onsite1(:) 1202 complex(gwpc),allocatable :: ur_ae2(:),ur_ae_onsite2(:),ur_ps_onsite2(:) 1203 complex(dpc), allocatable :: coeffW_BZ(:,:,:,:,:,:) 1204 logical,allocatable :: bbp_mask(:,:) 1205 type(pawcprj_type),allocatable :: Cprj1_kmq(:,:),Cprj2_k(:,:) 1206 type(pawpwij_t),allocatable :: Pwij(:),Pwij_fft(:) 1207 !************************************************************************ 1208 1209 DBG_ENTER("COLL") 1210 1211 call timab(331,1,tsec) ! cchi0 1212 call cwtime(cpu_time,wall_time,gflops,"start") 1213 1214 nsppol = Wfd%nsppol; nspinor = Wfd%nspinor 1215 is_metallic = ebands_has_metal_scheme(qp_ebands) 1216 1217 ucrpa_bands(1)=dtset%ucrpa_bands(1) 1218 ucrpa_bands(2)=dtset%ucrpa_bands(2) 1219 luwindow=.false. 1220 if(abs(dtset%ucrpa_window(1)+1_dp)>tol8.or.(abs(dtset%ucrpa_window(2)+1_dp)>tol8)) then 1221 luwindow=.true. 1222 endif 1223 !write(6,*)"ucrpa_bands",ucrpa_bands; write(6,*)"ucrpa_window",dtset%ucrpa_window; write(6,*)"luwindow",luwindow 1224 1225 ! For cRPA calculation of U: read forlb.ovlp 1226 if(dtset%ucrpa>=1 .AND. dtset%plowan_compute <10) then 1227 call read_plowannier(Cryst,bandinf,bandsup,coeffW_BZ,itypatcor,Kmesh,lcor,luwindow,& 1228 nspinor,nsppol,pawang,dtset%prtvol,ucrpa_bands) 1229 endif 1230 ! End of reading forlb.ovlp 1231 1232 if ( ANY(ngfft_gw(1:3) /= Wfd%ngfft(1:3)) ) call wfd%change_ngfft(Cryst,Psps,ngfft_gw) 1233 gw_mgfft = MAXVAL(ngfft_gw(1:3)) 1234 gw_fftalga = ngfft_gw(7)/100 !; gw_fftalgc=MOD(ngfft_gw(7),10) 1235 1236 if (Dtset%pawcross==1) mgfftf = MAXVAL(ngfftf(1:3)) 1237 1238 ! == Copy some values === 1239 comm = Wfd%comm 1240 mband = Wfd%mband 1241 nfft = Wfd%nfft 1242 ABI_CHECK(Wfd%nfftot==nfftot_gw,"Wrong nfftot_gw") 1243 1244 dim_rtwg=1 !; if (nspinor==2) dim_rtwg=2 ! can reduce size depending on Ep%nI and Ep%nj 1245 size_chi0 = Ep%npwe*Ep%nI*Ep%npwe*Ep%nJ*Ep%nomega 1246 1247 qp_eig => qp_ebands%eig; qp_occ => qp_ebands%occ 1248 1249 ! Initialize the completeness correction 1250 if (Ep%gwcomp==1) then 1251 en_high=MAXVAL(qp_eig(Ep%nbnds,:,:)) + Ep%gwencomp 1252 write(msg,'(a,f8.2,a)')' Using completeness correction with the energy ',en_high*Ha_eV,' [eV]' 1253 call wrtout(std_out, msg) 1254 1255 ! Allocation of wfwfg and green_enhigh_w moved inside openmp loop 1256 ! Init the largest G-sphere contained in the FFT box for the wavefunctions. 1257 call Gsph_FFT%in_fftbox(Cryst,Wfd%ngfft) 1258 1259 !call Gsph_FFT%print(unit=std_out,prtvol=10) 1260 1261 ABI_MALLOC(gspfft_igfft,(Gsph_FFT%ng)) 1262 ABI_MALLOC(dummy_gbound,(2*gw_mgfft+8,2)) 1263 1264 ! Mapping between G-sphere and FFT box. 1265 call Gsph_FFT%fft_tabs((/0,0,0/),Wfd%mgfft,Wfd%ngfft,dummy,dummy_gbound,gspfft_igfft) 1266 ABI_FREE(dummy_gbound) 1267 1268 if (Psps%usepaw==1) then ! * Prepare the onsite contributions on the GW FFT mesh. 1269 ABI_MALLOC(gw_gfft,(3,nfft)) 1270 q0=zero 1271 call get_gftt(ngfft_gw,q0,Cryst%gmet,gw_gsq,gw_gfft) ! Get the set of plane waves in the FFT Box. 1272 ABI_MALLOC(Pwij_fft,(Psps%ntypat)) 1273 call pawpwij_init(Pwij_fft,nfft,(/zero,zero,zero/),gw_gfft,Cryst%rprimd,Psps,Pawtab,Paw_pwff) 1274 end if 1275 end if 1276 1277 ! Setup weights (2 for spin unpolarized system, 1 for polarized). 1278 ! spin_fact is used to normalize the occupation factors to one. Consider also the AFM case. 1279 select case (nsppol) 1280 case (1) 1281 weight = two / Kmesh%nbz; spin_fact = half 1282 if (Wfd%nspden==2) then 1283 weight = one / Kmesh%nbz; spin_fact = half 1284 end if 1285 if (nspinor==2) then 1286 weight = one / Kmesh%nbz; spin_fact = one 1287 end if 1288 case (2) 1289 weight = one / Kmesh%nbz; spin_fact = one 1290 case default 1291 ABI_BUG("Wrong nsppol") 1292 end select 1293 1294 ! Weight for points in the IBZ_q. 1295 wtk_ltg(:) = 1 1296 if (Ep%symchi==1) then 1297 do ik_bz=1,Ltg_q%nbz 1298 wtk_ltg(ik_bz)=0 1299 if (Ltg_q%ibzq(ik_bz)/=1) CYCLE ! Only k-points in the IBZ_q. 1300 wtk_ltg(ik_bz)=SUM(Ltg_q%wtksym(:,:,ik_bz)) 1301 end do 1302 end if 1303 1304 write(msg,'(a,i2,2a,i2)')& 1305 ' Using spectral method for the imaginary part = ',Ep%spmeth,ch10,& 1306 ' Using symmetries to sum only over the IBZ_q = ',Ep%symchi 1307 call wrtout(std_out, msg) 1308 1309 if (use_tr) then 1310 ! Special care has to be taken in metals and/or spin dependent systems 1311 ! as Wfs_val might contain unoccupied states. 1312 call wrtout(std_out,' Using faster algorithm based on time reversal symmetry.') 1313 else 1314 call wrtout(std_out,' Using slow algorithm without time reversal symmetry.') 1315 end if 1316 1317 ! TODO this table can be calculated for each k-point 1318 my_nbbpks=0; allup="All"; got=0 1319 ABI_MALLOC(bbp_ks_distrb,(mband,mband,Kmesh%nbz,nsppol)) 1320 call wrtout(std_out, sjoin(' Memory needed for bbp_ks_distrb: ', & 1321 ftoa(four*mband**2*Kmesh%nbz*nsppol*b2Mb, fmt="f8.1"), ' [Mb] <<< MEM')) 1322 1323 ABI_MALLOC(bbp_mask,(mband, mband)) 1324 1325 do spin=1,nsppol 1326 do ik_bz=1,Kmesh%nbz 1327 if (Ep%symchi == 1) then 1328 if (Ltg_q%ibzq(ik_bz) /= 1) CYCLE ! Only IBZ_q 1329 end if 1330 1331 ! Get ik_ibz, non-symmorphic phase, ph_mkt, and symmetries from ik_bz. 1332 call kmesh%get_BZ_item(ik_bz,kbz,ik_ibz,isym_k,itim_k) 1333 1334 ! Get index of k-q in the BZ, stop if not found as the weight=one/nkbz is not correct. 1335 call kmesh%get_BZ_diff(kbz,qpoint,ikmq_bz,g0,nfound) 1336 ABI_CHECK(nfound==1,"Check kmesh") 1337 1338 ! Get ikmq_ibz, non-symmorphic phase, ph_mkmqt, and symmetries from ikmq_bz. 1339 call kmesh%get_BZ_item(ikmq_bz,kmq_bz,ikmq_ibz,isym_kmq,itim_kmq) 1340 1341 call chi0_bbp_mask(ikmq_ibz, ik_ibz, spin, spin_fact, use_tr, & 1342 ep%gwcomp, ep%spmeth, ep%nbnds, mband, qp_ebands, bbp_mask) 1343 1344 call wfd%distribute_kb_kpbp(ikmq_ibz,ik_ibz,spin,allup,my_nbbp,bbp_ks_distrb(:,:,ik_bz,spin),got=got,bbp_mask=bbp_mask) 1345 my_nbbpks = my_nbbpks + my_nbbp 1346 end do 1347 end do 1348 1349 ABI_FREE(bbp_mask) 1350 1351 write(msg,'(a,i0,a)')" Will sum ",my_nbbpks," (b,b',k,s) states in chi0." 1352 call wrtout(std_out, msg) 1353 1354 if (Psps%usepaw==1) then 1355 ABI_MALLOC(Pwij,(Psps%ntypat)) 1356 call pawpwij_init(Pwij,Ep%npwepG0,qpoint,Gsph_epsG0%gvec,Cryst%rprimd,Psps,Pawtab,Paw_pwff) 1357 ! Allocate statements moved to inside openmp loop 1358 end if 1359 1360 SELECT CASE (Ep%spmeth) 1361 CASE (0) 1362 call wrtout(std_out,' Calculating chi0(q,omega,G,G")') 1363 ! Allocation of green_w moved inside openmp loop 1364 1365 CASE (1, 2) 1366 call wrtout(std_out,' Calculating Im chi0(q,omega,G,G")') 1367 1368 ! Find Max and min resonant transitions for this q, report also treated by this proc. 1369 call make_transitions(Wfd,1,Ep%nbnds,nbvw,nsppol,Ep%symchi,Cryst%timrev,GW_TOL_DOCC,& 1370 & max_rest,min_rest,my_max_rest,my_min_rest,Kmesh,Ltg_q,qp_eig,qp_occ,qpoint,bbp_ks_distrb) 1371 ! 1372 ! Calculate frequency dependent weights for Hilbert transform. 1373 ABI_MALLOC(omegasf,(Ep%nomegasf)) 1374 ABI_MALLOC(kkweight,(Ep%nomegasf,Ep%nomega)) 1375 !my_wl=1; my_wr=Ep%nomegasf 1376 call setup_spectral(Ep%nomega,Ep%omega,Ep%nomegasf,omegasf,max_rest,min_rest,my_max_rest,my_min_rest,& 1377 & 0,Ep%zcut,zero,my_wl,my_wr,kkweight) 1378 1379 if (.not.use_tr) then 1380 ABI_BUG('spectral method requires time-reversal') 1381 end if 1382 1383 memreq = two*gwpc*Ep%npwe**2*(my_wr-my_wl+1)*b2Gb 1384 write(msg,'(a,f10.4,a)')' memory required per spectral point: ',two*gwpc*Ep%npwe**2*b2Mb,' [Mb]' 1385 call wrtout(std_out,msg) 1386 write(msg,'(a,f10.4,a)')' memory required by sf_chi0: ',memreq,' [Gb]' 1387 call wrtout(std_out,msg) 1388 if (memreq > two) then 1389 ABI_WARNING(' Memory required for sf_chi0 is larger than 2.0 Gb!') 1390 end if 1391 ABI_MALLOC_OR_DIE(sf_chi0,(Ep%npwe,Ep%npwe,my_wl:my_wr), ierr) 1392 sf_chi0=czero_gw 1393 1394 CASE DEFAULT 1395 ABI_BUG("Wrong spmeth") 1396 END SELECT 1397 1398 nkpt_summed=Kmesh%nbz 1399 if (Ep%symchi == 1) then 1400 nkpt_summed = Ltg_q%nibz_ltg 1401 call Ltg_q%print(std_out, Dtset%prtvol) 1402 end if 1403 1404 write(msg,'(a,i0,a)')' Calculation status: ',nkpt_summed,' k-points to be completed' 1405 call wrtout(std_out, msg) 1406 1407 ! ============================================ 1408 ! === Begin big fat loop over transitions === 1409 ! ============================================ 1410 chi0=czero_gw; chi0_sumrule=zero 1411 1412 ! === Loop on spin to calculate trace $\chi_{up,up}+\chi_{down,down}$ === 1413 ! Only $\chi_{up,up} for AFM. 1414 do spin=1,nsppol 1415 if (ALL(bbp_ks_distrb(:,:,:,spin) /= Wfd%my_rank)) CYCLE 1416 1417 ! Allocation of arrays that are private to loop 1418 if (Ep%gwcomp==1) then 1419 ABI_MALLOC(wfwfg,(nfft*nspinor**2)) 1420 end if 1421 if (Ep%gwcomp==1) then 1422 ABI_MALLOC(green_enhigh_w,(Ep%nomega)) 1423 end if 1424 if (Ep%spmeth==0) then 1425 ABI_MALLOC(green_w,(Ep%nomega)) 1426 end if 1427 if (Psps%usepaw==1) then 1428 ABI_MALLOC(Cprj2_k ,(Cryst%natom,nspinor)) 1429 call pawcprj_alloc(Cprj2_k, 0,Wfd%nlmn_atm) 1430 ABI_MALLOC(Cprj1_kmq,(Cryst%natom,nspinor)) 1431 call pawcprj_alloc(Cprj1_kmq,0,Wfd%nlmn_atm) 1432 if (Dtset%pawcross==1) then 1433 ABI_MALLOC(ur_ae1,(nfftf_tot*nspinor)) 1434 ABI_MALLOC(ur_ae_onsite1,(nfftf_tot*nspinor)) 1435 ABI_MALLOC(ur_ps_onsite1,(nfftf_tot*nspinor)) 1436 ABI_MALLOC(ur_ae2,(nfftf_tot*nspinor)) 1437 ABI_MALLOC(ur_ae_onsite2,(nfftf_tot*nspinor)) 1438 ABI_MALLOC(ur_ps_onsite2,(nfftf_tot*nspinor)) 1439 ABI_MALLOC(igfftepsG0f,(Ep%npwepG0)) 1440 ABI_MALLOC(tabrf_k,(nfftf_tot)) 1441 ABI_MALLOC(tabrf_kmq,(nfftf_tot)) 1442 end if 1443 end if 1444 1445 ABI_MALLOC(rhotwg,(Ep%npwepG0*dim_rtwg)) 1446 ABI_MALLOC(tabr_k,(nfft)) 1447 ABI_MALLOC(tabr_kmq,(nfft)) 1448 ABI_MALLOC(ur1_kmq_ibz,(nfft*nspinor)) 1449 ABI_MALLOC(ur2_k_ibz,(nfft*nspinor)) 1450 ABI_MALLOC(usr1_kmq,(nfft*nspinor)) 1451 ABI_MALLOC(ur2_k, (nfft*nspinor)) 1452 ABI_MALLOC(igfftepsG0,(Ep%npwepG0)) 1453 1454 ! Loop over k-points in the BZ. 1455 do ik_bz=1,Kmesh%nbz 1456 1457 if (Ep%symchi==1) then 1458 if (Ltg_q%ibzq(ik_bz)/=1) CYCLE ! Only IBZ_q 1459 end if 1460 1461 if (ALL(bbp_ks_distrb(:,:,ik_bz,spin) /= Wfd%my_rank)) CYCLE 1462 1463 write(msg,'(2(a,i4),a,i2,a,i3)')' ik= ',ik_bz,'/',Kmesh%nbz,' spin= ',spin,' done by mpi rank:',Wfd%my_rank 1464 call wrtout(std_out,msg) 1465 1466 ! Get ik_ibz, non-symmorphic phase, ph_mkt, and symmetries from ik_bz. 1467 call kmesh%get_BZ_item(ik_bz,kbz,ik_ibz,isym_k,itim_k,ph_mkt,umklp_k,isirred_k) 1468 1469 call kmesh%get_BZ_diff(kbz,qpoint,ikmq_bz,G0,nfound) 1470 if (nfound==0) then 1471 ABI_ERROR("Cannot find kbz - qpoint in Kmesh") 1472 end if 1473 1474 ! Get ikmq_ibz, non-symmorphic phase, ph_mkmqt, and symmetries from ikmq_bz. 1475 call kmesh%get_BZ_item(ikmq_bz,kmq_bz,ikmq_ibz,isym_kmq,itim_kmq,ph_mkmqt,umklp_kmq,isirred_kmq) 1476 1477 !BEGIN DEBUG 1478 !if (ANY(umklp_k /=0)) then 1479 ! write(msg,'(a,3i2)')" umklp_k /= 0 ",umklp_k 1480 ! ABI_ERROR(msg) 1481 !end if 1482 !if (ANY( g0 /= -umklp_kmq + umklp_k) ) then 1483 !if (ANY( g0 /= -umklp_kmq ) ) then 1484 ! write(msg,'(a,3(1x,3i2))')" g0 /= -umklp_kmq + umklp_k ",g0, umklp_kmq, umklp_k 1485 ! ABI_ERROR(msg) 1486 !end if 1487 !g0 = -umklp_k + umklp_kmq 1488 !g0 = +umklp_k - umklp_kmq 1489 !if (ANY (ABS(g0) > Ep%mg0) ) then 1490 ! write(msg,'(a,6(1x,i0))')" ABS(g0) > Ep%mg0 ",g0,Ep%mg0 1491 ! ABI_ERROR(msg) 1492 !end if 1493 !END DEBUG 1494 1495 ! Copy tables for rotated FFT points 1496 tabr_k(:) =ktabr(:,ik_bz) 1497 spinrot_k(:)=Cryst%spinrot(:,isym_k) 1498 1499 tabr_kmq(:)=ktabr(:,ikmq_bz) 1500 spinrot_kmq(:)=Cryst%spinrot(:,isym_kmq) 1501 1502 if (Dtset%pawcross==1) then 1503 tabrf_k(:) =ktabrf(:,ik_bz) 1504 tabrf_kmq(:)=ktabrf(:,ikmq_bz) 1505 end if 1506 ! 1507 ! Tables for the FFT of the oscillators. 1508 ! a) FFT index of G-G0. 1509 ! b) gw_gbound table for the zero-padded FFT performed in rhotwg. 1510 ABI_MALLOC(gw_gbound,(2*gw_mgfft+8,2)) 1511 call Gsph_epsG0%fft_tabs(g0,gw_mgfft,ngfft_gw,use_padfft,gw_gbound,igfftepsG0) 1512 if ( ANY(gw_fftalga == [2, 4]) ) use_padfft=0 ! Pad-FFT is not coded in rho_tw_g 1513 #ifdef FC_IBM 1514 ! XLF does not deserve this optimization (problem with [v67mbpt][t03]) 1515 use_padfft = 0 1516 #endif 1517 if (use_padfft==0) then 1518 ABI_FREE(gw_gbound) 1519 ABI_MALLOC(gw_gbound,(2*gw_mgfft+8,2*use_padfft)) 1520 end if 1521 1522 if (Dtset%pawcross==1) then 1523 ABI_MALLOC(gboundf,(2*mgfftf+8,2)) 1524 call Gsph_epsG0%fft_tabs(g0,mgfftf,ngfftf,use_padfftf,gboundf,igfftepsG0f) 1525 if (ANY(gw_fftalga == [2, 4])) use_padfftf=0 1526 if (use_padfftf==0) then 1527 ABI_FREE(gboundf) 1528 ABI_MALLOC(gboundf,(2*mgfftf+8,2*use_padfftf)) 1529 end if 1530 end if 1531 1532 nbmax=Ep%nbnds 1533 do band1=1,nbmax ! Loop over "conduction" states. 1534 if (ALL(bbp_ks_distrb(band1,:,ik_bz,spin) /= Wfd%my_rank)) CYCLE 1535 1536 call wfd%get_ur(band1,ikmq_ibz,spin,ur1_kmq_ibz) 1537 1538 if (Psps%usepaw==1) then 1539 call wfd%get_cprj(band1,ikmq_ibz,spin,Cryst,Cprj1_kmq,sorted=.FALSE.) 1540 call paw_symcprj(ikmq_bz,nspinor,1,Cryst,Kmesh,Pawtab,Pawang,Cprj1_kmq) 1541 if (Dtset%pawcross==1) then 1542 call wfdf%paw_get_aeur(band1,ikmq_ibz,spin,Cryst,Paw_onsite,Psps,Pawtab,Pawfgrtab,ur_ae1,ur_ae_onsite1,ur_ps_onsite1) 1543 end if 1544 end if 1545 1546 e_b1_kmq=qp_eig(band1,ikmq_ibz,spin) 1547 f_b1_kmq= qp_occ(band1,ikmq_ibz,spin) 1548 1549 do band2=1,nbmax ! Loop over "valence" states. 1550 if (luwindow.AND.dtset%ucrpa==1 & 1551 .AND.((qp_ebands%eig(band1,ik_ibz ,spin)-qp_ebands%fermie)<=dtset%ucrpa_window(2)) & 1552 .AND.((qp_ebands%eig(band1,ik_ibz ,spin)-qp_ebands%fermie)>=dtset%ucrpa_window(1)) & 1553 .AND.((qp_ebands%eig(band2,ikmq_ibz,spin)-qp_ebands%fermie)<=dtset%ucrpa_window(2)) & 1554 .AND.((qp_ebands%eig(band2,ikmq_ibz,spin)-qp_ebands%fermie)>=dtset%ucrpa_window(1))) CYCLE 1555 1556 if (bbp_ks_distrb(band1,band2,ik_bz,spin) /= Wfd%my_rank) CYCLE 1557 1558 deltaf_b1kmq_b2k=spin_fact*(f_b1_kmq-qp_occ(band2,ik_ibz,spin)) 1559 1560 if (Ep%gwcomp==0) then ! Skip negligible transitions. 1561 if (ABS(deltaf_b1kmq_b2k) < GW_TOL_DOCC) CYCLE 1562 else 1563 ! When the completeness correction is used, 1564 ! we need to also consider transitions with vanishing deltaf 1565 !if (qp_occ(band2,ik_ibz,spin) < GW_TOL_DOCC) CYCLE 1566 ! 1567 ! Rangel This is to compute chi correctly when using the extrapolar method 1568 if (qp_occ(band2,ik_ibz,spin) < GW_TOL_DOCC .and. (ABS(deltaf_b1kmq_b2k) < GW_TOL_DOCC .or. band1<band2)) CYCLE 1569 end if 1570 1571 deltaeGW_b1kmq_b2k=e_b1_kmq-qp_eig(band2,ik_ibz,spin) 1572 1573 call wfd%get_ur(band2,ik_ibz,spin,ur2_k_ibz) 1574 1575 if (Psps%usepaw==1) then 1576 call wfd%get_cprj(band2,ik_ibz,spin,Cryst,Cprj2_k,sorted=.FALSE.) 1577 call paw_symcprj(ik_bz,nspinor,1,Cryst,Kmesh,Pawtab,Pawang,Cprj2_k) 1578 if (Dtset%pawcross==1) then 1579 call wfdf%paw_get_aeur(band2,ik_ibz,spin,Cryst,Paw_onsite,Psps,Pawtab,Pawfgrtab,ur_ae2,ur_ae_onsite2,ur_ps_onsite2) 1580 end if 1581 end if 1582 1583 SELECT CASE (Ep%spmeth) 1584 CASE (0) 1585 ! Standard Adler-Wiser expression. 1586 ! Add the small imaginary of the Time-Ordered RF only for non-zero real omega ! FIXME What about metals? 1587 if (.not.use_tr) then 1588 ! Have to sum over all possible resonant and anti-resonant transitions. 1589 do io=1,Ep%nomega 1590 green_w(io) = g0g0w(Ep%omega(io),deltaf_b1kmq_b2k,deltaeGW_b1kmq_b2k,Ep%zcut,GW_TOL_W0,one_pole) 1591 end do 1592 1593 else 1594 if (Ep%gwcomp==0) then ! cannot be completely skipped in case of completeness correction 1595 if (band1<band2) CYCLE ! Here we GAIN a factor ~2 1596 end if 1597 1598 do io=1,Ep%nomega 1599 !Rangel: In metals, the intra-band transitions term does not contain the antiresonant part 1600 !green_w(io) = g0g0w(Ep%omega(io),deltaf_b1kmq_b2k,deltaeGW_b1kmq_b2k,Ep%zcut,GW_TOL_W0) 1601 if (band1==band2) then 1602 green_w(io) = g0g0w(Ep%omega(io),deltaf_b1kmq_b2k,deltaeGW_b1kmq_b2k,Ep%zcut,GW_TOL_W0,one_pole) 1603 else 1604 green_w(io) = g0g0w(Ep%omega(io),deltaf_b1kmq_b2k,deltaeGW_b1kmq_b2k,Ep%zcut,GW_TOL_W0,two_poles) 1605 end if 1606 1607 if (Ep%gwcomp==1) then ! Calculate the completeness correction 1608 numerator= -spin_fact*qp_occ(band2,ik_ibz,spin) 1609 deltaeGW_enhigh_b2k = en_high-qp_eig(band2,ik_ibz,spin) 1610 1611 if (REAL(Ep%omega(io))<GW_TOL_W0) then ! Completeness correction is NOT valid for real frequencies 1612 green_enhigh_w(io) = g0g0w(Ep%omega(io),numerator,deltaeGW_enhigh_b2k,Ep%zcut,GW_TOL_W0,two_poles) 1613 else 1614 green_enhigh_w(io) = local_czero_gw 1615 end if 1616 ! 1617 !Rangel Correction for metals 1618 !if (deltaf_b1kmq_b2k<0.d0) then 1619 if (band1>=band2 .and. ABS(deltaf_b1kmq_b2k) > GW_TOL_DOCC ) then 1620 green_w(io)= green_w(io) - green_enhigh_w(io) 1621 else ! Disregard green_w, since it is already accounted for through the time-reversal 1622 green_w(io)= - green_enhigh_w(io) 1623 end if 1624 end if !gwcomp==1 1625 end do !io 1626 1627 if (Ep%gwcomp==1.and.band1==band2) then 1628 ! Add the "delta part" of the extrapolar method. TODO doesnt work for spinor 1629 call calc_wfwfg(tabr_k,itim_k,spinrot_k,nfft,nspinor,ngfft_gw,ur2_k_ibz,ur2_k_ibz,wfwfg) 1630 1631 if (Psps%usepaw==1) then 1632 call paw_rho_tw_g(cryst,Pwij_fft, nfft,dim_rtwg,nspinor,gw_gfft,Cprj2_k,Cprj2_k,wfwfg) 1633 1634 ! Add PAW cross term 1635 if (Dtset%pawcross==1) then 1636 call paw_cross_rho_tw_g(nspinor,Ep%npwepG0,nfftf_tot,ngfftf,1,use_padfftf,igfftepsG0f,gboundf,& 1637 & ur_ae2,ur_ae_onsite2,ur_ps_onsite2,itim_kmq,tabrf_kmq,ph_mkmqt,spinrot_kmq,& 1638 & ur_ae2,ur_ae_onsite2,ur_ps_onsite2,itim_k ,tabrf_k ,ph_mkt ,spinrot_k,dim_rtwg,wfwfg) 1639 end if 1640 end if 1641 1642 qzero=.FALSE. 1643 call completechi0_deltapart(ik_bz,qzero,Ep%symchi,Ep%npwe,Gsph_FFT%ng,Ep%nomega,nspinor,& 1644 & nfft,ngfft_gw,gspfft_igfft,gsph_FFT,Ltg_q,green_enhigh_w,wfwfg,chi0) 1645 1646 end if 1647 end if ! use_tr 1648 1649 CASE (1, 2) 1650 ! Spectral method, WARNING time-reversal here is always assumed! 1651 if (deltaeGW_b1kmq_b2k<0) CYCLE 1652 call approxdelta(Ep%nomegasf,omegasf,deltaeGW_b1kmq_b2k,Ep%spsmear,iomegal,iomegar,wl,wr,Ep%spmeth) 1653 END SELECT 1654 1655 ! Form rho-twiddle(r)=u^*_{b1,kmq_bz}(r) u_{b2,kbz}(r) and its FFT transform. 1656 call rho_tw_g(nspinor,Ep%npwepG0,nfft,ndat1,ngfft_gw,1,use_padfft,igfftepsG0,gw_gbound,& 1657 & ur1_kmq_ibz,itim_kmq,tabr_kmq,ph_mkmqt,spinrot_kmq,& 1658 & ur2_k_ibz, itim_k ,tabr_k ,ph_mkt ,spinrot_k,dim_rtwg,rhotwg) 1659 1660 if (Psps%usepaw==1) then 1661 ! Add PAW on-site contribution, projectors are already in the BZ. 1662 call paw_rho_tw_g(cryst, Pwij, Ep%npwepG0,dim_rtwg,nspinor,Gsph_epsG0%gvec,Cprj1_kmq,Cprj2_k,rhotwg) 1663 1664 ! Add PAW cross term 1665 if (Dtset%pawcross==1) then 1666 call paw_cross_rho_tw_g(nspinor,Ep%npwepG0,nfftf_tot,ngfftf,1,use_padfftf,igfftepsG0f,gboundf,& 1667 & ur_ae1,ur_ae_onsite1,ur_ps_onsite1,itim_kmq,tabrf_kmq,ph_mkmqt,spinrot_kmq,& 1668 & ur_ae2,ur_ae_onsite2,ur_ps_onsite2,itim_k ,tabrf_k ,ph_mkt ,spinrot_k,dim_rtwg,rhotwg) 1669 end if 1670 end if 1671 1672 SELECT CASE (Ep%spmeth) 1673 1674 CASE (0) ! Adler-Wiser. 1675 !debug if(dtset%ucrpa==2) then 1676 if(dtset%ucrpa>=1.and..not.luwindow) then 1677 fac=one 1678 fac2=one 1679 fac3=one 1680 m1=-1 1681 m2=-1 1682 call flush_unit(std_out) 1683 if(dtset%ucrpa<=2) then 1684 if ( band1<=ucrpa_bands(2).AND.band1>=ucrpa_bands(1)& 1685 & .AND.band2<=ucrpa_bands(2).AND.band2>=ucrpa_bands(1)) then 1686 if (dtset%plowan_compute >=10) then 1687 band1c=band1-wan%bandi_wan+1 1688 band2c=band2-wan%bandi_wan+1 1689 do iat1=1, wan%natom_wan 1690 do iat2=1, wan%natom_wan 1691 do ispinor1=1,wan%nspinor 1692 do ispinor2=1,wan%nspinor 1693 do il1=1,wan%nbl_atom_wan(iat1) 1694 do il2=1,wan%nbl_atom_wan(iat2) 1695 do m1=1,2*wan%latom_wan(iat1)%lcalc(il1)+1 1696 do m2=1,2*wan%latom_wan(iat2)%lcalc(il2)+1 1697 fac=fac - real(wan%psichi(ik_bz,band1c,iat1)%atom(il1)%matl(m1,spin,ispinor1)*& 1698 &conjg(wan%psichi(ik_bz,band1c,iat1)%atom(il1)%matl(m1,spin,ispinor1))*& 1699 &wan%psichi(ikmq_bz,band2c,iat2)%atom(il2)%matl(m2,spin,ispinor2)*& 1700 &conjg(wan%psichi(ikmq_bz,band2c,iat2)%atom(il2)%matl(m2,spin,ispinor2))) 1701 enddo !m2 1702 enddo !m1 1703 enddo !il2 1704 enddo !il1 1705 enddo !ispinor2 1706 enddo !isspinor1 1707 enddo !iat2 1708 enddo !iat1 1709 else !plowan_compute>=10 1710 do iat=1, cryst%nattyp(itypatcor) 1711 do ispinor1=1,nspinor 1712 do ispinor2=1,nspinor 1713 do m1=1,2*lcor+1 1714 do m2=1,2*lcor+1 1715 fac=fac - real(coeffW_BZ(iat,spin,band1,ik_bz,ispinor1,m1)*& 1716 & conjg(coeffW_BZ(iat,spin,band1,ik_bz,ispinor1,m1))* & 1717 & coeffW_BZ(iat,spin,band2,ikmq_bz,ispinor2,m2)*& 1718 & conjg(coeffW_BZ(iat,spin,band2,ikmq_bz,ispinor2,m2))) 1719 enddo !m2 1720 enddo !m1 1721 enddo !ispinor2 1722 enddo !ispinor1 1723 enddo !iat 1724 endif !plowan_compute>=10 1725 if(dtset%ucrpa==1) fac=zero 1726 endif 1727 else if (dtset%ucrpa>=3) then 1728 if (band1<=ucrpa_bands(2).AND.band1>=ucrpa_bands(1)) then 1729 do iat=1, cryst%nattyp(itypatcor) 1730 do ispinor1=1,nspinor 1731 do m1=1,2*lcor+1 1732 fac2=fac2-real(coeffW_BZ(iat,spin,band1,ik_bz,ispinor1,m1)*& 1733 & conjg(coeffW_BZ(iat,spin,band1,ik_bz,ispinor1,m1))) 1734 enddo 1735 enddo 1736 enddo 1737 if(dtset%ucrpa==4) fac2=zero 1738 endif 1739 if (band2<=ucrpa_bands(2).AND.band2>=ucrpa_bands(1)) then 1740 do iat=1, cryst%nattyp(itypatcor) 1741 do ispinor1=1,nspinor 1742 do m1=1,2*lcor+1 1743 fac3=fac3-real(coeffW_BZ(iat,spin,band2,ikmq_bz,ispinor1,m1)*& 1744 & conjg(coeffW_BZ(iat,spin,band2,ikmq_bz,ispinor1,m1))) 1745 enddo 1746 enddo 1747 enddo 1748 if(dtset%ucrpa==4) fac3=zero 1749 endif 1750 fac=real(fac2*fac3) 1751 endif 1752 1753 ! if(dtset%prtvol>=10) write(6,'(6i3,e15.5,a)') ik_bz,ikmq_bz,band1,band2,m1,m2,fac," q=/0" 1754 ! if(dtset%prtvol>=10.and.abs(fac-one)>0.00001) & 1755 !& write(6,'(a,6i4,e15.5,a)') "*****FAC*********",ik_bz,ikmq_bz,band1,band2,m1,m2,fac," q/=0" 1756 ! if(dtset%prtvol>=10) write(6,'(a,6i4,e15.5,a)') "*****FAC*********",ik_bz,ikmq_bz,band1,band2,m1,m2,fac," q/=0" 1757 green_w=green_w*fac 1758 endif 1759 1760 call assemblychi0_sym(is_metallic,ik_bz,nspinor,Ep,Ltg_q,green_w,Ep%npwepG0,rhotwg,Gsph_epsG0,chi0) 1761 1762 CASE (1, 2) 1763 ! Spectral method (not yet adapted for nspinor=2) 1764 call assemblychi0sf(ik_bz,Ep%symchi,Ltg_q,Ep%npwepG0,Ep%npwe,rhotwg,Gsph_epsG0,& 1765 & deltaf_b1kmq_b2k,my_wl,iomegal,wl,my_wr,iomegar,wr,Ep%nomegasf,sf_chi0) 1766 1767 CASE DEFAULT 1768 ABI_BUG("Wrong spmeth") 1769 END SELECT 1770 1771 ! Accumulating the sum rule on chi0. Eq. (5.284) in G.D. Mahan Many-Particle Physics 3rd edition. [[cite:Mahan2000]] 1772 ! TODO Does not work with spinor 1773 factor=spin_fact*qp_occ(band2,ik_ibz,spin) 1774 call accumulate_chi0sumrule(ik_bz,Ep%symchi,Ep%npwe,factor,deltaeGW_b1kmq_b2k,& 1775 & Ltg_q,Gsph_epsG0,Ep%npwepG0,rhotwg,chi0_sumrule) 1776 1777 ! Include also the completeness correction in the sum rule 1778 if (Ep%gwcomp==1) then 1779 factor=-spin_fact*qp_occ(band2,ik_ibz,spin) 1780 call accumulate_chi0sumrule(ik_bz,Ep%symchi,Ep%npwe,factor,deltaeGW_enhigh_b2k,& 1781 & Ltg_q,Gsph_epsG0,Ep%npwepG0,rhotwg,chi0_sumrule) 1782 if (band1==Ep%nbnds) then 1783 chi0_sumrule(:)=chi0_sumrule(:) + wtk_ltg(ik_bz)*spin_fact*qp_occ(band2,ik_ibz,spin)*deltaeGW_enhigh_b2k 1784 end if 1785 end if 1786 1787 end do !band2 1788 end do !band1 1789 1790 ABI_FREE(gw_gbound) 1791 if (Dtset%pawcross==1) then 1792 ABI_FREE(gboundf) 1793 end if 1794 end do !ik_bz 1795 1796 ! Deallocation of arrays private to the spin loop. 1797 ABI_FREE(igfftepsG0) 1798 ABI_FREE(ur1_kmq_ibz) 1799 ABI_FREE(ur2_k_ibz) 1800 ABI_FREE(usr1_kmq) 1801 ABI_FREE(ur2_k) 1802 ABI_FREE(rhotwg) 1803 ABI_FREE(tabr_k) 1804 ABI_FREE(tabr_kmq) 1805 1806 ABI_SFREE(green_w) 1807 ABI_SFREE(wfwfg) 1808 ABI_SFREE(green_enhigh_w) 1809 if (Psps%usepaw==1) then 1810 call pawcprj_free(Cprj2_k) 1811 ABI_FREE(Cprj2_k) 1812 call pawcprj_free(Cprj1_kmq) 1813 ABI_FREE(Cprj1_kmq) 1814 if (Dtset%pawcross==1) then 1815 ABI_FREE(ur_ae1) 1816 ABI_FREE(ur_ae_onsite1) 1817 ABI_FREE(ur_ps_onsite1) 1818 ABI_FREE(ur_ae2) 1819 ABI_FREE(ur_ae_onsite2) 1820 ABI_FREE(ur_ps_onsite2) 1821 ABI_FREE(tabrf_k) 1822 ABI_FREE(tabrf_kmq) 1823 ABI_FREE(gboundf) 1824 ABI_FREE(igfftepsG0f) 1825 end if 1826 end if 1827 end do !spin 1828 1829 ! After big loop over transitions, now MPI 1830 ! Master took care of the contribution in case of metallic|spin polarized systems. 1831 SELECT CASE (Ep%spmeth) 1832 CASE (0) 1833 ! Adler-Wiser 1834 ! Collective sum of the contributions of each node. 1835 ! Looping on frequencies to avoid problems with the size of the MPI packet 1836 do io=1,Ep%nomega 1837 call xmpi_sum(chi0(:,:,io),comm,ierr) 1838 end do 1839 1840 CASE (1, 2) 1841 ! Spectral method. 1842 call hilbert_transform(Ep%npwe,Ep%nomega,Ep%nomegasf,my_wl,my_wr,kkweight,sf_chi0,chi0,Ep%spmeth) 1843 1844 ! Deallocate here before xmpi_sum 1845 ABI_SFREE(sf_chi0) 1846 1847 ! Collective sum of the contributions. 1848 ! Looping over frequencies to avoid problems with the size of the MPI packet 1849 do io=1,Ep%nomega 1850 call xmpi_sum(chi0(:,:,io),comm,ierr) 1851 end do 1852 1853 CASE DEFAULT 1854 ABI_BUG("Wrong spmeth") 1855 END SELECT 1856 1857 ! Divide by the volume 1858 !$OMP PARALLEL WORKSHARE 1859 chi0=chi0*weight/Cryst%ucvol 1860 !$OMP END PARALLEL WORKSHARE 1861 1862 ! === Collect the sum rule === 1863 ! * The pi factor comes from Im[1/(x-ieta)] = pi delta(x) 1864 call xmpi_sum(chi0_sumrule,comm,ierr) 1865 chi0_sumrule=chi0_sumrule*pi*weight/Cryst%ucvol 1866 ! 1867 ! ************************************************* 1868 ! **** Now each node has chi0(q,G,Gp,Ep%omega) **** 1869 ! ************************************************* 1870 1871 ! Impose Hermiticity (valid only for zero or purely imaginary frequencies) 1872 ! MG what about metals, where we have poles around zero? 1873 ! FB because of the intraband term, chi0 is never hermitian in case of metals 1874 ! FIXME: as of today, hermitianity is also enforced for metallic systems 1875 !if (.not. is_metallic) then 1876 do io=1,Ep%nomega 1877 if (ABS(REAL(Ep%omega(io))) <0.00001) then 1878 do ig2=1,Ep%npwe 1879 do ig1=1,ig2-1 1880 chi0(ig2,ig1,io) = GWPC_CONJG(chi0(ig1,ig2,io)) 1881 end do 1882 end do 1883 end if 1884 end do 1885 !endif 1886 1887 ! === Symmetrize chi0 in case of AFM system === 1888 ! Reconstruct $chi0{\down,\down}$ from $chi0{\up,\up}$. 1889 ! Works only in case of magnetic group Shubnikov type IV. 1890 if (Cryst%use_antiferro) then 1891 call symmetrize_afm_chi0(Cryst,Gsph_epsG0,Ltg_q,Ep%npwe,Ep%nomega,chi0=chi0) 1892 end if 1893 1894 ! ===================== 1895 ! ==== Free memory ==== 1896 ! ===================== 1897 ABI_FREE(bbp_ks_distrb) 1898 1899 ABI_SFREE(gw_gfft) 1900 ABI_SFREE(kkweight) 1901 ABI_SFREE(omegasf) 1902 ABI_SFREE(gspfft_igfft) 1903 1904 call Gsph_FFT%free() 1905 1906 ! deallocation for PAW. 1907 if (Psps%usepaw==1) then 1908 call pawpwij_free(Pwij) 1909 ABI_FREE(Pwij) 1910 if (allocated(Pwij_fft)) then 1911 call pawpwij_free(Pwij_fft) 1912 ABI_FREE(Pwij_fft) 1913 end if 1914 end if 1915 1916 if(dtset%ucrpa>=1 .AND. dtset%plowan_compute<10) then 1917 ABI_FREE(coeffW_BZ) 1918 endif 1919 1920 call timab(331,2,tsec) 1921 call cwtime(cpu_time,wall_time,gflops,"stop") 1922 write(std_out,'(2(a,f9.1))')" cpu_time = ",cpu_time,", wall_time = ",wall_time 1923 1924 DBG_EXIT("COLL") 1925 1926 end subroutine cchi0
ABINIT/cchi0q0 [ Functions ]
NAME
cchi0q0
FUNCTION
Calculate chi0 in the limit q --> 0
INPUTS
use_tr=If .TRUE. Wfs_val are allocate and only resonant transitions are evaluated (assumes time reversal symmetry) Dtset <type(dataset_type)>=all input variables in this dataset Ep= datatype gathering differening parameters related to the calculation of the inverse dielectric matrix Gsph_epsG0<gvectors_data_type>: Info on the G-sphere used to describe chi0/espilon (including umklapp) %ng=number of G vectors %rottbm1(ng,2,nsym)=contains the index (IS^{-1}) G in the array gvec %phmGt(ng,nsym)=phase factor e^{-iG.\tau} needed to symmetrize oscillator matrix elements and chi0 %gmet(3,3)=reciprocal space metric ($\textrm{bohr}^{-2}$). %gprimd(3,3)=dimensional reciprocal space primitive translations (b^-1) Ep%inclvkb=flag to include (or not) the grad of Vkb Ltg_q= little group datatype nbvw=number of bands in the arrays wfrv,wfgv Kmesh<kmesh_t> The k-point mesh %kbz(3,nbz)=k-point coordinates, full Brillouin zone %tab(nbz)= table giving for each k-point in the BZ (kBZ), the corresponding irreducible point (kIBZ), where kBZ= (IS) kIBZ and I is either the inversion or the identity %tabi(nbzx)= for each point in the BZ defines whether inversion has to be considered in the relation kBZ=(IS) kIBZ (1 => only S; -1 => -S) %tabo(nbzx)= the symmetry operation S that takes kIBZ to each kBZ %tabp(nbzx)= phase factor associated to tnons e^{-i 2 \pi k\cdot R{^-1}t} ktabr(nfftot_gw,Kmesh%nbz) index of R^-(r-t) in the FFT array, where k_BZ = (IS) k_IBZ and S = \transpose R^{-1} Ep%nbnds=number of bands ngfft_gw(18)= array containing all the information for 3D FFT for the oscillator strengths. Ep%nomega=number of frequencies Cryst<crystal_t>= data type gathering info on symmetries and unit cell %natom=number of atoms %nsym=number of symmetry operations %symrec(3,3,nsym)=symmetry operations in reciprocal space %typat(natom)=type of each atom %xred(3,natom)=reduced coordinated of atoms %rprimd(3,3)=dimensional primitive translations in real space (bohr) %timrev=2 if time-reversal symmetry can be used, 1 otherwise Ep%npwe=number of planewaves for sigma exchange (input variable) nfftot_gw=Total number of points in the GW FFT grid Ep%omega(Ep%nomega)=frequencies Psps <type(pseudopotential_type)>=variables related to pseudopotentials %mpsang=1+maximum angular momentum for nonlocal pseudopotential Pawang<pawang_type> angular mesh discretization and related data: Pawrad(ntypat*usepaw)<Pawrad_type>=paw radial mesh and related data Paw_ij(natom*usepaw)<Paw_ij_type)>=paw arrays given on (i,j) channels qp_ebands<ebands_t>=Quasiparticle energies and occupations (for the moment real quantities) %mband=MAX number of bands over k-points and spin (==Ep%nbnds) %occ(mband,nkpt,nsppol)=QP occupation numbers, for each k point in IBZ, and each band %eig(mband,nkpt,nsppol)=GW energies, for self-consistency purposes ks_ebands<ebands_t>=KS energies and occupations. %eig(mband,nkpt,nsppol)=KS energies Paw_pwff<pawpwff_t>=Form factor used to calculate the onsite mat. elements of a plane wave. Wfd<wfdgw_t>=Object used to access the wavefunctions
OUTPUT
chi0(Ep%npwe,Ep%npwe,Ep%nomega)=independent-particle susceptibility matrix for wavevector qq, and frequencies defined by Ep%omega chi0_lwing(Ep%npwe*Ep%nI,Ep%nomega,3)= Lower wings chi0_uwing(Ep%npwe*Ep%nJ,Ep%nomega,3)= Upper wings chi0_head(3,3,Ep%nomega)=Head of chi0.
NOTES
*) The terms "head", "wings" and "body" of chi(G,Gp) refer to G=Gp=0, either G or Gp=0, and neither=0 respectively *) Symmetry conventions: 1) symmetry in real space is defined as: R_t f(r) = f(R^-1(r-t)) 2) S=\transpose R^-1 3) kbz=S kibz The wavefunctions for the k-point in the BZ are (assuming nondegenerate states): u(G,b, Sk) = u ( S^-1G,b,k)* e^{-i(Sk+G)*t) u(G,b,-Sk) = u*(-S^-1G,b,k)* e^{ i(Sk-G)*t) u(r,b, Sk) = u (R^-1(r-t),b,k) e^{-iSk*t} u(r,b,-Sk) = u*(R^-1(r-t),b,k) e^{ iSK*t} The gradient of Vnl(K,Kp) for the k-point in the BZ should be: gradvnl(SG,SGp,Sk)=S gradvnl(G,Gp,kibz)
TODO
Check npwepG0 before activating umklapps
SOURCE
166 subroutine cchi0q0(use_tr,Dtset,Cryst,Ep,Psps,Kmesh,qp_ebands,ks_ebands,Gsph_epsG0,& 167 Pawang,Pawrad,Pawtab,Paw_ij,Paw_pwff,Pawfgrtab,Paw_onsite,ktabr,ktabrf,nbvw,ngfft_gw,& 168 nfftot_gw,ngfftf,nfftf_tot,chi0,chi0_head,chi0_lwing,chi0_uwing,Ltg_q,chi0_sumrule,Wfd,Wfdf,wan) 169 170 !Arguments ------------------------------------ 171 !scalars 172 integer,intent(in) :: nbvw,nfftot_gw,nfftf_tot 173 logical,intent(in) :: use_tr 174 type(ebands_t),target,intent(in) :: qp_ebands,ks_ebands 175 type(crystal_t),intent(in) :: Cryst 176 type(Dataset_type),intent(in) :: Dtset 177 type(littlegroup_t),intent(in) :: Ltg_q 178 type(em1params_t),intent(in) :: Ep 179 type(kmesh_t),intent(in) :: Kmesh 180 type(gsphere_t),intent(in) :: Gsph_epsG0 181 type(Pseudopotential_type),intent(in) :: Psps 182 type(Pawang_type),intent(in) :: Pawang 183 type(wfdgw_t),target,intent(inout) :: Wfd,Wfdf 184 !arrays 185 integer,intent(in) :: ktabr(nfftot_gw,Kmesh%nbz),ktabrf(nfftf_tot*Dtset%pawcross,Kmesh%nbz) 186 integer,intent(in) :: ngfft_gw(18),ngfftf(18) 187 real(dp),intent(out) :: chi0_sumrule(Ep%npwe) 188 complex(gwpc),intent(out) :: chi0(Ep%npwe,Ep%npwe,Ep%nomega) 189 complex(dpc),intent(out) :: chi0_lwing(Ep%npwe*Ep%nI,Ep%nomega,3) 190 complex(dpc),intent(out) :: chi0_uwing(Ep%npwe*Ep%nJ,Ep%nomega,3) 191 complex(dpc),intent(out) :: chi0_head(3,3,Ep%nomega) 192 type(Pawrad_type),intent(in) :: Pawrad(Psps%ntypat*Psps%usepaw) 193 type(Pawtab_type),intent(in) :: Pawtab(Psps%ntypat*Psps%usepaw) 194 type(Paw_ij_type),intent(in) :: Paw_ij(Cryst%natom*Psps%usepaw) 195 type(pawpwff_t),intent(in) :: Paw_pwff(Psps%ntypat*Psps%usepaw) 196 type(pawfgrtab_type),intent(inout) :: Pawfgrtab(Cryst%natom*Psps%usepaw) 197 type(plowannier_type),intent(inout) :: wan 198 type(paw_pwaves_lmn_t),intent(in) :: Paw_onsite(Cryst%natom) 199 200 !Local variables ------------------------------ 201 !scalars 202 integer,parameter :: tim_fourdp=1,enough=10,two_poles=2,one_pole=1,ndat1=1 203 integer :: bandinf,bandsup,lcor,nspinor,npw_k,istwf_k,mband,nfft,band1c,band2c 204 integer :: band1,band2,iat1,iat2,iat,ig,ig1,ig2,itim_k,ik_bz,ik_ibz,io,iqlwl,ispinor1,ispinor2,isym_k,il1,il2 205 integer :: itypatcor,m1,m2,nkpt_summed,dim_rtwg,use_padfft,gw_fftalga,use_padfftf,mgfftf 206 integer :: my_nbbp,my_nbbpks,spin,nsppol,iq,nq 207 integer :: comm,ierr,my_wl,my_wr,iomegal,iomegar,gw_mgfft,dummy 208 real(dp) :: cpu_time,wall_time,gflops 209 real(dp) :: fac,fac1,fac2,fac3,fac4,spin_fact,deltaf_b1b2,weight,factor 210 real(dp) :: max_rest,min_rest,my_max_rest,my_min_rest, qlen 211 real(dp) :: en_high,deltaeGW_enhigh_b2,wl,wr,numerator,deltaeGW_b1b2,gw_gsq,memreq 212 complex(dpc) :: deltaeKS_b1b2 213 logical :: qzero,luwindow,is_metallic 214 character(len=500) :: msg_tmp,msg,allup 215 type(gsphere_t) :: Gsph_FFT 216 type(wave_t),pointer :: wave1, wave2 217 !arrays 218 integer,ABI_CONTIGUOUS pointer :: kg_k(:,:) 219 integer :: ucrpa_bands(2), got(Wfd%nproc) 220 integer :: wtk_ltg(Kmesh%nbz) 221 integer,allocatable :: tabr_k(:),tabrf_k(:), igffteps0(:),gspfft_igfft(:),igfftepsG0f(:) 222 integer,allocatable :: gw_gfft(:,:),gw_gbound(:,:),dummy_gbound(:,:),gboundf(:,:), bbp_ks_distrb(:,:,:,:) 223 real(dp) :: kbz(3),spinrot_kbz(4),q0(3) 224 real(dp),ABI_CONTIGUOUS pointer :: ks_eig(:,:,:),qp_eig(:,:,:),qp_occ(:,:,:) 225 real(dp),allocatable :: omegasf(:), qdirs(:,:) 226 complex(gwpc) :: rhotwx(3,Wfd%nspinor**2) 227 complex(gwpc),allocatable :: rhotwg(:) 228 complex(dpc),allocatable :: green_w(:),green_enhigh_w(:) 229 complex(dpc),allocatable :: sf_lwing(:,:,:),sf_uwing(:,:,:),sf_head(:,:,:) 230 complex(dpc) :: chq(3), wng(3) 231 complex(dpc) :: ph_mkt 232 complex(dpc),allocatable :: kkweight(:,:) 233 complex(gwpc),allocatable :: ur1_kibz(:),ur2_kibz(:), usr1_k(:),ur2_k(:), wfwfg(:), sf_chi0(:,:,:) 234 complex(gwpc),allocatable :: ur_ae1(:),ur_ae_onsite1(:),ur_ps_onsite1(:) 235 complex(gwpc),allocatable :: ur_ae2(:),ur_ae_onsite2(:),ur_ps_onsite2(:) 236 complex(gwpc),ABI_CONTIGUOUS pointer :: ug1(:),ug2(:) 237 complex(dpc), allocatable :: coeffW_BZ(:,:,:,:,:,:), head_qvals(:) 238 logical :: gradk_not_done(Kmesh%nibz) 239 logical,allocatable :: bbp_mask(:,:) 240 type(pawcprj_type),allocatable :: Cprj1_bz(:,:),Cprj2_bz(:,:), Cprj1_ibz(:,:),Cprj2_ibz(:,:) 241 type(pawpwij_t),allocatable :: Pwij(:),Pwij_fft(:) 242 type(pawhur_t),allocatable :: Hur(:) 243 type(vkbr_t),allocatable :: vkbr(:) 244 245 !************************************************************************ 246 247 DBG_ENTER("COLL") 248 249 call cwtime(cpu_time, wall_time, gflops, "start") 250 251 ! Change FFT mesh if needed 252 if (ANY(ngfft_gw(1:3) /= Wfd%ngfft(1:3))) call wfd%change_ngfft(Cryst,Psps,ngfft_gw) 253 254 gw_mgfft = MAXVAL(ngfft_gw(1:3)); gw_fftalga = ngfft_gw(7)/100 !; gw_fftalgc=MOD(ngfft_gw(7),10) 255 if (Dtset%pawcross==1) mgfftf = MAXVAL(ngfftf(1:3)) 256 257 ! Copy important variables. 258 comm = Wfd%comm; nsppol = Wfd%nsppol; nspinor = Wfd%nspinor; mband = Wfd%mband; nfft = Wfd%nfft 259 ABI_CHECK(Wfd%nfftot == nfftot_gw, "Wrong nfftot_gw") 260 dim_rtwg = 1 !; if (nspinor==2) dim_rtwg=2 ! Can reduce size depending on Ep%nI and Ep%nj 261 262 is_metallic = ebands_has_metal_scheme(qp_ebands) 263 ucrpa_bands(1)=dtset%ucrpa_bands(1) 264 ucrpa_bands(2)=dtset%ucrpa_bands(2) 265 luwindow=.false. 266 if (abs(dtset%ucrpa_window(1)+1_dp)>tol8.or.(abs(dtset%ucrpa_window(2)+1_dp)>tol8)) luwindow=.true. 267 268 ! For cRPA calculation of U: read forlb.ovlp 269 if(dtset%ucrpa>=1 .AND. dtset%plowan_compute<10) then 270 call read_plowannier(Cryst,bandinf,bandsup,coeffW_BZ,itypatcor,Kmesh,lcor,luwindow,& 271 nspinor,nsppol,pawang,dtset%prtvol,ucrpa_bands) 272 endif 273 274 ks_eig => ks_ebands%eig 275 qp_eig => qp_ebands%eig; qp_occ => qp_ebands%occ 276 277 chi0_lwing = czero; chi0_uwing = czero; chi0_head = czero 278 279 if (Psps%usepaw == 0) then 280 if (Ep%inclvkb /= 0) then 281 ! Include the term <n,k|[Vnl,iqr]|n"k>' for q -> 0. 282 ABI_CHECK(nspinor == 1, "nspinor with inclvkb not coded") 283 else 284 ABI_WARNING('Neglecting <n,k|[Vnl,iqr]|m,k>') 285 end if 286 287 else 288 ! For PAW+DFT+U, precalculate <\phi_i|[Hu,r]|phi_j\> 289 ABI_MALLOC(HUr, (Cryst%natom)) 290 if (Dtset%usepawu /= 0) then 291 call pawhur_init(hur,nsppol,Dtset%pawprtvol,Cryst,Pawtab,Pawang,Pawrad,Paw_ij) 292 end if 293 end if 294 295 ! Initialize the completeness correction. 296 ABI_MALLOC(green_enhigh_w, (Ep%nomega)) 297 green_enhigh_w = czero 298 299 if (Ep%gwcomp == 1) then 300 en_high = MAXVAL(qp_eig(Ep%nbnds,:,:))+Ep%gwencomp 301 write(msg,'(a,f8.2,a)')' Using completeness correction with energy ',en_high*Ha_eV,' [eV] ' 302 call wrtout(std_out, msg) 303 ABI_MALLOC(wfwfg,(nfft*nspinor**2)) 304 305 ! Init the largest G-sphere contained in the FFT box for the wavefunctions. 306 call Gsph_FFT%in_fftbox(Cryst,Wfd%ngfft) 307 call Gsph_FFT%print(unit=std_out, prtvol=10) 308 309 ABI_MALLOC(gspfft_igfft,(Gsph_FFT%ng)) 310 ABI_MALLOC(dummy_gbound,(2*gw_mgfft+8,2)) 311 312 ! Mapping between G-sphere and FFT box. 313 call Gsph_FFT%fft_tabs([0, 0, 0],Wfd%mgfft,Wfd%ngfft,dummy,dummy_gbound,gspfft_igfft) 314 ABI_FREE(dummy_gbound) 315 316 if (Psps%usepaw==1) then 317 ! Prepare the onsite contributions on the GW FFT mesh. 318 ABI_MALLOC(gw_gfft,(3,nfft)) 319 q0=zero 320 call get_gftt(ngfft_gw,q0,Cryst%gmet,gw_gsq,gw_gfft) ! The set of plane waves in the FFT Box. 321 ABI_MALLOC(Pwij_fft,(Psps%ntypat)) 322 call pawpwij_init(Pwij_fft,nfft,(/zero,zero,zero/),gw_gfft,Cryst%rprimd,Psps,Pawtab,Paw_pwff) 323 end if 324 end if 325 326 ! Setup weight (2 for spin unpolarized systems, 1 for polarized). 327 ! spin_fact is used to normalize the occupation factors to one. 328 ! Consider also the AFM case. 329 select case (nsppol) 330 case (1) 331 weight = two / Kmesh%nbz; spin_fact = half 332 if (Wfd%nspden == 2) then 333 weight = one / Kmesh%nbz; spin_fact = half 334 end if 335 if (nspinor == 2) then 336 weight = one / Kmesh%nbz; spin_fact = one 337 end if 338 case (2) 339 weight = one / Kmesh%nbz; spin_fact = one 340 case default 341 ABI_BUG(sjoin("Wrong nsppol:", itoa(nsppol))) 342 end select 343 344 ! k-weights for points in the IBZ_q 345 wtk_ltg(:) = 1 346 if (Ep%symchi == 1) then 347 do ik_bz=1,Ltg_q%nbz 348 wtk_ltg(ik_bz) = 0 349 if (Ltg_q%ibzq(ik_bz) /= 1) CYCLE ! Only k in IBZ_q 350 wtk_ltg(ik_bz) = sum(Ltg_q%wtksym(:,:,ik_bz)) 351 end do 352 end if 353 354 write(msg,'(a,i3,a)')' Q-points for long wave-length limit. # ',Ep%nqlwl,ch10 355 do iqlwl=1,Ep%nqlwl 356 write(msg_tmp,'(1x,i5,a,2x,3f12.6,a)') iqlwl,')',Ep%qlwl(:,iqlwl),ch10 357 msg=TRIM(msg)//msg_tmp 358 end do 359 call wrtout(std_out, msg) 360 361 write(msg,'(a,i2,2a,i2)')& 362 ' Using spectral method for the imaginary part = ',Ep%spmeth,ch10,& 363 ' Using symmetries to sum only over the IBZ_q = ',Ep%symchi 364 call wrtout(std_out, msg) 365 366 if (use_tr) then 367 call wrtout(std_out,' Using faster algorithm based on time reversal symmetry.') 368 else 369 call wrtout(std_out,' Using slow algorithm without time reversal symmetry.') 370 end if 371 372 ! Evaluate oscillator matrix elements btw partial waves. Note q=Gamma 373 if (Psps%usepaw == 1) then 374 ABI_MALLOC(Pwij,(Psps%ntypat)) 375 call pawpwij_init(Pwij,Ep%npwepG0, [zero, zero, zero],Gsph_epsG0%gvec,Cryst%rprimd,Psps,Pawtab,Paw_pwff) 376 377 ABI_MALLOC(Cprj1_bz,(Cryst%natom,nspinor)) 378 call pawcprj_alloc(Cprj1_bz,0,Wfd%nlmn_atm) 379 ABI_MALLOC(Cprj2_bz,(Cryst%natom,nspinor)) 380 call pawcprj_alloc(Cprj2_bz,0,Wfd%nlmn_atm) 381 382 ABI_MALLOC(Cprj1_ibz,(Cryst%natom,nspinor)) 383 call pawcprj_alloc(Cprj1_ibz,0,Wfd%nlmn_atm) 384 ABI_MALLOC(Cprj2_ibz,(Cryst%natom,nspinor)) 385 call pawcprj_alloc(Cprj2_ibz,0,Wfd%nlmn_atm) 386 if (Dtset%pawcross==1) then 387 ABI_MALLOC(ur_ae1,(nfftf_tot*nspinor)) 388 ABI_MALLOC(ur_ae_onsite1,(nfftf_tot*nspinor)) 389 ABI_MALLOC(ur_ps_onsite1,(nfftf_tot*nspinor)) 390 ABI_MALLOC(ur_ae2,(nfftf_tot*nspinor)) 391 ABI_MALLOC(ur_ae_onsite2,(nfftf_tot*nspinor)) 392 ABI_MALLOC(ur_ps_onsite2,(nfftf_tot*nspinor)) 393 ABI_MALLOC(igfftepsG0f,(Ep%npwepG0)) 394 ABI_MALLOC(tabrf_k,(nfftf_tot)) 395 end if 396 end if 397 398 ABI_MALLOC(rhotwg,(Ep%npwe*dim_rtwg)) 399 ABI_MALLOC(tabr_k,(nfft)) 400 ABI_MALLOC(ur1_kibz,(nfft*nspinor)) 401 ABI_MALLOC(ur2_kibz,(nfft*nspinor)) 402 403 ABI_MALLOC(usr1_k,(nfft*nspinor)) 404 ABI_MALLOC(ur2_k,(nfft*nspinor)) 405 ! 406 ! Tables for the FFT of the oscillators. 407 ! a) FFT index of the G sphere (only vertical transitions, unlike cchi0, no need to shift the sphere). 408 ! b) gw_gbound table for the zero-padded FFT performed in rhotwg. 409 410 ABI_MALLOC(igffteps0,(Gsph_epsG0%ng)) 411 ABI_MALLOC(gw_gbound,(2*gw_mgfft+8,2)) 412 call Gsph_epsG0%fft_tabs([0, 0, 0], gw_mgfft,ngfft_gw,use_padfft,gw_gbound,igffteps0) 413 if ( ANY(gw_fftalga == [2, 4]) ) use_padfft=0 ! Pad-FFT is not coded in rho_tw_g 414 #ifdef FC_IBM 415 ! XLF does not deserve this optimization (problem with [v67mbpt][t03]) 416 use_padfft = 0 417 #endif 418 if (use_padfft==0) then 419 ABI_FREE(gw_gbound) 420 ABI_MALLOC(gw_gbound,(2*gw_mgfft+8,2*use_padfft)) 421 end if 422 if (Dtset%pawcross==1) then 423 ABI_MALLOC(gboundf,(2*mgfftf+8,2)) 424 call Gsph_epsG0%fft_tabs((/0,0,0/),mgfftf,ngfftf,use_padfftf,gboundf,igfftepsG0f) 425 if ( ANY(gw_fftalga == (/2,4/)) ) use_padfftf=0 426 if (use_padfftf==0) then 427 ABI_FREE(gboundf) 428 ABI_MALLOC(gboundf,(2*mgfftf+8,2*use_padfftf)) 429 end if 430 end if 431 432 ! TODO this table can be calculated for each k-point 433 my_nbbpks=0; allup="All"; got=0 434 ABI_MALLOC(bbp_ks_distrb,(mband,mband,Kmesh%nbz,nsppol)) 435 call wrtout(std_out, sjoin(' Memory needed for bbp_ks_distrb: ', & 436 ftoa(four*mband**2*Kmesh%nbz*nsppol*b2Mb, fmt="f8.1"), ' [Mb] <<< MEM')) 437 438 ABI_MALLOC(bbp_mask, (mband, mband)) 439 440 do spin=1,nsppol 441 do ik_bz=1,Kmesh%nbz 442 443 if (Ep%symchi==1) then 444 if (Ltg_q%ibzq(ik_bz)/=1) CYCLE ! Only IBZ_q 445 end if 446 ik_ibz=Kmesh%tab(ik_bz) 447 448 call chi0_bbp_mask(ik_ibz, ik_ibz, spin, spin_fact, use_tr, & 449 ep%gwcomp, ep%spmeth, ep%nbnds, mband, qp_ebands, bbp_mask) 450 451 call wfd%distribute_bbp(ik_ibz,spin,allup,my_nbbp,bbp_ks_distrb(:,:,ik_bz,spin),got=got,bbp_mask=bbp_mask) 452 my_nbbpks = my_nbbpks + my_nbbp 453 end do 454 end do 455 456 ABI_FREE(bbp_mask) 457 458 write(msg,'(a,i0,a)')" Will sum ",my_nbbpks," (b,b',k,s) states in chi0q0." 459 call wrtout(std_out, msg) 460 461 SELECT CASE (Ep%spmeth) 462 CASE (0) 463 call wrtout(std_out,' Calculating chi0(q=(0,0,0),omega,G,G")') 464 ABI_MALLOC(green_w, (Ep%nomega)) 465 466 CASE (1, 2) 467 call wrtout(std_out,' Calculating Im chi0(q=(0,0,0),omega,G,G")') 468 ! 469 ! === Find max and min resonant transitions for this q, report values for this processor === 470 call make_transitions(Wfd,1,Ep%nbnds,nbvw,nsppol,Ep%symchi,Cryst%timrev,GW_TOL_DOCC,& 471 max_rest,min_rest,my_max_rest,my_min_rest,Kmesh,Ltg_q,qp_eig,qp_occ, [zero,zero,zero], bbp_ks_distrb) 472 473 ! === Calculate frequency dependent weights for Kramers Kronig transform === 474 ABI_MALLOC(omegasf, (Ep%nomegasf)) 475 ABI_MALLOC(kkweight, (Ep%nomegasf,Ep%nomega)) 476 !my_wl=1; my_wr=Ep%nomegasf 477 call setup_spectral(Ep%nomega,Ep%omega,Ep%nomegasf,omegasf,max_rest,min_rest,my_max_rest,my_min_rest,& 478 0,Ep%zcut,zero,my_wl,my_wr,kkweight) 479 480 ABI_CHECK(use_tr, 'Hilbert transform requires time-reversal') 481 482 ! allocate head and wings of the spectral function. 483 ABI_MALLOC(sf_head,(3,3,my_wl:my_wr)) 484 ABI_MALLOC(sf_lwing,(Ep%npwe,my_wl:my_wr,3)) 485 ABI_MALLOC(sf_uwing,(Ep%npwe,my_wl:my_wr,3)) 486 sf_head=czero; sf_lwing=czero; sf_uwing=czero 487 488 memreq = two*gwpc*Ep%npwe**2*(my_wr-my_wl+1)*b2Gb 489 write(msg,'(a,f10.4,a)')' memory required per spectral point: ',two*gwpc*Ep%npwe**2*b2Mb,' [Mb]' 490 call wrtout(std_out, msg) 491 write(msg,'(a,f10.4,a)')' memory required by sf_chi0q0: ',memreq,' [Gb]' 492 call wrtout(std_out, msg) 493 if (memreq > two) then 494 ABI_WARNING(' Memory required for sf_chi0q0 is larger than 2.0 Gb!') 495 end if 496 ABI_MALLOC_OR_DIE(sf_chi0,(Ep%npwe,Ep%npwe,my_wl:my_wr), ierr) 497 sf_chi0=czero_gw 498 499 CASE DEFAULT 500 ABI_BUG("Wrong spmeth") 501 END SELECT 502 503 nkpt_summed = Kmesh%nbz 504 if (Ep%symchi /= 0) then 505 nkpt_summed = Ltg_q%nibz_ltg 506 call Ltg_q%print(std_out, Dtset%prtvol) 507 end if 508 call wrtout(std_out, sjoin(' Calculation status: ', itoa(nkpt_summed), ' k-points to be completed')) 509 510 ABI_MALLOC(vkbr, (Kmesh%nibz)) 511 gradk_not_done = .TRUE. 512 513 ! ============================================ 514 ! === Begin big fat loop over transitions ==== 515 ! ============================================ 516 chi0 = czero_gw; chi0_sumrule = zero 517 518 ! Loop on spin to calculate $\chi_{\up,\up} + \chi_{\down,\down}$ 519 do spin=1,nsppol 520 if (ALL(bbp_ks_distrb(:,:,:,spin) /= Wfd%my_rank)) CYCLE 521 522 ! Loop over k-points in the BZ. 523 do ik_bz=1,Kmesh%nbz 524 if (Ep%symchi == 1) then 525 if (Ltg_q%ibzq(ik_bz) /= 1) CYCLE ! Only IBZ_q 526 end if 527 528 if (ALL(bbp_ks_distrb(:,:,ik_bz,spin) /= Wfd%my_rank)) CYCLE 529 530 write(msg,'(2(a,i4),a,i2,a,i3)')' ik= ',ik_bz,'/',Kmesh%nbz,' spin=',spin,' done by mpi rank:',Wfd%my_rank 531 call wrtout(std_out, msg) 532 533 ! Get ik_ibz, non-symmorphic phase and symmetries from ik_bz. 534 call kmesh%get_BZ_item(ik_bz,kbz,ik_ibz,isym_k,itim_k,ph_mkt) 535 tabr_k=ktabr(:,ik_bz) ! Table for rotated FFT points 536 spinrot_kbz(:)=Cryst%spinrot(:,isym_k) 537 if (Dtset%pawcross==1) tabrf_k(:) = ktabrf(:,ik_bz) 538 539 istwf_k = Wfd%istwfk(ik_ibz) 540 npw_k = Wfd%npwarr(ik_ibz) 541 kg_k => Wfd%Kdata(ik_ibz)%kg_k 542 543 if (psps%usepaw == 0 .and. Ep%inclvkb /= 0 .and. gradk_not_done(ik_ibz)) then 544 ! Include term <n,k|[Vnl,iqr]|n"k>' for q -> 0. 545 call vkbr_init(vkbr(ik_ibz), Cryst, Psps, Ep%inclvkb, istwf_k, npw_k, Kmesh%ibz(:,ik_ibz), kg_k) 546 gradk_not_done(ik_ibz) = .FALSE. 547 end if 548 549 ! Loop over "conduction" states. 550 do band1=1,Ep%nbnds 551 if (ALL(bbp_ks_distrb(band1,:,ik_bz,spin) /= Wfd%my_rank)) CYCLE 552 553 ABI_CHECK(wfd%get_wave_ptr(band1, ik_ibz, spin, wave1, msg) == 0, msg) 554 ug1 => wave1%ug 555 call wfd%get_ur(band1,ik_ibz,spin,ur1_kibz) 556 557 if (Psps%usepaw==1) then 558 call wfd%get_cprj(band1,ik_ibz,spin,Cryst,Cprj1_ibz,sorted=.FALSE.) 559 call pawcprj_copy(Cprj1_ibz,Cprj1_bz) 560 call paw_symcprj(ik_bz,nspinor,1,Cryst,Kmesh,Pawtab,Pawang,Cprj1_bz) 561 if (Dtset%pawcross==1) then 562 call wfdf%paw_get_aeur(band1,ik_ibz,spin,Cryst,Paw_onsite,Psps,Pawtab,Pawfgrtab,ur_ae1,ur_ae_onsite1,ur_ps_onsite1) 563 end if 564 end if 565 566 ! Loop over "valence" states. 567 do band2=1,Ep%nbnds 568 569 ! write(std_out,*) "ik,band1,band2",ik_bz,band1,band2 570 ! ----------------- cRPA for U 571 !debug if (.not.luwindow.AND.dtset%ucrpa==1.AND.band1<=ucrpa_bands(2).AND.band1>=ucrpa_bands(1)& 572 !debug& .AND.band2<=ucrpa_bands(2).AND.band2>=ucrpa_bands(1)) CYCLE 573 if (luwindow.AND.dtset%ucrpa==1 & 574 .AND.((ks_ebands%eig(band1,ik_ibz,spin)-ks_ebands%fermie)<=dtset%ucrpa_window(2)) & 575 .AND.((ks_ebands%eig(band1,ik_ibz,spin)-ks_ebands%fermie)>=dtset%ucrpa_window(1)) & 576 .AND.((ks_ebands%eig(band2,ik_ibz,spin)-ks_ebands%fermie)<=dtset%ucrpa_window(2)) & 577 .AND.((ks_ebands%eig(band2,ik_ibz,spin)-ks_ebands%fermie)>=dtset%ucrpa_window(1))) CYCLE 578 !----------------- cRPA for U 579 580 if (bbp_ks_distrb(band1,band2,ik_bz,spin) /= Wfd%my_rank) CYCLE 581 582 deltaeKS_b1b2 = ks_eig(band1, ik_ibz, spin) - ks_eig(band2, ik_ibz, spin) 583 deltaf_b1b2 = spin_fact * (qp_occ(band1, ik_ibz, spin) - qp_occ(band2, ik_ibz, spin)) 584 deltaeGW_b1b2 = qp_eig(band1, ik_ibz, spin) - qp_eig(band2, ik_ibz, spin) 585 586 if (Ep%gwcomp == 0) then 587 ! Skip negligible transitions. 588 if (abs(deltaf_b1b2) < GW_TOL_DOCC) CYCLE 589 else 590 ! when the completeness trick is used, we need to also consider transitions with vanishing deltaf 591 ! Rangel Correction for metals 592 if (qp_occ(band2,ik_ibz,spin) < GW_TOL_DOCC .and. ( ABS(deltaf_b1b2)< GW_TOL_DOCC .or. band1 < band2)) CYCLE 593 end if 594 595 ABI_CHECK(wfd%get_wave_ptr(band2, ik_ibz, spin, wave2, msg) == 0, msg) 596 ug2 => wave2%ug 597 call wfd%get_ur(band2,ik_ibz,spin,ur2_kibz) 598 599 if (Psps%usepaw==1) then 600 call wfd%get_cprj(band2,ik_ibz,spin,Cryst,Cprj2_ibz,sorted=.FALSE.) 601 call pawcprj_copy(Cprj2_ibz,Cprj2_bz) 602 call paw_symcprj(ik_bz,nspinor,1,Cryst,Kmesh,Pawtab,Pawang,Cprj2_bz) 603 if (Dtset%pawcross==1) then 604 call wfdf%paw_get_aeur(band2,ik_ibz,spin,Cryst,Paw_onsite,Psps,Pawtab,Pawfgrtab,& 605 ur_ae2,ur_ae_onsite2,ur_ps_onsite2) 606 end if 607 end if 608 609 SELECT CASE (Ep%spmeth) 610 CASE (0) 611 ! Adler-Wiser expression. 612 ! Add small imaginary of the Time-Ordered response function but only for non-zero real omega 613 ! FIXME What about metals? 614 615 if (.not. use_tr) then 616 ! Adler-Wiser without time-reversal. 617 do io=1,Ep%nomega 618 green_w(io) = g0g0w(Ep%omega(io), deltaf_b1b2, deltaeGW_b1b2, Ep%zcut, GW_TOL_W0, one_pole) 619 end do 620 621 else 622 if (Ep%gwcomp == 0) then ! cannot be completely skipped in case of completeness correction 623 if (band1 < band2) CYCLE ! Here we GAIN a factor ~2 624 end if 625 626 do io=1,Ep%nomega 627 ! Rangel: In metals, the intra-band transitions term does not contain the antiresonant part 628 ! if(abs(deltaeGW_b1b2)>GW_TOL_W0) green_w(io) = g0g0w(Ep%omega(io),deltaf_b1b2,deltaeGW_b1b2,Ep%zcut,GW_TOL_W0) 629 if (band1 == band2) green_w(io) = g0g0w(Ep%omega(io), deltaf_b1b2, deltaeGW_b1b2, Ep%zcut, GW_TOL_W0, one_pole) 630 if (band1 /= band2) green_w(io) = g0g0w(Ep%omega(io), deltaf_b1b2, deltaeGW_b1b2, Ep%zcut, GW_TOL_W0, two_poles) 631 632 if (Ep%gwcomp == 1) then 633 ! Calculate the completeness correction 634 numerator= -spin_fact * qp_occ(band2,ik_ibz,spin) 635 deltaeGW_enhigh_b2 = en_high - qp_eig(band2,ik_ibz,spin) 636 ! Completeness correction is NOT valid for real frequencies 637 if (REAL(Ep%omega(io)) < GW_TOL_W0) then 638 green_enhigh_w(io) = g0g0w(Ep%omega(io),numerator,deltaeGW_enhigh_b2,Ep%zcut,GW_TOL_W0,two_poles) 639 else 640 green_enhigh_w(io) = czero_gw 641 endif 642 ! 643 ! Rangel Correction for metals 644 if (band1 >= band2 .and. abs(deltaf_b1b2) > GW_TOL_DOCC) then 645 green_w(io)= green_w(io) - green_enhigh_w(io) 646 else 647 ! Disregard green_w, since it is already accounted for through the time-reversal 648 green_w(io)= - green_enhigh_w(io) 649 end if 650 end if !gwcomp == 1 651 end do !io 652 653 if (Ep%gwcomp == 1 .and. band1 == band2) then 654 ! Add the "delta part", symmetrization is done inside the routine. 655 call calc_wfwfg(tabr_k,itim_k,spinrot_kbz,nfft,nspinor,ngfft_gw,ur2_kibz,ur2_kibz,wfwfg) 656 657 if (Psps%usepaw==1) then 658 call paw_rho_tw_g(cryst, Pwij_fft, nfft,dim_rtwg,nspinor,gw_gfft,Cprj2_bz,Cprj2_bz,wfwfg) 659 660 ! Add PAW cross term 661 if (Dtset%pawcross==1) then 662 call paw_cross_rho_tw_g(nspinor,Ep%npwepG0,nfftf_tot,ngfftf,1,use_padfftf,igfftepsG0f,gboundf,& 663 ur_ae2,ur_ae_onsite2,ur_ps_onsite2,itim_k,tabrf_k,ph_mkt,spinrot_kbz,& 664 ur_ae2,ur_ae_onsite2,ur_ps_onsite2,itim_k,tabrf_k,ph_mkt,spinrot_kbz,& 665 dim_rtwg,wfwfg) 666 end if 667 end if 668 669 qzero = .TRUE. 670 call completechi0_deltapart(ik_bz,qzero,Ep%symchi,Ep%npwe,Gsph_FFT%ng,Ep%nomega,nspinor,& 671 nfft,ngfft_gw,gspfft_igfft,Gsph_FFT,Ltg_q,green_enhigh_w,wfwfg,chi0) 672 end if 673 end if ! use_tr 674 675 CASE (1, 2) 676 ! Spectral method, here time-reversal is always assumed. 677 if (deltaeGW_b1b2 < 0) CYCLE 678 call approxdelta(Ep%nomegasf,omegasf,deltaeGW_b1b2,Ep%spsmear,iomegal,iomegar,wl,wr,Ep%spmeth) 679 END SELECT 680 681 ! FFT of u^*_{b1,k}(r) u_{b2,k}(r) and (q,G=0) limit using small q and k.p perturbation theory 682 call rho_tw_g(nspinor,Ep%npwe,nfft,ndat1,ngfft_gw,1,use_padfft,igffteps0,gw_gbound,& 683 ur1_kibz,itim_k,tabr_k,ph_mkt,spinrot_kbz,& 684 ur2_kibz,itim_k,tabr_k,ph_mkt,spinrot_kbz,& 685 dim_rtwg,rhotwg) 686 687 if (psps%usepaw == 0) then 688 ! Matrix elements of i[H,r] for NC pseudopotentials. 689 rhotwx = nc_ihr_comm(vkbr(ik_ibz), cryst, psps, npw_k, nspinor, istwf_k, Ep%inclvkb, & 690 Kmesh%ibz(:,ik_ibz), ug1, ug2, kg_k) 691 692 else 693 ! 1) Add PAW onsite contribution, projectors are already in the BZ. 694 call paw_rho_tw_g(cryst, Pwij, Ep%npwe,dim_rtwg,nspinor,Gsph_epsG0%gvec,Cprj1_bz,Cprj2_bz,rhotwg) 695 696 ! 2) Matrix elements of i[H,r] for PAW. 697 rhotwx = paw_ihr(spin,nspinor,npw_k,istwf_k,Kmesh%ibz(:,ik_ibz),Cryst,Pawtab,ug1,ug2,kg_k,Cprj1_ibz,Cprj2_ibz,HUr) 698 699 ! Add PAW cross term 700 if (Dtset%pawcross==1) then 701 call paw_cross_rho_tw_g(nspinor,Ep%npwepG0,nfftf_tot,ngfftf,1,use_padfftf,igfftepsG0f,gboundf,& 702 ur_ae1,ur_ae_onsite1,ur_ps_onsite1,itim_k,tabrf_k,ph_mkt,spinrot_kbz,& 703 ur_ae2,ur_ae_onsite2,ur_ps_onsite2,itim_k,tabrf_k,ph_mkt,spinrot_kbz,& 704 dim_rtwg,rhotwg) 705 706 ! Add cross-term contribution to the commutator 707 if (Dtset%userib/=111) then 708 call paw_cross_ihr_comm(rhotwx,nspinor,nfftf_tot,Cryst,Pawfgrtab,Paw_onsite,& 709 ur_ae1,ur_ae2,ur_ae_onsite1,ur_ae_onsite2,Cprj1_ibz,Cprj2_ibz) 710 end if 711 end if 712 end if 713 714 ! Treat a possible degeneracy between v and c. 715 if (abs(deltaeKS_b1b2) > GW_TOL_W0) then 716 rhotwx = -rhotwx / deltaeKS_b1b2 717 else 718 rhotwx = czero_gw 719 end if 720 721 SELECT CASE (Ep%spmeth) 722 CASE (0) 723 ! ---------------- Ucrpa (begin) 724 if(dtset%ucrpa>=1.and..not.luwindow) then 725 fac=one 726 fac1=zero 727 fac2=zero 728 fac3=zero 729 fac4=one 730 m1=-1 731 m2=-1 732 if(dtset%ucrpa<=2) then 733 call flush_unit(std_out) 734 if ( band1<=ucrpa_bands(2).AND.band1>=ucrpa_bands(1)& 735 & .AND.band2<=ucrpa_bands(2).AND.band2>=ucrpa_bands(1)) then 736 ! if(dtset%prtvol>=10)write(6,*)"calculation is in progress",band1,band2,ucrpa_bands(1),ucrpa_bands(2) 737 if (dtset%plowan_compute>=10) then 738 band1c=band1-wan%bandi_wan+1 739 band2c=band2-wan%bandi_wan+1 740 do iat1=1,wan%natom_wan 741 do ispinor1=1,wan%nspinor 742 do il1=1,wan%nbl_atom_wan(iat1) 743 do m1=1,2*(wan%latom_wan(iat1)%lcalc(il1))+1 744 fac1=fac1 + real(wan%psichi(ik_bz,band1c,iat1)%atom(il1)%matl(m1,spin,ispinor1))*& 745 &conjg(wan%psichi(ik_bz,band1c,iat1)%atom(il1)%matl(m1,spin,ispinor1)) 746 fac2=fac2 + real(wan%psichi(ik_bz,band2c,iat1)%atom(il1)%matl(m1,spin,ispinor1))*& 747 &conjg(wan%psichi(ik_bz,band2c,iat1)%atom(il1)%matl(m1,spin,ispinor1)) 748 do iat2=1,wan%natom_wan 749 do ispinor2=1,wan%nspinor 750 do il2=1,wan%nbl_atom_wan(iat2) 751 do m2=1,2*(wan%latom_wan(iat2)%lcalc(il2))+1 752 fac=fac - real(wan%psichi(ik_bz,band1c,iat1)%atom(il1)%matl(m1,spin,ispinor1)*& 753 & conjg(wan%psichi(ik_bz,band1c,iat1)%atom(il1)%matl(m1,spin,ispinor1))*& 754 & wan%psichi(ik_bz,band2c,iat2)%atom(il2)%matl(m2,spin,ispinor2)*& 755 & conjg(wan%psichi(ik_bz,band2c,iat2)%atom(il2)%matl(m2,spin,ispinor2))) 756 fac3=fac3+ real(wan%psichi(ik_bz,band1c,iat1)%atom(il1)%matl(m1,spin,ispinor1))*& 757 &conjg(wan%psichi(ik_bz,band1c,iat1)%atom(il1)%matl(m1,spin,ispinor1))*& 758 &wan%psichi(ik_bz,band2c,iat2)%atom(il2)%matl(m2,spin,ispinor2)*& 759 &conjg(wan%psichi(ik_bz,band2c,iat2)%atom(il2)%matl(m2,spin,ispinor2)) 760 enddo !m2 761 enddo !il2 762 enddo !ispinor2 763 enddo !iat2 764 enddo !m1 765 enddo !il1 766 enddo !ispinor1 767 enddo !iat 768 else !plowan_compute 769 do iat=1, cryst%nattyp(itypatcor) 770 do ispinor1=1,nspinor 771 do m1=1,2*lcor+1 772 fac1=fac1+ real(coeffW_BZ(iat,spin,band1,ik_bz,ispinor1,m1)*& 773 &conjg(coeffW_BZ(iat,spin,band1,ik_bz,ispinor1,m1))) 774 fac2=fac2+ real(coeffW_BZ(iat,spin,band2,ik_bz,ispinor1,m1)*& 775 & conjg(coeffW_BZ(iat,spin,band2,ik_bz,ispinor1,m1))) 776 do ispinor2=1,nspinor 777 do m2=1,2*lcor+1 778 fac=fac - real(coeffW_BZ(iat,spin,band1,ik_bz,ispinor1,m1)*& 779 & conjg(coeffW_BZ(iat,spin,band1,ik_bz,ispinor1,m1))*& 780 & coeffW_BZ(iat,spin,band2,ik_bz,ispinor2,m2)*& 781 & conjg(coeffW_BZ(iat,spin,band2,ik_bz,ispinor2,m2))) 782 fac3=fac3 + real(coeffW_BZ(iat,spin,band1,ik_bz,ispinor1,m1)*& 783 & conjg(coeffW_BZ(iat,spin,band1,ik_bz,ispinor1,m1))* & 784 & coeffW_BZ(iat,spin,band2,ik_bz,ispinor2,m2)*& 785 & conjg(coeffW_BZ(iat,spin,band2,ik_bz,ispinor2,m2))) 786 ! if(dtset%prtvol>=10)write(6,*) fac,fac3 787 enddo !m2 788 enddo !ispinor2 789 ! if(dtset%prtvol>=10)write(6,*) fac,fac3,fac1,fac2,fac1*fac2 790 enddo !m1 791 enddo !ispinor1 792 enddo !iat 793 endif !plowan_compute>=10 794 fac4=fac 795 ! fac=zero 796 if(dtset%ucrpa==1) fac=zero 797 ! write(6,'(a,6i4,e15.5,a)') "*****FAC*********",ik_bz,ik_bz,band1,band2,m1,m2,fac," q==0" 798 endif 799 else if (dtset%ucrpa==3) then 800 if (band1<=ucrpa_bands(2).AND.band1>=ucrpa_bands(1)) then 801 do iat=1, cryst%nattyp(itypatcor) 802 do ispinor1=1,nspinor 803 do m1=1,2*lcor+1 804 fac2=fac2-real(coeffW_BZ(iat,spin,band1,ik_bz,ispinor1,m1)*& 805 & conjg(coeffW_BZ(iat,spin,band1,ik_bz,ispinor1,m1))) 806 enddo 807 enddo 808 enddo 809 if(dtset%ucrpa==4) fac2=zero 810 endif 811 if (band2<=ucrpa_bands(2).AND.band2>=ucrpa_bands(1)) then 812 do iat=1, cryst%nattyp(itypatcor) 813 do ispinor1=1,nspinor 814 do m1=1,2*lcor+1 815 fac3=fac3-real(coeffW_BZ(iat,spin,band2,ik_bz,ispinor1,m1)*& 816 & conjg(coeffW_BZ(iat,spin,band2,ik_bz,ispinor1,m1))) 817 enddo 818 enddo 819 enddo 820 if(dtset%ucrpa==4) fac3=zero 821 endif 822 fac=real(fac2*fac3) 823 endif 824 !if(dtset%prtvol>=10) write(6,'(6i4,e15.5,a)') ik_bz,ik_bz,band1,band2,m1,m2,fac," q==0" 825 !if(abs(fac-one)>0.00001) write(6,'(a,6i4,e15.5,a)') "*****FAC*********",ik_bz,ik_bz,band1,band2,m1,m2,fac," q==0" 826 ! if(dtset%prtvol>=10) write(6,'(a,6i4,e15.5,a)') "*****FAC*********",ik_bz,ik_bz,band1,band2,m1,m2,fac4," q==0" 827 green_w=green_w*fac 828 endif 829 ! ---------------- Ucrpa (end) 830 831 ! Adler-Wiser expression, to be consistent here we use the KS eigenvalues (?) 832 call accumulate_chi0_q0(is_metallic,ik_bz,isym_k,itim_k,Ep%gwcomp,nspinor,Ep%npwepG0,Ep,& 833 Cryst,Ltg_q,Gsph_epsG0,chi0,rhotwx,rhotwg,green_w,green_enhigh_w,deltaf_b1b2,chi0_head,chi0_lwing,chi0_uwing) 834 835 CASE (1, 2) 836 ! Spectral method, to be consistent here we use the KS eigenvalues. 837 call accumulate_sfchi0_q0(ik_bz,isym_k,itim_k,nspinor,Ep%symchi,Ep%npwepG0,Ep%npwe,Cryst,Ltg_q,& 838 Gsph_epsG0,deltaf_b1b2,my_wl,iomegal,wl,my_wr,iomegar,wr,rhotwx,rhotwg,Ep%nomegasf,& 839 sf_chi0,sf_head,sf_lwing,sf_uwing) 840 841 CASE DEFAULT 842 ABI_BUG("Wrong spmeth") 843 END SELECT 844 845 ! Accumulating the sum rule on chi0. Eq. (5.284) in G.D. Mahan Many-Particle Physics 3rd edition. [[cite:Mahan2000]] 846 factor = spin_fact * qp_occ(band2,ik_ibz,spin) 847 848 call accumulate_chi0sumrule(ik_bz,Ep%symchi,Ep%npwe,factor,deltaeGW_b1b2,& 849 Ltg_q,Gsph_epsG0,Ep%npwepG0,rhotwg,chi0_sumrule) 850 851 if (Ep%gwcomp == 1) then 852 ! Include also the completeness correction in the sum rule. 853 factor=-spin_fact*qp_occ(band2,ik_ibz,spin) 854 call accumulate_chi0sumrule(ik_bz,Ep%symchi,Ep%npwe,factor,deltaeGW_enhigh_b2,& 855 Ltg_q,Gsph_epsG0,Ep%npwepG0,rhotwg,chi0_sumrule) 856 if (band1 == Ep%nbnds) then 857 chi0_sumrule(:) = chi0_sumrule(:) + wtk_ltg(ik_bz)*spin_fact*qp_occ(band2,ik_ibz,spin)*deltaeGW_enhigh_b2 858 end if 859 end if 860 861 end do ! band2 862 end do ! band1 863 864 if (Psps%usepaw == 0 .and. Ep%inclvkb /= 0 .and. Ep%symchi == 1) then 865 call vkbr_free(vkbr(ik_ibz)) ! Not need anymore as we loop only over IBZ. 866 end if 867 end do !ik_bz 868 end do !spin 869 870 ABI_FREE(igffteps0) 871 872 call vkbr_free(vkbr) 873 ABI_FREE(vkbr) 874 875 ! === After big fat loop over transitions, now MPI === 876 ! * Master took care of the contribution in case of (metallic|spin) polarized systems. 877 select case (Ep%spmeth) 878 case (0) 879 ! Adler-Wiser expression 880 ! Sum contributions from each proc. 881 ! Looping on frequencies to avoid problems with the size of the MPI packet. 882 do io=1,Ep%nomega 883 call xmpi_sum(chi0(:,:,io),comm,ierr) 884 end do 885 886 case (1, 2) 887 ! Spectral method. 888 call hilbert_transform(Ep%npwe,Ep%nomega,Ep%nomegasf,my_wl,my_wr,kkweight,sf_chi0,chi0,Ep%spmeth) 889 890 ABI_SFREE(sf_chi0) 891 892 ! Sum contributions from each proc 893 ! Looping on frequencies to avoid problems with the size of the MPI packet 894 do io=1,Ep%nomega 895 call xmpi_sum(chi0(:,:,io),comm,ierr) 896 end do 897 898 call hilbert_transform_headwings(Ep%npwe,Ep%nomega,Ep%nomegasf,& 899 my_wl,my_wr,kkweight,sf_lwing,sf_uwing,sf_head,chi0_lwing,& 900 chi0_uwing,chi0_head,Ep%spmeth) 901 902 case default 903 ABI_BUG(sjoin("Wrong spmeth:", itoa(ep%spmeth))) 904 end select 905 906 ! Divide by the volume 907 !$OMP PARALLEL WORKSHARE 908 chi0 = chi0 * weight / Cryst%ucvol 909 !$OMP END PARALLEL WORKSHARE 910 911 ! Collect sum rule. pi comes from Im[1/(x-ieta)] = pi delta(x) 912 call xmpi_sum(chi0_sumrule, comm, ierr) 913 chi0_sumrule = chi0_sumrule * pi * weight / Cryst%ucvol 914 915 ! Collect head and wings. 916 call xmpi_sum(chi0_head, comm, ierr) 917 call xmpi_sum(chi0_lwing, comm, ierr) 918 call xmpi_sum(chi0_uwing, comm, ierr) 919 920 chi0_head = chi0_head * weight / cryst%ucvol 921 ! Tensor in terms of reciprocal lattice vectors. 922 do io=1,Ep%nomega 923 chi0_head(:,:,io) = matmul(chi0_head(:,:,io), cryst%gmet) * (two_pi**2) 924 end do 925 chi0_lwing = chi0_lwing * weight / cryst%ucvol 926 chi0_uwing = chi0_uwing * weight / cryst%ucvol 927 928 ! =============================================== 929 ! ==== Symmetrize chi0 in case of AFM system ==== 930 ! =============================================== 931 ! Reconstruct $chi0{\down,\down}$ from $chi0{\up,\up}$. 932 ! Works only in the case of magnetic group Shubnikov type IV. 933 if (cryst%use_antiferro) then 934 call symmetrize_afm_chi0(Cryst, Gsph_epsG0, Ltg_q, Ep%npwe, Ep%nomega, chi0=chi0, & 935 chi0_head=chi0_head, chi0_lwing=chi0_lwing, chi0_uwing=chi0_uwing) 936 end if 937 938 ! =================================================== 939 ! ==== Construct head and wings from the tensor ===== 940 ! =================================================== 941 do io=1,Ep%nomega 942 do ig=2,Ep%npwe 943 wng = chi0_uwing(ig,io,:) 944 chi0(1,ig,io) = vdotw(Ep%qlwl(:,1), wng, Cryst%gmet,"G") 945 wng = chi0_lwing(ig,io,:) 946 chi0(ig,1,io) = vdotw(Ep%qlwl(:,1), wng, Cryst%gmet,"G") 947 end do 948 chq = matmul(chi0_head(:,:,io), Ep%qlwl(:,1)) 949 chi0(1,1,io) = vdotw(Ep%qlwl(:,1), chq, Cryst%gmet,"G") ! Use user-defined small q 950 end do 951 952 if (wfd%my_rank == 0 .and. dtset%prtvol > 20) then 953 qlen = tol3 954 call cryst%get_redcart_qdirs(nq, qdirs, qlen=qlen) 955 ABI_MALLOC(head_qvals, (nq)) 956 call wrtout([std_out, ab_out], "Head of the irreducible polarizability for q --> 0", pre_newlines=1) 957 call wrtout([std_out, ab_out], sjoin(" q0_len:", ftoa(qlen), "(Bohr^-1)")) 958 write(msg, "(*(a14))") "omega_re (eV)", "omega_im (eV)", "[100]", "[010]", "[001]", "x", "y", "z" 959 call wrtout([std_out, ab_out], msg) 960 do io=1,Ep%nomega 961 do iq=1,nq 962 chq = matmul(chi0_head(:,:,io), qdirs(:,iq)) 963 head_qvals(iq) = vdotw(qdirs(:, iq), chq, cryst%gmet, "G") 964 end do 965 write(msg, "(*(es12.5,2x))") ep%omega(io) * Ha_eV, real(head_qvals(:)) 966 call wrtout([std_out, ab_out], msg) 967 ! Write imag part to std_out 968 write(msg, "(*(es12.5,2x))") ep%omega(io) * Ha_eV, aimag(head_qvals(:)) 969 call wrtout(std_out, msg) 970 end do 971 call wrtout([std_out, ab_out], " ") 972 ABI_FREE(qdirs) 973 ABI_FREE(head_qvals) 974 end if 975 !stop 976 977 ! Impose Hermiticity (valid only for zero or purely imaginary frequencies) 978 ! MG: what about metals, where we have poles around zero? 979 ! FB: because of the intraband term, chi0 is never hermitian in case of metals 980 if (.not. is_metallic) then 981 do io=1,Ep%nomega 982 if (ABS(REAL(Ep%omega(io))) < 0.00001) then 983 do ig2=1,Ep%npwe 984 do ig1=1,ig2-1 985 chi0(ig2,ig1,io)=GWPC_CONJG(chi0(ig1,ig2,io)) 986 end do 987 end do 988 end if 989 end do 990 end if 991 992 ! ===================== 993 ! ==== Free memory ==== 994 ! ===================== 995 ABI_FREE(bbp_ks_distrb) 996 ABI_FREE(rhotwg) 997 ABI_FREE(tabr_k) 998 ABI_FREE(ur1_kibz) 999 ABI_FREE(ur2_kibz) 1000 ABI_FREE(usr1_k) 1001 ABI_FREE(ur2_k) 1002 ABI_FREE(gw_gbound) 1003 1004 if (Dtset%pawcross==1) then 1005 ABI_FREE(gboundf) 1006 end if 1007 1008 ABI_SFREE(green_enhigh_w) 1009 ABI_SFREE(gw_gfft) 1010 ABI_SFREE(wfwfg) 1011 ABI_SFREE(kkweight) 1012 ABI_SFREE(omegasf) 1013 ABI_SFREE(green_w) 1014 ABI_SFREE(sf_head) 1015 ABI_SFREE(sf_lwing) 1016 ABI_SFREE(sf_uwing) 1017 ABI_SFREE(gspfft_igfft) 1018 1019 call Gsph_FFT%free() 1020 1021 if (Psps%usepaw==1) then 1022 ! deallocation for PAW. 1023 call pawcprj_free(Cprj1_bz) 1024 ABI_FREE(Cprj1_bz) 1025 call pawcprj_free(Cprj2_bz) 1026 ABI_FREE(Cprj2_bz) 1027 call pawcprj_free(Cprj1_ibz) 1028 ABI_FREE(Cprj1_ibz) 1029 call pawcprj_free(Cprj2_ibz) 1030 ABI_FREE(Cprj2_ibz) 1031 call pawpwij_free(Pwij) 1032 ABI_FREE(Pwij) 1033 if (allocated(Pwij_fft)) then 1034 call pawpwij_free(Pwij_fft) 1035 ABI_FREE(Pwij_fft) 1036 end if 1037 call pawhur_free(Hur) 1038 ABI_FREE(Hur) 1039 if (Dtset%pawcross==1) then 1040 ABI_FREE(ur_ae1) 1041 ABI_FREE(ur_ae_onsite1) 1042 ABI_FREE(ur_ps_onsite1) 1043 ABI_FREE(ur_ae2) 1044 ABI_FREE(ur_ae_onsite2) 1045 ABI_FREE(ur_ps_onsite2) 1046 ABI_FREE(tabrf_k) 1047 ABI_FREE(gboundf) 1048 ABI_FREE(igfftepsG0f) 1049 end if 1050 end if 1051 1052 if(dtset%ucrpa >= 1 .AND. dtset%plowan_compute < 10) then 1053 ABI_FREE(coeffW_BZ) 1054 endif 1055 1056 call cwtime(cpu_time, wall_time, gflops, "stop") 1057 write(std_out,'(2(a,f9.1))')" cpu_time = ",cpu_time,", wall_time = ",wall_time 1058 1059 DBG_EXIT("COLL") 1060 1061 end subroutine cchi0q0
ABINIT/chi0q0_intraband [ Functions ]
NAME
chi0q0_intraband
FUNCTION
Calculate chi0 in the limit q-->0
INPUTS
use_tr=If .TRUE. Wfs_val are allocate and only resonant transitions are evaluated (assumes time reversal symmetry) Ep= datatype gathering differening parameters related to the calculation of the inverse dielectric matrix Gsph_epsG0<gvectors_data_type>: Info on the G-sphere used to describe chi0/espilon (including umklapp) %ng=number of G vectors %rottbm1(ng,2,nsym)=contains the index (IS^{-1}) G in the array gvec %phmGt(ng,nsym)=phase factor e^{-iG.\tau} needed to symmetrize oscillator matrix elements and chi0 %gmet(3,3)=reciprocal space metric ($\textrm{bohr}^{-2}$). %gprimd(3,3)=dimensional reciprocal space primitive translations (b^-1) Ep%nbnds=number of bands ngfft_gw(18)= array containing all the information for 3D FFT for the oscillator strengths. Ep%nomega=number of frequencies Cryst<crystal_t>= data type gathering info on symmetries and unit cell %natom=number of atoms %nsym=number of symmetry operations %symrec(3,3,nsym)=symmetry operations in reciprocal space %typat(natom)=type of each atom %xred(3,natom)=reduced coordinated of atoms %rprimd(3,3)=dimensional primitive translations in real space (bohr) %timrev=2 if time-reversal symmetry can be used, 1 otherwise Ep%npwe=number of planewaves for sigma exchange (input variable) Ep%nsppol=1 for unpolarized, 2 for spin-polarized Ep%omega(Ep%nomega)=frequencies Psps <type(pseudopotential_type)>=variables related to pseudopotentials %mpsang=1+maximum angular momentum for nonlocal pseudopotential Pawang<pawang_type> angular mesh discretization and related data: Pawrad(ntypat*usepaw)<Pawrad_type>=paw radial mesh and related data Paw_ij(natom*usepaw)<Paw_ij_type)>=paw arrays given on (i,j) channels BSt<ebands_t>=Quasiparticle energies and occupations (for the moment real quantities) %mband=MAX number of bands over k-points and spin (==Ep%nbnds) %occ(mband,nkpt,nsppol)=QP occupation numbers, for each k point in IBZ, and each band %eig(mband,nkpt,nsppol)=GW energies, for self-consistency purposes Paw_pwff<pawpwff_t>=Form factor used to calculate the onsite mat. elements of a plane wave.
OUTPUT
chi0(Ep%npwe,Ep%npwe,Ep%nomega)=independent-particle susceptibility matrix for wavevector qq, and frequencies defined by Ep%omega
NOTES
*) The terms "head", "wings" and "body" of chi(G,Gp) refer to G=Gp=0, either G or Gp=0, and neither=0 respectively
TODO
Check npwepG0 before Switching on umklapp
SOURCE
1984 subroutine chi0q0_intraband(Wfd,Cryst,Ep,Psps,BSt,Gsph_epsG0,Pawang,Pawrad,Pawtab,Paw_ij,Paw_pwff,use_tr,usepawu,& 1985 & ngfft_gw,chi0,chi0_head,chi0_lwing,chi0_uwing) 1986 1987 !Arguments ------------------------------------ 1988 !scalars 1989 integer,intent(in) :: usepawu 1990 logical,intent(in) :: use_tr 1991 type(ebands_t),intent(in) :: BSt 1992 type(crystal_t),intent(in) :: Cryst 1993 type(em1params_t),intent(in) :: Ep 1994 type(gsphere_t),intent(in) :: Gsph_epsG0 1995 type(Pseudopotential_type),intent(in) :: Psps 1996 type(Pawang_type),intent(in) :: Pawang 1997 type(wfdgw_t),target,intent(inout) :: Wfd 1998 !arrays 1999 integer,intent(in) :: ngfft_gw(18) 2000 complex(gwpc),intent(out) :: chi0(Ep%npwe*Ep%nI,Ep%npwe*Ep%nJ,Ep%nomega) 2001 complex(dpc),intent(out) :: chi0_lwing(Ep%npwe*Ep%nI,Ep%nomega,3) 2002 complex(dpc),intent(out) :: chi0_uwing(Ep%npwe*Ep%nJ,Ep%nomega,3) 2003 complex(dpc),intent(out) :: chi0_head(3,3,Ep%nomega) 2004 type(Pawrad_type),intent(in) :: Pawrad(Psps%ntypat*Psps%usepaw) 2005 type(Pawtab_type),intent(in) :: Pawtab(Psps%ntypat*Psps%usepaw) 2006 type(Paw_ij_type),intent(in) :: Paw_ij(Cryst%natom*Psps%usepaw) 2007 type(pawpwff_t),intent(in) :: Paw_pwff(Psps%ntypat*Psps%usepaw) 2008 2009 !Local variables ------------------------------ 2010 !scalars 2011 integer,parameter :: tim_fourdp1=1,two_poles=2,one_pole=1,ndat1=1 2012 integer,parameter :: unitdos0=0,option1=1,NOMEGA_PRINTED=15 2013 integer :: nqlwl,nband_k,iomega,istwf_k,npw_k,my_nband,lbidx 2014 integer :: band,itim_k,ik_bz,ik_ibz,io,isym_k,spin,iqlwl !ig,ig1,ig2,my_nbbp,my_nbbpks 2015 integer :: nkpt_summed,dim_rtwg,use_padfft,gw_fftalga,ifft 2016 integer :: kptopt,isym,nsppol,nspinor 2017 integer :: comm,ierr,gw_mgfft,use_umklp,inclvkb 2018 real(dp) :: spin_fact,deltaf_b1b2,weight 2019 real(dp) :: deltaeGW_b1b2,zcut 2020 real(dp),parameter :: dummy_dosdeltae=HUGE(zero) 2021 real(dp) :: o_entropy,o_nelect,maxocc 2022 complex(dpc) :: ph_mkt 2023 logical :: iscompatibleFFT,is_metallic !,ltest 2024 character(len=500) :: msg,msg_tmp !,allup 2025 type(kmesh_t) :: Kmesh 2026 type(littlegroup_t) :: Ltg_q 2027 type(vkbr_t) :: vkbr 2028 type(wave_t),pointer :: wave 2029 !arrays 2030 integer :: my_band_list(Wfd%mband) 2031 integer,ABI_CONTIGUOUS pointer :: kg_k(:,:) 2032 integer,allocatable :: ktabr(:,:),irottb(:,:) 2033 !integer :: got(Wfd%nproc) 2034 integer,allocatable :: tabr_k(:),igffteps0(:),gw_gbound(:,:) 2035 real(dp),parameter :: q0(3)=(/zero,zero,zero/) 2036 real(dp) :: kpt(3),dedk(3),kbz(3),spinrot_kbz(4) 2037 !real(dp),ABI_CONTIGUOUS pointer :: ks_eig(:,:,:),qp_eig(:,:,:),qp_occ(:,:,:) 2038 real(dp) :: shift_ene(BSt%mband,BSt%nkpt,BSt%nsppol) 2039 real(dp) :: delta_occ(BSt%mband,BSt%nkpt,BSt%nsppol) 2040 !real(dp) :: eigen_vec(BSt%bantot) 2041 real(dp) :: o_doccde(BSt%bantot) 2042 real(dp) :: eigen_pdelta_vec(BSt%bantot),eigen_mdelta_vec(BSt%bantot) 2043 real(dp) :: o_occ_pdelta(BSt%bantot),o_occ_mdelta(BSt%bantot) 2044 real(dp) :: delta_ene(BSt%mband,BSt%nkpt,BSt%nsppol) 2045 real(dp) :: test_docc(BSt%mband,BSt%nkpt,BSt%nsppol) 2046 real(dp),allocatable :: qlwl(:,:) 2047 complex(gwpc) :: comm_kbbs(3,Wfd%nspinor**2) 2048 complex(dpc),allocatable :: ihr_comm(:,:,:,:,:) 2049 complex(gwpc),allocatable :: rhotwg(:) 2050 complex(dpc) :: green_w(Ep%nomega) 2051 complex(gwpc),allocatable :: ur1(:) 2052 complex(gwpc),ABI_CONTIGUOUS pointer :: ug(:) 2053 logical :: bmask(Wfd%mband) 2054 type(pawcprj_type),allocatable :: Cprj1_bz(:,:),Cprj1_ibz(:,:),Cp_bks(:,:) 2055 type(pawpwij_t),allocatable :: Pwij(:) 2056 type(pawhur_t),allocatable :: Hur(:) 2057 2058 !************************************************************************ 2059 2060 DBG_ENTER("COLL") 2061 2062 nsppol = Wfd%nsppol 2063 nspinor = Wfd%nspinor 2064 is_metallic = ebands_has_metal_scheme(BSt) 2065 2066 gw_mgfft = MAXVAL(ngfft_gw(1:3)) 2067 gw_fftalga = ngfft_gw(7)/100 !; gw_fftalgc=MOD(ngfft_gw(7),10) 2068 2069 ! Calculate <k,b1|i[H,r]|k',b2>. 2070 inclvkb=2; if (Wfd%usepaw==1) inclvkb=0 2071 ABI_MALLOC(ihr_comm,(3,nspinor**2,Wfd%mband,Wfd%nkibz,nsppol)) 2072 ihr_comm = czero 2073 2074 if (Wfd%usepaw==1) then 2075 ABI_MALLOC(Cp_bks,(Cryst%natom,nspinor)) 2076 call pawcprj_alloc(Cp_bks,0,Wfd%nlmn_atm) 2077 ABI_MALLOC(HUr,(Cryst%natom)) 2078 if (usepawu/=0) then ! For PAW+DFT+U, precalculate <\phi_i|[Hu,r]|phi_j\>. 2079 call pawhur_init(hur,nsppol,Wfd%pawprtvol,Cryst,Pawtab,Pawang,Pawrad,Paw_ij) 2080 end if 2081 end if 2082 2083 do spin=1,nsppol 2084 do ik_ibz=1,Wfd%nkibz 2085 npw_k = Wfd%npwarr(ik_ibz) 2086 nband_k= Wfd%nband(ik_ibz,spin) 2087 kpt = Wfd%kibz(:,ik_ibz) 2088 kg_k => Wfd%Kdata(ik_ibz)%kg_k 2089 istwf_k = Wfd%istwfk(ik_ibz) 2090 ABI_CHECK(istwf_k==1,"istwf_k/=1 not coded") 2091 2092 ! Distribute bands. 2093 bmask=.FALSE.; bmask(1:nband_k)=.TRUE. ! TODO only bands around EF should be included. 2094 call wfd%distribute_bands(ik_ibz,spin,my_nband,my_band_list,bmask=bmask) 2095 if (my_nband==0) CYCLE 2096 2097 if (Wfd%usepaw==0.and.inclvkb/=0) then ! Include term <n,k|[Vnl,iqr]|n"k>' for q->0. 2098 call vkbr_init(vkbr,Cryst,Psps,inclvkb,istwf_k,npw_k,kpt,kg_k) 2099 end if 2100 2101 do lbidx=1,my_nband 2102 band=my_band_list(lbidx) 2103 2104 ABI_CHECK(wfd%get_wave_ptr(band, ik_ibz, spin, wave, msg) == 0, msg) 2105 ug => wave%ug 2106 2107 if (Wfd%usepaw==0) then 2108 ! Matrix elements of i[H,r] for NC pseudopotentials. 2109 comm_kbbs = nc_ihr_comm(vkbr,cryst,psps,npw_k,nspinor,istwf_k,inclvkb,Kmesh%ibz(:,ik_ibz),ug,ug,kg_k) 2110 else 2111 ! Matrix elements of i[H,r] for PAW. 2112 call wfd%get_cprj(band,ik_ibz,spin,Cryst,Cp_bks,sorted=.FALSE.) 2113 comm_kbbs = paw_ihr(spin,nspinor,npw_k,istwf_k,Kmesh%ibz(:,ik_ibz),Cryst,Pawtab,ug,ug,kg_k,Cp_bks,Cp_bks,HUr) 2114 end if 2115 2116 ihr_comm(:,:,band,ik_ibz,spin) = comm_kbbs 2117 end do 2118 2119 call vkbr_free(vkbr) ! Not need anymore as we loop only over IBZ. 2120 end do 2121 end do 2122 ! 2123 ! Gather the commutator on each node. 2124 call xmpi_sum(ihr_comm,Wfd%comm,ierr) 2125 2126 if (Wfd%usepaw==1) then 2127 call pawcprj_free(Cp_bks) 2128 ABI_FREE(Cp_bks) 2129 call pawhur_free(Hur) 2130 ABI_FREE(Hur) 2131 end if 2132 2133 nqlwl=1 2134 ABI_MALLOC(qlwl,(3,nqlwl)) 2135 !qlwl = GW_Q0_DEFAULT(3) 2136 qlwl(:,1) = (/0.00001_dp, 0.00002_dp, 0.00003_dp/) 2137 ! 2138 write(msg,'(a,i3,a)')' Q-points for long wave-length limit in chi0q_intraband. # ',nqlwl,ch10 2139 do iqlwl=1,nqlwl 2140 write(msg_tmp,'(1x,i5,a,2x,3f12.6,a)') iqlwl,')',qlwl(:,iqlwl),ch10 2141 msg=TRIM(msg)//msg_tmp 2142 end do 2143 call wrtout(std_out, msg) 2144 ! 2145 ! delta_ene = e_{b,k-q} - e_{b,k} = -q. <b,k| i[H,r] |b,k> + O(q^2). 2146 delta_ene = zero 2147 do spin=1,nsppol 2148 do ik_ibz=1,Wfd%nkibz 2149 do band=1,Wfd%nband(ik_ibz,spin) 2150 dedk = REAL(ihr_comm(:,1,band,ik_ibz,spin)) 2151 delta_ene(band,ik_ibz,spin) = -vdotw(qlwl(:,1),dedk,Cryst%gmet,"G") 2152 end do 2153 end do 2154 end do 2155 2156 maxocc=two/(nsppol*nspinor) 2157 2158 ! Calculate the occupations at f(e+delta/2). 2159 shift_ene = BSt%eig + half*delta_ene 2160 2161 call pack_eneocc(BSt%nkpt,BSt%nsppol,BSt%mband,BSt%nband,BSt%bantot,shift_ene,eigen_pdelta_vec) 2162 2163 if (BSt%occopt < 9) then 2164 call getnel(o_doccde,dummy_dosdeltae,eigen_pdelta_vec,o_entropy,BSt%fermie,BSt%fermie,maxocc,BSt%mband,BSt%nband,& 2165 & o_nelect,BSt%nkpt,BSt%nsppol,o_occ_pdelta,BSt%occopt,option1,BSt%tphysel,BSt%tsmear,unitdos0,BSt%wtk,1,BSt%nband(1)) 2166 else 2167 ABI_ERROR('occopt 9 not implemented for GW calculations') 2168 end if 2169 ! 2170 ! Calculate the occupations at f(e-delta/2). 2171 shift_ene = BSt%eig - half*delta_ene 2172 2173 call pack_eneocc(BSt%nkpt,BSt%nsppol,BSt%mband,BSt%nband,BSt%bantot,shift_ene,eigen_mdelta_vec) 2174 2175 if (BSt%occopt < 9) then 2176 call getnel(o_doccde,dummy_dosdeltae,eigen_mdelta_vec,o_entropy,BSt%fermie,BSt%fermie,maxocc,BSt%mband,BSt%nband,& 2177 & o_nelect,BSt%nkpt,BSt%nsppol,o_occ_mdelta,BSt%occopt,option1,BSt%tphysel,BSt%tsmear,unitdos0,BSt%wtk,1,BSt%nband(1)) 2178 write(std_out,*)"nelect2: ",o_nelect 2179 else 2180 ABI_ERROR("occopt 9 not implemented for GW calculations") 2181 end if 2182 2183 ! f(e-delta/2) - f(e+delta/2). 2184 o_occ_pdelta = o_occ_mdelta - o_occ_pdelta 2185 2186 call unpack_eneocc(BSt%nkpt,BSt%nsppol,BSt%mband,BSt%nband,o_occ_pdelta,delta_occ) 2187 ! 2188 ! Expand f(e-delta/2) - f(e+delta/2) up to the first order in the small q. 2189 do spin=1,nsppol 2190 do ik_ibz=1,Wfd%nkibz 2191 do band=1,Wfd%nband(ik_ibz,spin) 2192 dedk = REAL(ihr_comm(:,1,band,ik_ibz,spin)) 2193 test_docc(band,ik_ibz,spin) = +vdotw(qlwl(:,1),dedk,Cryst%gmet,"G") * BSt%doccde(band,ik_ibz,spin) 2194 write(std_out,'(a,3(i0,1x),1x,3es16.8)')" spin,ik_ibz,band, delta_occ: ",& 2195 & spin,ik_ibz,band,delta_occ(band,ik_ibz,spin),& 2196 & test_docc(band,ik_ibz,spin),delta_occ(band,ik_ibz,spin)-test_docc(band,ik_ibz,spin) 2197 end do 2198 end do 2199 end do 2200 2201 ! ABI_ERROR("DONE") 2202 ! do spin=1,nsppol 2203 ! do ik_ibz=1,Wfd%nkibz 2204 ! nband_k = Wfd%nband(ik_ibz,spin) 2205 ! do band=1,nband_k 2206 ! write(std_out,'(a,3i3,2es14.6)')" spin, band, ik_ibz, delta_ene, delta_occ ",& 2207 !& spin,band,ik_ibz,delta_ene(band,ik_ibz,spin),delta_occ(band,ik_ibz,spin) 2208 ! end do 2209 ! end do 2210 ! end do 2211 2212 ABI_FREE(ihr_comm) 2213 ABI_FREE(qlwl) 2214 2215 if ( ANY(ngfft_gw(1:3) /= Wfd%ngfft(1:3)) ) call wfd%change_ngfft(Cryst,Psps,ngfft_gw) 2216 2217 ! TODO take into account the case of random k-meshes. 2218 kptopt=3 2219 call Kmesh%init(Cryst,Wfd%nkibz,Wfd%kibz,kptopt) 2220 ! 2221 !=== Get the FFT index of $ (R^{-1}(r-\tau)) $ === 2222 !* S= $\transpose R^{-1}$ and k_BZ = S k_IBZ 2223 !* irottb is the FFT index of $ R^{-1} (r-\tau) $ used to symmetrize u_Sk. 2224 ABI_MALLOC(irottb,(Wfd%nfftot,Cryst%nsym)) 2225 2226 call rotate_FFT_mesh(Cryst%nsym,Cryst%symrel,Cryst%tnons,Wfd%ngfft,irottb,iscompatibleFFT) 2227 ABI_CHECK(iscompatibleFFT,"FFT mesh not compatible with symmetries") 2228 2229 ABI_MALLOC(ktabr,(Wfd%nfftot,Kmesh%nbz)) 2230 do ik_bz=1,Kmesh%nbz 2231 isym=Kmesh%tabo(ik_bz) 2232 do ifft=1,Wfd%nfftot 2233 ktabr(ifft,ik_bz)=irottb(ifft,isym) 2234 end do 2235 end do 2236 ABI_FREE(irottb) 2237 ! 2238 ! === Setup weight (2 for spin unpolarized systems, 1 for polarized) === 2239 ! spin_fact is used to normalize the occupation factors to one. 2240 ! Consider also the AFM case. 2241 SELECT CASE (nsppol) 2242 CASE (1) 2243 weight=two/Kmesh%nbz; spin_fact=half 2244 if (Wfd%nspden==2) then 2245 weight=one/Kmesh%nbz; spin_fact=half 2246 end if 2247 if (nspinor==2) then 2248 weight=one/Kmesh%nbz; spin_fact=one 2249 end if 2250 CASE (2) 2251 weight=one/Kmesh%nbz; spin_fact=one 2252 CASE DEFAULT 2253 ABI_BUG("Wrong nsppol") 2254 END SELECT 2255 2256 use_umklp = 0 2257 call Ltg_q%init(q0, Kmesh%nbz, Kmesh%bz, Cryst, use_umklp, Ep%npwepG0, gvec=Gsph_epsG0%gvec) 2258 2259 write(msg,'(a,i2)')' Using symmetries to sum only over the IBZ_q = ',Ep%symchi 2260 call wrtout(std_out, msg) 2261 ! 2262 ! Evaluate oscillator matrix elements btw partial waves. Note that q=Gamma is used. 2263 if (Psps%usepaw==1) then 2264 ABI_MALLOC(Pwij,(Psps%ntypat)) 2265 call pawpwij_init(Pwij,Ep%npwepG0, [zero,zero,zero], Gsph_epsG0%gvec,Cryst%rprimd,Psps,Pawtab,Paw_pwff) 2266 2267 ABI_MALLOC(Cprj1_bz ,(Cryst%natom,nspinor)) 2268 call pawcprj_alloc(Cprj1_bz, 0,Wfd%nlmn_atm) 2269 ABI_MALLOC(Cprj1_ibz,(Cryst%natom,nspinor)) 2270 call pawcprj_alloc(Cprj1_ibz,0,Wfd%nlmn_atm) 2271 end if 2272 2273 ABI_MALLOC(rhotwg,(Ep%npwe*nspinor**2)) 2274 ABI_MALLOC(tabr_k,(Wfd%nfftot)) 2275 ABI_MALLOC(ur1,(Wfd%nfft*nspinor)) 2276 ! 2277 ! Tables for the FFT of the oscillators. 2278 ! a) FFT index of the G sphere (only vertical transitions, unlike cchi0, no need to shift the sphere). 2279 ! b) gw_gbound table for the zero-padded FFT performed in rhotwg. 2280 ABI_MALLOC(gw_gbound,(2*gw_mgfft+8,2)) 2281 ABI_MALLOC(igffteps0,(Gsph_epsG0%ng)) 2282 2283 call Gsph_epsG0%fft_tabs([0, 0, 0], gw_mgfft,ngfft_gw,use_padfft,gw_gbound,igffteps0) 2284 if ( ANY(gw_fftalga == [2, 4]) ) use_padfft=0 ! Pad-FFT is not coded in rho_tw_g 2285 if (use_padfft==0) then 2286 ABI_FREE(gw_gbound) 2287 ABI_MALLOC(gw_gbound,(2*gw_mgfft+8,2*use_padfft)) 2288 end if 2289 2290 nkpt_summed=Kmesh%nbz 2291 if (Ep%symchi/=0) then 2292 nkpt_summed=Ltg_q%nibz_ltg 2293 call Ltg_q%print(std_out, Wfd%prtvol) 2294 end if 2295 ! 2296 ! ============================================ 2297 ! === Begin big fat loop over transitions ==== 2298 ! ============================================ 2299 chi0 = czero_gw 2300 chi0_head = czero_gw; chi0_lwing = czero_gw; chi0_uwing = czero_gw 2301 dim_rtwg=1; if (nspinor==2) dim_rtwg=2 !can reduce size depending on Ep%nI and Ep%nj 2302 2303 zcut = Ep%zcut 2304 zcut = 0.1/Ha_eV 2305 write(std_out,*)" using zcut ",zcut*Ha_eV," [eV]" 2306 2307 ! Loop on spin to calculate $ \chi_{\up,\up} + \chi_{\down,\down} 2308 do spin=1,nsppol 2309 ! Loop over k-points in the BZ. 2310 do ik_bz=1,Kmesh%nbz 2311 if (Ep%symchi==1) then 2312 if (Ltg_q%ibzq(ik_bz)/=1) CYCLE ! Only IBZ_q 2313 end if 2314 2315 ! Get ik_ibz, non-symmorphic phase and symmetries from ik_bz. 2316 call kmesh%get_BZ_item(ik_bz,kbz,ik_ibz,isym_k,itim_k,ph_mkt) 2317 tabr_k=ktabr(:,ik_bz) ! Table for rotated FFT points 2318 spinrot_kbz(:)=Cryst%spinrot(:,isym_k) 2319 nband_k=Wfd%nband(ik_ibz,spin) 2320 2321 ! Distribute bands. 2322 bmask=.FALSE.; bmask(1:nband_k)=.TRUE. ! TODO only bands around EF should be included. 2323 call wfd%distribute_bands(ik_ibz,spin,my_nband,my_band_list,bmask=bmask) 2324 if (my_nband==0) CYCLE 2325 2326 write(msg,'(2(a,i4),a,i2,a,i3)')' ik = ',ik_bz,' / ',Kmesh%nbz,' spin = ',spin,' done by processor ',Wfd%my_rank 2327 call wrtout(std_out, msg) 2328 2329 do lbidx=1,my_nband 2330 ! Loop over bands treated by this node. 2331 band=my_band_list(lbidx) 2332 call wfd%get_ur(band,ik_ibz,spin,ur1) 2333 2334 if (Psps%usepaw==1) then 2335 call wfd%get_cprj(band,ik_ibz,spin,Cryst,Cprj1_ibz,sorted=.FALSE.) 2336 call pawcprj_copy(Cprj1_ibz,Cprj1_bz) 2337 call paw_symcprj(ik_bz,nspinor,1,Cryst,Kmesh,Pawtab,Pawang,Cprj1_bz) 2338 end if 2339 2340 deltaf_b1b2 = spin_fact*delta_occ(band,ik_ibz,spin) 2341 deltaeGW_b1b2= delta_ene(band,ik_ibz,spin) 2342 2343 ! Add small imaginary of the Time-Ordered resp function but only for non-zero real omega FIXME What about metals? 2344 if (.not.use_tr) then 2345 do io=1,Ep%nomega 2346 !green_w(io) = g0g0w(Ep%omega(io),deltaf_b1b2,deltaeGW_b1b2,zcut,-one,one_pole) 2347 green_w(io) = g0g0w(Ep%omega(io),deltaf_b1b2,deltaeGW_b1b2,zcut,GW_TOL_W0,one_pole) 2348 end do 2349 else 2350 do io=1,Ep%nomega ! This expression implements time-reversal even when the input k-mesh breaks it. 2351 !green_w(io) = half * g0g0w(Ep%omega(io),deltaf_b1b2,deltaeGW_b1b2,zcut,-one,two_poles) 2352 green_w(io) = half * g0g0w(Ep%omega(io),deltaf_b1b2,deltaeGW_b1b2,zcut,GW_TOL_W0,two_poles) 2353 end do !io 2354 end if ! use_tr 2355 2356 ! FFT of u^*_{b1,k}(r) u_{b2,k}(r). 2357 call rho_tw_g(nspinor,Ep%npwe,Wfd%nfft,ndat1,ngfft_gw,1,use_padfft,igffteps0,gw_gbound,& 2358 & ur1,itim_k,tabr_k,ph_mkt,spinrot_kbz,& 2359 & ur1,itim_k,tabr_k,ph_mkt,spinrot_kbz,& 2360 & dim_rtwg,rhotwg) 2361 2362 if (Psps%usepaw==1) then 2363 ! Add PAW onsite contribution, projectors are already in the BZ. 2364 call paw_rho_tw_g(cryst,Pwij,Ep%npwe,dim_rtwg,nspinor,Gsph_epsG0%gvec,Cprj1_bz,Cprj1_bz,rhotwg) 2365 end if 2366 2367 ! ==== Adler-Wiser expression, to be consistent here we use the KS eigenvalues (?) ==== 2368 call assemblychi0_sym(is_metallic,ik_bz,nspinor,Ep,Ltg_q,green_w,Ep%npwepG0,rhotwg,Gsph_epsG0,chi0) 2369 end do !band 2370 end do !ik_bz 2371 end do !spin 2372 2373 ! Collect body, head and wings within comm 2374 comm=Wfd%comm 2375 do io=1,Ep%nomega 2376 call xmpi_sum(chi0(:,:,io),comm,ierr) 2377 end do 2378 call xmpi_sum(chi0_head,comm,ierr) 2379 call xmpi_sum(chi0_lwing,comm,ierr) 2380 call xmpi_sum(chi0_uwing,comm,ierr) 2381 2382 ! Divide by the volume 2383 chi0 = chi0 * weight/Cryst%ucvol 2384 chi0_head = chi0_head * weight/Cryst%ucvol 2385 do io=1,Ep%nomega ! Tensor in the basis of the reciprocal lattice vectors. 2386 chi0_head(:,:,io) = MATMUL(chi0_head(:,:,io),Cryst%gmet) * (two_pi**2) 2387 end do 2388 chi0_lwing = chi0_lwing * weight/Cryst%ucvol 2389 chi0_uwing = chi0_uwing * weight/Cryst%ucvol 2390 ! 2391 ! =============================================== 2392 ! ==== Symmetrize chi0 in case of AFM system ==== 2393 ! =============================================== 2394 ! * Reconstruct $chi0{\down,\down}$ from $chi0{\up,\up}$. 2395 ! * Works only in the case of magnetic group Shubnikov type IV. 2396 if (Cryst%use_antiferro) then 2397 call symmetrize_afm_chi0(Cryst, Gsph_epsG0, Ltg_q, Ep%npwe, Ep%nomega, chi0=chi0, & 2398 chi0_head=chi0_head, chi0_lwing=chi0_lwing, chi0_uwing=chi0_uwing) 2399 end if 2400 ! 2401 ! ================================================== 2402 ! ==== Construct head and wings from the tensor ==== 2403 ! ================================================== 2404 !do io=1,Ep%nomega 2405 ! do ig=2,Ep%npwe 2406 ! wng = chi0_uwing(ig,io,:) 2407 ! chi0(1,ig,io) = vdotw(Ep%qlwl(:,1),wng,Cryst%gmet,"G") 2408 ! wng = chi0_lwing(ig,io,:) 2409 ! chi0(ig,1,io) = vdotw(Ep%qlwl(:,1),wng,Cryst%gmet,"G") 2410 ! end do 2411 ! chq = MATMUL(chi0_head(:,:,io), Ep%qlwl(:,1)) 2412 ! chi0(1,1,io) = vdotw(Ep%qlwl(:,1),chq,Cryst%gmet,"G") ! Use user-defined small q 2413 !end do 2414 !call wfd_barrier(Wfd) 2415 2416 ! Impose Hermiticity (valid only for zero or purely imaginary frequencies) 2417 ! MG what about metals, where we have poles around zero? 2418 !do io=1,Ep%nomega 2419 ! if (ABS(REAL(Ep%omega(io)))<0.00001) then 2420 ! do ig2=1,Ep%npwe 2421 ! do ig1=1,ig2-1 2422 ! chi0(ig2,ig1,io)=CONJG(chi0(ig1,ig2,io)) 2423 ! end do 2424 ! end do 2425 ! end if 2426 !end do 2427 2428 do iomega=1,MIN(Ep%nomega,NOMEGA_PRINTED) 2429 write(msg,'(1x,a,i4,a,2f9.4,a)')' chi0_intra(G,G'') at the ',iomega,' th omega',Ep%omega(iomega)*Ha_eV,' [eV]' 2430 call wrtout(std_out, msg) 2431 call print_arr(chi0(:,:,iomega),unit=std_out) 2432 end do 2433 2434 ! ===================== 2435 ! ==== Free memory ==== 2436 ! ===================== 2437 ABI_FREE(rhotwg) 2438 ABI_FREE(tabr_k) 2439 ABI_FREE(ur1) 2440 ABI_FREE(gw_gbound) 2441 ABI_FREE(ktabr) 2442 ABI_FREE(igffteps0) 2443 2444 ! deallocation for PAW. 2445 if (Psps%usepaw==1) then 2446 call pawcprj_free(Cprj1_bz) 2447 ABI_FREE(Cprj1_bz) 2448 call pawcprj_free(Cprj1_ibz) 2449 ABI_FREE(Cprj1_ibz) 2450 call pawpwij_free(Pwij) 2451 ABI_FREE(Pwij) 2452 end if 2453 2454 call Ltg_q%free() 2455 call Kmesh%free() 2456 2457 DBG_EXIT("COLL") 2458 2459 end subroutine chi0q0_intraband
ABINIT/m_chi0 [ Modules ]
NAME
m_chi0
FUNCTION
COPYRIGHT
Copyright (C) 1999-2024 ABINIT group (GMR, VO, LR, RWG, MG, RShaltaf) This file is distributed under the terms of the GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt .
SOURCE
16 #if defined HAVE_CONFIG_H 17 #include "config.h" 18 #endif 19 20 #include "abi_common.h" 21 22 module m_chi0 23 24 use defs_basis 25 use m_abicore 26 use m_xmpi 27 use m_errors 28 use m_hide_blas 29 use m_time 30 use m_wfd 31 use m_dtset 32 33 use defs_datatypes, only : pseudopotential_type, ebands_t 34 use m_fstrings, only : ftoa, sjoin, itoa 35 use m_gwdefs, only : GW_TOL_DOCC, GW_TOL_W0, czero_gw, em1params_t, g0g0w 36 use m_numeric_tools, only : imin_loc, print_arr 37 use m_geometry, only : normv, vdotw 38 use m_crystal, only : crystal_t 39 use m_fft_mesh, only : rotate_FFT_mesh, get_gftt 40 use m_occ, only : getnel 41 use m_ebands, only : pack_eneocc, unpack_eneocc, ebands_has_metal_scheme 42 use m_bz_mesh, only : kmesh_t, littlegroup_t 43 use m_gsphere, only : gsphere_t 44 use m_io_tools, only : flush_unit 45 use m_oscillators, only : rho_tw_g, calc_wfwfg 46 use m_vkbr, only : vkbr_t, vkbr_free, vkbr_init, nc_ihr_comm 47 use m_chi0tk, only : hilbert_transform, setup_spectral, assemblychi0_sym, assemblychi0sf, symmetrize_afm_chi0, & 48 approxdelta, completechi0_deltapart, accumulate_chi0sumrule, make_transitions, & 49 chi0_bbp_mask, accumulate_chi0_q0, accumulate_sfchi0_q0, hilbert_transform_headwings 50 use m_pawang, only : pawang_type 51 use m_pawrad, only : pawrad_type 52 use m_pawtab, only : pawtab_type 53 use m_paw_ij, only : paw_ij_type 54 use m_pawfgrtab, only : pawfgrtab_type 55 use m_pawcprj, only : pawcprj_type, pawcprj_alloc, pawcprj_free, pawcprj_copy 56 use m_pawpwij, only : pawpwff_t, pawpwij_t, pawpwij_init, pawpwij_free, paw_rho_tw_g, paw_cross_rho_tw_g 57 use m_paw_sym, only : paw_symcprj 58 use m_paw_pwaves_lmn, only : paw_pwaves_lmn_t 59 use m_paw_hr, only : pawhur_t, pawhur_free, pawhur_init, paw_ihr, paw_cross_ihr_comm 60 use m_read_plowannier, only : read_plowannier 61 use m_plowannier, only : plowannier_type 62 63 implicit none 64 65 private