The Elk Code
writew90spn.f90
Go to the documentation of this file.
1 
2 ! Copyright (C) 2015 Manh Duc Le, 2017-18 Arsenii Gerasimov, Yaroslav Kvashnin
3 ! and Lars Nordstrom. This file is distributed under the terms of the GNU
4 ! General Public License. See the file COPYING for license details.
5 
6 subroutine writew90spn
7 use modmain
8 use modw90
9 implicit none
10 ! local variables
11 integer ik,ist,jst,i,j
12 complex(8) z1
13 character(256) fname
14 ! allocatable arrays
15 complex(8), allocatable :: evecsv(:,:),smat(:,:,:,:)
16 if (.not.spinpol) return
17 fname=trim(seedname)//'.spn'
18 open(50,file=trim(fname),action='WRITE',form='FORMATTED')
19 write(50,'("Generated by Elk version ",I0,".",I0,".",I0)') version
20 write(50,'(3I8)') num_bands,nkptnr
21 allocate(evecsv(nstsv,nstsv),smat(nstsv,nstsv,2,2))
22 do ik=1,nkptnr
23 ! generate the spin operator matrix elements
24 ! (note that this is *not* a density matrix)
25  call getevecsv(filext,0,vkl(:,ik),evecsv)
26  call gensmatk(evecsv,smat)
27  do j=1,num_bands
28  jst=idxw90(j)
29  do i=1,j
30  ist=idxw90(i)
31  z1=smat(ist,jst,1,2)+smat(ist,jst,2,1)
32  write(50,'(2G18.10)') z1
33  z1=smat(ist,jst,1,2)-smat(ist,jst,2,1)
34  z1=cmplx(z1%im,-z1%re,8)
35  write(50,'(2G18.10)') z1
36  z1=smat(ist,jst,1,1)-smat(ist,jst,2,2)
37  write(50,'(2G18.10)') z1
38  end do
39  end do
40 end do
41 deallocate(evecsv,smat)
42 close(50)
43 write(*,*)
44 write(*,'("Info(writew90spn): created file ",A)') trim(fname)
45 end subroutine
46 
integer num_bands
Definition: modw90.f90:20
subroutine getevecsv(fext, ikp, vpl, evecsv)
Definition: getevecsv.f90:7
character(256) filext
Definition: modmain.f90:1301
logical spinpol
Definition: modmain.f90:228
integer nkptnr
Definition: modmain.f90:463
subroutine gensmatk(evecsv, smat)
Definition: gensmatk.f90:7
character(256) seedname
Definition: modw90.f90:12
integer nstsv
Definition: modmain.f90:889
Definition: modw90.f90:6
integer, dimension(:), allocatable idxw90
Definition: modw90.f90:22
real(8), dimension(:,:), allocatable vkl
Definition: modmain.f90:471
integer, dimension(3), parameter version
Definition: modmain.f90:1289
subroutine writew90spn
Definition: writew90spn.f90:7