TABLE OF CONTENTS
ABINIT/fermi_green [ Functions ]
NAME
fermi_green
FUNCTION
Compute Fermi level for DMFT or LDA.
COPYRIGHT
Copyright (C) 2006-2018 ABINIT group (BAmadon) This file is distributed under the terms of the GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt .
INPUTS
fermie: current value f_precision: precision of f required ITRLST: =1 if last iteration of DMFT opt_noninter : if one wants the LDA fermi level max_iter : max number of iterations.
OUTPUT
fermie: output value
CHILDREN
wrtout, spline_c,invfourier, nfourier
PARENTS
dmft_solve
CHILDREN
compute_green,integrate_green,newton,wrtout
SOURCE
36 #if defined HAVE_CONFIG_H 37 #include "config.h" 38 #endif 39 40 41 #include "abi_common.h" 42 43 subroutine fermi_green(cryst_struc,green,paw_dmft,pawang,& 44 & self) 45 46 use m_profiling_abi 47 48 use defs_basis 49 use defs_abitypes 50 use m_pawang, only : pawang_type 51 use m_crystal, only : crystal_t 52 use m_green, only : green_type,integrate_green,compute_green 53 use m_paw_dmft, only: paw_dmft_type 54 use m_self, only : self_type 55 56 !This section has been created automatically by the script Abilint (TD). 57 !Do not modify the following lines by hand. 58 #undef ABI_FUNC 59 #define ABI_FUNC 'fermi_green' 60 use interfaces_14_hidewrite 61 use interfaces_68_dmft, except_this_one => fermi_green 62 !End of the abilint section 63 64 implicit none 65 66 !Arguments ------------------------------------ 67 !scalars 68 type(crystal_t),intent(in) :: cryst_struc 69 type(green_type),intent(inout) :: green 70 type(paw_dmft_type), intent(inout) :: paw_dmft 71 !type(MPI_type), intent(in) :: mpi_enreg 72 type(pawang_type),intent(in) :: pawang 73 type(self_type), intent(inout) :: self 74 75 !Local variables------------------------------- 76 integer :: ierr_hh,opt_noninter,max_iter 77 real(dp):: x_precision,f_precision,fermi_old 78 ! real(dp) :: hx 79 character(len=500) :: message 80 !************************************************************************ 81 ! 82 write(message,'(a,8x,a)') ch10," == Compute Fermi level" 83 call wrtout(std_out,message,'COLL') 84 85 !============= 86 !headers 87 !============= 88 write(message,'(2a)') ch10, " |---Newton method to search Fermi level ------------|" 89 call wrtout(std_out,message,'COLL') 90 write(message,'(2a,f13.6)') ch10, " |--- Initial value for Fermi level",paw_dmft%fermie 91 call wrtout(std_out,message,'COLL') 92 93 !======================================== 94 !Define precision and nb of iterations 95 !========================================= 96 fermi_old=paw_dmft%fermie 97 ierr_hh=0 98 f_precision=paw_dmft%dmft_chpr 99 !f_precision=0.01 100 x_precision=tol5 101 !if(option==1) then 102 !f_precision=(erreursurlacharge)/100d0 103 !else 104 !f_precision=tol11 105 !endif 106 max_iter=1 ! for tests only 107 !write(6,*) "for tests max_iter=1" 108 max_iter=50 109 opt_noninter=4 110 111 !===================== 112 !Call newton method 113 !===================== 114 write(message,'(a,4x,a,e13.6)') ch10," Precision required :",f_precision 115 call wrtout(std_out,message,'COLL') 116 if (f_precision<10_dp) then 117 call newton(cryst_struc,green,paw_dmft,pawang,self,& 118 & paw_dmft%fermie,x_precision,max_iter,f_precision,ierr_hh,opt_noninter) 119 end if 120 121 !=========================== 122 !Deals with errors signals 123 !=========================== 124 if(ierr_hh==-314) then 125 write(message,'(a)') "Warning, check Fermi level" 126 call wrtout(std_out,message,'COLL') 127 ! call leave_new('COLL') 128 write(message,'(2a,f13.6)') ch10, " |--- Final value for Fermi level (check)",paw_dmft%fermie 129 call wrtout(std_out,message,'COLL') 130 else if (ierr_hh==-123) then 131 write(message,'(a,f13.6)') " Fermi level is put to",fermi_old 132 paw_dmft%fermie=fermi_old 133 call wrtout(std_out,message,'COLL') 134 135 ! ===================================== 136 ! If fermi level search was successful 137 ! ===================================== 138 else 139 write(message,'(a,4x,a,e13.6)') ch10," Precision achieved on Fermi Level :",x_precision 140 call wrtout(std_out,message,'COLL') 141 write(message,'(4x,a,e13.6)') " Precision achieved on number of electrons :",f_precision 142 call wrtout(std_out,message,'COLL') 143 write(message,'(2a,f13.6)') ch10, " |--- Final value for Fermi level",paw_dmft%fermie 144 call wrtout(std_out,message,'COLL') 145 end if 146 147 !======================================================== 148 !Check convergence of fermi level during DMFT iterations 149 !======================================================== 150 if(paw_dmft%idmftloop>=2) then 151 if(abs(paw_dmft%fermie-fermi_old).le.paw_dmft%dmft_fepr) then 152 ! write(message,'(a,8x,a,e9.2,a,8x,a,e12.5)') ch10,"|fermie(n)-fermie(n-1)|=<",paw_dmft%dmft_fepr,ch10,& 153 write(message,'(a,8x,a,e9.2,a,e9.2,a,8x,a,e12.5)') ch10,"|fermie(n)-fermie(n-1)|=",& 154 & abs(paw_dmft%fermie-fermi_old),"<",paw_dmft%dmft_fepr,ch10,& 155 & "=> DMFT Loop: Fermi level is converged to:",paw_dmft%fermie 156 call wrtout(std_out,message,'COLL') 157 green%ifermie_cv=1 158 else 159 write(message,'(a,8x,a,2f12.5)') ch10,"DMFT Loop: Fermi level is not converged:",& 160 & paw_dmft%fermie 161 call wrtout(std_out,message,'COLL') 162 green%ifermie_cv=0 163 end if 164 end if 165 write(message,'(2a)') ch10, " |---------------------------------------------------|" 166 call wrtout(std_out,message,'COLL') 167 ! 168 169 !========================================================== 170 !Recompute full green function including non diag elements 171 !========================================================== 172 call compute_green(cryst_struc,green,paw_dmft,pawang,0,self,opt_self=1,opt_nonxsum=1) 173 call integrate_green(cryst_struc,green,paw_dmft,pawang,prtopt=0,opt_ksloc=3) !,opt_nonxsum=1) 174 175 return 176 end subroutine fermi_green