File: linengy.f90

package info (click to toggle)
elkcode 5.4.24-2
  • links: PTS, VCS
  • area: main
  • in suites: buster, sid
  • size: 12,840 kB
  • sloc: f90: 48,415; fortran: 22,457; perl: 965; makefile: 384; sh: 369; python: 105; ansic: 67
file content (133 lines) | stat: -rw-r--r-- 3,531 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133

! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
! This file is distributed under the terms of the GNU General Public License.
! See the file COPYING for license details.

!BOP
! !ROUTINE: linengy
! !INTERFACE:
subroutine linengy
! !USES:
use modmain
use modmpi
! !DESCRIPTION:
!   Calculates the new linearisation energies for both the APW and local-orbital
!   radial functions. See the routine {\tt findband}.
!
! !REVISION HISTORY:
!   Created May 2003 (JKD)
!EOP
!BOC
implicit none
! local variables
logical fnd
integer is,ia,ja,ias,jas
integer nr,nri,iro,ir,l,i
integer ilo,io,jo,nnf
! automatic arrays
logical done(natmmax)
real(8) vr(nrmtmax)
nnf=0
! begin loops over atoms and species
do is=1,nspecies
  nr=nrmt(is)
  nri=nrmti(is)
  iro=nri+1
  done(:)=.false.
  do ia=1,natoms(is)
    if (done(ia)) cycle
    ias=idxas(ia,is)
    i=1
    do ir=1,nri
      vr(ir)=vsmt(i,ias)*y00
      i=i+lmmaxi
    end do
    do ir=iro,nr
      vr(ir)=vsmt(i,ias)*y00
      i=i+lmmaxo
    end do
!-----------------------!
!     APW functions     !
!-----------------------!
    do l=0,lmaxapw
      do io=1,apword(l,is)
        if (apwve(io,l,is)) then
! check if previous radial functions have same default energies
          do jo=1,io-1
            if (apwve(jo,l,is)) then
              if (abs(apwe0(io,l,is)-apwe0(jo,l,is)).lt.1.d-4) then
                apwe(io,l,ias)=apwe(jo,l,ias)
                goto 10
              end if
            end if
          end do
! find the band energy starting from default
          apwe(io,l,ias)=apwe0(io,l,is)
          call findband(solsc,l,nr,rsp(:,is),vr,epsband,demaxbnd, &
           apwe(io,l,ias),fnd)
          if (.not.fnd) nnf=nnf+1
        else
! set linearisation energy automatically
          if (autolinengy) apwe(io,l,ias)=efermi+dlefe
        end if
10 continue
      end do
    end do
!---------------------------------!
!     local-orbital functions     !
!---------------------------------!
    do ilo=1,nlorb(is)
      do io=1,lorbord(ilo,is)
        if (lorbve(io,ilo,is)) then
! check if previous radial functions have same default energies
          do jo=1,io-1
            if (lorbve(jo,ilo,is)) then
              if (abs(lorbe0(io,ilo,is)-lorbe0(jo,ilo,is)).lt.1.d-4) then
                lorbe(io,ilo,ias)=lorbe(jo,ilo,ias)
                goto 20
              end if
            end if
          end do
          l=lorbl(ilo,is)
! find the band energy starting from default
          lorbe(io,ilo,ias)=lorbe0(io,ilo,is)
          call findband(solsc,l,nr,rsp(:,is),vr,epsband,demaxbnd, &
           lorbe(io,ilo,ias),fnd)
          if (.not.fnd) nnf=nnf+1
        else
! set linearisation energy automatically
          if (autolinengy) lorbe(io,ilo,ias)=efermi+dlefe
        end if
20 continue
      end do
    end do
    done(ia)=.true.
! copy to equivalent atoms
    do ja=1,natoms(is)
      if ((.not.done(ja)).and.(eqatoms(ia,ja,is))) then
        jas=idxas(ja,is)
        do l=0,lmaxapw
          do io=1,apword(l,is)
            apwe(io,l,jas)=apwe(io,l,ias)
          end do
        end do
        do ilo=1,nlorb(is)
          do io=1,lorbord(ilo,is)
            lorbe(io,ilo,jas)=lorbe(io,ilo,ias)
          end do
        end do
        done(ja)=.true.
      end if
    end do
! end loops over atoms and species
  end do
end do
if (mp_mpi.and.(nnf.gt.0)) then
  write(*,*)
  write(*,'("Warning(linengy): could not find ",I3," linearisation energies &
   &in s.c. loop ",I5)') nnf,iscl
end if
return
end subroutine
!EOC