The Elk Code
 
Loading...
Searching...
No Matches
genfdufr.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: genfdufr
8! !INTERFACE:
9subroutine genfdufr(idu)
10! !USES:
11use modmain
12use moddftu
13! !INPUT/OUTPUT PARAMETERS:
14! idu : DFT+U entry (in,integer)
15! !DESCRIPTION:
16! Generates the radial functions used to calculate the Slater integrals
17! through a Yukawa potential.
18!
19! !REVISION HISTORY:
20! Created April 2008 from genapwfr (Francesco Cricchio)
21!EOP
22!BOC
23implicit none
24integer, intent(in) :: idu
25! local variables
26integer is,ia,ias
27integer nr,nri,nn,l
28real(8) t1
29! automatic arrays
30real(8) vr(nrmtmax),p0(nrmtmax),p1(nrmtmax),q0(nrmtmax),q1(nrmtmax)
31is=isldu(1,idu)
32l=isldu(2,idu)
33nr=nrmt(is)
34nri=nrmti(is)
35do ia=1,natoms(is)
36 ias=idxas(ia,is)
37 call rfmtlm(1,nr,nri,vsmt(:,ias),vr)
38 vr(1:nr)=vr(1:nr)*y00
39! integrate the radial Schrodinger equation
40 call rschrodint(solsc,l,efdu(l,ias),nr,rlmt(:,1,is),vr,nn,p0,p1,q0,q1)
41! divide the radial wavefunction by r
42 p0(1:nr)=p0(1:nr)*rlmt(1:nr,-1,is)
43! normalise radial functions
44 t1=sum(wr2mt(1:nr,is)*p0(1:nr)**2)
45 if (t1 < 1.d-20) then
46 write(*,*)
47 write(*,'("Error(genfdufr): degenerate radial functions")')
48 write(*,'(" for species ",I4)') is
49 write(*,'(" atom ",I4)') ia
50 write(*,'(" and angular momentum ",I4)') l
51 write(*,*)
52 stop
53 end if
54 t1=1.d0/sqrt(abs(t1))
55 p0(1:nr)=t1*p0(1:nr)
56! store in global array
57 fdufr(1:nr,l,ias)=p0(1:nr)
58end do
59end subroutine
60!EOC
61
subroutine genfdufr(idu)
Definition genfdufr.f90:10
integer, dimension(2, maxdftu) isldu
Definition moddftu.f90:40
real(8), dimension(:,:,:), allocatable fdufr
Definition moddftu.f90:58
real(8), dimension(:,:), allocatable efdu
Definition moddftu.f90:56
integer, dimension(maxspecies) nrmti
Definition modmain.f90:211
real(8), parameter y00
Definition modmain.f90:1233
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
integer, dimension(maxspecies) natoms
Definition modmain.f90:36
integer, dimension(maxatoms, maxspecies) idxas
Definition modmain.f90:42
real(8), dimension(:,:), pointer, contiguous vsmt
Definition modmain.f90:649
real(8) solsc
Definition modmain.f90:1252
real(8), dimension(:,:), allocatable wr2mt
Definition modmain.f90:183
real(8), dimension(:,:,:), allocatable rlmt
Definition modmain.f90:179
pure subroutine rfmtlm(lm, nr, nri, rfmt, fr)
Definition rfmtlm.f90:7
pure subroutine rschrodint(sol, l, e, nr, r, vr, nn, p0, p1, q0, q1)