The Elk Code
writespecies.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 writespecies(symb,name,zn,mass,rmin,rm,rmax,nrm,nst,n,l,k,occ,eval)
7 use modmain
8 use modmpi
9 implicit none
10 ! arguments
11 character(*), intent(in) :: symb,name
12 real(8), intent(in) :: zn,mass
13 real(8), intent(in) :: rmin,rm,rmax
14 integer, intent(in) :: nrm,nst
15 integer, intent(in) :: n(nst),l(nst),k(nst)
16 real(8), intent(in) :: occ(nst)
17 real(8), intent(in) :: eval(nst)
18 ! local variables
19 integer lmax,nlo
20 integer ist,jst,i
21 logical core(maxstsp),lorb(maxstsp)
22 ! default APW band energy
23 real(8), parameter :: e0=0.15d0
24 ! find which states belong to core
25 core(1:nst)=(eval(1:nst) < ecvcut)
26 ! check that the state for same n and l but different k is also core
27 do ist=1,nst
28  if (core(ist)) then
29  do jst=1,nst
30  if ((n(ist) == n(jst)).and.(l(ist) == l(jst))) core(jst)=.true.
31  end do
32  end if
33 end do
34 lmax=1
35 do ist=1,nst
36  if (.not.core(ist)) lmax=max(lmax,l(ist))
37 end do
38 ! determine the local orbitals
39 nlo=lmax+1
40 lorb(:)=.false.
41 do ist=1,nst
42  if (.not.core(ist)) then
43  if ((l(ist) == 0).or.(l(ist) < k(ist))) then
44  if ((eval(ist) < esccut).or.(l(ist) >= 2)) then
45  lorb(ist)=.true.
46  nlo=nlo+1
47  end if
48  end if
49  end if
50 end do
51 if (mp_mpi) then
52  open(55,file=trim(symb)//'.in',form='FORMATTED')
53  write(55,'(" ''",A,"''",T45,": spsymb")') trim(symb)
54  write(55,'(" ''",A,"''",T45,": spname")') trim(name)
55  write(55,'(G14.6,T45,": spzn")') zn
56  write(55,'(G18.10,T45,": spmass")') mass
57  write(55,'(G14.6,2F10.4,I6,T45,": rminsp, rmt, rmaxsp, nrmt")') rmin,rm, &
58  rmax,nrm
59  write(55,'(I4,T45,": nstsp")') nst
60  write(55,'(3I4,F10.5," ",L1,T45,": nsp, lsp, ksp, occsp, spcore")') n(1), &
61  l(1),k(1),occ(1),core(1)
62  do ist=2,nst
63  write(55,'(3I4,F10.5," ",L1)') n(ist),l(ist),k(ist),occ(ist),core(ist)
64  end do
65  write(55,'(I4,T45,": apword")') 1
66  write(55,'(F10.4,I4," ",L1,T45,": apwe0, apwdm, apwve")') e0,0,.false.
67  write(55,'(I4,T45,": nlx")') 0
68  write(55,'(I4,T45,": nlorb")') nlo
69  do i=0,lmax
70  write(55,'(2I4,T45,": lorbl, lorbord")') i,2
71  write(55,'(F10.4,I4," ",L1,T45,": lorbe0, lorbdm, lorbve")') e0,0,.false.
72  write(55,'(F10.4,I4," ",L1)') e0,1,.false.
73  end do
74  do ist=1,nst
75  if (lorb(ist)) then
76  write(55,'(2I4,T45,": lorbl, lorbord")') l(ist),2
77  write(55,'(F10.4,I4," ",L1,T45,": lorbe0, lorbdm, lorbve")') e0,0,.false.
78  write(55,'(F10.4,I4," ",L1)') eval(ist),0,.true.
79  end if
80  end do
81  close(55)
82 end if
83 ! synchronise MPI processes
84 call mpi_barrier(mpicom,ierror)
85 end subroutine
86 
real(8) ecvcut
Definition: modmain.f90:117
subroutine writespecies(symb, name, zn, mass, rmin, rm, rmax, nrm, nst, n, l, k, occ, eval)
Definition: writespecies.f90:7
logical mp_mpi
Definition: modmpi.f90:17
real(8) esccut
Definition: modmain.f90:119
Definition: modmpi.f90:6
integer mpicom
Definition: modmpi.f90:11
integer ierror
Definition: modmpi.f90:19