The Elk Code
engyfdu.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: engyfdu
8 ! !INTERFACE:
9 subroutine engyfdu(idu)
10 ! !USES:
11 use modmain
12 use moddftu
13 use modmpi
14 ! !INPUT/OUTPUT PARAMETERS:
15 ! idu : DFT+U entry (in,integer)
16 ! !DESCRIPTION:
17 ! Calculates the energies of radial functions to be used to calculate the
18 ! Slater integrals. By convention those energies are chosen to be the ones at
19 ! the center of the band.
20 !
21 ! !REVISION HISTORY:
22 ! Created April 2008 (F. Cricchio)
23 !EOP
24 !BOC
25 implicit none
26 integer, intent(in) :: idu
27 ! local variables
28 integer is,ia,ja,ias,jas
29 integer nr,nri,nnf,l
30 logical fnd
31 ! automatic arrays
32 real(8) vr(nrmtmax)
33 nnf=0
34 is=isldu(1,idu)
35 l=isldu(2,idu)
36 nr=nrmt(is)
37 nri=nrmti(is)
38 do ia=1,natoms(is)
39 ! perform calculation for only the first equivalent atom
40  if (.not.tfeqat(ia,is)) cycle
41  ias=idxas(ia,is)
42  call rfmtlm(1,nr,nri,vsmt(:,ias),vr)
43  vr(1:nr)=vr(1:nr)*y00
44 ! find the center of the band starting from -0.5 Ha
45  efdu(l,ias)=-0.5d0
46  call findband(solsc,l,nrmt(is),rsp(1,is),vr,epsband,demaxbnd,efdu(l,ias),fnd)
47  if (.not.fnd) nnf=nnf+1
48 ! copy to equivalent atoms
49  do ja=1,natoms(is)
50  if (eqatoms(ia,ja,is).and.(ia /= ja)) then
51  jas=idxas(ja,is)
52  efdu(l,jas)=efdu(l,ias)
53  end if
54  end do
55 ! end loops over atoms and species
56 end do
57 if (mp_mpi.and.(nnf > 0)) then
58  write(*,*)
59  write(*,'("Warning(engyfdu): could not find ",I3," energies")') nnf
60 end if
61 end subroutine
62 !EOC
63 
real(8) epsband
Definition: modmain.f90:820
real(8), dimension(:,:), allocatable efdu
Definition: moddftu.f90:56
logical mp_mpi
Definition: modmpi.f90:17
logical, dimension(:,:), allocatable tfeqat
Definition: modmain.f90:372
integer, dimension(maxatoms, maxspecies) idxas
Definition: modmain.f90:42
logical, dimension(:,:,:), allocatable eqatoms
Definition: modmain.f90:370
integer, dimension(2, maxdftu) isldu
Definition: moddftu.f90:40
pure subroutine rfmtlm(lm, nr, nri, rfmt, fr)
Definition: rfmtlm.f90:10
real(8) solsc
Definition: modmain.f90:1253
integer, dimension(maxspecies) natoms
Definition: modmain.f90:36
Definition: modmpi.f90:6
subroutine engyfdu(idu)
Definition: engyfdu.f90:10
real(8), dimension(:,:), allocatable rsp
Definition: modmain.f90:135
real(8) demaxbnd
Definition: modmain.f90:823
real(8), dimension(:,:), pointer, contiguous vsmt
Definition: modmain.f90:649
integer, dimension(maxspecies) nrmti
Definition: modmain.f90:211
real(8), parameter y00
Definition: modmain.f90:1236
subroutine findband(sol, l, nr, r, vr, eps, demax, e, fnd)
Definition: findband.f90:10
integer, dimension(maxspecies) nrmt
Definition: modmain.f90:150