The Elk Code
bandstr.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: bandstr
8 ! !INTERFACE:
9 subroutine bandstr
10 ! !USES:
11 use modmain
12 use modomp
13 ! !DESCRIPTION:
14 ! Produces a band structure along the path in reciprocal space which connects
15 ! the vertices in the array {\tt vvlp1d}. The band structure is obtained from
16 ! the second-variational eigenvalues and is written to the file {\tt BAND.OUT}
17 ! with the Fermi energy set to zero. If required, band structures are plotted
18 ! to files {\tt BAND\_Sss\_Aaaaa.OUT} for atom {\tt aaaa} of species {\tt ss},
19 ! which include the band characters for each $l$ component of that atom in
20 ! columns 4 onwards. Column 3 contains the sum over $l$ of the characters.
21 ! Vertex location lines are written to {\tt BANDLINES.OUT}.
22 !
23 ! !REVISION HISTORY:
24 ! Created June 2003 (JKD)
25 !EOP
26 !BOC
27 implicit none
28 ! local variables
29 logical tspndg,tlmdg
30 integer ik,ist,ispn,idm
31 integer is,ia,ias,nthd
32 integer l,m,lm,iv
33 real(8) emin,emax,sm
34 real(8) v(ndmag),t1,t2
35 character(256) fname
36 complex(8) su2(2,2),z1
37 ! allocatable arrays
38 real(8), allocatable :: evalfv(:,:)
39 ! low precision for band character array saves memory
40 real(4), allocatable :: bc(:,:,:,:),elm(:,:)
41 complex(8), allocatable :: ulm(:,:,:),dmat(:,:,:,:,:)
42 complex(8), allocatable :: apwalm(:,:,:,:,:),evecfv(:,:,:),evecsv(:,:)
43 ! initialise universal variables
44 call init0
45 call init1
46 ! calculate only the diagonal parts of the density matrices by default
47 tspndg=.true.
48 tlmdg=.true.
49 select case(task)
50 case(21)
51  allocate(bc(0:lmaxdb,natmtot,nstsv,nkpt))
52 case(22)
53  allocate(bc(lmmaxdb,natmtot,nstsv,nkpt))
54 ! generate unitary matrices which convert the Yₗₘ basis into the irreducible
55 ! representation basis of the symmetry group at each atomic site
56  if (lmirep) then
57  allocate(elm(lmmaxdb,natmtot),ulm(lmmaxdb,lmmaxdb,natmtot))
58  call genlmirep(elm,ulm)
59 ! write the eigenvalues of a pseudorandom matrix symmetrised by the site
60 ! symmetries in the Yₗₘ basis
61  call writeelmirep('.OUT',elm)
62 ! require full density matrix in (l,m) degrees of freedom
63  tlmdg=.false.
64  end if
65 case(23)
66  if (.not.spinpol) goto 10
67  allocate(bc(nspinor,natmtot,nstsv,nkpt))
68 ! compute the SU(2) operator used for rotating the density matrix to the
69 ! desired spin-quantisation axis; if axis is not z then the full matrix in
70 ! spin degrees of freedom is required
71  call sqasu2(sqaxis,tspndg,su2)
72 case(24)
73  if (.not.spinpol) goto 10
74  allocate(bc(ndmag,natmtot,nstsv,nkpt))
75 end select
76 ! read density and potentials from file
77 call readstate
78 ! Fourier transform Kohn-Sham potential to G-space
79 call genvsig
80 ! find the new linearisation energies
81 call linengy
82 ! generate the APW and local-orbital radial functions and integrals
83 call genapwlofr
84 ! generate the spin-orbit coupling radial functions
85 call gensocfr
86 ! begin parallel loop over k-points
87 call holdthd(nkpt,nthd)
88 !$OMP PARALLEL DEFAULT(SHARED) &
89 !$OMP PRIVATE(evalfv,evecfv,evecsv) &
90 !$OMP PRIVATE(dmat,apwalm,ispn,ias,ist) &
91 !$OMP PRIVATE(l,m,lm,sm,v,t1,t2,z1) &
92 !$OMP NUM_THREADS(nthd)
93 allocate(evalfv(nstfv,nspnfv))
94 allocate(evecfv(nmatmax,nstfv,nspnfv),evecsv(nstsv,nstsv))
95 if (task >= 21) then
96  allocate(dmat(lmmaxdb,nspinor,lmmaxdb,nspinor,nstsv))
97  allocate(apwalm(ngkmax,apwordmax,lmmaxapw,natmtot,nspnfv))
98 end if
99 !$OMP DO SCHEDULE(DYNAMIC)
100 do ik=1,nkpt
101 !$OMP CRITICAL(bandstr_)
102  write(*,'("Info(bandstr): ",I6," of ",I6," k-points")') ik,nkpt
103 !$OMP END CRITICAL(bandstr_)
104 ! solve the first- and second-variational eigenvalue equations
105  call eveqn(ik,evalfv,evecfv,evecsv)
106 ! compute the band characters if required
107  if (task >= 21) then
108 ! find the matching coefficients
109  do ispn=1,nspnfv
110  call match(ngk(ispn,ik),vgkc(:,:,ispn,ik),gkc(:,ispn,ik), &
111  sfacgk(:,:,ispn,ik),apwalm(:,:,:,:,ispn))
112  end do
113 ! average band character over spin and m for all atoms
114  do ias=1,natmtot
115 ! generate the density matrix
116  call gendmatk(tspndg,tlmdg,0,lmaxdb,ias,nstsv,[0],ngk(:,ik),apwalm, &
117  evecfv,evecsv,lmmaxdb,dmat)
118 ! convert (l,m) part to an irreducible representation if required
119  if (.not.tlmdg) call dmatulm(ulm(:,:,ias),dmat)
120 ! spin rotate the density matrices to desired spin-quantisation axis
121  if (.not.tspndg) call dmatsu2(lmmaxdb,su2,dmat)
122  do ist=1,nstsv
123  select case(task)
124  case(21)
125 ! l character of band
126  lm=0
127  do l=0,lmaxdb
128  sm=0.d0
129  do m=-l,l
130  lm=lm+1
131  do ispn=1,nspinor
132  sm=sm+dble(dmat(lm,ispn,lm,ispn,ist))
133  end do
134  end do
135  bc(l,ias,ist,ik)=sm
136  end do
137  case(22)
138 ! (l,m) character of band
139  do lm=1,lmmaxdb
140  sm=0.d0
141  do ispn=1,nspinor
142  sm=sm+dble(dmat(lm,ispn,lm,ispn,ist))
143  end do
144  bc(lm,ias,ist,ik)=sm
145  end do
146  case(23)
147 ! spin character of band
148  do ispn=1,nspinor
149  sm=0.d0
150  do lm=1,lmmaxdb
151  sm=sm+dble(dmat(lm,ispn,lm,ispn,ist))
152  end do
153  bc(ispn,ias,ist,ik)=sm
154  end do
155  case(24)
156 ! magnetic moment character of band
157  v(1:ndmag)=0.d0
158  do lm=1,lmmaxdb
159  t1=dble(dmat(lm,1,lm,1,ist))
160  t2=dble(dmat(lm,2,lm,2,ist))
161  v(ndmag)=v(ndmag)+t1-t2
162  if (ncmag) then
163  z1=2.d0*dmat(lm,1,lm,2,ist)
164  v(1)=v(1)+z1%re
165  v(2)=v(2)-z1%im
166  end if
167  end do
168  bc(1:ndmag,ias,ist,ik)=v(1:ndmag)
169  end select
170  end do
171  end do
172  end if
173 ! end loop over k-points
174 end do
175 !$OMP END DO
176 deallocate(evalfv,evecfv,evecsv)
177 if (task >= 21) deallocate(dmat,apwalm)
178 !$OMP END PARALLEL
179 call freethd(nthd)
180 ! subtract the Fermi energy
181 evalsv(:,:)=evalsv(:,:)-efermi
182 ! find the minimum and maximum eigenvalues
183 emin=minval(evalsv(:,:))
184 emax=maxval(evalsv(:,:))
185 t1=(emax-emin)*0.5d0
186 emin=emin-t1
187 emax=emax+t1
188 ! output the band structure
189 if (task == 20) then
190  open(50,file='BAND.OUT',form='FORMATTED',action='WRITE')
191  do ist=1,nstsv
192  do ik=1,nkpt
193  write(50,'(2G18.10)') dpp1d(ik),evalsv(ist,ik)
194  end do
195  write(50,*)
196  end do
197  close(50)
198  write(*,*)
199  write(*,'("Info(bandstr):")')
200  write(*,'(" Band structure plot written to BAND.OUT")')
201 else
202  do ias=1,natmtot
203  is=idxis(ias)
204  ia=idxia(ias)
205  write(fname,'("BAND_S",I2.2,"_A",I4.4,".OUT")') is,ia
206  open(50,file=trim(fname),form='FORMATTED',action='WRITE')
207  do ist=1,nstsv
208  do ik=1,nkpt
209  select case(task)
210  case(21)
211 ! sum band character over l to find total atomic character
212  sm=sum(bc(0:lmaxdb,ias,ist,ik))
213  write(50,'(2G18.10,F12.6)',advance='NO') dpp1d(ik),evalsv(ist,ik),sm
214  do l=0,lmaxdb
215  write(50,'(F12.6)',advance='NO') bc(l,ias,ist,ik)
216  end do
217  write(50,*)
218  case(22)
219  write(50,'(2G18.10)',advance='NO') dpp1d(ik),evalsv(ist,ik)
220  do lm=1,lmmaxdb
221  write(50,'(F12.6)',advance='NO') bc(lm,ias,ist,ik)
222  end do
223  write(50,*)
224  case(23)
225  write(50,'(2G18.10,2F12.6)') dpp1d(ik),evalsv(ist,ik), &
226  (bc(ispn,ias,ist,ik),ispn=1,nspinor)
227  case(24)
228  write(50,'(2G18.10,3F12.6)') dpp1d(ik),evalsv(ist,ik), &
229  (bc(idm,ias,ist,ik),idm=1,ndmag)
230  end select
231  end do
232  write(50,*)
233  end do
234  close(50)
235  end do
236  write(*,*)
237  write(*,'("Info(bandstr):")')
238  write(*,'(" Band structure plot written to BAND_Sss_Aaaaa.OUT")')
239  write(*,'(" for all species and atoms")')
240  write(*,*)
241  write(*,'(" Columns in the file are :")')
242  select case(task)
243  case(21)
244  write(*,'(" distance, eigenvalue, total atomic character, l character &
245  &(l = 0...",I1,")")') lmaxdb
246  case(22)
247  write(*,'(" distance, eigenvalue, (l,m) character &
248  &(l = 0...",I1,", m = -l...l)")') lmaxdb
249  if (lmirep) then
250  write(*,*)
251  write(*,'(" Eigenvalues of a random matrix symmetrised with the site")')
252  write(*,'(" symmetries in the Yₗₘ basis written to ELMIREP.OUT for all")')
253  write(*,'(" species and atoms. Degenerate eigenvalues correspond to")')
254  write(*,'(" irreducible representations of each site symmetry group")')
255  end if
256  case(23)
257  write(*,'(" distance, eigenvalue, spin-up and spin-down characters")')
258  write(*,*)
259  write(*,'(" Spin-quantisation axis : ",3G18.10)') sqaxis(:)
260  case(24)
261  write(*,'(" distance, eigenvalue, moment character")')
262  end select
263 end if
264 write(*,*)
265 write(*,'(" Fermi energy is at zero in plot")')
266 ! output the vertex location lines
267 open(50,file='BANDLINES.OUT',form='FORMATTED',action='WRITE')
268 do iv=1,nvp1d
269  write(50,'(2G18.10)') dvp1d(iv),emin
270  write(50,'(2G18.10)') dvp1d(iv),emax
271  write(50,*)
272 end do
273 close(50)
274 write(*,*)
275 write(*,'(" Vertex location lines written to BANDLINES.OUT")')
276 if (task >= 21) deallocate(bc)
277 if ((task == 22).and.lmirep) deallocate(elm,ulm)
278 return
279 10 continue
280 write(*,*)
281 write(*,'("Error(bandstr): spin-unpolarised calculation")')
282 write(*,*)
283 stop
284 end subroutine
285 !EOC
286 
integer nmatmax
Definition: modmain.f90:853
real(8) efermi
Definition: modmain.f90:907
real(8), dimension(:,:), allocatable evalsv
Definition: modmain.f90:921
integer task
Definition: modmain.f90:1299
logical spinpol
Definition: modmain.f90:228
integer lmmaxapw
Definition: modmain.f90:199
integer nkpt
Definition: modmain.f90:461
integer ngkmax
Definition: modmain.f90:499
real(8), dimension(:), allocatable dpp1d
Definition: modmain.f90:1126
subroutine match(ngp, vgpc, gpc, sfacgp, apwalm)
Definition: match.f90:10
Definition: modomp.f90:6
subroutine dmatulm(ulm, dmat)
Definition: dmatulm.f90:10
subroutine sqasu2(sqaxis, tsqaz, su2)
Definition: sqasu2.f90:7
subroutine bandstr
Definition: bandstr.f90:10
subroutine genvsig
Definition: genvsig.f90:10
subroutine genapwlofr
Definition: genapwlofr.f90:7
subroutine gensocfr
Definition: gensocfr.f90:7
complex(8), dimension(:,:,:,:), allocatable sfacgk
Definition: modmain.f90:509
subroutine linengy
Definition: linengy.f90:10
subroutine genlmirep(elm, ulm)
Definition: genlmirep.f90:10
integer lmmaxdb
Definition: modmain.f90:1081
integer nstsv
Definition: modmain.f90:889
integer, dimension(:,:), allocatable ngk
Definition: modmain.f90:497
real(8), dimension(:), allocatable dvp1d
Definition: modmain.f90:1122
logical lmirep
Definition: modmain.f90:1098
real(8), dimension(3) sqaxis
Definition: modmain.f90:1101
subroutine eveqn(ik, evalfv, evecfv, evecsv)
Definition: eveqn.f90:10
integer lmaxdb
Definition: modmain.f90:1081
integer nspinor
Definition: modmain.f90:267
real(8), dimension(:,:,:,:), allocatable vgkc
Definition: modmain.f90:505
subroutine init1
Definition: init1.f90:10
subroutine gendmatk(tspndg, tlmdg, lmin, lmax, ias, nst, idx, ngp, apwalm, evecfv, evecsv, ld, dmat)
Definition: gendmatk.f90:8
integer apwordmax
Definition: modmain.f90:760
integer, dimension(maxatoms *maxspecies) idxis
Definition: modmain.f90:44
real(8), dimension(:,:,:), allocatable gkc
Definition: modmain.f90:507
subroutine writeelmirep(fext, elm)
Definition: writeelmirep.f90:7
subroutine readstate
Definition: readstate.f90:10
subroutine freethd(nthd)
Definition: modomp.f90:106
subroutine holdthd(nloop, nthd)
Definition: modomp.f90:78
integer, dimension(maxatoms *maxspecies) idxia
Definition: modmain.f90:45
subroutine init0
Definition: init0.f90:10
integer natmtot
Definition: modmain.f90:40
logical ncmag
Definition: modmain.f90:240
subroutine dmatsu2(lmmax, su2, dmat)
Definition: dmatsu2.f90:10
integer nvp1d
Definition: modmain.f90:1114
integer nstfv
Definition: modmain.f90:887
integer nspnfv
Definition: modmain.f90:289