The Elk Code
linengy.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2002-2005 J. K. Dewhurst, S. Sharma and C. Ambrosch-Draxl.
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: linengy
8 ! !INTERFACE:
9 subroutine linengy
10 ! !USES:
11 use modmain
12 use modmpi
13 use modomp
14 ! !DESCRIPTION:
15 ! Calculates the new linearisation energies for both the APW and local-orbital
16 ! radial functions. See the routine {\tt findband}.
17 !
18 ! !REVISION HISTORY:
19 ! Created May 2003 (JKD)
20 !EOP
21 !BOC
22 implicit none
23 ! local variables
24 logical fnd
25 integer is,ia,ja,ias,jas
26 integer nr,nri,iro,i0,i1
27 integer l,io,jo,ilo,nnf,nthd
28 ! automatic arrays
29 real(8) vr(nrmtmax)
30 nnf=0
31 ! begin loops over atoms and species
32 call holdthd(natmtot,nthd)
33 !$OMP PARALLEL DO DEFAULT(SHARED) &
34 !$OMP PRIVATE(vr,is,ia,nr,nri,iro,i0,i1) &
35 !$OMP PRIVATE(l,io,jo,fnd,ilo,ja,jas) &
36 !$OMP REDUCTION(+:nnf) &
37 !$OMP SCHEDULE(DYNAMIC) &
38 !$OMP NUM_THREADS(nthd)
39 do ias=1,natmtot
40  is=idxis(ias)
41  ia=idxia(ias)
42 ! perform calculation for only the first equivalent atom
43  if (.not.tfeqat(ia,is)) cycle
44  nr=nrmt(is)
45  nri=nrmti(is)
46  iro=nri+1
47  i1=lmmaxi*(nri-1)+1
48  vr(1:nri)=vsmt(1:i1:lmmaxi,ias)*y00
49  i0=i1+lmmaxi
50  i1=lmmaxo*(nr-iro)+i0
51  vr(iro:nr)=vsmt(i0:i1:lmmaxo,ias)*y00
52 !-----------------------!
53 ! APW functions !
54 !-----------------------!
55  do l=0,lmaxapw
56  do io=1,apword(l,is)
57  if (apwve(io,l,is)) then
58 ! check if previous radial functions have same default energies
59  do jo=1,io-1
60  if (apwve(jo,l,is)) then
61  if (abs(apwe0(io,l,is)-apwe0(jo,l,is)) < 1.d-4) then
62  apwe(io,l,ias)=apwe(jo,l,ias)
63  goto 10
64  end if
65  end if
66  end do
67 ! find the band energy starting from default
68  apwe(io,l,ias)=apwe0(io,l,is)
69  call findband(solsc,l,nr,rlmt(:,1,is),vr,epsband,demaxbnd, &
70  apwe(io,l,ias),fnd)
71  if (.not.fnd) nnf=nnf+1
72  else
73 ! set linearisation energy automatically
74  if (autolinengy) apwe(io,l,ias)=efermi+dlefe
75  end if
76 10 continue
77  end do
78  end do
79 !---------------------------------!
80 ! local-orbital functions !
81 !---------------------------------!
82  do ilo=1,nlorb(is)
83  do io=1,lorbord(ilo,is)
84  if (lorbve(io,ilo,is)) then
85 ! check if previous radial functions have same default energies
86  do jo=1,io-1
87  if (lorbve(jo,ilo,is)) then
88  if (abs(lorbe0(io,ilo,is)-lorbe0(jo,ilo,is)) < 1.d-4) then
89  lorbe(io,ilo,ias)=lorbe(jo,ilo,ias)
90  goto 20
91  end if
92  end if
93  end do
94  l=lorbl(ilo,is)
95 ! find the band energy starting from default
96  lorbe(io,ilo,ias)=lorbe0(io,ilo,is)
97  call findband(solsc,l,nr,rlmt(:,1,is),vr,epsband,demaxbnd, &
98  lorbe(io,ilo,ias),fnd)
99  if (.not.fnd) nnf=nnf+1
100  else
101 ! set linearisation energy automatically
102  if (autolinengy) lorbe(io,ilo,ias)=efermi+dlefe
103  end if
104 20 continue
105  end do
106  end do
107 ! copy to equivalent atoms
108  do ja=1,natoms(is)
109  if (eqatoms(ia,ja,is).and.(ia /= ja)) then
110  jas=idxas(ja,is)
111  do l=0,lmaxapw
112  do io=1,apword(l,is)
113  apwe(io,l,jas)=apwe(io,l,ias)
114  end do
115  end do
116  do ilo=1,nlorb(is)
117  do io=1,lorbord(ilo,is)
118  lorbe(io,ilo,jas)=lorbe(io,ilo,ias)
119  end do
120  end do
121  end if
122  end do
123 ! end loop over species and atoms
124 end do
125 !$OMP END PARALLEL DO
126 call freethd(nthd)
127 if (mp_mpi.and.(nnf > 0)) then
128  write(*,*)
129  write(*,'("Warning(linengy): could not find ",I3," linearisation energies &
130  &in s.c. loop ",I5)') nnf,iscl
131 end if
132 end subroutine
133 !EOC
134 
real(8) epsband
Definition: modmain.f90:820
real(8) efermi
Definition: modmain.f90:907
integer, dimension(maxspecies) nlorb
Definition: modmain.f90:786
logical mp_mpi
Definition: modmpi.f90:17
integer lmmaxo
Definition: modmain.f90:203
logical, dimension(:,:), allocatable tfeqat
Definition: modmain.f90:372
integer, dimension(maxatoms, maxspecies) idxas
Definition: modmain.f90:42
real(8), dimension(:,:,:), allocatable rlmt
Definition: modmain.f90:179
real(8), dimension(maxlorbord, maxlorb, maxspecies) lorbe0
Definition: modmain.f90:804
integer iscl
Definition: modmain.f90:1051
Definition: modomp.f90:6
logical, dimension(:,:,:), allocatable eqatoms
Definition: modmain.f90:370
integer lmaxapw
Definition: modmain.f90:197
subroutine linengy
Definition: linengy.f90:10
logical, dimension(maxlorbord, maxlorb, maxspecies) lorbve
Definition: modmain.f90:812
real(8), dimension(maxapword, 0:maxlapw, maxspecies) apwe0
Definition: modmain.f90:766
real(8), dimension(:,:,:), allocatable lorbe
Definition: modmain.f90:808
real(8), dimension(:,:,:), allocatable apwe
Definition: modmain.f90:768
integer, dimension(0:maxlapw, maxspecies) apword
Definition: modmain.f90:758
real(8) solsc
Definition: modmain.f90:1253
integer, dimension(maxspecies) natoms
Definition: modmain.f90:36
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
integer lmmaxi
Definition: modmain.f90:207
logical, dimension(maxapword, 0:maxlapw, maxspecies) apwve
Definition: modmain.f90:772
Definition: modmpi.f90:6
integer, dimension(maxlorb, maxspecies) lorbord
Definition: modmain.f90:792
real(8) dlefe
Definition: modmain.f90:830
logical autolinengy
Definition: modmain.f90:828
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
integer, dimension(maxatoms *maxspecies) idxia
Definition: modmain.f90:45
integer natmtot
Definition: modmain.f90:40
integer, dimension(maxlorb, maxspecies) lorbl
Definition: modmain.f90:796
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