The Elk Code
 
Loading...
Searching...
No Matches
genfdu.f90
Go to the documentation of this file.
1
2! Copyright (C) 2008 F. Bultmark, F. Cricchio, L. Nordstrom and J. K. Dewhurst.
3! This file is distributed under the terms of the GNU General Public License.
4! See the file COPYING for license details.
5
6!BOP
7! !ROUTINE: genfdu
8! !INTERFACE:
9subroutine genfdu(idu,u,j,f)
10! !USES:
11use moddftu
12use modmpi
13! !INPUT/OUTPUT PARAMETERS:
14! idu : DFT+U entry (in,integer)
15! u : parameter U (inout,real)
16! j : parameter J (inout,real)
17! f : Slater parameters (inout,real)
18! !DESCRIPTION:
19! Calculate the Slater parameters for DFT+$U$ calculation with different
20! approaches, see {\it Phys. Rev. B} {\bf 80}, 035121 (2009). The relations
21! among Slater and Racah parameters are from E. U. Condon and G. H. Shortley,
22! {\it The Theory of Atomic Spectra}, The University Press, Cambridge (1935).
23!
24! !REVISION HISTORY:
25! Created July 2008 (Francesco Cricchio)
26!EOP
27!BOC
28implicit none
29! arguments
30integer, intent(in) :: idu
31real(8), intent(inout) :: u,j
32real(8), intent(inout) :: f(0:2*lmaxdm)
33! local variables
34integer is,l,k,q
35real(8) ufix,r1,r2
36real(8) a(3,3),v1(3),v2(3)
37! automatic arrays
38real(8) e(0:lmaxdm)
39! external functions
40real(8), external :: fyukawa,fyukawa0
41is=isldu(1,idu)
42l=isldu(2,idu)
43if (l > 3) return
44! load input parameters to calculate Slater integrals
45u=ujdu(1,idu)
46j=ujdu(2,idu)
47f(:)=fdu(:,idu)
48e(:)=edu(:,idu)
49if (inpdftu < 4) then
50! F(0) = U for any l-shell
51 if (inpdftu == 1) f(0)=u
52 select case(l)
53 case(0)
54! s electrons only f(0)=u
55 if (inpdftu == 3) then
56 f(0)=e(0)
57 u=f(0)
58 end if
59 case(1)
60! p electrons
61 if (inpdftu == 1) then
62! F(2) = 5.0 * J
63 f(2)=5.d0*j
64 else if (inpdftu == 3) then
65! F(0) = E(0) + J= E(0) + 5/3 * E(1)
66 f(0)=e(0)+(5.d0/3.d0)*e(1)
67! F(2) = 5 * J = 25/3 * E1, Eq. 101
68 f(2)=(25.d0/3.d0)*e(1)
69 end if
70 case(2)
71! d electrons
72 if (inpdftu == 1) then
73! r1 = F(4)/F(2), see PRB 52, R5467 (1995)
74 r1=0.625d0
75 f(2)=(14.d0*j)/(1.d0+r1)
76 f(4)=f(2)*r1
77 else if (inpdftu == 3) then
78! copy Racah parameters
79 v1(1:3)=e(0:2)
80! transformation matrix from Racah to Slater parameters
81! obtained from inversion of Eq. 110-112, LN Notes 29-12-08
82 a(1,1)=1.d0
83 a(1,2)=1.4d0
84 a(1,3)=0.d0
85 a(2,1)=0.d0
86 a(2,2)=0.1428571428571428d0
87 a(2,3)=1.285714285714286d0
88 a(3,1)=0.d0
89 a(3,2)=2.8571428571428571d-2
90 a(3,3)=-0.1428571428571428d0
91! multiply transformation matrix by Racah parameters
92 call r3mv(a,v1,v2)
93! Slater parameters, Eq. 104-105, LN Notes 29-12-08
94 f(0)=v2(1)
95 f(2)=49.d0*v2(2)
96 f(4)=441.d0*v2(3)
97 end if
98 case(3)
99! f electrons
100 if (inpdftu == 1) then
101! r2 = F(6)/F(2), r1 = F(4)/F(2), see PRB 50, 16861 (1994)
102 r1=451.d0/675.d0
103 r2=1001.d0/2025.d0
104 f(2)=6435.d0*j/(286.d0+195.d0*r1+250.d0*r2)
105 f(4)=f(2)*r1
106 f(6)=f(2)*r2
107 else if (inpdftu == 3) then
108! F(0) = E(0) + 9/7 * E(1) , Eq. 119, LN Notes 29-12-08
109 f(0)=e(0)+(9.d0/7.d0)*e(1)
110! copy Racah parameters
111 v1(1:3)=e(1:3)
112! transformation matrix from Racah to Slater parameters
113! obtained from inversion of Eq. 120-122, LN Notes 29-12-08
114 a(1,1)=2.3809523809523808d-2
115 a(1,2)=3.404761904761904d0
116 a(1,3)=0.2619047619047619d0
117 a(2,1)=1.2987012987012984d-2
118 a(2,2)=-1.688311688311688d0
119 a(2,3)=5.1948051948051951d-2
120 a(3,1)=2.1645021645021645d-3
121 a(3,2)=7.5757575757575760d-2
122 a(3,3)=-1.5151515151515152d-2
123! multiply transformation matrix by Racah parameters
124 call r3mv(a,v1,v2)
125! Slater parameters, Eq. 115-117, LN Notes 29-12-08
126 f(2)=225.d0*v2(1)
127 f(4)=1089.d0*v2(2)
128 f(6)=(184041.d0/25.d0)*v2(3)
129 end if
130 end select
131else if (inpdftu >= 4) then
132! define energies for Slater parameters
133 call engyfdu(idu)
134! calculate radial functions for Slater parameters
135 call genfdufr(idu)
136 if (inpdftu == 5) then
137 ufix=udufix(idu)
138! calculate the screening length λ corresponding to udufix
139! lamdu0 is in/out and is initialized to 0 in init0
140 call findlambda(is,l,ufix,lamdu0(idu),lamdu(idu))
141 end if
142 do q=0,l
143 k=2*q
144 if (lamdu(idu) < 1.d-2) then
145! unscreened Slater parameters
146 f(k)=fyukawa0(is,l,k)
147 else
148! screened Slater parameters
149 f(k)=fyukawa(is,l,k,lamdu(idu))
150 end if
151 end do
152end if
153! calculate U and J from Slater integrals
154if (inpdftu /= 1) then
155 u=f(0)
156 select case(l)
157 case(0)
158 j=0.d0
159 case(1)
160! J = 1/5 * F(2)
161 j=(1.d0/5.d0)*f(2)
162 case(2)
163! J = 1/14 * ( F(2) + F(4) ), Eq. 106, LN Notes 29-12-08
164 j=(1.d0/14.d0)*(f(2)+f(4))
165 case(3)
166! J= 2/45 * F(2) + 1/33 * F(4) + 50/1287 * F(6), Eq. 118, LN Notes 29-12-08
167 j=(2.d0/45.d0)*f(2)+(1.d0/33.d0)*f(4)+(50.d0/1287.d0)*f(6)
168 end select
169end if
170! save calculated parameters
171! (except Racah parameters that are provided only as input)
172ujdu(1,idu)=u
173ujdu(2,idu)=j
174fdu(:,idu)=f(:)
175end subroutine
176!EOC
177
subroutine engyfdu(idu)
Definition engyfdu.f90:10
subroutine findlambda(is, l, ufix, lambda0, lambda)
real(8) function fyukawa0(is, l, k)
Definition fyukawa0.f90:10
real(8) function fyukawa(is, l, k, lambda)
Definition fyukawa.f90:10
subroutine genfdu(idu, u, j, f)
Definition genfdu.f90:10
subroutine genfdufr(idu)
Definition genfdufr.f90:10
real(8), dimension(0:lmaxdm, maxdftu) edu
Definition moddftu.f90:50
real(8), dimension(0:2 *lmaxdm, maxdftu) fdu
Definition moddftu.f90:48
real(8), dimension(maxdftu) lamdu0
Definition moddftu.f90:54
real(8), dimension(maxdftu) udufix
Definition moddftu.f90:60
real(8), dimension(2, maxdftu) ujdu
Definition moddftu.f90:42
integer, dimension(2, maxdftu) isldu
Definition moddftu.f90:40
integer inpdftu
Definition moddftu.f90:34
real(8), dimension(maxdftu) lamdu
Definition moddftu.f90:52
pure subroutine r3mv(a, x, y)
Definition r3mv.f90:10