The Elk Code
genfspecies.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2014 J. K. Dewhurst, S. Sharma and E. K. U. Gross.
3 ! This file is distributed under the terms of the GNU General Public License.
4 ! See the file COPYING for license details.
5 
6 subroutine genfspecies(zn,symb)
7 use modmain
8 use modmpi
9 implicit none
10 ! arguments
11 real(8), intent(in) :: zn
12 character(*), intent(in) :: symb
13 ! local variables
14 integer, parameter :: nit=4
15 integer nst,ist,jst
16 integer nmax,in,il,ik
17 integer nrm,nr,ir,it
18 integer n(maxstsp),l(maxstsp),k(maxstsp)
19 integer idx(maxstsp),iv(maxstsp)
20 real(8) rm,rmin,rmax
21 real(8) mass,t1,t2,t3
22 real(8) occ(maxstsp),eval(maxstsp),rv(maxstsp)
23 character(64) name
24 ! allocatable arrays
25 real(8), allocatable :: r(:),rho(:),vr(:),rwf(:,:,:)
26 ! external functions
27 real(8), external :: massnucl
28 name='Fractional species'
29 ! set up the initial occupation numbers
30 occ(:)=0.d0
31 t1=abs(zn)
32 nmax=1
33 ist=0
34 do in=1,maxstsp
35  do il=0,in-1
36  do ik=max(il,1),il+1
37  t2=dble(2*ik)
38  t2=min(t2,t1)
39  ist=ist+1
40  n(ist)=in
41  l(ist)=il
42  k(ist)=ik
43  occ(ist)=t2
44  if (t2 > epsocc) nmax=in
45  t1=t1-t2
46  if (ist == maxstsp) then
47  if (t1 > epsocc) then
48  write(*,*)
49  write(*,'("Error(genfspecies): too many states for fractional &
50  &species ",A)') trim(symb)
51  write(*,*)
52  stop
53  else
54  goto 10
55  end if
56  end if
57  end do
58  end do
59 end do
60 10 continue
61 ! minimum radius
62 rmin=2.d-6/sqrt(abs(zn))
63 ! initial maximum radius
64 rmax=100.d0
65 ! initial muffin-tin radius
66 rm=2.d0
67 ! number of points to muffin-tin radius
68 nrm=100*(nmax+1)
69 ! iterate the solution but not to self-consistency
70 do it=1,nit
71 ! number of points to effective infinity
72  t1=log(rm/rmin)
73  t2=log(rmax/rmin)
74  t3=dble(nrm)*t2/t1
75  nr=int(t3)
76  allocate(r(nr),rho(nr),vr(nr),rwf(nr,2,maxstsp))
77 ! generate logarithmic radial mesh
78  t2=t1/dble(nrm-1)
79  do ir=1,nr
80  r(ir)=rmin*exp(dble(ir-1)*t2)
81  end do
82 ! solve the Kohn-Sham-Dirac equation for the atom
83  call atom(sol,.true.,zn,maxstsp,n,l,k,occ,3,0,nr,r,eval,rho,vr,rwf)
84 ! check for spurious eigenvalues
85  do ist=2,maxstsp
86  if (eval(ist) < eval(1)) eval(ist)=1.d6
87  end do
88 ! recompute the effective infinity
89  do ir=nr,1,-1
90  if (rho(ir) > 1.d-20) then
91  rmax=1.75d0*r(ir)
92  exit
93  end if
94  end do
95 ! estimate the muffin-tin radius
96  do ir=nr,1,-1
97  if (rho(ir) > 2.d-2) then
98  rm=r(ir)
99  exit
100  end if
101  end do
102  if (rm < 1.d0) rm=1.d0
103  if (rm > 3.2d0) rm=3.2d0
104 ! sort the eigenvalues
105  call sortidx(maxstsp,eval,idx)
106 ! recompute the occupation numbers
107  occ(:)=0.d0
108  t1=abs(zn)
109  do ist=1,maxstsp
110  jst=idx(ist)
111  ik=k(jst)
112  t2=dble(2*ik)
113  t2=min(t2,t1)
114  occ(jst)=t2
115  t1=t1-t2
116  end do
117  deallocate(r,rho,vr,rwf)
118 end do
119 ! rearrange the arrays
120 iv(:)=n(:)
121 n(:)=iv(idx(:))
122 iv(:)=l(:)
123 l(:)=iv(idx(:))
124 iv(:)=k(:)
125 k(:)=iv(idx(:))
126 rv(:)=occ(:)
127 occ(:)=rv(idx(:))
128 rv(:)=eval(:)
129 eval(:)=rv(idx(:))
130 ! find the number of occupied states
131 nst=0
132 do ist=1,maxstsp
133  if (occ(ist) < epsocc) then
134  nst=ist
135  exit
136  end if
137 end do
138 ! estimate the nuclear mass
139 mass=massnucl(zn)
140 ! convert from 'atomic mass units' to atomic units
141 mass=mass*amu
142 ! write the species file
143 call writespecies(symb,name,zn,mass,rmin,rm,rmax,nrm,nst,n,l,k,occ,eval)
144 if (mp_mpi) then
145  write(*,'("Info(genfspecies): wrote fractional species file ",A,".in")') &
146  trim(symb)
147 end if
148 end subroutine
149 
subroutine genfspecies(zn, symb)
Definition: genfspecies.f90:7
subroutine writespecies(symb, name, zn, mass, rmin, rm, rmax, nrm, nst, n, l, k, occ, eval)
Definition: writespecies.f90:7
subroutine atom(sol, ptnucl, zn, nst, n, l, k, occ, xctype, xcgrad, nr, r, eval, rho, vr, rwf)
Definition: atom.f90:10
logical mp_mpi
Definition: modmpi.f90:17
real(8) epsocc
Definition: modmain.f90:903
real(8), parameter amu
Definition: modmain.f90:1283
real(8), parameter sol
Definition: modmain.f90:1251
pure subroutine sortidx(n, x, idx)
Definition: sortidx.f90:10
Definition: modmpi.f90:6