The Elk Code
 
Loading...
Searching...
No Matches
readspecies.f90
Go to the documentation of this file.
1
2! Copyright (C) 2002-2008 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
6subroutine readspecies
7use modmain
8implicit none
9! local variables
10integer is,ist,ios
11integer nlx,ilx,lx,ilo
12integer io,jo,ko,l,i,j
13e0min=0.d0
14do is=1,nspecies
15 open(50,file=trim(sppath)//trim(spfname(is)),status='OLD',form='FORMATTED', &
16 action='READ',iostat=ios)
17 if (ios /= 0) then
18 write(*,*)
19 write(*,'("Error(readspecies): error opening species file ",A)') &
20 trim(sppath)//trim(spfname(is))
21 write(*,*)
22 stop
23 end if
24 read(50,*) spsymb(is)
25 read(50,*) spname(is)
26 read(50,*) spzn(is)
27 read(50,*) spmass(is)
28 read(50,*) rminsp(is),rmt(is),rmaxsp(is),nrmt(is)
29 if (rminsp(is) <= 0.d0) then
30 write(*,*)
31 write(*,'("Error(readspecies): rminsp <= 0 : ",G18.10)') rminsp(is)
32 write(*,'(" for species ",I4)') is
33 write(*,*)
34 stop
35 end if
36 if (rmt(is) <= rminsp(is)) then
37 write(*,*)
38 write(*,'("Error(readspecies): rmt <= rminsp : ",2G18.10)') rmt(is), &
39 rminsp(is)
40 write(*,'(" for species ",I4)') is
41 write(*,*)
42 stop
43 end if
44 if (rmaxsp(is) < rmt(is)) then
45 write(*,*)
46 write(*,'("Error(readspecies): rmaxsp < rmt : ",2G18.10)') rmaxsp(is), &
47 rmt(is)
48 write(*,*)
49 stop
50 end if
51 if (nrmt(is) < 20) then
52 write(*,*)
53 write(*,'("Error(readspecies): nrmt too small : ",I8)') nrmt(is)
54 write(*,'(" for species ",I4)') is
55 write(*,*)
56 stop
57 end if
58! multiply nrmt by the scale factor
59 nrmt(is)=nint(dble(nrmt(is))*nrmtscf)
60! reduce the minimum radial mesh point by the same factor
61 rminsp(is)=rminsp(is)/nrmtscf
62 read(50,*) nstsp(is)
63 if ((nstsp(is) < 1).or.(nstsp(is) > maxstsp)) then
64 write(*,*)
65 write(*,'("Error(readspecies): nstsp out of range : ",I8)') nstsp(is)
66 write(*,'(" for species ",I4)') is
67 write(*,*)
68 stop
69 end if
70 do ist=1,nstsp(is)
71 read(50,*) nsp(ist,is),lsp(ist,is),ksp(ist,is),occsp(ist,is),spcore(ist,is)
72 if (nsp(ist,is) < 1) then
73 write(*,*)
74 write(*,'("Error(readspecies): nsp < 1 : ",I8)') nsp(ist,is)
75 write(*,'(" for species ",I4)') is
76 write(*,'(" and state ",I4)') ist
77 write(*,*)
78 stop
79 end if
80 if (lsp(ist,is) < 0) then
81 write(*,*)
82 write(*,'("Error(readspecies): lsp < 0 : ",I8)') lsp(ist,is)
83 write(*,'(" for species ",I4)') is
84 write(*,'(" and state ",I4)') ist
85 write(*,*)
86 stop
87 end if
88 if (ksp(ist,is) < 1) then
89 write(*,*)
90 write(*,'("Error(readspecies): ksp < 1 : ",I8)') ksp(ist,is)
91 write(*,'(" for species ",I4)') is
92 write(*,'(" and state ",I4)') ist
93 write(*,*)
94 stop
95 end if
96 if (occsp(ist,is) < 0.d0) then
97 write(*,*)
98 write(*,'("Error(readspecies): occsp < 0 : ",G18.10)') occsp(ist,is)
99 write(*,'(" for species ",I4)') is
100 write(*,'(" and state ",I4)') ist
101 write(*,*)
102 stop
103 end if
104 end do
105 read(50,*) apword(0,is)
106 if (apword(0,is) < 1) then
107 write(*,*)
108 write(*,'("Error(readspecies): apword < 1 : ",I8)') apword(0,is)
109 write(*,'(" for species ",I4)') is
110 write(*,*)
111 stop
112 end if
113 if (apword(0,is) > maxapword) then
114 write(*,*)
115 write(*,'("Error(readspecies): apword too large : ",I8)') apword(0,is)
116 write(*,'(" for species ",I4)') is
117 write(*,'("Adjust maxapword in modmain and recompile code")')
118 write(*,*)
119 stop
120 end if
121! set the APW orders for l>0
122 apword(1:lmaxapw,is)=apword(0,is)
123 do io=1,apword(0,is)
124 read(50,*) apwe0(io,0,is),apwdm(io,0,is),apwve(io,0,is)
125 if (apwdm(io,0,is) < 0) then
126 write(*,*)
127 write(*,'("Error(readspecies): apwdm < 0 : ",I8)') apwdm(io,0,is)
128 write(*,'(" for species ",I4)') is
129 write(*,'(" and order ",I4)') io
130 write(*,*)
131 stop
132 end if
133! set the APW linearisation energies, derivative orders and variability for l>0
134 apwe0(io,1:lmaxapw,is)=apwe0(io,0,is)
135 apwdm(io,1:lmaxapw,is)=apwdm(io,0,is)
136 apwve(io,1:lmaxapw,is)=apwve(io,0,is)
137 e0min=min(e0min,apwe0(io,0,is))
138 end do
139 read(50,*) nlx
140 if (nlx < 0) then
141 write(*,*)
142 write(*,'("Error(readspecies): nlx < 0 : ",I8)') nlx
143 write(*,'(" for species ",I4)') is
144 write(*,*)
145 stop
146 end if
147 do ilx=1,nlx
148 read(50,*) lx,io
149 if (lx < 0) then
150 write(*,*)
151 write(*,'("Error(readspecies): lx < 0 : ",I8)') lx
152 write(*,'(" for species ",I4)') is
153 write(*,'(" and exception number ",I4)') ilx
154 write(*,*)
155 stop
156 end if
157 if (lx > lmaxapw) then
158 write(*,*)
159 write(*,'("Error(readspecies): lx > lmaxapw : ",I8)') lx
160 write(*,'(" for species ",I4)') is
161 write(*,'(" and exception number ",I4)') ilx
162 write(*,*)
163 stop
164 end if
165 apword(lx,is)=io
166 if (apword(lx,is) < 1) then
167 write(*,*)
168 write(*,'("Error(readspecies): apword < 1 : ",I8)') apword(lx,is)
169 write(*,'(" for species ",I4)') is
170 write(*,'(" and exception number ",I4)') ilx
171 write(*,*)
172 stop
173 end if
174 if (apword(lx,is) > maxapword) then
175 write(*,*)
176 write(*,'("Error(readspecies): apword too large : ",I8)') apword(lx,is)
177 write(*,'(" for species ",I4)') is
178 write(*,'(" and exception number ",I4)') ilx
179 write(*,'("Adjust maxapword in modmain and recompile code")')
180 write(*,*)
181 stop
182 end if
183 do io=1,apword(lx,is)
184 read(50,*) apwe0(io,lx,is),apwdm(io,lx,is),apwve(io,lx,is)
185 if (apwdm(io,lx,is) < 0) then
186 write(*,*)
187 write(*,'("Error(readspecies): apwdm < 0 : ",I8)') apwdm(io,lx,is)
188 write(*,'(" for species ",I4)') is
189 write(*,'(" exception number ",I4)') ilx
190 write(*,'(" and order ",I4)') io
191 write(*,*)
192 stop
193 end if
194 e0min=min(e0min,apwe0(io,lx,is))
195 end do
196 end do
197! add excess order to APW functions if required
198 if (nxoapwlo > 0) then
199 do l=0,lmaxapw
200 jo=apword(l,is)
201 ko=jo+nxoapwlo
202 if (ko > maxapword) ko=maxapword
203 i=0
204 do io=jo+1,ko
205 i=i+1
206 apwe0(io,l,is)=apwe0(jo,l,is)
207 apwdm(io,l,is)=apwdm(jo,l,is)+i
208 apwve(io,l,is)=apwve(jo,l,is)
209 end do
210 apword(l,is)=ko
211 end do
212 end if
213 read(50,*) nlorb(is)
214 if (nlorb(is) < 0) then
215 write(*,*)
216 write(*,'("Error(readspecies): nlorb < 0 : ",I8)') nlorb(is)
217 write(*,'(" for species ",I4)') is
218 write(*,*)
219 stop
220 end if
221 if (nlorb(is) > maxlorb) then
222 write(*,*)
223 write(*,'("Error(readspecies): nlorb too large : ",I8)') nlorb(is)
224 write(*,'(" for species ",I4)') is
225 write(*,'("Adjust maxlorb in modmain and recompile code")')
226 write(*,*)
227 stop
228 end if
229 do ilo=1,nlorb(is)
230 read(50,*) lorbl(ilo,is),lorbord(ilo,is)
231 if (lorbl(ilo,is) < 0) then
232 write(*,*)
233 write(*,'("Error(readspecies): lorbl < 0 : ",I8)') lorbl(ilo,is)
234 write(*,'(" for species ",I4)') is
235 write(*,'(" and local-orbital ",I4)') ilo
236 write(*,*)
237 stop
238 end if
239 if (lorbl(ilo,is) > lmaxo) then
240 write(*,*)
241 write(*,'("Error(readspecies): lorbl > lmaxo : ",2I8)') lorbl(ilo,is), &
242 lmaxo
243 write(*,'(" for species ",I4)') is
244 write(*,'(" and local-orbital ",I4)') ilo
245 write(*,*)
246 stop
247 end if
248 if (lorbord(ilo,is) < 2) then
249 write(*,*)
250 write(*,'("Error(readspecies): lorbord < 2 : ",I8)') lorbord(ilo,is)
251 write(*,'(" for species ",I4)') is
252 write(*,'(" and local-orbital ",I4)') ilo
253 write(*,*)
254 stop
255 end if
256 if (lorbord(ilo,is) > maxlorbord) then
257 write(*,*)
258 write(*,'("Error(readspecies): lorbord too large : ",I8)') lorbord(ilo,is)
259 write(*,'(" for species ",I4)') is
260 write(*,'(" and local-orbital ",I4)') ilo
261 write(*,'("Adjust maxlorbord in modmain and recompile code")')
262 write(*,*)
263 stop
264 end if
265 do io=1,lorbord(ilo,is)
266 read(50,*) lorbe0(io,ilo,is),lorbdm(io,ilo,is),lorbve(io,ilo,is)
267 if (lorbdm(io,ilo,is) < 0) then
268 write(*,*)
269 write(*,'("Error(readspecies): lorbdm < 0 : ",I8)') lorbdm(io,ilo,is)
270 write(*,'(" for species ",I4)') is
271 write(*,'(" local-orbital ",I4)') ilo
272 write(*,'(" and order ",I4)') io
273 write(*,*)
274 stop
275 end if
276 e0min=min(e0min,lorbe0(io,ilo,is))
277 end do
278 end do
279! add excess local-orbitals if required
280 if (nxlo > 0) then
281 lx=-1
282 do ilo=1,nlorb(is)
283 do io=1,lorbord(ilo,is)
284 if (lorbe0(io,ilo,is) < 0.d0) goto 10
285 end do
286 if (lorbl(ilo,is) > lx) lx=lorbl(ilo,is)
28710 continue
288 end do
289 ilo=nlorb(is)
290 do i=1,nxlo
291 if (ilo == maxlorb) exit
292 l=lx+i
293 if (l > lmaxo) exit
294 ilo=ilo+1
295 lorbl(ilo,is)=l
296 lorbord(ilo,is)=apword(l,is)+1
297 do io=1,lorbord(ilo,is)
298 lorbe0(io,ilo,is)=apwe0(1,l,is)
299 lorbdm(io,ilo,is)=io-1
300 lorbve(io,ilo,is)=apwve(1,l,is)
301 end do
302 end do
303 nlorb(is)=ilo
304 end if
305! add excess order to local-orbitals if required
306 if (nxoapwlo > 0) then
307 do ilo=1,nlorb(is)
308! find the maximum energy derivative
309 jo=1
310 j=lorbdm(jo,ilo,is)
311 do io=1,lorbord(ilo,is)
312 i=lorbdm(io,ilo,is)
313 if (i > j) then
314 jo=io
315 j=i
316 end if
317 end do
318 ko=lorbord(ilo,is)+nxoapwlo
319 if (ko > maxlorbord) ko=maxlorbord
320 i=0
321 do io=lorbord(ilo,is)+1,ko
322 i=i+1
323 lorbe0(io,ilo,is)=lorbe0(jo,ilo,is)
324 lorbdm(io,ilo,is)=lorbdm(jo,ilo,is)+i
325 lorbve(io,ilo,is)=lorbve(jo,ilo,is)
326 end do
327 lorbord(ilo,is)=ko
328 end do
329 end if
330 close(50)
331end do
332if (rmtall > 0.d0) then
333! set all muffin-tin radii to single value if required
335else
336! scale the muffin-tin radii if required
337 if (rmtscf /= 1.d0) rmt(1:nspecies)=rmtscf*rmt(1:nspecies)
338! apply averaging scheme to the muffin-tin radii
339 if (mrmtav > 0) call rmtavrg
340end if
341! make a copy of the muffin-tin radii
343! add conduction state local-orbitals if required
344if (lorbcnd) call addlorbcnd
345! maximum number of local-orbitals over all species
346nlomax=maxval(nlorb(1:nspecies))
347! generate the index which arranges the local-orbitals in ascending energy
348call genidxelo
349! subtract 2 Hartree from the minimum energy
350e0min=e0min-2.d0
351end subroutine
352
subroutine addlorbcnd
Definition addlorbcnd.f90:7
subroutine genidxelo
Definition genidxelo.f90:7
integer, dimension(maxstsp, maxspecies) lsp
Definition modmain.f90:123
real(8), dimension(maxapword, 0:maxlapw, maxspecies) apwe0
Definition modmain.f90:766
integer nxlo
Definition modmain.f90:838
integer, dimension(maxstsp, maxspecies) nsp
Definition modmain.f90:121
integer, parameter maxstsp
Definition modmain.f90:111
real(8) nrmtscf
Definition modmain.f90:148
character(256), dimension(maxspecies) spfname
Definition modmain.f90:74
integer, dimension(maxspecies) nrmt
Definition modmain.f90:150
logical, dimension(maxstsp, maxspecies) spcore
Definition modmain.f90:127
real(8), dimension(maxspecies) rmt
Definition modmain.f90:162
integer, parameter maxapword
Definition modmain.f90:754
integer nxoapwlo
Definition modmain.f90:836
integer, dimension(0:maxlapw, maxspecies) apword
Definition modmain.f90:758
real(8), dimension(maxlorbord, maxlorb, maxspecies) lorbe0
Definition modmain.f90:804
integer, dimension(maxlorbord, maxlorb, maxspecies) lorbdm
Definition modmain.f90:810
real(8), dimension(maxstsp, maxspecies) occsp
Definition modmain.f90:133
integer, dimension(maxspecies) nstsp
Definition modmain.f90:113
logical lorbcnd
Definition modmain.f90:832
integer nlomax
Definition modmain.f90:788
real(8) e0min
Definition modmain.f90:825
integer, dimension(maxapword, 0:maxlapw, maxspecies) apwdm
Definition modmain.f90:770
integer lmaxapw
Definition modmain.f90:197
character(64), dimension(maxspecies) spsymb
Definition modmain.f90:78
integer, dimension(maxlorb, maxspecies) lorbord
Definition modmain.f90:792
integer, parameter maxlorbord
Definition modmain.f90:782
real(8) rmtall
Definition modmain.f90:158
integer, dimension(maxstsp, maxspecies) ksp
Definition modmain.f90:125
integer mrmtav
Definition modmain.f90:156
logical, dimension(maxapword, 0:maxlapw, maxspecies) apwve
Definition modmain.f90:772
character(256) sppath
Definition modmain.f90:72
integer lmaxo
Definition modmain.f90:201
real(8), dimension(maxspecies) rmaxsp
Definition modmain.f90:105
real(8), dimension(maxspecies) rminsp
Definition modmain.f90:103
real(8) rmtscf
Definition modmain.f90:154
real(8), dimension(maxspecies) rmt0
Definition modmain.f90:162
integer, dimension(maxspecies) nlorb
Definition modmain.f90:786
integer nspecies
Definition modmain.f90:34
logical, dimension(maxlorbord, maxlorb, maxspecies) lorbve
Definition modmain.f90:812
character(64), dimension(maxspecies) spname
Definition modmain.f90:76
integer, parameter maxlorb
Definition modmain.f90:780
real(8), dimension(maxspecies) spmass
Definition modmain.f90:101
real(8), dimension(maxspecies) spzn
Definition modmain.f90:80
integer, dimension(maxlorb, maxspecies) lorbl
Definition modmain.f90:796
subroutine readspecies
subroutine rmtavrg
Definition rmtavrg.f90:7